flow_ebos FIP: fix the field totals code

this unifies the code paths of the code that calculates the FIP field
totals for the parallel and the sequential cases and makes the code
more robust because it does not hard-code the presence of an intensive
quantities cache anymore. also, rock compressibility is now also
included in the field totals instead of just the FIP regions. this was
forgotten in the last FIP PR because the region values are calculated
in a different class using completely different code. (i.e., regions
are done by the model, field totals by the simulator. that design
should win an award, IMO.)

with this patch, the field totals for Norne and SPE1 seem to match
those produced by E100 _very_ closely and also parallel and sequential
runs for Norne and SPE1 of flow_ebos produce exactly the same
numbers. (This is probably the case for all decks, but I haven't
tested anything else.)
This commit is contained in:
Andreas Lauser 2017-03-13 10:04:35 +01:00
parent 5194cc5122
commit aa9966d7b8

View File

@ -652,56 +652,54 @@ protected:
totals[i] += fip[reg][i];
}
}
const int numCells = Opm::AutoDiffGrid::numCells(grid());
const auto& pv = geo_.poreVolume();
const auto& gridView = ebosSimulator_.gridManager().gridView();
const auto& comm = gridView.comm();
double pv_hydrocarbon_sum = 0.0;
double p_pv_hydrocarbon_sum = 0.0;
if ( ! is_parallel_run_ )
ElementContext elemCtx(ebosSimulator_);
const auto& elemEndIt = gridView.template end</*codim=*/0>();
for (auto elemIt = gridView.template begin</*codim=*/0>();
elemIt != elemEndIt;
++elemIt)
{
for (int cellIdx = 0; cellIdx < numCells; ++cellIdx) {
const auto& intQuants = *ebosSimulator_.model().cachedIntensiveQuantities(cellIdx, /*timeIdx=*/0);
const auto& fs = intQuants.fluidState();
const double& p = fs.pressure(FluidSystem::oilPhaseIdx).value();
const double hydrocarbon = fs.saturation(FluidSystem::oilPhaseIdx).value() + fs.saturation(FluidSystem::gasPhaseIdx).value();
totals[5] += pv[cellIdx];
pv_hydrocarbon_sum += pv[cellIdx] * hydrocarbon;
p_pv_hydrocarbon_sum += p * pv[cellIdx] * hydrocarbon;
}
}
else
{
#if HAVE_MPI
const auto & pinfo =
boost::any_cast<const ParallelISTLInformation&>(solver_.parallelInformation());
// Mask with 1 for owned cell and 0 otherwise
const auto& mask = pinfo.updateOwnerMask(pv);
for (int cellIdx = 0; cellIdx < numCells; ++cellIdx) {
const auto& intQuants = *ebosSimulator_.model().cachedIntensiveQuantities(cellIdx, /*timeIdx=*/0);
const auto& fs = intQuants.fluidState();
const double& p = fs.pressure(FluidSystem::oilPhaseIdx).value();
const double hydrocarbon = fs.saturation(FluidSystem::oilPhaseIdx).value() + fs.saturation(FluidSystem::gasPhaseIdx).value();
if( mask[cellIdx] )
{
totals[5] += pv[cellIdx];
pv_hydrocarbon_sum += pv[cellIdx] * hydrocarbon;
p_pv_hydrocarbon_sum += p * pv[cellIdx] * hydrocarbon;
}
const auto& elem = *elemIt;
if (elem.partitionType() != Dune::InteriorEntity) {
continue;
}
totals[5] = pinfo.communicator().sum(totals[5]);
pv_hydrocarbon_sum = pinfo.communicator().sum(pv_hydrocarbon_sum);
p_pv_hydrocarbon_sum= pinfo.communicator().sum(p_pv_hydrocarbon_sum);
#else
OPM_THROW(std::logic_error, "Requested a parallel run without MPI available!");
#endif
elemCtx.updatePrimaryStencil(elem);
elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
const unsigned cellIdx = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
const auto& intQuants = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0);
const auto& fs = intQuants.fluidState();
const double p = fs.pressure(FluidSystem::oilPhaseIdx).value();
const double hydrocarbon = fs.saturation(FluidSystem::oilPhaseIdx).value() + fs.saturation(FluidSystem::gasPhaseIdx).value();
// calculate the pore volume of the current cell. Note that the
// porosity returned by the intensive quantities is defined as the
// ratio of pore space to total cell volume and includes all pressure
// dependent (-> rock compressibility) and static modifiers (MULTPV,
// MULTREGP, NTG, PORV, MINPV and friends). Also note that because of
// this, the porosity returned by the intensive quantities can be
// outside of the physical range [0, 1] in pathetic cases.
const double pv =
ebosSimulator_.model().dofTotalVolume(cellIdx)
* intQuants.porosity().value();
totals[5] += pv;
pv_hydrocarbon_sum += pv*hydrocarbon;
p_pv_hydrocarbon_sum += p*pv*hydrocarbon;
}
pv_hydrocarbon_sum = comm.sum(pv_hydrocarbon_sum);
p_pv_hydrocarbon_sum = comm.sum(p_pv_hydrocarbon_sum);
totals[5] = comm.sum(totals[5]);
totals[6] = (p_pv_hydrocarbon_sum / pv_hydrocarbon_sum);
return totals;
}
@ -747,7 +745,7 @@ protected:
<< std::fixed << std::setprecision(0)
<< " : PORV =" << std::setw(14) << cip[5] << " RB :\n";
if (!reg) {
ss << " : Pressure is weighted by hydrocarbon pore voulme :\n"
ss << " : Pressure is weighted by hydrocarbon pore volume :\n"
<< " : Pore volumes are taken at reference conditions :\n";
}
ss << " :--------------- Oil STB ---------------:-- Wat STB --:--------------- Gas MSCF ---------------:\n";