|
void | computeTotalMobility (const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &totmob) |
| Computes total mobility for a set of s/c values. More...
|
|
void | computeTotalMobilityOmega (const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &totmob, std::vector< double > &omega) |
| Computes total mobility and omega for a set of s/c values. More...
|
|
void | computeFractionalFlow (const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &fractional_flows) |
|
void | computeFractionalFlow (const Opm::BlackoilPropertiesInterface &props, const Opm::PolymerProperties &polyprops, 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, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &fractional_flows) |
|
void | computeInjectedProduced (const IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const PolymerState &state, const std::vector< double > &transport_src, const std::vector< double > &inj_c, const double dt, double *injected, double *produced, double &polyinj, double &polyprod) |
| Computes injected and produced volumes of all phases, and injected and produced polymer mass. 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 | computeInjectedProduced (const BlackoilPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const PolymerBlackoilState &state, const std::vector< double > &transport_src, const std::vector< double > &inj_c, const double dt, double *injected, double *produced, double &polyinj, double &polyprod) |
| Computes injected and produced volumes of all phases, and injected and produced polymer mass - in the compressible case. 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...
|
|
double | computePolymerMass (const std::vector< double > &pv, const std::vector< double > &s, const std::vector< double > &c, const double dps) |
| Computes total (free) polymer mass over all grid cells. More...
|
|
double | computePolymerAdsorbed (const IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< double > &pv, const std::vector< double > &cmax) |
| Computes total absorbed polymer mass over all grid cells. More...
|
|
double | computePolymerAdsorbed (const UnstructuredGrid &grid, const BlackoilPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const PolymerBlackoilState &state, const RockCompressibility *rock_comp) |
| Computes total absorbed polymer mass over all grid cells. With compressibility. More...
|
|
◆ computeFractionalFlow() [1/2]
void Opm::computeFractionalFlow |
( |
const Opm::BlackoilPropertiesInterface & |
props, |
|
|
const Opm::PolymerProperties & |
polyprops, |
|
|
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, |
|
|
const std::vector< double > & |
c, |
|
|
const std::vector< double > & |
cmax, |
|
|
std::vector< double > & |
fractional_flows |
|
) |
| |
Computes the fractional flow for each cell in the cells argument - Parameters
-
[in] | props | rock and fluid properties |
[in] | polyprops | polymer 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) |
[in] | c | concentration values |
[in] | cmax | max polymer concentration experienced by cell |
[out] | fractional_flow | the fractional flow for each phase for each cell. |
◆ computeFractionalFlow() [2/2]
void Opm::computeFractionalFlow |
( |
const Opm::IncompPropertiesInterface & |
props, |
|
|
const Opm::PolymerProperties & |
polyprops, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
s, |
|
|
const std::vector< double > & |
c, |
|
|
const std::vector< double > & |
cmax, |
|
|
std::vector< double > & |
fractional_flows |
|
) |
| |
Computes the fractional flow for each cell in the cells argument - Parameters
-
[in] | props | rock and fluid properties |
[in] | polyprops | polymer properties |
[in] | cells | cells with which the saturation values are associated |
[in] | s | saturation values (for all phases) |
[in] | c | concentration values |
[in] | cmax | max polymer concentration experienced by cell |
[out] | fractional_flow | the fractional flow for each phase for each cell. |
◆ computeInjectedProduced() [1/2]
void Opm::computeInjectedProduced |
( |
const BlackoilPropertiesInterface & |
props, |
|
|
const Opm::PolymerProperties & |
polyprops, |
|
|
const PolymerBlackoilState & |
state, |
|
|
const std::vector< double > & |
transport_src, |
|
|
const std::vector< double > & |
inj_c, |
|
|
const double |
dt, |
|
|
double * |
injected, |
|
|
double * |
produced, |
|
|
double & |
polyinj, |
|
|
double & |
polyprod |
|
) |
| |
Computes injected and produced volumes of all phases, and injected and produced polymer mass - in the compressible case. 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] | polyprops | polymer properties |
[in] | state | state variables (pressure, fluxes etc.) |
[in] | src | if < 0: total reservoir volume outflow, if > 0: first phase surface volume inflow. |
[in] | inj_c | injected concentration by cell |
[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. |
[out] | polyinj | injected mass of polymer |
[out] | polyprod | produced mass of polymer |
◆ computeInjectedProduced() [2/2]
void Opm::computeInjectedProduced |
( |
const IncompPropertiesInterface & |
props, |
|
|
const Opm::PolymerProperties & |
polyprops, |
|
|
const PolymerState & |
state, |
|
|
const std::vector< double > & |
transport_src, |
|
|
const std::vector< double > & |
inj_c, |
|
|
const double |
dt, |
|
|
double * |
injected, |
|
|
double * |
produced, |
|
|
double & |
polyinj, |
|
|
double & |
polyprod |
|
) |
| |
Computes injected and produced volumes of all phases, and injected and produced polymer mass. 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] | polyprops | polymer properties |
[in] | state | state variables (pressure, fluxes etc.) |
[in] | src | if < 0: total reservoir volume outflow, if > 0: first phase reservoir volume inflow. |
[in] | inj_c | injected concentration by cell |
[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. |
[out] | polyinj | injected mass of polymer |
[out] | polyprod | produced mass of polymer |
◆ computePolymerAdsorbed() [1/2]
double Opm::computePolymerAdsorbed |
( |
const IncompPropertiesInterface & |
props, |
|
|
const Opm::PolymerProperties & |
polyprops, |
|
|
const std::vector< double > & |
pv, |
|
|
const std::vector< double > & |
cmax |
|
) |
| |
Computes total absorbed polymer mass over all grid cells.
- Parameters
-
[in] | props | fluid and rock properties. |
[in] | polyprops | polymer properties |
[in] | pv | the pore volume by cell. |
[in] | cmax | max polymer concentration for cell |
- Returns
- total absorbed polymer mass.
◆ computePolymerAdsorbed() [2/2]
double Opm::computePolymerAdsorbed |
( |
const UnstructuredGrid & |
grid, |
|
|
const BlackoilPropertiesInterface & |
props, |
|
|
const Opm::PolymerProperties & |
polyprops, |
|
|
const PolymerBlackoilState & |
state, |
|
|
const RockCompressibility * |
rock_comp |
|
) |
| |
Computes total absorbed polymer mass over all grid cells. With compressibility.
- Parameters
-
[in] | grid | grid |
[in] | props | fluid and rock properties. |
[in] | polyprops | polymer properties |
[in] | state | State variables |
[in] | rock_comp | Rock compressibility (optional) |
- Returns
- total absorbed polymer mass.
◆ computePolymerMass()
double Opm::computePolymerMass |
( |
const std::vector< double > & |
pv, |
|
|
const std::vector< double > & |
s, |
|
|
const std::vector< double > & |
c, |
|
|
const double |
dps |
|
) |
| |
Computes total (free) polymer mass over all grid cells.
- Parameters
-
[in] | pv | the pore volume by cell. |
[in] | s | saturation values (for all P phases) |
[in] | c | polymer concentration |
[in] | dps | dead pore space |
- Returns
- total polymer mass in grid.
◆ computeTotalMobility()
void Opm::computeTotalMobility |
( |
const Opm::IncompPropertiesInterface & |
props, |
|
|
const Opm::PolymerProperties & |
polyprops, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
s, |
|
|
const std::vector< double > & |
c, |
|
|
const std::vector< double > & |
cmax, |
|
|
std::vector< double > & |
totmob |
|
) |
| |
Computes total mobility for a set of s/c values.
- Parameters
-
[in] | props | rock and fluid properties |
[in] | polyprops | polymer properties |
[in] | cells | cells with which the saturation values are associated |
[in] | s | saturation values (for all phases) |
[in] | c | polymer concentration |
[out] | totmob | total mobilities. |
◆ computeTotalMobilityOmega()
void Opm::computeTotalMobilityOmega |
( |
const Opm::IncompPropertiesInterface & |
props, |
|
|
const Opm::PolymerProperties & |
polyprops, |
|
|
const std::vector< int > & |
cells, |
|
|
const std::vector< double > & |
s, |
|
|
const std::vector< double > & |
c, |
|
|
const std::vector< double > & |
cmax, |
|
|
std::vector< double > & |
totmob, |
|
|
std::vector< double > & |
omega |
|
) |
| |
Computes total mobility and omega for a set of s/c values.
- Parameters
-
[in] | props | rock and fluid properties |
[in] | polyprops | polymer properties |
[in] | cells | cells with which the saturation values are associated |
[in] | s | saturation values (for all phases) |
[in] | c | polymer concentration |
[in] | cmax | max polymer concentration experienced by cell |
[out] | totmob | total mobility |
[out] | omega | mobility-weighted (or fractional-flow weighted) fluid densities. |
|