|
opm-upscaling
|
Helper class for relperm upscaling applications. More...
#include <RelPermUtils.hpp>
Public Member Functions | |
| RelPermUpscaleHelper (int mpi_rank, std::map< std::string, std::string > &options_) | |
| Default constructor. | |
| void | collectResults () |
| Collect results from all MPI nodes. | |
| std::vector< std::vector< double > > | getRelPerm (int phase) const |
| Calculate relperm values from phase permeabilities. | |
| void | upscaleSinglePhasePermeability () |
| Upscale single phase permeability. | |
| void | sanityCheckInput (const Opm::Deck &deck, const double minPerm, const double maxPerm, const double minPoro) |
| Do sanity checks for input file. | |
| void | checkCriticalSaturations () |
| Check that input relperm curevs specify critical saturations. | |
| void | setupBoundaryConditions () |
| Setup requested boundary conditions. | |
| double | tesselateGrid (const Opm::Deck &deck) |
| Tesselate grid. | |
| void | calculateCellPressureGradients () |
| Find cell center pressure gradient for every cell. | |
| void | calculateMinMaxCapillaryPressure () |
| Calculate minimum and maximum capillary pressures. | |
| void | upscaleCapillaryPressure () |
| Upscale capillary pressure. | |
| std::tuple< double, double > | upscalePermeability (int mpi_rank) |
| Upscale permeabilities. | |
Static Public Member Functions | |
| template<class String> | |
| static Deck | parseEclipseFile (const String &eclipseFileName) |
| Form Deck object from input file (ECLIPSE format). | |
Public Attributes | |
| bool | isMaster |
| Whether this is the master MPI node or not. | |
| bool | doEclipseCheck |
| Whether to check that input relperm curves include relperm at critical saturation points. | |
| bool | anisotropic_input |
| Whether input eclipse file has diagonal anisotrophy. | |
| int | points |
| Number of saturation points to upscale for. | |
| int | tensorElementCount |
| Number of independent elements in resulting tensor. | |
| bool | upscaleBothPhases |
| Whether to upscale both phases. | |
| std::vector< double > | WaterSaturation |
| Re-upscaled water saturation for the computed pressure points. | |
| SinglePhaseUpscaler::permtensor_t | permTensor |
| Tensor of upscaled results. | |
| std::vector< int > | satnums |
| Cell satnums. | |
| std::vector< MonotCubicInterpolator > | InvJfunctions |
| Inverse of the loaded J-functions. | |
| std::array< std::array< std::vector< MonotCubicInterpolator >, 2 >, 3 > | Krfunctions |
| Relperm-curves for each (component->phase->stone type). | |
| SinglePhaseUpscaler::BoundaryConditionType | boundaryCondition |
| Boundary conditions to use. | |
| std::vector< MonotCubicInterpolator > | SwPcfunctions |
| Holds Sw(Pc) for each rocktype. | |
| std::string | saturationstring |
| Fluid system type. | |
| size_t | tesselatedCells |
| Number of "active" cells (Sintef interpretation of "active"). | |
| double | volume |
| Total volume. | |
| double | poreVolume |
| Total pore volume. | |
| std::vector< double > | pressurePoints |
| Vector of capillary pressure points between Swor and Swir. | |
Helper class for relperm upscaling applications.
| Opm::RelPermUpscaleHelper::RelPermUpscaleHelper | ( | int | mpi_rank, |
| std::map< std::string, std::string > & | options_ ) |
Default constructor.
| [in] | mpi_rank | Rank of this process (for parallel simulations). |
| [in] | options_ | Options structure. |
Uses the following options: fluids
| void Opm::RelPermUpscaleHelper::calculateCellPressureGradients | ( | ) |
Find cell center pressure gradient for every cell.
Uses the following options: gravity, waterDensity, oilDensity
| void Opm::RelPermUpscaleHelper::calculateMinMaxCapillaryPressure | ( | ) |
Calculate minimum and maximum capillary pressures.
Uses the following options: maxPermContrast, minPerm, gravity, linsolver_tolerance
| std::vector< std::vector< double > > Opm::RelPermUpscaleHelper::getRelPerm | ( | int | phase | ) | const |
Calculate relperm values from phase permeabilities.
| [in] | phase | The phase to calculate values for (0-indexed). |
First index is voigt index, second index is pressure point.
|
static |
Form Deck object from input file (ECLIPSE format).
| String | String type for representing pathnames. Typically std::string
const char*
|
| [in] | eclipseFileName | Name of input ECLIPSE file. |
| void Opm::RelPermUpscaleHelper::sanityCheckInput | ( | const Opm::Deck & | deck, |
| const double | minPerm, | ||
| const double | maxPerm, | ||
| const double | minPoro ) |
Do sanity checks for input file.
| [in] | deck | The deck to sanity check. |
| [in] | minPerm | Minimum permeability. |
| [in] | maxPerm | Maximum permeability. |
| [in] | minPoro | Minimum porosity. |
Throws error string.
| void Opm::RelPermUpscaleHelper::setupBoundaryConditions | ( | ) |
Setup requested boundary conditions.
Uses the following options: bc
| double Opm::RelPermUpscaleHelper::tesselateGrid | ( | const Opm::Deck & | deck | ) |
Tesselate grid.
| [in] | deck | The grid to tesselate. |
| [in] | options | Option structure. |
Uses the following options: linsolver_tolerance, linsolver_verbosity, linsolver_type, linsolver_max_iterations, linsolver_smooth_steps, linsolver_prolongate_factor, minPerm
| void Opm::RelPermUpscaleHelper::upscaleCapillaryPressure | ( | ) |
Upscale capillary pressure.
Uses the following options: saturationThreshold
| std::tuple< double, double > Opm::RelPermUpscaleHelper::upscalePermeability | ( | int | mpi_rank | ) |
Upscale permeabilities.
| [in] | mpi_rank | MPI rank of this process. |
Uses the following options: minPerm, maxPermContrast