|
std::pair< std::vector< double >, std::vector< double > > | computeFandPhi (const std::vector< double > &pv, const std::vector< double > &ftof, const std::vector< double > &rtof) |
| Compute flow-capacity/storage-capacity based on time-of-flight. More...
|
|
double | computeLorenz (const std::vector< double > &flowcap, const std::vector< double > &storagecap) |
| Compute the Lorenz coefficient based on the F-Phi curve. More...
|
|
std::pair< std::vector< double >, std::vector< double > > | computeSweep (const std::vector< double > &flowcap, const std::vector< double > &storagecap) |
| Compute sweep efficiency versus dimensionless time (PVI). More...
|
|
std::vector< std::tuple< int, int, double > > | computeWellPairs (const Wells &wells, const std::vector< double > &porevol, const std::vector< double > &ftracer, const std::vector< double > &btracer) |
| Compute volumes associated with injector-producer pairs. More...
|
|
void | extractParallelGridInformationToISTL (boost::any &anyComm, const UnstructuredGrid &grid) |
| Extracts the information about the data decomposition from the grid for dune-istl. More...
|
|
template<class T1 > |
auto | accumulateMaskedValues (const T1 &container, const std::vector< double > *maskContainer) -> decltype(container[0] *(*maskContainer)[0]) |
| Accumulates entries masked with 1. More...
|
|
PhaseUsage | phaseUsageFromDeck (const Opm::EclipseState &eclipseState) |
|
PhaseUsage | phaseUsageFromDeck (const Opm::Deck &deck) |
|
template<class Props > |
static void | initSaturation (const std::vector< int > &cells, const Props &props, SimulationDataContainer &state, ExtremalSat satType) |
|
template<class State > |
void | initStateBasic (const UnstructuredGrid &grid, const IncompPropertiesInterface &props, const ParameterGroup ¶m, const double gravity, State &state) |
| Initialize a twophase state from parameters. More...
|
|
template<class FaceCells , class CCI , class FCI , class State > |
void | initStateBasic (int number_of_cells, const int *global_cell, const int *cartdims, int number_of_faces, FaceCells face_cells, FCI begin_face_centroids, CCI begin_cell_centroids, int dimensions, const IncompPropertiesInterface &props, const ParameterGroup ¶m, const double gravity, State &state) |
|
template<class State > |
void | initStateBasic (const UnstructuredGrid &grid, const BlackoilPropertiesInterface &props, const ParameterGroup ¶m, const double gravity, State &state) |
| Initialize a blackoil state from parameters. More...
|
|
template<class FaceCells , class FCI , class CCI , class State > |
void | initStateBasic (int number_of_cells, const int *global_cell, const int *cartdims, int number_of_faces, FaceCells face_cells, FCI begin_face_centroids, CCI begin_cell_centroids, int dimensions, const BlackoilPropertiesInterface &props, const ParameterGroup ¶m, const double gravity, State &state) |
|
template<class Props , class State > |
void | initStateFromDeck (const UnstructuredGrid &grid, const Props &props, const EclipseState &es, const double gravity, State &state) |
|
template<class Props , class State > |
void | initBlackoilStateFromDeck (const UnstructuredGrid &grid, const Props &props, const Opm::EclipseState &es, const double gravity, State &state) |
| Initialize a blackoil state from input deck. More...
|
|
template<class FaceCells , class FCI , class CCI , class Props , class State > |
void | initBlackoilStateFromDeck (int number_of_cells, const int *global_cell, int number_of_faces, FaceCells face_cells, FCI begin_face_centroids, CCI begin_cell_centroids, int dimensions, const Props &props, const Opm::EclipseState &es, const double gravity, State &state) |
| Initialize a blackoil state from input deck. More...
|
|
template<class Props > |
static void | initSaturation (const std::vector< int > &cells, const Props &props, SimulationDataContainer &state, ExtremalSat satType) |
|
template<class Props , class State > |
void | initStateFromDeck (const UnstructuredGrid &grid, const Props &props, const Opm::EclipseState &es, const double gravity, State &state) |
| Initialize a state from input deck. More...
|
|
template<class FaceCells , class FCI , class CCI , class Props , class State > |
void | initStateFromDeck (int number_of_cells, const int *global_cell, int number_of_faces, FaceCells face_cells, FCI begin_face_centroids, CCI begin_cell_centroids, int dimensions, const Props &props, const Opm::EclipseState &es, const double gravity, State &state) |
| Initialize a state from input deck. More...
|
|
template<class Props , class State > |
void | initBlackoilSurfvol (const UnstructuredGrid &grid, const Props &props, State &state) |
|
template<class Props , class State > |
void | initBlackoilSurfvol (int number_of_cells, const Props &props, State &state) |
|
template<class Props , class State > |
void | initBlackoilSurfvolUsingRSorRV (const UnstructuredGrid &grid, const Props &props, State &state) |
|
template<class Props , class State > |
void | initBlackoilSurfvolUsingRSorRV (int number_of_cells, const Props &props, State &state) |
|
void | initHydroCarbonState (BlackoilState &state, const PhaseUsage &pu, const int num_cells, const bool has_disgas, const bool has_vapoil) |
|
void | computePorevolume (const UnstructuredGrid &grid, const double *porosity, std::vector< double > &porevol) |
| Computes pore volume of all cells in a grid. More...
|
|
template<class T > |
void | computePorevolume (int number_of_cells, T begin_cell_volume, const double *porosity, std::vector< double > &porevol) |
| Computes pore volume of all cells in a grid. More...
|
|
void | computePorevolume (const UnstructuredGrid &grid, const double *porosity, const RockCompressibility &rock_comp, const std::vector< double > &pressure, std::vector< double > &porevol) |
| Computes pore volume of all cells in a grid, with rock compressibility effects. More...
|
|
template<class T > |
void | computePorevolume (int number_of_cells, T begin_cell_volume, const double *porosity, const RockCompressibility &rock_comp, const std::vector< double > &pressure, std::vector< double > &porevol) |
| Computes pore volume of all cells in a grid, with rock compressibility effects. More...
|
|
void | computePorosity (const UnstructuredGrid &grid, const double *porosity_standard, const RockCompressibility &rock_comp, const std::vector< double > &pressure, std::vector< double > &porosity) |
| Computes porosity of all cells in a grid, with rock compressibility effects. More...
|
|
void | computeSaturatedVol (const std::vector< double > &pv, const std::vector< double > &s, double *sat_vol) |
| Computes total saturated volumes over all grid cells. More...
|
|
void | computeAverageSat (const std::vector< double > &pv, const std::vector< double > &s, double *aver_sat) |
| Computes average saturations over all grid cells. More...
|
|
void | computeInjectedProduced (const IncompPropertiesInterface &props, const std::vector< double > &s, const std::vector< double > &src, const double dt, double *injected, double *produced) |
| Computes injected and produced volumes of all phases. Note 1: assumes that only the first phase is injected. Note 2: assumes that transport has been done with an implicit method, i.e. that the current state gives the mobilities used for the preceding timestep. More...
|
|
void | computeTotalMobility (const Opm::IncompPropertiesInterface &props, const std::vector< int > &cells, const std::vector< double > &s, std::vector< double > &totmob) |
| Computes total mobility for a set of saturation values. More...
|
|
void | computeTotalMobilityOmega (const Opm::IncompPropertiesInterface &props, const std::vector< int > &cells, const std::vector< double > &s, std::vector< double > &totmob, std::vector< double > &omega) |
| Computes total mobility and omega for a set of saturation values. More...
|
|
void | computePhaseMobilities (const Opm::IncompPropertiesInterface &props, const std::vector< int > &cells, const std::vector< double > &s, std::vector< double > &pmobc) |
| Computes phase mobilities for a set of saturation values. More...
|
|
void | computeFractionalFlow (const Opm::IncompPropertiesInterface &props, const std::vector< int > &cells, const std::vector< double > &saturations, std::vector< double > &fractional_flows) |
|
void | computeTransportSource (const UnstructuredGrid &grid, const std::vector< double > &src, const std::vector< double > &faceflux, const double inflow_frac, const Wells *wells, const std::vector< double > &well_perfrates, std::vector< double > &transport_src) |
|
void | estimateCellVelocity (const UnstructuredGrid &grid, const std::vector< double > &face_flux, std::vector< double > &cell_velocity) |
| Estimates a scalar cell velocity from face fluxes. More...
|
|
template<class CC , class FC , class FC1 , class CV > |
void | estimateCellVelocity (int number_of_cells, int number_of_faces, FC begin_face_centroids, FC1 face_cells, CC begin_cell_centroids, CV begin_cell_volumes, int dimension, const std::vector< double > &face_flux, std::vector< double > &cell_velocity) |
| Estimates a scalar cell velocity from face fluxes. More...
|
|
void | toWaterSat (const std::vector< double > &sboth, std::vector< double > &sw) |
|
void | toBothSat (const std::vector< double > &sw, std::vector< double > &sboth) |
|
void | wellsToSrc (const Wells &wells, const int num_cells, std::vector< double > &src) |
|
void | computeWDP (const Wells &wells, const UnstructuredGrid &grid, const std::vector< double > &saturations, const double *densities, const double gravity, const bool per_grid_cell, std::vector< double > &wdp) |
|
template<class T > |
void | computeWDP (const Wells &wells, int number_of_cells, T begin_cell_centroids, const std::vector< double > &saturations, const double *densities, const double gravity, const bool per_grid_cell, std::vector< double > &wdp) |
|
void | computeFlowRatePerWell (const Wells &wells, const std::vector< double > &flow_rates_per_cell, std::vector< double > &flow_rates_per_well) |
|
void | computePhaseFlowRatesPerWell (const Wells &wells, const std::vector< double > &flow_rates_per_well_cell, const std::vector< double > &fractional_flows, std::vector< double > &phase_flow_per_well) |
|
void | computeInjectedProduced (const BlackoilPropertiesInterface &props, const BlackoilState &state, const std::vector< double > &transport_src, const double dt, double *injected, double *produced) |
| Computes injected and produced surface volumes of all phases. Note 1: assumes that only the first phase is injected. Note 2: assumes that transport has been done with an implicit method, i.e. that the current state gives the mobilities used for the preceding timestep. Note 3: Gives surface volume values, not reservoir volumes (as the incompressible version of the function does). Also, assumes that transport_src is given in surface volumes for injector terms! More...
|
|
void | computeTotalMobility (const Opm::BlackoilPropertiesInterface &props, const std::vector< int > &cells, const std::vector< double > &p, const std::vector< double > &z, const std::vector< double > &s, std::vector< double > &totmob) |
| Computes total mobility for a set of saturation values. More...
|
|
void | computeTotalMobilityOmega (const Opm::BlackoilPropertiesInterface &props, const std::vector< int > &cells, const std::vector< double > &p, const std::vector< double > &z, const std::vector< double > &s, std::vector< double > &totmob, std::vector< double > &omega) |
| Computes total mobility and omega for a set of saturation values. More...
|
|
void | computePhaseMobilities (const Opm::BlackoilPropertiesInterface &props, const std::vector< int > &cells, const std::vector< double > &p, const std::vector< double > &T, const std::vector< double > &z, const std::vector< double > &s, std::vector< double > &pmobc) |
| Computes phase mobilities for a set of saturation values. More...
|
|
void | computeFractionalFlow (const Opm::BlackoilPropertiesInterface &props, const std::vector< int > &cells, const std::vector< double > &p, const std::vector< double > &T, const std::vector< double > &z, const std::vector< double > &s, std::vector< double > &fractional_flows) |
|
void | computeSurfacevol (const int n, const int np, const double *A, const double *saturation, double *surfacevol) |
|
void | computeSaturation (const BlackoilPropertiesInterface &props, BlackoilState &state) |
| Computes saturation from surface volume densities. More...
|
|
void | computeTransportSource (const BlackoilPropertiesInterface &props, const Wells *wells, const WellState &well_state, std::vector< double > &transport_src) |
|
std::shared_ptr< WellsGroupInterface > | createWellWellsGroup (const Well *well, size_t timeStep, const PhaseUsage &phase_usage) |
|
std::shared_ptr< WellsGroupInterface > | createGroupWellsGroup (const Group &group, size_t timeStep, const PhaseUsage &phase_usage) |
|
◆ ExtremalSat
Will initialize the first and second component of the SATURATION field in all the cells in the set @cells. The @props object will be queried, and depending on the value @satType either the minimum or the maximum saturation is applied to thee first component in the SATURATION field. For the second component (1 - first_sat) is used.
◆ HydroCarbonState
Enumerator |
---|
GasOnly | |
GasAndOil | |
OilOnly | |
◆ accumulateMaskedValues()
template<class T1 >
auto Opm::accumulateMaskedValues |
( |
const T1 & |
container, |
|
|
const std::vector< double > * |
maskContainer |
|
) |
| -> decltype(container[0]*(*maskContainer)[0])
|
Accumulates entries masked with 1.
- Parameters
-
container | The container whose values to accumulate. |
maskContainer | null pointer or a pointer to a container with entries 0 and 1. Only values at indices with a 1 stored will be accumulated. If null then all values will be accumulated |
- Returns
- the summ of all entries that should be represented.
◆ computeAverageSat()
void Opm::computeAverageSat |
( |
const std::vector< double > & |
pv, |
|
|
const std::vector< double > & |
s, |
|
|
double * |
aver_sat |
|
) |
| |
Computes average saturations over all grid cells.
- Parameters
-
[in] | pv | the pore volume by cell. |
[in] | s | saturation values (for all P phases) |
[out] | aver_sat | must point to a valid array with P elements, where P = s.size()/pv.size(). For each phase p, we compute aver_sat_p = (sum_i s_p_i pv_i) / (sum_i pv_i). |
◆ computeFandPhi()
std::pair< std::vector< double >, std::vector< double > > Opm::computeFandPhi |
( |
const std::vector< double > & |
pv, |
|
|
const std::vector< double > & |
ftof, |
|
|
const std::vector< double > & |
rtof |
|
) |
| |
Compute flow-capacity/storage-capacity based on time-of-flight.
The F-Phi curve is an analogue to the fractional flow curve in a 1D displacement. It can be used to compute other interesting diagnostic quantities such as the Lorenz coefficient. For a technical description see Shavali et al. (SPE 146446), Shook and Mitchell (SPE 124625).
- Parameters
-
[in] | pv | pore volumes of each cell |
[in] | ftof | forward (time from injector) time-of-flight values for each cell |
[in] | rtof | reverse (time to producer) time-of-flight values for each cell |
- Returns
- a pair of vectors, the first containing F (flow capacity) the second containing Phi (storage capacity).
◆ computeFlowRatePerWell()
void Opm::computeFlowRatePerWell |
( |
const Wells & |
wells, |
|
|
const std::vector< double > & |
flow_rates_per_cell, |
|
|
std::vector< double > & |
flow_rates_per_well |
|
) |
| |
Computes (sums) the flow rate for each well. - Parameters
-
[in] | wells | The wells for which the flow rate should be computed. |
[in] | flow_rates_per_cell | Flow rates per well cells. Should ordered the same way as wells. |
[out] | flow_rates_per_well | Will contain the summed up flow_rates for each well. |
◆ computeFractionalFlow() [1/2]
void Opm::computeFractionalFlow |
( |
const Opm::BlackoilPropertiesInterface & |
props, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
p, |
|
|
const std::vector< double > & |
T, |
|
|
const std::vector< double > & |
z, |
|
|
const std::vector< double > & |
s, |
|
|
std::vector< double > & |
fractional_flows |
|
) |
| |
Computes the fractional flow for each cell in the cells argument - Parameters
-
[in] | props | rock and fluid properties |
[in] | cells | cells with which the saturation values are associated |
[in] | p | pressure (one value per cell) |
[in] | T | temperature(one value per cell) |
[in] | z | surface-volume values (for all P phases) |
[in] | s | saturation values (for all phases) |
[out] | fractional_flow | the fractional flow for each phase for each cell. |
◆ computeFractionalFlow() [2/2]
void Opm::computeFractionalFlow |
( |
const Opm::IncompPropertiesInterface & |
props, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
saturations, |
|
|
std::vector< double > & |
fractional_flows |
|
) |
| |
Computes the fractional flow for each cell in the cells argument - Parameters
-
[in] | props | rock and fluid properties |
[in] | cells | cells with which the saturation values are associated |
[in] | saturations | saturation values (for all phases) |
[out] | fractional_flow | the fractional flow for each phase for each cell. |
◆ computeInjectedProduced() [1/2]
void Opm::computeInjectedProduced |
( |
const BlackoilPropertiesInterface & |
props, |
|
|
const BlackoilState & |
state, |
|
|
const std::vector< double > & |
transport_src, |
|
|
const double |
dt, |
|
|
double * |
injected, |
|
|
double * |
produced |
|
) |
| |
Computes injected and produced surface volumes of all phases. Note 1: assumes that only the first phase is injected. Note 2: assumes that transport has been done with an implicit method, i.e. that the current state gives the mobilities used for the preceding timestep. Note 3: Gives surface volume values, not reservoir volumes (as the incompressible version of the function does). Also, assumes that transport_src is given in surface volumes for injector terms!
- Parameters
-
[in] | props | fluid and rock properties. |
[in] | state | state variables (pressure, sat, surfvol) |
[in] | transport_src | if < 0: total resv outflow, if > 0: first phase surfv inflow |
[in] | dt | timestep used |
[out] | injected | must point to a valid array with P elements, where P = s.size()/src.size(). |
[out] | produced | must also point to a valid array with P elements. |
◆ computeInjectedProduced() [2/2]
void Opm::computeInjectedProduced |
( |
const IncompPropertiesInterface & |
props, |
|
|
const std::vector< double > & |
s, |
|
|
const std::vector< double > & |
src, |
|
|
const double |
dt, |
|
|
double * |
injected, |
|
|
double * |
produced |
|
) |
| |
Computes injected and produced volumes of all phases. Note 1: assumes that only the first phase is injected. Note 2: assumes that transport has been done with an implicit method, i.e. that the current state gives the mobilities used for the preceding timestep.
- Parameters
-
[in] | props | fluid and rock properties. |
[in] | s | saturation values (for all P phases) |
[in] | src | if < 0: total outflow, if > 0: first phase inflow. |
[in] | dt | timestep used |
[out] | injected | must point to a valid array with P elements, where P = s.size()/src.size(). |
[out] | produced | must also point to a valid array with P elements. |
◆ computeLorenz()
double Opm::computeLorenz |
( |
const std::vector< double > & |
flowcap, |
|
|
const std::vector< double > & |
storagecap |
|
) |
| |
Compute the Lorenz coefficient based on the F-Phi curve.
The Lorenz coefficient is a measure of heterogeneity. It is equal to twice the area between the F-Phi curve and the F = Phi line. The coefficient can vary from zero to one. If the coefficient is zero (so the F-Phi curve is a straight line) we have perfect piston-like displacement while a coefficient of one indicates infinitely heterogenous displacement (essentially no sweep).
Note: The coefficient is analogous to the Gini coefficient of economic theory, where the name Lorenz curve is applied to what we call the F-Phi curve.
- Parameters
-
- Returns
- the Lorenz coefficient
◆ computePhaseFlowRatesPerWell()
void Opm::computePhaseFlowRatesPerWell |
( |
const Wells & |
wells, |
|
|
const std::vector< double > & |
flow_rates_per_well_cell, |
|
|
const std::vector< double > & |
fractional_flows, |
|
|
std::vector< double > & |
phase_flow_per_well |
|
) |
| |
Computes the phase flow rate per well - Parameters
-
[in] | wells | The wells for which the flow rate should be computed |
[in] | flow_rates_per_well_cell | The total flow rate for each cell (ordered the same way as the wells struct |
[in] | fractional_flows | the fractional flow for each cell in each well |
[out] | phase_flow_per_well | Will contain the phase flow per well |
◆ computePhaseMobilities() [1/2]
void Opm::computePhaseMobilities |
( |
const Opm::BlackoilPropertiesInterface & |
props, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
p, |
|
|
const std::vector< double > & |
T, |
|
|
const std::vector< double > & |
z, |
|
|
const std::vector< double > & |
s, |
|
|
std::vector< double > & |
pmobc |
|
) |
| |
Computes phase mobilities for a set of saturation values.
- Parameters
-
[in] | props | rock and fluid properties |
[in] | cells | cells with which the saturation values are associated |
[in] | p | pressure (one value per cell) |
[in] | T | temperature (one value per cell) |
[in] | z | surface-volume values (for all P phases) |
[in] | s | saturation values (for all phases) |
[out] | pmobc | phase mobilities (for all phases). |
◆ computePhaseMobilities() [2/2]
void Opm::computePhaseMobilities |
( |
const Opm::IncompPropertiesInterface & |
props, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
s, |
|
|
std::vector< double > & |
pmobc |
|
) |
| |
Computes phase mobilities for a set of saturation values.
- Parameters
-
[in] | props | rock and fluid properties |
[in] | cells | cells with which the saturation values are associated |
[in] | s | saturation values (for all phases) |
[out] | pmobc | phase mobilities (for all phases). |
◆ computePorevolume() [1/4]
void Opm::computePorevolume |
( |
const UnstructuredGrid & |
grid, |
|
|
const double * |
porosity, |
|
|
const RockCompressibility & |
rock_comp, |
|
|
const std::vector< double > & |
pressure, |
|
|
std::vector< double > & |
porevol |
|
) |
| |
Computes pore volume of all cells in a grid, with rock compressibility effects.
- Parameters
-
[in] | grid | a grid |
[in] | porosity | array of grid.number_of_cells porosity values (at reference pressure) |
[in] | rock_comp | rock compressibility properties |
[in] | pressure | pressure by cell |
[out] | porevol | the pore volume by cell. |
◆ computePorevolume() [2/4]
void Opm::computePorevolume |
( |
const UnstructuredGrid & |
grid, |
|
|
const double * |
porosity, |
|
|
std::vector< double > & |
porevol |
|
) |
| |
Computes pore volume of all cells in a grid.
- Parameters
-
[in] | grid | a grid |
[in] | porosity | array of grid.number_of_cells porosity values |
[out] | porevol | the pore volume by cell. |
◆ computePorevolume() [3/4]
template<class T >
void Opm::computePorevolume |
( |
int |
number_of_cells, |
|
|
T |
begin_cell_volumes, |
|
|
const double * |
porosity, |
|
|
const RockCompressibility & |
rock_comp, |
|
|
const std::vector< double > & |
pressure, |
|
|
std::vector< double > & |
porevol |
|
) |
| |
Computes pore volume of all cells in a grid, with rock compressibility effects.
- Parameters
-
[in] | number_of_cells | The number of cells of the grid. |
[in] | Pointer | to/ Iterator at the first cell volume. |
[in] | porosity | array of grid.number_of_cells porosity values |
[in] | rock_comp | rock compressibility properties |
[in] | pressure | pressure by cell |
[out] | porevol | the pore volume by cell. |
[in] | number_of_cells | The number of cells of the grid. |
[in] | porosity | array of grid.number_of_cells porosity values |
[in] | rock_comp | rock compressibility properties |
[in] | pressure | pressure by cell |
[out] | porevol | the pore volume by cell. |
References Opm::RockCompressibility::poroMult().
◆ computePorevolume() [4/4]
template<class T >
void Opm::computePorevolume |
( |
int |
number_of_cells, |
|
|
T |
begin_cell_volume, |
|
|
const double * |
porosity, |
|
|
std::vector< double > & |
porevol |
|
) |
| |
Computes pore volume of all cells in a grid.
- Parameters
-
[in] | number_of_cells | The number of cells of the grid. |
[in] | begin_cell_volume | Iterator to the volume of the first cell. |
[in] | porosity | array of grid.number_of_cells porosity values |
[out] | porevol | the pore volume by cell. |
◆ computePorosity()
void Opm::computePorosity |
( |
const UnstructuredGrid & |
grid, |
|
|
const double * |
porosity_standard, |
|
|
const RockCompressibility & |
rock_comp, |
|
|
const std::vector< double > & |
pressure, |
|
|
std::vector< double > & |
porosity |
|
) |
| |
Computes porosity of all cells in a grid, with rock compressibility effects.
- Parameters
-
[in] | grid | a grid |
[in] | porosity_standard | array of grid.number_of_cells porosity values (at reference presure) |
[in] | rock_comp | rock compressibility properties |
[in] | pressure | pressure by cell |
[out] | porosity | porosity (at reservoir condition) |
◆ computeSaturatedVol()
void Opm::computeSaturatedVol |
( |
const std::vector< double > & |
pv, |
|
|
const std::vector< double > & |
s, |
|
|
double * |
sat_vol |
|
) |
| |
Computes total saturated volumes over all grid cells.
- Parameters
-
[in] | pv | the pore volume by cell. |
[in] | s | saturation values (for all P phases) |
[out] | sat_vol | must point to a valid array with P elements, where P = s.size()/pv.size(). For each phase p, we compute sat_vol_p = sum_i s_p_i pv_i |
◆ computeSaturation()
◆ computeSurfacevol()
void Opm::computeSurfacevol |
( |
const int |
n, |
|
|
const int |
np, |
|
|
const double * |
A, |
|
|
const double * |
saturation, |
|
|
double * |
surfacevol |
|
) |
| |
Computes the surface volume densities from saturations by the formula z = A s for a number of data points, where z is the surface volume density, s is the saturation (both as column vectors) and A is the phase-to-component relation matrix. - Parameters
-
[in] | n | number of data points |
[in] | np | number of phases, must be 2 or 3 |
[in] | A | array containing n square matrices of size num_phases, in Fortran ordering, typically the output of a call to the matrix() method of a BlackoilProperties* class. |
[in] | saturation | concatenated saturation values (for all P phases) |
[out] | surfacevol | concatenated surface-volume values (for all P phases) |
◆ computeSweep()
std::pair< std::vector< double >, std::vector< double > > Opm::computeSweep |
( |
const std::vector< double > & |
flowcap, |
|
|
const std::vector< double > & |
storagecap |
|
) |
| |
Compute sweep efficiency versus dimensionless time (PVI).
The sweep efficiency is analogue to 1D displacement using the F-Phi curve as flux function.
- Parameters
-
- Returns
- a pair of vectors, the first containing Ev (sweep efficiency) the second containing tD (dimensionless time).
◆ computeTotalMobility() [1/2]
void Opm::computeTotalMobility |
( |
const Opm::BlackoilPropertiesInterface & |
props, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
p, |
|
|
const std::vector< double > & |
z, |
|
|
const std::vector< double > & |
s, |
|
|
std::vector< double > & |
totmob |
|
) |
| |
Computes total mobility for a set of saturation values.
- Parameters
-
[in] | props | rock and fluid properties |
[in] | cells | cells with which the saturation values are associated |
[in] | p | pressure (one value per cell) |
[in] | z | surface-volume values (for all P phases) |
[in] | s | saturation values (for all phases) |
[out] | totmob | total mobilities. |
◆ computeTotalMobility() [2/2]
void Opm::computeTotalMobility |
( |
const Opm::IncompPropertiesInterface & |
props, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
s, |
|
|
std::vector< double > & |
totmob |
|
) |
| |
Computes total mobility for a set of saturation values.
- Parameters
-
[in] | props | rock and fluid properties |
[in] | cells | cells with which the saturation values are associated |
[in] | s | saturation values (for all phases) |
[out] | totmob | total mobilities. |
◆ computeTotalMobilityOmega() [1/2]
void Opm::computeTotalMobilityOmega |
( |
const Opm::BlackoilPropertiesInterface & |
props, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
p, |
|
|
const std::vector< double > & |
z, |
|
|
const std::vector< double > & |
s, |
|
|
std::vector< double > & |
totmob, |
|
|
std::vector< double > & |
omega |
|
) |
| |
Computes total mobility and omega for a set of saturation values.
- Parameters
-
[in] | props | rock and fluid properties |
[in] | cells | cells with which the saturation values are associated |
[in] | p | pressure (one value per cell) |
[in] | z | surface-volume values (for all P phases) |
[in] | s | saturation values (for all phases) |
[out] | totmob | total mobility |
[out] | omega | fractional-flow weighted fluid densities. |
◆ computeTotalMobilityOmega() [2/2]
void Opm::computeTotalMobilityOmega |
( |
const Opm::IncompPropertiesInterface & |
props, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
s, |
|
|
std::vector< double > & |
totmob, |
|
|
std::vector< double > & |
omega |
|
) |
| |
Computes total mobility and omega for a set of saturation values.
- Parameters
-
[in] | props | rock and fluid properties |
[in] | cells | cells with which the saturation values are associated |
[in] | s | saturation values (for all phases) |
[out] | totmob | total mobility |
[out] | omega | fractional-flow weighted fluid densities. |
◆ computeTransportSource() [1/2]
Compute two-phase transport source terms from well terms. Note: Unlike the incompressible version of this function, this version computes surface volume injection rates, production rates are still total reservoir volumes. - Parameters
-
[in] | props | Fluid and rock properties. |
[in] | wells | Wells data structure. |
[in] | well_state | Well pressures and fluxes. |
[out] | transport_src | The transport source terms. They are to be interpreted depending on sign: (+) positive inflow of first (water) phase (surface volume), (-) negative total outflow of both phases (reservoir volume). |
◆ computeTransportSource() [2/2]
void Opm::computeTransportSource |
( |
const UnstructuredGrid & |
grid, |
|
|
const std::vector< double > & |
src, |
|
|
const std::vector< double > & |
faceflux, |
|
|
const double |
inflow_frac, |
|
|
const Wells * |
wells, |
|
|
const std::vector< double > & |
well_perfrates, |
|
|
std::vector< double > & |
transport_src |
|
) |
| |
Compute two-phase transport source terms from face fluxes, and pressure equation source terms. This puts boundary flows into the source terms for the transport equation. - Parameters
-
[in] | grid | The grid used. |
[in] | src | Pressure eq. source terms. The sign convention is: (+) positive total inflow (positive velocity divergence) (-) negative total outflow |
[in] | faceflux | Signed face fluxes, typically the result from a flow solver. |
[in] | inflow_frac | Fraction of inflow (boundary and source terms) that consists of first phase. Example: if only water is injected, inflow_frac == 1.0. Note: it is not possible (with this method) to use different fractions for different inflow sources, be they source terms of boundary flows. |
[in] | wells | Wells data structure, or null if no wells. |
[in] | well_perfrates | Volumetric flow rates per well perforation. |
[out] | transport_src | The transport source terms. They are to be interpreted depending on sign: (+) positive inflow of first phase (water) (-) negative total outflow of both phases |
◆ computeWDP() [1/2]
void Opm::computeWDP |
( |
const Wells & |
wells, |
|
|
const UnstructuredGrid & |
grid, |
|
|
const std::vector< double > & |
saturations, |
|
|
const double * |
densities, |
|
|
const double |
gravity, |
|
|
const bool |
per_grid_cell, |
|
|
std::vector< double > & |
wdp |
|
) |
| |
Computes the WDP for each well. - Parameters
-
[in] | wells | Wells that need their wdp calculated. |
[in] | grid | The associated grid to make cell lookups. |
[in] | saturations | A vector of weights for each cell for each phase in the grid (or well, see per_grid_cell parameter). So for cell i, saturations[i*densities.size() + p] should give the weight of phase p in cell i. |
[in] | densities | Density for each phase. |
[out] | wdp | Will contain, for each well, the wdp of the well. |
[in] | per_grid_cell | Whether or not the saturations are per grid cell or per well cell. |
◆ computeWDP() [2/2]
template<class T >
void Opm::computeWDP |
( |
const Wells & |
wells, |
|
|
int |
number_of_cells, |
|
|
T |
begin_cell_centroids, |
|
|
const std::vector< double > & |
saturations, |
|
|
const double * |
densities, |
|
|
const double |
gravity, |
|
|
const bool |
per_grid_cell, |
|
|
std::vector< double > & |
wdp |
|
) |
| |
Computes the WDP for each well. - Parameters
-
[in] | wells | Wells that need their wdp calculated. |
[in] | number_of_cells | The number of cells in the grid. |
[in] | begin_cell_centroids | Pointer/Iterator to the first cell centroid. |
[in] | saturations | A vector of weights for each cell for each phase in the grid (or well, see per_grid_cell parameter). So for cell i, saturations[i*densities.size() + p] should give the weight of phase p in cell i. |
[in] | densities | Density for each phase. |
[out] | wdp | Will contain, for each well, the wdp of the well. |
[in] | per_grid_cell | Whether or not the saturations are per grid cell or per well cell. |
References Wells::depth_ref, Wells::number_of_wells, Wells::well_cells, and Wells::well_connpos.
◆ computeWellPairs()
std::vector< std::tuple< int, int, double > > Opm::computeWellPairs |
( |
const Wells & |
wells, |
|
|
const std::vector< double > & |
porevol, |
|
|
const std::vector< double > & |
ftracer, |
|
|
const std::vector< double > & |
btracer |
|
) |
| |
Compute volumes associated with injector-producer pairs.
- Parameters
-
[in] | wells | wells structure, containing NI injector wells and NP producer wells. |
[in] | porevol | pore volume of each grid cell |
[in] | ftracer | array of forward (injector) tracer values, NI per cell |
[in] | btracer | array of backward (producer) tracer values, NP per cell |
- Returns
- a vector of tuples, one tuple for each injector-producer pair, where the first and second elements are well indices for the injector and producer, and the third element is the pore volume associated with that pair.
◆ createGroupWellsGroup()
Creates the WellsGroupInterface for the given Group - Parameters
-
[in] | group | the Group to construct object for |
[in] | timeStep | the time step in question |
[in] | the | phase usage |
◆ createWellWellsGroup()
Creates the WellsGroupInterface for the given well - Parameters
-
[in] | well | the Well to construct object for |
[in] | timeStep | the time step in question |
[in] | the | phase usage |
◆ estimateCellVelocity() [1/2]
void Opm::estimateCellVelocity |
( |
const UnstructuredGrid & |
grid, |
|
|
const std::vector< double > & |
face_flux, |
|
|
std::vector< double > & |
cell_velocity |
|
) |
| |
Estimates a scalar cell velocity from face fluxes.
- Parameters
-
[in] | grid | a grid |
[in] | face_flux | signed per-face fluxes |
[out] | cell_velocity | the estimated velocities. |
◆ estimateCellVelocity() [2/2]
template<class CC , class FC , class FC1 , class CV >
void Opm::estimateCellVelocity |
( |
int |
number_of_cells, |
|
|
int |
number_of_faces, |
|
|
FC |
begin_face_centroids, |
|
|
FC1 |
face_cells, |
|
|
CC |
begin_cell_centroids, |
|
|
CV |
begin_cell_volumes, |
|
|
int |
dimension, |
|
|
const std::vector< double > & |
face_flux, |
|
|
std::vector< double > & |
cell_velocity |
|
) |
| |
Estimates a scalar cell velocity from face fluxes.
- Parameters
-
[in] | number_of_cells | The number of cells of the grid |
[in] | number_of_faces | The number of cells of the grid |
[in] | begin_face_centroids | Iterator pointing to first face centroid. |
[in] | face_cells | Mapping from faces to connected cells. |
[in] | dimensions | The dimensions of the grid. |
[in] | begin_cell_centroids | Iterator pointing to first cell centroid. |
[in] | face_flux | signed per-face fluxes |
[out] | cell_velocity | the estimated velocities. |
[in] | number_of_cells | The number of cells of the grid |
[in] | begin_face_centroids | Iterator pointing to first face centroid. |
[in] | face_cells | Mapping from faces to connected cells. |
[in] | dimensions | The dimensions of the grid. |
[in] | begin_cell_centroids | Iterator pointing to first cell centroid. |
[in] | face_flux | signed per-face fluxes |
[out] | cell_velocity | the estimated velocities. |
◆ extractParallelGridInformationToISTL()
void Opm::extractParallelGridInformationToISTL |
( |
boost::any & |
anyComm, |
|
|
const UnstructuredGrid & |
grid |
|
) |
| |
|
inline |
Extracts the information about the data decomposition from the grid for dune-istl.
In the case that grid is a parallel grid this method will query it to get the information about the data decompoisition and convert it to the format expected by the linear algebra of dune-istl. \warn for UnstructuredGrid this function doesn't do anything. - Parameters
-
anyComm | The handle to store the information in. If grid is a parallel grid then this will ecapsulate an instance of ParallelISTLInformation. |
grid | The grid to inspect. |
◆ initBlackoilStateFromDeck() [1/2]
template<class Props , class State >
void Opm::initBlackoilStateFromDeck |
( |
const UnstructuredGrid & |
grid, |
|
|
const Props & |
props, |
|
|
const Opm::EclipseState & |
es, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
Initialize a blackoil state from input deck.
Initialize a two-phase water-oil blackoil state from input deck. If EQUIL is present:
- saturation is set according to the water-oil contact,
- pressure is set to hydrostatic equilibrium. Otherwise:
- saturation is set according to SWAT,
- pressure is set according to PRESSURE. In addition, this function sets surfacevol.
References initBlackoilStateFromDeck().
Referenced by initBlackoilStateFromDeck().
◆ initBlackoilStateFromDeck() [2/2]
template<class FaceCells , class FCI , class CCI , class Props , class State >
void Opm::initBlackoilStateFromDeck |
( |
int |
number_of_cells, |
|
|
const int * |
global_cell, |
|
|
int |
number_of_faces, |
|
|
FaceCells |
face_cells, |
|
|
FCI |
begin_face_centroids, |
|
|
CCI |
begin_cell_centroids, |
|
|
int |
dimensions, |
|
|
const Props & |
props, |
|
|
const Opm::EclipseState & |
es, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
◆ initBlackoilSurfvol() [1/2]
template<class Props , class State >
void Opm::initBlackoilSurfvol |
( |
const UnstructuredGrid & |
grid, |
|
|
const Props & |
props, |
|
|
State & |
state |
|
) |
| |
Initialize surface volume from pressure and saturation by z = As. Here saturation is used as an initial guess for z in the computation of A.
References initBlackoilSurfvol().
Referenced by initBlackoilSurfvol().
◆ initBlackoilSurfvol() [2/2]
template<class Props , class State >
void Opm::initBlackoilSurfvol |
( |
int |
number_of_cells, |
|
|
const Props & |
props, |
|
|
State & |
state |
|
) |
| |
◆ initBlackoilSurfvolUsingRSorRV() [1/2]
template<class Props , class State >
void Opm::initBlackoilSurfvolUsingRSorRV |
( |
const UnstructuredGrid & |
grid, |
|
|
const Props & |
props, |
|
|
State & |
state |
|
) |
| |
◆ initBlackoilSurfvolUsingRSorRV() [2/2]
template<class Props , class State >
void Opm::initBlackoilSurfvolUsingRSorRV |
( |
int |
number_of_cells, |
|
|
const Props & |
props, |
|
|
State & |
state |
|
) |
| |
◆ initHydroCarbonState()
void Opm::initHydroCarbonState |
( |
BlackoilState & |
state, |
|
|
const PhaseUsage & |
pu, |
|
|
const int |
num_cells, |
|
|
const bool |
has_disgas, |
|
|
const bool |
has_vapoil |
|
) |
| |
|
inline |
◆ initSaturation() [1/2]
template<class Props >
static void Opm::initSaturation |
( |
const std::vector< int > & |
cells, |
|
|
const Props & |
props, |
|
|
SimulationDataContainer & |
state, |
|
|
ExtremalSat |
satType |
|
) |
| |
|
static |
◆ initSaturation() [2/2]
template<class Props >
static void Opm::initSaturation |
( |
const std::vector< int > & |
cells, |
|
|
const Props & |
props, |
|
|
SimulationDataContainer & |
state, |
|
|
ExtremalSat |
satType |
|
) |
| |
|
static |
◆ initStateBasic() [1/4]
template<class State >
void Opm::initStateBasic |
( |
const UnstructuredGrid & |
grid, |
|
|
const BlackoilPropertiesInterface & |
props, |
|
|
const ParameterGroup & |
param, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
Initialize a blackoil state from parameters.
Initialize a blackoil state from parameters. The following parameters are accepted (defaults):
- convection_testcase (false) – Water in the 'left' part of the grid.
- ref_pressure (100) – Initial pressure in bar for all cells (if convection_testcase is true), or pressure at woc depth.
- water_oil_contact (none) – Depth of water-oil contact (woc). If convection_testcase is true, the saturation is initialised as indicated, and pressure is initialised to a constant value ('ref_pressure'). Otherwise we have 2 cases:
If 'water_oil_contact' is given, saturation is initialised accordingly.
- Water saturation is set to minimum. In both cases, pressure is initialised hydrostatically. In case 2., the depth of the first cell is used as reference depth.
References initStateBasic().
◆ initStateBasic() [2/4]
template<class State >
void Opm::initStateBasic |
( |
const UnstructuredGrid & |
grid, |
|
|
const IncompPropertiesInterface & |
props, |
|
|
const ParameterGroup & |
param, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
Initialize a twophase state from parameters.
Initialize a two-phase state from parameters. The following parameters are accepted (defaults):
- convection_testcase (false) – Water in the 'left' part of the grid.
- ref_pressure (100) – Initial pressure in bar for all cells (if convection_testcase is true), or pressure at woc depth.
- segregation_testcase (false) – Water above the woc instead of below.
- water_oil_contact (none) – Depth of water-oil contact (woc).
- init_saturation (none) – Initial water saturation for all cells.
If convection_testcase is true, the saturation is initialised as indicated, and pressure is initialised to a constant value ('ref_pressure'). If segregation_testcase is true, the saturation is initialised as indicated, and pressure is initialised hydrostatically. Otherwise we have 3 cases:
- If 'water_oil_contact' is given, saturation is initialised accordingly.
- If 'water_oil_contact' is not given, but 'init_saturation' is given, water saturation is set to that value everywhere.
- If neither are given, water saturation is set to minimum.
In all three cases, pressure is initialised hydrostatically. In case 2) and 3), the depth of the first cell is used as reference depth.
References initStateBasic().
Referenced by initStateBasic().
◆ initStateBasic() [3/4]
template<class FaceCells , class FCI , class CCI , class State >
void Opm::initStateBasic |
( |
int |
number_of_cells, |
|
|
const int * |
global_cell, |
|
|
const int * |
cartdims, |
|
|
int |
number_of_faces, |
|
|
FaceCells |
face_cells, |
|
|
FCI |
begin_face_centroids, |
|
|
CCI |
begin_cell_centroids, |
|
|
int |
dimensions, |
|
|
const BlackoilPropertiesInterface & |
props, |
|
|
const ParameterGroup & |
param, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
Initialize a blackoil state from parameters. The following parameters are accepted (defaults):
- convection_testcase (false) – Water in the 'left' part of the grid.
- ref_pressure (100) – Initial pressure in bar for all cells (if convection_testcase is true), or pressure at woc depth.
- water_oil_contact (none) – Depth of water-oil contact (woc). If convection_testcase is true, the saturation is initialised as indicated, and pressure is initialised to a constant value ('ref_pressure'). Otherwise we have 2 cases:
If 'water_oil_contact' is given, saturation is initialised accordingly.
- Water saturation is set to minimum. In both cases, pressure is initialised hydrostatically. In case 2., the depth of the first cell is used as reference depth.
References initSaturation(), MaxSat, MinSat, Opm::BlackoilPropertiesInterface::numCells(), and Opm::BlackoilPropertiesInterface::numPhases().
◆ initStateBasic() [4/4]
template<class FaceCells , class CCI , class FCI , class State >
void Opm::initStateBasic |
( |
int |
number_of_cells, |
|
|
const int * |
global_cell, |
|
|
const int * |
cartdims, |
|
|
int |
number_of_faces, |
|
|
FaceCells |
face_cells, |
|
|
FCI |
begin_face_centroids, |
|
|
CCI |
begin_cell_centroids, |
|
|
int |
dimensions, |
|
|
const IncompPropertiesInterface & |
props, |
|
|
const ParameterGroup & |
param, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
Initialize a two-phase state from parameters. The following parameters are accepted (defaults):
- convection_testcase (false) – Water in the 'left' part of the grid.
- ref_pressure (100) – Initial pressure in bar for all cells (if convection_testcase is true), or pressure at woc depth.
- segregation_testcase (false) – Water above the woc instead of below.
- water_oil_contact (none) – Depth of water-oil contact (woc).
- init_saturation (none) – Initial water saturation for all cells.
If convection_testcase is true, the saturation is initialised as indicated, and pressure is initialised to a constant value ('ref_pressure'). If segregation_testcase is true, the saturation is initialised as indicated, and pressure is initialised hydrostatically. Otherwise we have 3 cases:
- If 'water_oil_contact' is given, saturation is initialised accordingly.
- If 'water_oil_contact' is not given, but 'init_saturation' is given, water saturation is set to that value everywhere.
- If neither are given, water saturation is set to minimum.
In all three cases, pressure is initialised hydrostatically. In case 2) and 3), the depth of the first cell is used as reference depth.
References Opm::IncompPropertiesInterface::density(), initSaturation(), MaxSat, MinSat, Opm::IncompPropertiesInterface::numCells(), and Opm::IncompPropertiesInterface::numPhases().
◆ initStateFromDeck() [1/3]
template<class Props , class State >
void Opm::initStateFromDeck |
( |
const UnstructuredGrid & |
grid, |
|
|
const Props & |
props, |
|
|
const EclipseState & |
es, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
Initialize a two-phase state from input deck. If EQUIL is present:
- saturation is set according to the water-oil contact,
- pressure is set to hydrostatic equilibrium. Otherwise:
- saturation is set according to SWAT,
- pressure is set according to PRESSURE.
Referenced by initBlackoilStateFromDeck(), and initStateFromDeck().
◆ initStateFromDeck() [2/3]
template<class Props , class State >
void Opm::initStateFromDeck |
( |
const UnstructuredGrid & |
grid, |
|
|
const Props & |
props, |
|
|
const Opm::EclipseState & |
es, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
◆ initStateFromDeck() [3/3]
template<class FaceCells , class FCI , class CCI , class Props , class State >
void Opm::initStateFromDeck |
( |
int |
number_of_cells, |
|
|
const int * |
global_cell, |
|
|
int |
number_of_faces, |
|
|
FaceCells |
face_cells, |
|
|
FCI |
begin_face_centroids, |
|
|
CCI |
begin_cell_centroids, |
|
|
int |
dimensions, |
|
|
const Props & |
props, |
|
|
const Opm::EclipseState & |
es, |
|
|
const double |
gravity, |
|
|
State & |
state |
|
) |
| |
◆ phaseUsageFromDeck() [1/2]
PhaseUsage Opm::phaseUsageFromDeck |
( |
const Opm::Deck & |
deck | ) |
|
|
inline |
Looks at presence of WATER, OIL and GAS keywords in deck to determine active phases.
References Opm::BlackoilPhases::Aqua, Opm::BlackoilPhases::Energy, GAS, Opm::PhaseUsage::has_energy, Opm::PhaseUsage::has_polymer, Opm::PhaseUsage::has_solvent, Opm::BlackoilPhases::Liquid, Opm::BlackoilPhases::MaxNumPhases, Opm::PhaseUsage::num_phases, Opm::BlackoilPhases::NumCryptoPhases, OIL, Opm::PhaseUsage::phase_pos, Opm::PhaseUsage::phase_used, Opm::BlackoilPhases::Polymer, Opm::BlackoilPhases::Solvent, Opm::BlackoilPhases::Vapour, and WATER.
◆ phaseUsageFromDeck() [2/2]
PhaseUsage Opm::phaseUsageFromDeck |
( |
const Opm::EclipseState & |
eclipseState | ) |
|
|
inline |
Looks at presence of WATER, OIL and GAS keywords in state object to determine active phases.
References Opm::BlackoilPhases::Aqua, Opm::BlackoilPhases::Energy, GAS, Opm::PhaseUsage::has_energy, Opm::PhaseUsage::has_polymer, Opm::PhaseUsage::has_solvent, Opm::BlackoilPhases::Liquid, Opm::BlackoilPhases::MaxNumPhases, Opm::PhaseUsage::num_phases, Opm::BlackoilPhases::NumCryptoPhases, OIL, Opm::PhaseUsage::phase_pos, Opm::PhaseUsage::phase_used, Opm::BlackoilPhases::Polymer, Opm::BlackoilPhases::Solvent, Opm::BlackoilPhases::Vapour, and WATER.
Referenced by Opm::SaturationPropsFromDeck::init(), and initStateFromDeck().
◆ toBothSat()
void Opm::toBothSat |
( |
const std::vector< double > & |
sw, |
|
|
std::vector< double > & |
sboth |
|
) |
| |
Make a vector of interleaved water and oil saturations from a vector of water saturations.
◆ toWaterSat()
void Opm::toWaterSat |
( |
const std::vector< double > & |
sboth, |
|
|
std::vector< double > & |
sw |
|
) |
| |
Extract a vector of water saturations from a vector of interleaved water and oil saturations.
◆ wellsToSrc()
void Opm::wellsToSrc |
( |
const Wells & |
wells, |
|
|
const int |
num_cells, |
|
|
std::vector< double > & |
src |
|
) |
| |
Create a src vector equivalent to a wells structure. For this to be valid, the wells must be all rate-controlled and single-perforation.
|