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/core/utility/MonotCubicInterpolator.hpp>
29 #include <array>
30 #include <tuple>
31 #include <vector>
32 
33 namespace Opm {
38  double getVoigtValue(const SinglePhaseUpscaler::permtensor_t& K, int voigt_idx);
39 
44  void setVoigtValue(SinglePhaseUpscaler::permtensor_t& K, int voigt_idx, double val);
45 
48  public:
49  bool isMaster;
52  int points;
55  std::vector<double> WaterSaturation;
57  std::vector<int> satnums;
58  std::vector<MonotCubicInterpolator> InvJfunctions;
59  std::array<std::array<std::vector<MonotCubicInterpolator>,2>,3> Krfunctions;
62  std::vector<MonotCubicInterpolator> SwPcfunctions;
63  std::string saturationstring;
64  size_t tesselatedCells;
65  double volume;
66  double poreVolume;
67  std::vector<double> pressurePoints;
68 
73  RelPermUpscaleHelper(int mpi_rank, std::map<std::string,std::string>& options_);
74 
76  void collectResults();
77 
82  std::vector<std::vector<double>> getRelPerm(int phase) const;
83 
86 
93  void sanityCheckInput(Opm::DeckConstPtr deck,
94  double minPerm, double maxPerm, double minPoro);
95 
98 
102 
110  double tesselateGrid(Opm::DeckConstPtr deck);
111 
115  void calculateCellPressureGradients(const std::array<int,3>& res);
116 
121 
125 
130  std::tuple<double,double> upscalePermeability(int mpi_rank);
131  private:
137  bool checkCurve(MonotCubicInterpolator& func);
138 
139  double Swir;
140  double Swor;
141  double Pcmin;
142  double Pcmax;
143  double critRelpThresh;
144  std::vector<int> node_vs_pressurepoint;
145  std::array<std::vector<std::vector<double>>,2> PhasePerm;
146  SinglePhaseUpscaler::permtensor_t permTensorInv;
147  SinglePhaseUpscaler upscaler;
148  std::array<std::vector<double>,3> perms;
149  std::vector<double> poros;
150  std::vector<double> zcorns;
151  double minSinglePhasePerm;
152  std::vector<double> cellPoreVolumes;
153  MonotCubicInterpolator WaterSaturationVsCapPressure;
154 
155  // Reference: http://www.spe.org/spe-site/spe/spe/papers/authors/Metric_Standard.pdf
156  const double milliDarcyToSqMetre;
157 
158  std::vector<double> dP;
159  std::map<std::string,std::string>& options;
160  };
161 }
162 
163 #endif
void upscaleCapillaryPressure()
Upscale capillary pressure.
SinglePhaseUpscaler::BoundaryConditionType boundaryCondition
Boundary conditions to use.
Definition: RelPermUtils.hpp:61
std::vector< double > pressurePoints
Vector of capillary pressure points between Swor and Swir.
Definition: RelPermUtils.hpp:67
Definition: applier.hpp:18
void setVoigtValue(SinglePhaseUpscaler::permtensor_t &K, int voigt_idx, double val)
Set value in tensor.
std::vector< int > satnums
Cell satnums.
Definition: RelPermUtils.hpp:57
ResProp::MutablePermTensor permtensor_t
A type for the upscaled permeability.
Definition: UpscalerBase.hpp:63
bool doEclipseCheck
Whether to check that input relperm curves include relperm at critical saturation points...
Definition: RelPermUtils.hpp:50
double poreVolume
Total pore volume.
Definition: RelPermUtils.hpp:66
bool upscaleBothPhases
Whether to upscale both phases.
Definition: RelPermUtils.hpp:54
int tensorElementCount
Number of independent elements in resulting tensor.
Definition: RelPermUtils.hpp:53
RelPermUpscaleHelper(int mpi_rank, std::map< std::string, std::string > &options_)
Default constructor.
double volume
Total volume.
Definition: RelPermUtils.hpp:65
bool isMaster
Whether this is the master MPI node or not.
Definition: RelPermUtils.hpp:49
void setupBoundaryConditions()
Setup requested boundary conditions.
std::string saturationstring
Fluid system type.
Definition: RelPermUtils.hpp:63
void collectResults()
Collect results from all MPI nodes.
void upscaleSinglePhasePermeability()
Upscale single phase permeability.
void calculateCellPressureGradients(const std::array< int, 3 > &res)
Find cell center pressure gradient for every cell.
void calculateMinMaxCapillaryPressure()
Calculate minimum and maximum capillary pressures.
Helper class for relperm upscaling applications.
Definition: RelPermUtils.hpp:47
std::vector< MonotCubicInterpolator > InvJfunctions
Definition: RelPermUtils.hpp:58
size_t tesselatedCells
Number of "active" cells (Sintef interpretation of "active")
Definition: RelPermUtils.hpp:64
double getVoigtValue(const SinglePhaseUpscaler::permtensor_t &K, int voigt_idx)
Get value from tensor.
void checkCriticalSaturations()
Check that input relperm curevs specify critical saturations.
bool anisotropic_input
Whether input eclipse file has diagonal anisotrophy.
Definition: RelPermUtils.hpp:51
std::array< std::array< std::vector< MonotCubicInterpolator >, 2 >, 3 > Krfunctions
Relperm-curves for each (component->phase->stone type)
Definition: RelPermUtils.hpp:60
std::vector< std::vector< double > > getRelPerm(int phase) const
Calculate relperm values from phase permeabilities.
A class for doing single phase (permeability) upscaling.
Definition: SinglePhaseUpscaler.hpp:50
int points
Number of saturation points to upscale for.
Definition: RelPermUtils.hpp:52
Opm::DeckConstPtr deck(parser->parseFile(file, parseMode))
std::vector< MonotCubicInterpolator > SwPcfunctions
Holds Sw(Pc) for each rocktype.
Definition: RelPermUtils.hpp:62
SinglePhaseUpscaler::permtensor_t permTensor
Tensor of upscaled results.
Definition: RelPermUtils.hpp:56
BoundaryConditionType
Definition: UpscalerBase.hpp:65
void sanityCheckInput(Opm::DeckConstPtr deck, double minPerm, double maxPerm, double minPoro)
Do sanity checks for input file.
std::vector< double > WaterSaturation
Re-upscaled water saturation for the computed pressure points.
Definition: RelPermUtils.hpp:55
std::tuple< double, double > upscalePermeability(int mpi_rank)
Upscale permeabilities.
double tesselateGrid(Opm::DeckConstPtr deck)
Tesselate grid.