Bugfix in calculation of adsorbed polymer for balance reports.

This commit is contained in:
Atgeirr Flø Rasmussen 2012-03-26 15:33:38 +02:00
parent c070f0ad4d
commit 294c8b03be
3 changed files with 9 additions and 4 deletions

View File

@ -940,7 +940,7 @@ main(int argc, char** argv)
// Report volume balances.
Opm::computeSaturatedVol(porevol, state.saturation(), satvol);
polymass = Opm::computePolymerMass(porevol, state.saturation(), state.concentration(), polyprop.deadPoreVol());
polymass_adsorbed = Opm::computePolymerAdsorbed(polyprop, porevol, state.cmax());
polymass_adsorbed = Opm::computePolymerAdsorbed(*props, polyprop, porevol, state.cmax());
Opm::computeInjectedProduced(*props, polyprop, state.saturation(), state.concentration(),
src, simtimer.currentStepLength(), inflow_c,
injected, produced, polyinj, polyprod);

View File

@ -189,19 +189,22 @@ namespace Opm
/// @brief Computes total absorbed polymer mass over all grid cells.
/// @param[in] props fluid and rock properties.
/// @param[in] polyprops polymer properties
/// @param[in] pv the pore volume by cell.
/// @param[in] cmax max polymer concentration for cell
/// @return total absorbed polymer mass.
double computePolymerAdsorbed(const Opm::PolymerProperties& polyprops,
double computePolymerAdsorbed(const IncompPropertiesInterface& props,
const Opm::PolymerProperties& polyprops,
const std::vector<double>& pv,
const std::vector<double>& cmax)
{
const int num_cells = pv.size();
const double rhor = polyprops.rockDensity();
const double* poro = props.porosity();
double abs_mass = 0.0;
for (int cell = 0; cell < num_cells; ++cell) {
abs_mass += polyprops.adsorbtion(cmax[cell])*pv[cell]*rhor;
abs_mass += polyprops.adsorbtion(cmax[cell])*pv[cell]*((1.0 - poro[cell])/poro[cell])*rhor;
}
return abs_mass;
}

View File

@ -103,11 +103,13 @@ namespace Opm
const double dps);
/// @brief Computes total absorbed polymer mass over all grid cells.
/// @param[in] props fluid and rock properties.
/// @param[in] polyprops polymer properties
/// @param[in] pv the pore volume by cell.
/// @param[in] cmax max polymer concentration for cell
/// @return total absorbed polymer mass.
double computePolymerAdsorbed(const Opm::PolymerProperties& polyprops,
double computePolymerAdsorbed(const IncompPropertiesInterface& props,
const Opm::PolymerProperties& polyprops,
const std::vector<double>& pv,
const std::vector<double>& cmax);