mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-01-16 01:51:55 -06:00
Fix for case with incompressible rock.
Made rock comp argument optional in computePolymerAdsorbed(). Check inside func for active rock comp.
This commit is contained in:
parent
dbd4c9e7b5
commit
c6b76e715d
@ -365,7 +365,7 @@ namespace Opm
|
||||
Opm::computeSaturatedVol(porevol, state.saturation(), satvol);
|
||||
polymass = Opm::computePolymerMass(porevol, state.saturation(), state.concentration(), poly_props_.deadPoreVol());
|
||||
polymass_adsorbed = Opm::computePolymerAdsorbed(grid_, props_, poly_props_,
|
||||
state, *rock_comp_props_);
|
||||
state, rock_comp_props_);
|
||||
tot_injected[0] += injected[0];
|
||||
tot_injected[1] += injected[1];
|
||||
tot_produced[0] += produced[0];
|
||||
|
@ -294,13 +294,17 @@ namespace Opm
|
||||
const BlackoilPropertiesInterface& props,
|
||||
const Opm::PolymerProperties& polyprops,
|
||||
const PolymerBlackoilState& state,
|
||||
const RockCompressibility& rock_comp
|
||||
const RockCompressibility* rock_comp
|
||||
)
|
||||
{
|
||||
const int num_cells = props.numCells();
|
||||
const double rhor = polyprops.rockDensity();
|
||||
std::vector<double> porosity;
|
||||
computePorosity(grid, props.porosity(), rock_comp, state.pressure(), porosity);
|
||||
if (rock_comp && rock_comp->isActive()) {
|
||||
computePorosity(grid, props.porosity(), *rock_comp, state.pressure(), porosity);
|
||||
} else {
|
||||
porosity.assign(props.porosity(), props.porosity() + num_cells);
|
||||
}
|
||||
double abs_mass = 0.0;
|
||||
const std::vector<double>& cmax = state.maxconcentration();
|
||||
for (int cell = 0; cell < num_cells; ++cell) {
|
||||
|
@ -163,13 +163,13 @@ namespace Opm
|
||||
/// @param[in] props fluid and rock properties.
|
||||
/// @param[in] polyprops polymer properties
|
||||
/// @param[in] state State variables
|
||||
/// @param[in] rock_comp Rock compressibility
|
||||
/// @param[in] rock_comp Rock compressibility (optional)
|
||||
/// @return total absorbed polymer mass.
|
||||
double computePolymerAdsorbed(const UnstructuredGrid& grid,
|
||||
const BlackoilPropertiesInterface& props,
|
||||
const Opm::PolymerProperties& polyprops,
|
||||
const PolymerBlackoilState& state,
|
||||
const RockCompressibility& rock_comp);
|
||||
const RockCompressibility* rock_comp);
|
||||
|
||||
|
||||
/// @brief Functor giving the injected amount of polymer as a function of time.
|
||||
|
Loading…
Reference in New Issue
Block a user