RelPermUtils.hpp
Go to the documentation of this file.
1/*
2 Copyright 2010 Statoil ASA.
3
4 This file is part of The Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
24#ifndef OPM_UPSCALING_RELPERM_UTILS_HPP
25#define OPM_UPSCALING_RELPERM_UTILS_HPP
26
27#include <opm/input/eclipse/Parser/Parser.hpp>
28#include <opm/input/eclipse/Deck/Deck.hpp>
29
30#include <opm/common/utility/numeric/MonotCubicInterpolator.hpp>
31
34
35#include <array>
36#include <map>
37#include <memory>
38#include <tuple>
39#include <vector>
40
41namespace Opm {
46 double getVoigtValue(const SinglePhaseUpscaler::permtensor_t& K, int voigt_idx);
47
52 void setVoigtValue(SinglePhaseUpscaler::permtensor_t& K, int voigt_idx, double val);
53
56 public:
57 bool isMaster;
60 int points;
63 std::vector<double> WaterSaturation;
65 std::vector<int> satnums;
66 std::vector<MonotCubicInterpolator> InvJfunctions;
68 std::array<std::array<std::vector<MonotCubicInterpolator>,2>,3> Krfunctions;
70 std::vector<MonotCubicInterpolator> SwPcfunctions;
71 std::string saturationstring;
73 double volume;
74 double poreVolume;
75 std::vector<double> pressurePoints;
76
85 template <class String>
86 static Deck
87 parseEclipseFile(const String& eclipseFileName);
88
93 RelPermUpscaleHelper(int mpi_rank, std::map<std::string,std::string>& options_);
94
97
102 std::vector<std::vector<double>> getRelPerm(int phase) const;
103
106
113 void sanityCheckInput(const Opm::Deck& deck,
114 const double minPerm,
115 const double maxPerm,
116 const double minPoro);
117
120
124
132 double tesselateGrid(const Opm::Deck& deck);
133
137
142
146
151 std::tuple<double,double> upscalePermeability(int mpi_rank);
152
153 private:
159 bool checkCurve(MonotCubicInterpolator& func);
160
161 double Swir;
162 double Swor;
163 double Pcmin;
164 double Pcmax;
165 double critRelpThresh;
166 std::vector<int> node_vs_pressurepoint;
167 std::array<std::vector<std::vector<double>>,2> PhasePerm;
169 SinglePhaseUpscaler upscaler;
170 std::array<std::vector<double>,3> perms;
171 std::vector<double> poros;
172 std::vector<double> zcorns;
173 double minSinglePhasePerm;
174 std::vector<double> cellPoreVolumes;
175 MonotCubicInterpolator WaterSaturationVsCapPressure;
176
177#if defined(UNITTEST_TRESPASS_PRIVATE_PROPERTY_DP)
178 public: // Intrusive unit testing of calculcateCellPressureGradients()
179#endif // UNITTEST_TRESPASS_PRIVATE_PROPERTY_DP
180 std::vector<double> dP;
181
182#if defined(UNITTEST_TRESPASS_PRIVATE_PROPERTY_DP)
183 private:
184#endif // UNITTEST_TRESPASS_PRIVATE_PROPERTY_DP
185
186 std::map<std::string,std::string>& options;
187 };
188
189 template <class String>
190 Deck
191 RelPermUpscaleHelper::parseEclipseFile(const String& eclipseFileName)
192 {
193 auto parser = Parser{};
195
196 return parser.parseFile(eclipseFileName);
197 }
198}
199
200#endif // OPM_UPSCALING_RELPERM_UTILS_HPP
Helper class for relperm upscaling applications.
Definition: RelPermUtils.hpp:55
std::tuple< double, double > upscalePermeability(int mpi_rank)
Upscale permeabilities.
SinglePhaseUpscaler::BoundaryConditionType boundaryCondition
Boundary conditions to use.
Definition: RelPermUtils.hpp:69
std::vector< int > satnums
Cell satnums.
Definition: RelPermUtils.hpp:65
double volume
Total volume.
Definition: RelPermUtils.hpp:73
void upscaleSinglePhasePermeability()
Upscale single phase permeability.
RelPermUpscaleHelper(int mpi_rank, std::map< std::string, std::string > &options_)
Default constructor.
void upscaleCapillaryPressure()
Upscale capillary pressure.
SinglePhaseUpscaler::permtensor_t permTensor
Tensor of upscaled results.
Definition: RelPermUtils.hpp:64
static Deck parseEclipseFile(const String &eclipseFileName)
Form Deck object from input file (ECLIPSE format)
Definition: RelPermUtils.hpp:191
std::vector< MonotCubicInterpolator > SwPcfunctions
Holds Sw(Pc) for each rocktype.
Definition: RelPermUtils.hpp:70
void checkCriticalSaturations()
Check that input relperm curevs specify critical saturations.
bool isMaster
Whether this is the master MPI node or not.
Definition: RelPermUtils.hpp:57
std::vector< MonotCubicInterpolator > InvJfunctions
Definition: RelPermUtils.hpp:66
std::vector< double > pressurePoints
Vector of capillary pressure points between Swor and Swir.
Definition: RelPermUtils.hpp:75
void sanityCheckInput(const Opm::Deck &deck, const double minPerm, const double maxPerm, const double minPoro)
Do sanity checks for input file.
size_t tesselatedCells
Number of "active" cells (Sintef interpretation of "active")
Definition: RelPermUtils.hpp:72
void calculateMinMaxCapillaryPressure()
Calculate minimum and maximum capillary pressures.
int tensorElementCount
Number of independent elements in resulting tensor.
Definition: RelPermUtils.hpp:61
std::array< std::array< std::vector< MonotCubicInterpolator >, 2 >, 3 > Krfunctions
Relperm-curves for each (component->phase->stone type)
Definition: RelPermUtils.hpp:68
int points
Number of saturation points to upscale for.
Definition: RelPermUtils.hpp:60
bool doEclipseCheck
Whether to check that input relperm curves include relperm at critical saturation points.
Definition: RelPermUtils.hpp:58
bool anisotropic_input
Whether input eclipse file has diagonal anisotrophy.
Definition: RelPermUtils.hpp:59
void collectResults()
Collect results from all MPI nodes.
std::vector< double > WaterSaturation
Re-upscaled water saturation for the computed pressure points.
Definition: RelPermUtils.hpp:63
double tesselateGrid(const Opm::Deck &deck)
Tesselate grid.
bool upscaleBothPhases
Whether to upscale both phases.
Definition: RelPermUtils.hpp:62
std::string saturationstring
Fluid system type.
Definition: RelPermUtils.hpp:71
std::vector< std::vector< double > > getRelPerm(int phase) const
Calculate relperm values from phase permeabilities.
double poreVolume
Total pore volume.
Definition: RelPermUtils.hpp:74
void setupBoundaryConditions()
Setup requested boundary conditions.
void calculateCellPressureGradients()
Find cell center pressure gradient for every cell.
A class for doing single phase (permeability) upscaling.
Definition: SinglePhaseUpscaler.hpp:51
ResProp::MutablePermTensor permtensor_t
A type for the upscaled permeability.
Definition: UpscalerBase.hpp:66
BoundaryConditionType
Definition: UpscalerBase.hpp:68
auto deck
Definition: elasticity_upscale_impl.hpp:590
Definition: ImplicitAssembly.hpp:43
double getVoigtValue(const SinglePhaseUpscaler::permtensor_t &K, int voigt_idx)
Get value from tensor.
void addNonStandardUpscalingKeywords(Parser &parser)
This function registers non-standard keywords used by opm-upscaling in a parser object.
void setVoigtValue(SinglePhaseUpscaler::permtensor_t &K, int voigt_idx, double val)
Set value in tensor.