2019-04-23 01:56:06 -05:00
|
|
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
|
2015-06-18 06:43:59 -05:00
|
|
|
// vi: set et ts=4 sw=4 sts=4:
|
2014-11-28 05:58:48 -06:00
|
|
|
/*
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
2019-10-30 11:51:49 -05:00
|
|
|
|
2014-11-28 05:58:48 -06:00
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 2 of the License, or
|
|
|
|
(at your option) any later version.
|
2019-10-30 11:51:49 -05:00
|
|
|
|
2014-11-28 05:58:48 -06:00
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
2019-10-30 11:51:49 -05:00
|
|
|
|
2014-11-28 05:58:48 -06:00
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
2019-10-30 11:51:49 -05:00
|
|
|
|
2016-03-14 07:21:47 -05:00
|
|
|
Consult the COPYING file in the top-level source directory of this
|
|
|
|
module for the precise wording of the license and the list of
|
|
|
|
copyright holders.
|
2014-11-28 05:58:48 -06:00
|
|
|
*/
|
|
|
|
/*!
|
|
|
|
* \file
|
2019-09-05 10:04:39 -05:00
|
|
|
* \copydoc Opm::EclOutputBlackOilModule
|
2014-11-28 05:58:48 -06:00
|
|
|
*/
|
|
|
|
#ifndef EWOMS_ECL_OUTPUT_BLACK_OIL_MODULE_HH
|
|
|
|
#define EWOMS_ECL_OUTPUT_BLACK_OIL_MODULE_HH
|
|
|
|
|
2022-12-13 05:54:27 -06:00
|
|
|
#include <opm/common/Exceptions.hpp>
|
2023-03-24 08:56:39 -05:00
|
|
|
#include <opm/common/TimingMacros.hpp>
|
2022-12-13 05:54:27 -06:00
|
|
|
#include <opm/common/OpmLog/OpmLog.hpp>
|
|
|
|
|
2023-03-24 08:56:39 -05:00
|
|
|
#include <opm/material/common/Valgrind.hpp>
|
|
|
|
#include <opm/material/fluidmatrixinteractions/EclEpsScalingPoints.hpp>
|
|
|
|
#include <opm/material/fluidstates/BlackOilFluidState.hpp>
|
|
|
|
#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
|
2018-03-12 08:31:36 -05:00
|
|
|
|
2023-03-24 08:56:39 -05:00
|
|
|
#include <opm/models/blackoil/blackoilproperties.hh>
|
|
|
|
#include <opm/models/discretization/common/fvbaseproperties.hh>
|
2019-09-16 03:58:20 -05:00
|
|
|
#include <opm/models/utils/parametersystem.hh>
|
2023-01-16 01:53:57 -06:00
|
|
|
#include <opm/models/utils/propertysystem.hh>
|
2014-11-28 05:58:48 -06:00
|
|
|
|
2017-12-07 05:38:48 -06:00
|
|
|
#include <opm/output/data/Cells.hpp>
|
|
|
|
#include <opm/output/eclipse/EclipseIO.hpp>
|
2020-11-17 01:01:13 -06:00
|
|
|
#include <opm/output/eclipse/Inplace.hpp>
|
2016-11-07 08:14:07 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
#include <ebos/eclgenericoutputblackoilmodule.hh>
|
|
|
|
|
2014-11-28 05:58:48 -06:00
|
|
|
#include <dune/common/fvector.hh>
|
|
|
|
|
2022-04-24 19:08:04 -05:00
|
|
|
#include <algorithm>
|
|
|
|
#include <array>
|
|
|
|
#include <cstddef>
|
|
|
|
#include <initializer_list>
|
|
|
|
#include <numeric>
|
|
|
|
#include <optional>
|
|
|
|
#include <stdexcept>
|
|
|
|
#include <tuple>
|
2014-11-28 05:58:48 -06:00
|
|
|
#include <type_traits>
|
2022-04-24 19:08:04 -05:00
|
|
|
#include <utility>
|
2014-11-28 05:58:48 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
namespace Opm::Properties
|
|
|
|
{
|
2018-06-14 09:06:05 -05:00
|
|
|
|
2017-12-07 05:38:48 -06:00
|
|
|
// create new type tag for the Ecl-output
|
2023-01-16 01:53:57 -06:00
|
|
|
namespace TTag
|
|
|
|
{
|
|
|
|
struct EclOutputBlackOil {
|
|
|
|
};
|
|
|
|
} // namespace TTag
|
2014-11-28 05:58:48 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
template <class TypeTag, class MyTypeTag>
|
2020-08-27 04:38:38 -05:00
|
|
|
struct ForceDisableFluidInPlaceOutput {
|
|
|
|
using type = UndefinedProperty;
|
|
|
|
};
|
2018-07-26 05:31:32 -05:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
template <class TypeTag>
|
2020-08-27 04:38:38 -05:00
|
|
|
struct ForceDisableFluidInPlaceOutput<TypeTag, TTag::EclOutputBlackOil> {
|
|
|
|
static constexpr bool value = false;
|
|
|
|
};
|
2018-06-14 09:06:05 -05:00
|
|
|
|
2021-08-25 11:24:05 -05:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
template <class TypeTag, class MyTypeTag>
|
2021-08-25 11:24:05 -05:00
|
|
|
struct ForceDisableResvFluidInPlaceOutput {
|
|
|
|
using type = UndefinedProperty;
|
|
|
|
};
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
template <class TypeTag>
|
2021-08-25 11:24:05 -05:00
|
|
|
struct ForceDisableResvFluidInPlaceOutput<TypeTag, TTag::EclOutputBlackOil> {
|
|
|
|
static constexpr bool value = false;
|
|
|
|
};
|
|
|
|
|
2020-08-21 06:42:08 -05:00
|
|
|
} // namespace Opm::Properties
|
2018-06-14 09:06:05 -05:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
namespace Opm
|
|
|
|
{
|
2014-11-28 05:58:48 -06:00
|
|
|
|
|
|
|
// forward declaration
|
|
|
|
template <class TypeTag>
|
|
|
|
class EcfvDiscretization;
|
|
|
|
|
|
|
|
/*!
|
2014-12-22 12:19:03 -06:00
|
|
|
* \ingroup EclBlackOilSimulator
|
2014-11-28 05:58:48 -06:00
|
|
|
*
|
|
|
|
* \brief Output module for the results black oil model writing in
|
|
|
|
* ECL binary format.
|
|
|
|
*/
|
|
|
|
template <class TypeTag>
|
2021-05-17 16:18:36 -05:00
|
|
|
class EclOutputBlackOilModule : public EclGenericOutputBlackoilModule<GetPropType<TypeTag, Properties::FluidSystem>,
|
|
|
|
GetPropType<TypeTag, Properties::Scalar>>
|
2014-11-28 05:58:48 -06:00
|
|
|
{
|
2020-08-26 03:49:52 -05:00
|
|
|
using Simulator = GetPropType<TypeTag, Properties::Simulator>;
|
|
|
|
using Discretization = GetPropType<TypeTag, Properties::Discretization>;
|
|
|
|
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
|
|
|
|
using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
|
|
|
|
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
|
|
|
|
using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
|
|
|
|
using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
|
2023-03-23 07:44:46 -05:00
|
|
|
using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
|
2020-08-26 03:49:52 -05:00
|
|
|
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
|
|
|
|
using GridView = GetPropType<TypeTag, Properties::GridView>;
|
|
|
|
using Element = typename GridView::template Codim<0>::Entity;
|
|
|
|
using ElementIterator = typename GridView::template Codim<0>::Iterator;
|
2021-05-17 16:18:36 -05:00
|
|
|
using BaseType = EclGenericOutputBlackoilModule<FluidSystem, Scalar>;
|
2022-11-17 08:01:14 -06:00
|
|
|
using Indices = GetPropType<TypeTag, Properties::Indices>;
|
2014-11-28 05:58:48 -06:00
|
|
|
|
2023-02-06 17:41:41 -06:00
|
|
|
enum { conti0EqIdx = Indices::conti0EqIdx };
|
2014-11-28 05:58:48 -06:00
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
|
|
|
enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
|
|
|
|
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
|
|
|
|
enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
|
|
|
|
enum { gasCompIdx = FluidSystem::gasCompIdx };
|
2015-02-03 06:58:18 -06:00
|
|
|
enum { oilCompIdx = FluidSystem::oilCompIdx };
|
2022-11-17 08:01:14 -06:00
|
|
|
enum { waterCompIdx = FluidSystem::waterCompIdx };
|
2020-08-27 02:13:30 -05:00
|
|
|
enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
|
2014-11-28 05:58:48 -06:00
|
|
|
|
|
|
|
public:
|
2022-04-24 19:08:04 -05:00
|
|
|
template <class CollectDataToIORankType>
|
2023-01-16 01:53:57 -06:00
|
|
|
EclOutputBlackOilModule(const Simulator& simulator,
|
2021-05-17 16:18:36 -05:00
|
|
|
const std::vector<std::size_t>& wbp_index_list,
|
2023-01-16 01:53:57 -06:00
|
|
|
const CollectDataToIORankType& collectToIORank)
|
2021-05-17 16:18:36 -05:00
|
|
|
: BaseType(simulator.vanguard().eclState(),
|
|
|
|
simulator.vanguard().schedule(),
|
|
|
|
simulator.vanguard().summaryConfig(),
|
|
|
|
simulator.vanguard().summaryState(),
|
|
|
|
getPropValue<TypeTag, Properties::EnableEnergy>(),
|
2021-06-01 04:55:54 -05:00
|
|
|
getPropValue<TypeTag, Properties::EnableTemperature>(),
|
2021-05-17 16:18:36 -05:00
|
|
|
getPropValue<TypeTag, Properties::EnableSolvent>(),
|
|
|
|
getPropValue<TypeTag, Properties::EnablePolymer>(),
|
|
|
|
getPropValue<TypeTag, Properties::EnableFoam>(),
|
|
|
|
getPropValue<TypeTag, Properties::EnableBrine>(),
|
2022-01-06 05:32:38 -06:00
|
|
|
getPropValue<TypeTag, Properties::EnableSaltPrecipitation>(),
|
2021-10-06 12:32:35 -05:00
|
|
|
getPropValue<TypeTag, Properties::EnableExtbo>(),
|
|
|
|
getPropValue<TypeTag, Properties::EnableMICP>())
|
2023-01-16 01:53:57 -06:00
|
|
|
, simulator_(simulator)
|
2018-01-12 08:32:43 -06:00
|
|
|
{
|
2022-04-24 19:08:04 -05:00
|
|
|
for (auto& region_pair : this->regions_) {
|
|
|
|
this->createLocalRegion_(region_pair.second);
|
|
|
|
}
|
2020-09-07 03:22:42 -05:00
|
|
|
|
2023-01-11 08:38:03 -06:00
|
|
|
auto isCartIdxOnThisRank = [&collectToIORank](int idx)
|
|
|
|
{
|
|
|
|
return collectToIORank.isCartIdxOnThisRank(idx);
|
|
|
|
};
|
|
|
|
|
|
|
|
this->setupBlockData(isCartIdxOnThisRank);
|
2018-07-26 05:31:32 -05:00
|
|
|
|
2020-12-01 02:09:06 -06:00
|
|
|
for (const auto& global_index : wbp_index_list) {
|
2022-04-24 19:08:04 -05:00
|
|
|
if (collectToIORank.isCartIdxOnThisRank(global_index - 1)) {
|
|
|
|
this->wbpData_.emplace(global_index, 0.0);
|
|
|
|
}
|
2020-12-01 02:09:06 -06:00
|
|
|
}
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
this->forceDisableFipOutput_ = EWOMS_GET_PARAM(TypeTag, bool, ForceDisableFluidInPlaceOutput);
|
2022-04-24 19:08:04 -05:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
this->forceDisableFipresvOutput_ = EWOMS_GET_PARAM(TypeTag, bool, ForceDisableResvFluidInPlaceOutput);
|
2018-07-26 05:31:32 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Register all run-time parameters for the Vtk output module.
|
|
|
|
*/
|
|
|
|
static void registerParameters()
|
|
|
|
{
|
2023-01-16 01:53:57 -06:00
|
|
|
EWOMS_REGISTER_PARAM(
|
|
|
|
TypeTag,
|
|
|
|
bool,
|
|
|
|
ForceDisableFluidInPlaceOutput,
|
|
|
|
"Do not print fluid-in-place values after each report step even if requested by the deck.");
|
|
|
|
EWOMS_REGISTER_PARAM(
|
|
|
|
TypeTag,
|
|
|
|
bool,
|
|
|
|
ForceDisableResvFluidInPlaceOutput,
|
|
|
|
"Do not print reservoir volumes values after each report step even if requested by the deck.");
|
2018-01-12 08:32:43 -06:00
|
|
|
}
|
2014-11-28 05:58:48 -06:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Allocate memory for the scalar fields we would like to
|
2017-12-07 05:38:48 -06:00
|
|
|
* write to ECL output files
|
2014-11-28 05:58:48 -06:00
|
|
|
*/
|
2023-01-16 01:53:57 -06:00
|
|
|
void
|
|
|
|
allocBuffers(unsigned bufferSize, unsigned reportStepNum, const bool substep, const bool log, const bool isRestart)
|
2014-11-28 05:58:48 -06:00
|
|
|
{
|
2023-01-16 01:53:57 -06:00
|
|
|
if (!std::is_same<Discretization, EcfvDiscretization<TypeTag>>::value)
|
2018-01-09 07:01:30 -06:00
|
|
|
return;
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
this->doAllocBuffers(bufferSize,
|
|
|
|
reportStepNum,
|
|
|
|
substep,
|
|
|
|
log,
|
|
|
|
isRestart,
|
2021-05-25 08:49:14 -05:00
|
|
|
simulator_.problem().vapparsActive(std::max(simulator_.episodeIndex(), 0)),
|
2021-05-17 16:18:36 -05:00
|
|
|
simulator_.problem().materialLawManager()->enableHysteresis(),
|
2022-11-17 08:01:14 -06:00
|
|
|
simulator_.problem().tracerModel().numTracers(),
|
|
|
|
simulator_.problem().eclWriter()->getOutputNnc().size());
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Modify the internal buffers according to the intensive quanties relevant
|
|
|
|
* for an element
|
|
|
|
*/
|
2016-11-07 08:14:07 -06:00
|
|
|
void processElement(const ElementContext& elemCtx)
|
2014-11-28 05:58:48 -06:00
|
|
|
{
|
2023-02-15 04:05:45 -06:00
|
|
|
OPM_TIMEBLOCK_LOCAL(processElement);
|
2023-01-16 01:53:57 -06:00
|
|
|
if (!std::is_same<Discretization, EcfvDiscretization<TypeTag>>::value)
|
2014-11-28 05:58:48 -06:00
|
|
|
return;
|
|
|
|
|
2019-03-12 09:51:41 -05:00
|
|
|
const auto& problem = elemCtx.simulator().problem();
|
2015-11-18 04:54:35 -06:00
|
|
|
for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
|
2017-12-07 05:38:48 -06:00
|
|
|
const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
|
|
|
|
const auto& fs = intQuants.fluidState();
|
2018-01-12 08:32:43 -06:00
|
|
|
|
2016-01-04 08:32:55 -06:00
|
|
|
typedef typename std::remove_const<typename std::remove_reference<decltype(fs)>::type>::type FluidState;
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
2015-12-31 06:20:41 -06:00
|
|
|
unsigned pvtRegionIdx = elemCtx.primaryVars(dofIdx, /*timeIdx=*/0).pvtRegionIndex();
|
2014-11-28 05:58:48 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->saturation_[phaseIdx].empty())
|
2018-01-09 07:01:30 -06:00
|
|
|
continue;
|
2016-11-15 11:00:48 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
this->saturation_[phaseIdx][globalDofIdx] = getValue(fs.saturation(phaseIdx));
|
|
|
|
Valgrind::CheckDefined(this->saturation_[phaseIdx][globalDofIdx]);
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
2018-01-09 07:01:30 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
if (!this->fluidPressure_.empty()) {
|
2019-10-09 08:24:23 -05:00
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
|
2023-01-16 01:53:57 -06:00
|
|
|
// Output oil pressure as default
|
|
|
|
this->fluidPressure_[globalDofIdx] = getValue(fs.pressure(oilPhaseIdx));
|
|
|
|
} else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
|
|
|
|
// Output gas if oil is not present
|
|
|
|
this->fluidPressure_[globalDofIdx] = getValue(fs.pressure(gasPhaseIdx));
|
|
|
|
} else {
|
|
|
|
// Output water if neither oil nor gas is present
|
|
|
|
this->fluidPressure_[globalDofIdx] = getValue(fs.pressure(waterPhaseIdx));
|
2019-10-09 08:24:23 -05:00
|
|
|
}
|
2023-01-16 01:53:57 -06:00
|
|
|
Valgrind::CheckDefined(this->fluidPressure_[globalDofIdx]);
|
2018-01-31 04:10:37 -06:00
|
|
|
}
|
2016-11-15 11:00:48 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->temperature_.empty()) {
|
|
|
|
this->temperature_[globalDofIdx] = getValue(fs.temperature(oilPhaseIdx));
|
|
|
|
Valgrind::CheckDefined(this->temperature_[globalDofIdx]);
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->gasDissolutionFactor_.empty()) {
|
2017-12-19 11:00:14 -06:00
|
|
|
Scalar SoMax = elemCtx.problem().maxOilSaturation(globalDofIdx);
|
2023-01-16 01:53:57 -06:00
|
|
|
this->gasDissolutionFactor_[globalDofIdx]
|
|
|
|
= FluidSystem::template saturatedDissolutionFactor<FluidState, Scalar>(
|
|
|
|
fs, oilPhaseIdx, pvtRegionIdx, SoMax);
|
2021-05-17 16:18:36 -05:00
|
|
|
Valgrind::CheckDefined(this->gasDissolutionFactor_[globalDofIdx]);
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->oilVaporizationFactor_.empty()) {
|
2017-12-19 11:00:14 -06:00
|
|
|
Scalar SoMax = elemCtx.problem().maxOilSaturation(globalDofIdx);
|
2023-01-16 01:53:57 -06:00
|
|
|
this->oilVaporizationFactor_[globalDofIdx]
|
|
|
|
= FluidSystem::template saturatedDissolutionFactor<FluidState, Scalar>(
|
|
|
|
fs, gasPhaseIdx, pvtRegionIdx, SoMax);
|
2021-05-17 16:18:36 -05:00
|
|
|
Valgrind::CheckDefined(this->oilVaporizationFactor_[globalDofIdx]);
|
2016-11-15 11:00:48 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->gasFormationVolumeFactor_.empty()) {
|
2023-01-16 01:53:57 -06:00
|
|
|
this->gasFormationVolumeFactor_[globalDofIdx] = 1.0
|
|
|
|
/ FluidSystem::template inverseFormationVolumeFactor<FluidState, Scalar>(
|
|
|
|
fs, gasPhaseIdx, pvtRegionIdx);
|
2021-05-17 16:18:36 -05:00
|
|
|
Valgrind::CheckDefined(this->gasFormationVolumeFactor_[globalDofIdx]);
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->saturatedOilFormationVolumeFactor_.empty()) {
|
2023-01-16 01:53:57 -06:00
|
|
|
this->saturatedOilFormationVolumeFactor_[globalDofIdx] = 1.0
|
|
|
|
/ FluidSystem::template saturatedInverseFormationVolumeFactor<FluidState, Scalar>(
|
|
|
|
fs, oilPhaseIdx, pvtRegionIdx);
|
2021-05-17 16:18:36 -05:00
|
|
|
Valgrind::CheckDefined(this->saturatedOilFormationVolumeFactor_[globalDofIdx]);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->oilSaturationPressure_.empty()) {
|
2023-01-16 01:53:57 -06:00
|
|
|
this->oilSaturationPressure_[globalDofIdx]
|
|
|
|
= FluidSystem::template saturationPressure<FluidState, Scalar>(fs, oilPhaseIdx, pvtRegionIdx);
|
2021-05-17 16:18:36 -05:00
|
|
|
Valgrind::CheckDefined(this->oilSaturationPressure_[globalDofIdx]);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->rs_.empty()) {
|
|
|
|
this->rs_[globalDofIdx] = getValue(fs.Rs());
|
|
|
|
Valgrind::CheckDefined(this->rs_[globalDofIdx]);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
2022-11-24 07:32:41 -06:00
|
|
|
if (!this->rsw_.empty()) {
|
|
|
|
this->rsw_[globalDofIdx] = getValue(fs.Rsw());
|
|
|
|
Valgrind::CheckDefined(this->rsw_[globalDofIdx]);
|
|
|
|
}
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->rv_.empty()) {
|
|
|
|
this->rv_[globalDofIdx] = getValue(fs.Rv());
|
|
|
|
Valgrind::CheckDefined(this->rv_[globalDofIdx]);
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
2022-04-29 08:35:16 -05:00
|
|
|
if (!this->pcow_.empty()) {
|
|
|
|
this->pcow_[globalDofIdx] = getValue(fs.pressure(oilPhaseIdx)) - getValue(fs.pressure(waterPhaseIdx));
|
|
|
|
Valgrind::CheckDefined(this->pcow_[globalDofIdx]);
|
|
|
|
}
|
|
|
|
if (!this->pcog_.empty()) {
|
|
|
|
this->pcog_[globalDofIdx] = getValue(fs.pressure(gasPhaseIdx)) - getValue(fs.pressure(oilPhaseIdx));
|
|
|
|
Valgrind::CheckDefined(this->pcog_[globalDofIdx]);
|
|
|
|
}
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2022-03-01 10:54:22 -06:00
|
|
|
if (!this->rvw_.empty()) {
|
|
|
|
this->rvw_[globalDofIdx] = getValue(fs.Rvw());
|
|
|
|
Valgrind::CheckDefined(this->rvw_[globalDofIdx]);
|
|
|
|
}
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->invB_[phaseIdx].empty())
|
2018-01-09 07:01:30 -06:00
|
|
|
continue;
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
this->invB_[phaseIdx][globalDofIdx] = getValue(fs.invB(phaseIdx));
|
|
|
|
Valgrind::CheckDefined(this->invB_[phaseIdx][globalDofIdx]);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->density_[phaseIdx].empty())
|
2018-01-09 07:01:30 -06:00
|
|
|
continue;
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
this->density_[phaseIdx][globalDofIdx] = getValue(fs.density(phaseIdx));
|
|
|
|
Valgrind::CheckDefined(this->density_[phaseIdx][globalDofIdx]);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->viscosity_[phaseIdx].empty())
|
2018-01-09 07:01:30 -06:00
|
|
|
continue;
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
if (!this->extboX_.empty() && phaseIdx == oilPhaseIdx)
|
2021-05-17 16:18:36 -05:00
|
|
|
this->viscosity_[phaseIdx][globalDofIdx] = getValue(intQuants.oilViscosity());
|
2023-01-16 01:53:57 -06:00
|
|
|
else if (!this->extboX_.empty() && phaseIdx == gasPhaseIdx)
|
2021-05-17 16:18:36 -05:00
|
|
|
this->viscosity_[phaseIdx][globalDofIdx] = getValue(intQuants.gasViscosity());
|
2020-11-04 09:08:17 -06:00
|
|
|
else
|
2021-05-17 16:18:36 -05:00
|
|
|
this->viscosity_[phaseIdx][globalDofIdx] = getValue(fs.viscosity(phaseIdx));
|
|
|
|
Valgrind::CheckDefined(this->viscosity_[phaseIdx][globalDofIdx]);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->relativePermeability_[phaseIdx].empty())
|
2018-01-09 07:01:30 -06:00
|
|
|
continue;
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
this->relativePermeability_[phaseIdx][globalDofIdx]
|
|
|
|
= getValue(intQuants.relativePermeability(phaseIdx));
|
2021-05-17 16:18:36 -05:00
|
|
|
Valgrind::CheckDefined(this->relativePermeability_[phaseIdx][globalDofIdx]);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2022-07-21 05:51:30 -05:00
|
|
|
if (!this->drsdtcon_.empty()) {
|
|
|
|
this->drsdtcon_[globalDofIdx] = problem.drsdtcon(globalDofIdx, elemCtx.simulator().episodeIndex());
|
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->sSol_.empty()) {
|
|
|
|
this->sSol_[globalDofIdx] = intQuants.solventSaturation().value();
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->cPolymer_.empty()) {
|
|
|
|
this->cPolymer_[globalDofIdx] = intQuants.polymerConcentration().value();
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->cFoam_.empty()) {
|
|
|
|
this->cFoam_[globalDofIdx] = intQuants.foamConcentration().value();
|
2019-04-26 03:38:37 -05:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->cSalt_.empty()) {
|
|
|
|
this->cSalt_[globalDofIdx] = fs.saltConcentration().value();
|
2019-11-07 02:39:42 -06:00
|
|
|
}
|
|
|
|
|
2022-01-06 05:32:38 -06:00
|
|
|
if (!this->pSalt_.empty()) {
|
|
|
|
this->pSalt_[globalDofIdx] = intQuants.saltSaturation().value();
|
|
|
|
}
|
|
|
|
|
|
|
|
if (!this->permFact_.empty()) {
|
|
|
|
this->permFact_[globalDofIdx] = intQuants.permFactor().value();
|
|
|
|
}
|
2023-01-16 01:53:57 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->extboX_.empty()) {
|
|
|
|
this->extboX_[globalDofIdx] = intQuants.xVolume().value();
|
2020-10-10 09:20:25 -05:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->extboY_.empty()) {
|
|
|
|
this->extboY_[globalDofIdx] = intQuants.yVolume().value();
|
2020-10-10 09:20:25 -05:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->extboZ_.empty()) {
|
|
|
|
this->extboZ_[globalDofIdx] = intQuants.zFraction().value();
|
2020-10-10 09:20:25 -05:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->mFracCo2_.empty()) {
|
2023-01-16 01:53:57 -06:00
|
|
|
const Scalar stdVolOil = getValue(fs.saturation(oilPhaseIdx)) * getValue(fs.invB(oilPhaseIdx))
|
|
|
|
+ getValue(fs.saturation(gasPhaseIdx)) * getValue(fs.invB(gasPhaseIdx)) * getValue(fs.Rv());
|
|
|
|
const Scalar stdVolGas = getValue(fs.saturation(gasPhaseIdx)) * getValue(fs.invB(gasPhaseIdx))
|
|
|
|
* (1.0 - intQuants.yVolume().value())
|
|
|
|
+ getValue(fs.saturation(oilPhaseIdx)) * getValue(fs.invB(oilPhaseIdx)) * getValue(fs.Rs())
|
|
|
|
* (1.0 - intQuants.xVolume().value());
|
|
|
|
const Scalar stdVolCo2 = getValue(fs.saturation(gasPhaseIdx)) * getValue(fs.invB(gasPhaseIdx))
|
|
|
|
* intQuants.yVolume().value()
|
|
|
|
+ getValue(fs.saturation(oilPhaseIdx)) * getValue(fs.invB(oilPhaseIdx)) * getValue(fs.Rs())
|
|
|
|
* intQuants.xVolume().value();
|
|
|
|
const Scalar rhoO = FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
|
|
|
|
const Scalar rhoG = FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
|
|
|
|
const Scalar rhoCO2 = intQuants.zRefDensity();
|
|
|
|
const Scalar stdMassTotal = 1.0e-10 + stdVolOil * rhoO + stdVolGas * rhoG + stdVolCo2 * rhoCO2;
|
|
|
|
this->mFracOil_[globalDofIdx] = stdVolOil * rhoO / stdMassTotal;
|
|
|
|
this->mFracGas_[globalDofIdx] = stdVolGas * rhoG / stdMassTotal;
|
|
|
|
this->mFracCo2_[globalDofIdx] = stdVolCo2 * rhoCO2 / stdMassTotal;
|
2020-10-10 09:20:25 -05:00
|
|
|
}
|
|
|
|
|
2021-10-06 12:32:35 -05:00
|
|
|
if (!this->cMicrobes_.empty()) {
|
|
|
|
this->cMicrobes_[globalDofIdx] = intQuants.microbialConcentration().value();
|
|
|
|
}
|
|
|
|
|
|
|
|
if (!this->cOxygen_.empty()) {
|
|
|
|
this->cOxygen_[globalDofIdx] = intQuants.oxygenConcentration().value();
|
|
|
|
}
|
|
|
|
|
|
|
|
if (!this->cUrea_.empty()) {
|
2023-01-16 01:53:57 -06:00
|
|
|
this->cUrea_[globalDofIdx] = 10
|
|
|
|
* intQuants.ureaConcentration()
|
|
|
|
.value(); // Reescaling back the urea concentration (see WellInterface_impl.hpp)
|
2021-10-06 12:32:35 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
if (!this->cBiofilm_.empty()) {
|
|
|
|
this->cBiofilm_[globalDofIdx] = intQuants.biofilmConcentration().value();
|
|
|
|
}
|
|
|
|
|
|
|
|
if (!this->cCalcite_.empty()) {
|
|
|
|
this->cCalcite_[globalDofIdx] = intQuants.calciteConcentration().value();
|
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->bubblePointPressure_.empty()) {
|
2017-12-07 05:38:48 -06:00
|
|
|
try {
|
2023-01-16 01:53:57 -06:00
|
|
|
this->bubblePointPressure_[globalDofIdx]
|
|
|
|
= getValue(FluidSystem::bubblePointPressure(fs, intQuants.pvtRegionIndex()));
|
|
|
|
} catch (const NumericalProblem&) {
|
2021-12-01 07:00:21 -06:00
|
|
|
const auto cartesianIdx = elemCtx.simulator().vanguard().cartesianIndex(globalDofIdx);
|
2021-05-17 16:18:36 -05:00
|
|
|
this->failedCellsPb_.push_back(cartesianIdx);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
}
|
2023-03-15 07:32:44 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->dewPointPressure_.empty()) {
|
2017-12-07 05:38:48 -06:00
|
|
|
try {
|
2023-01-16 01:53:57 -06:00
|
|
|
this->dewPointPressure_[globalDofIdx]
|
|
|
|
= getValue(FluidSystem::dewPointPressure(fs, intQuants.pvtRegionIndex()));
|
|
|
|
} catch (const NumericalProblem&) {
|
2021-12-01 07:00:21 -06:00
|
|
|
const auto cartesianIdx = elemCtx.simulator().vanguard().cartesianIndex(globalDofIdx);
|
2021-05-17 16:18:36 -05:00
|
|
|
this->failedCellsPd_.push_back(cartesianIdx);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->soMax_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->soMax_[globalDofIdx]
|
|
|
|
= std::max(getValue(fs.saturation(oilPhaseIdx)), problem.maxOilSaturation(globalDofIdx));
|
2019-03-12 09:51:41 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->swMax_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->swMax_[globalDofIdx]
|
|
|
|
= std::max(getValue(fs.saturation(waterPhaseIdx)), problem.maxWaterSaturation(globalDofIdx));
|
2019-03-12 09:51:41 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->minimumOilPressure_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->minimumOilPressure_[globalDofIdx]
|
|
|
|
= std::min(getValue(fs.pressure(oilPhaseIdx)), problem.minOilPressure(globalDofIdx));
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->overburdenPressure_.empty())
|
|
|
|
this->overburdenPressure_[globalDofIdx] = problem.overburdenPressure(globalDofIdx);
|
2019-03-12 09:51:41 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->rockCompPorvMultiplier_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->rockCompPorvMultiplier_[globalDofIdx]
|
|
|
|
= problem.template rockCompPoroMultiplier<Scalar>(intQuants, globalDofIdx);
|
2019-03-12 09:51:41 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->rockCompTransMultiplier_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->rockCompTransMultiplier_[globalDofIdx]
|
|
|
|
= problem.template rockCompTransMultiplier<Scalar>(intQuants, globalDofIdx);
|
2019-03-12 09:51:41 -05:00
|
|
|
|
|
|
|
const auto& matLawManager = problem.materialLawManager();
|
2018-01-09 07:01:30 -06:00
|
|
|
if (matLawManager->enableHysteresis()) {
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->pcSwMdcOw_.empty() && !this->krnSwMdcOw_.empty()) {
|
2017-12-07 05:38:48 -06:00
|
|
|
matLawManager->oilWaterHysteresisParams(
|
2023-01-16 01:53:57 -06:00
|
|
|
this->pcSwMdcOw_[globalDofIdx], this->krnSwMdcOw_[globalDofIdx], globalDofIdx);
|
2018-01-09 07:01:30 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->pcSwMdcGo_.empty() && !this->krnSwMdcGo_.empty()) {
|
2017-12-07 05:38:48 -06:00
|
|
|
matLawManager->gasOilHysteresisParams(
|
2023-01-16 01:53:57 -06:00
|
|
|
this->pcSwMdcGo_[globalDofIdx], this->krnSwMdcGo_[globalDofIdx], globalDofIdx);
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->ppcw_.empty()) {
|
|
|
|
this->ppcw_[globalDofIdx] = matLawManager->oilWaterScaledEpsInfoDrainage(globalDofIdx).maxPcow;
|
2023-01-16 01:53:57 -06:00
|
|
|
// printf("ppcw_[%d] = %lg\n", globalDofIdx, ppcw_[globalDofIdx]);
|
2020-01-13 08:46:50 -06:00
|
|
|
}
|
2023-03-15 07:32:44 -05:00
|
|
|
|
2017-12-07 05:38:48 -06:00
|
|
|
// hack to make the intial output of rs and rv Ecl compatible.
|
|
|
|
// For cells with swat == 1 Ecl outputs; rs = rsSat and rv=rvSat, in all but the initial step
|
|
|
|
// where it outputs rs and rv values calculated by the initialization. To be compatible we overwrite
|
|
|
|
// rs and rv with the values computed in the initially.
|
|
|
|
// Volume factors, densities and viscosities need to be recalculated with the updated rs and rv values.
|
|
|
|
// This can be removed when ebos has 100% controll over output
|
2023-03-15 07:32:44 -05:00
|
|
|
if ((elemCtx.simulator().episodeIndex() < 0) &&
|
|
|
|
FluidSystem::phaseIsActive(oilPhaseIdx) &&
|
|
|
|
FluidSystem::phaseIsActive(gasPhaseIdx))
|
|
|
|
{
|
2019-03-12 09:51:41 -05:00
|
|
|
const auto& fsInitial = problem.initialFluidState(globalDofIdx);
|
2017-12-07 05:38:48 -06:00
|
|
|
|
|
|
|
// use initial rs and rv values
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->rv_.empty())
|
|
|
|
this->rv_[globalDofIdx] = fsInitial.Rv();
|
2018-01-09 07:01:30 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->rs_.empty())
|
|
|
|
this->rs_[globalDofIdx] = fsInitial.Rs();
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2022-11-24 07:32:41 -06:00
|
|
|
if (!this->rsw_.empty())
|
|
|
|
this->rsw_[globalDofIdx] = fsInitial.Rsw();
|
|
|
|
|
2022-03-01 10:54:22 -06:00
|
|
|
if (!this->rvw_.empty())
|
|
|
|
this->rvw_[globalDofIdx] = fsInitial.Rvw();
|
|
|
|
|
2017-12-07 05:38:48 -06:00
|
|
|
// re-compute the volume factors, viscosities and densities if asked for
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->density_[oilPhaseIdx].empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->density_[oilPhaseIdx][globalDofIdx]
|
|
|
|
= FluidSystem::density(fsInitial, oilPhaseIdx, intQuants.pvtRegionIndex());
|
2023-03-15 07:32:44 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->density_[gasPhaseIdx].empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->density_[gasPhaseIdx][globalDofIdx]
|
|
|
|
= FluidSystem::density(fsInitial, gasPhaseIdx, intQuants.pvtRegionIndex());
|
2021-05-17 16:18:36 -05:00
|
|
|
|
|
|
|
if (!this->invB_[oilPhaseIdx].empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->invB_[oilPhaseIdx][globalDofIdx]
|
|
|
|
= FluidSystem::inverseFormationVolumeFactor(fsInitial, oilPhaseIdx, intQuants.pvtRegionIndex());
|
2023-03-15 07:32:44 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->invB_[gasPhaseIdx].empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->invB_[gasPhaseIdx][globalDofIdx]
|
|
|
|
= FluidSystem::inverseFormationVolumeFactor(fsInitial, gasPhaseIdx, intQuants.pvtRegionIndex());
|
2023-03-15 07:32:44 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->viscosity_[oilPhaseIdx].empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->viscosity_[oilPhaseIdx][globalDofIdx]
|
|
|
|
= FluidSystem::viscosity(fsInitial, oilPhaseIdx, intQuants.pvtRegionIndex());
|
2023-03-15 07:32:44 -05:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->viscosity_[gasPhaseIdx].empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
this->viscosity_[gasPhaseIdx][globalDofIdx]
|
|
|
|
= FluidSystem::viscosity(fsInitial, gasPhaseIdx, intQuants.pvtRegionIndex());
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
2023-03-24 08:56:23 -05:00
|
|
|
|
2018-02-07 07:30:43 -06:00
|
|
|
|
|
|
|
// Adding Well RFT data
|
2023-03-23 07:44:46 -05:00
|
|
|
const auto cartesianIdx = elemCtx.simulator().vanguard().cartesianIndex(globalDofIdx);
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->oilConnectionPressures_.count(cartesianIdx) > 0) {
|
|
|
|
this->oilConnectionPressures_[cartesianIdx] = getValue(fs.pressure(oilPhaseIdx));
|
2018-02-07 07:30:43 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->waterConnectionSaturations_.count(cartesianIdx) > 0) {
|
|
|
|
this->waterConnectionSaturations_[cartesianIdx] = getValue(fs.saturation(waterPhaseIdx));
|
2018-02-07 07:30:43 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->gasConnectionSaturations_.count(cartesianIdx) > 0) {
|
|
|
|
this->gasConnectionSaturations_[cartesianIdx] = getValue(fs.saturation(gasPhaseIdx));
|
2018-02-07 07:30:43 -06:00
|
|
|
}
|
2022-07-06 04:14:45 -05:00
|
|
|
if (this->wbpData_.count(cartesianIdx) > 0) {
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
|
|
|
|
this->wbpData_[cartesianIdx] = getValue(fs.pressure(oilPhaseIdx));
|
|
|
|
} else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
|
|
|
|
this->wbpData_[cartesianIdx] = getValue(fs.pressure(gasPhaseIdx));
|
|
|
|
} else if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
|
|
|
|
this->wbpData_[cartesianIdx] = getValue(fs.pressure(waterPhaseIdx));
|
|
|
|
}
|
|
|
|
}
|
2018-11-14 06:10:01 -06:00
|
|
|
|
|
|
|
// tracers
|
|
|
|
const auto& tracerModel = simulator_.problem().tracerModel();
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->tracerConcentrations_.empty()) {
|
2023-01-16 01:53:57 -06:00
|
|
|
for (int tracerIdx = 0; tracerIdx < tracerModel.numTracers(); tracerIdx++) {
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->tracerConcentrations_[tracerIdx].empty())
|
2018-11-14 06:10:01 -06:00
|
|
|
continue;
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
this->tracerConcentrations_[tracerIdx][globalDofIdx]
|
|
|
|
= tracerModel.tracerConcentration(tracerIdx, globalDofIdx);
|
2018-02-07 07:30:43 -06:00
|
|
|
}
|
2018-01-12 08:32:43 -06:00
|
|
|
}
|
2023-03-23 07:44:46 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
void processElementFlows(const ElementContext& elemCtx)
|
|
|
|
{
|
|
|
|
OPM_TIMEBLOCK_LOCAL(processElementBlockData);
|
2023-03-24 08:56:49 -05:00
|
|
|
if (!std::is_same_v<Discretization, EcfvDiscretization<TypeTag>>)
|
2023-03-23 07:44:46 -05:00
|
|
|
return;
|
2022-11-17 08:01:14 -06:00
|
|
|
|
2023-03-23 07:44:46 -05:00
|
|
|
const auto& problem = elemCtx.simulator().problem();
|
|
|
|
for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
|
|
|
|
|
|
|
|
unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
2023-03-24 08:56:23 -05:00
|
|
|
if (!problem.model().linearizer().getFlowsInfo().empty()) {
|
2022-11-17 08:01:14 -06:00
|
|
|
const auto& flowsInf = problem.model().linearizer().getFlowsInfo();
|
|
|
|
auto flowsInfos = flowsInf[globalDofIdx];
|
|
|
|
for (auto& flowsInfo : flowsInfos) {
|
|
|
|
if (flowsInfo.faceId == 1) {
|
|
|
|
if (!this->flowsi_[gasCompIdx].empty()) {
|
|
|
|
this->flowsi_[gasCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(gasCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->flowsi_[oilCompIdx].empty()) {
|
|
|
|
this->flowsi_[oilCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(oilCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->flowsi_[waterCompIdx].empty()) {
|
|
|
|
this->flowsi_[waterCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(waterCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
if (flowsInfo.faceId == 3) {
|
|
|
|
if (!this->flowsj_[gasCompIdx].empty()) {
|
|
|
|
this->flowsj_[gasCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(gasCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->flowsj_[oilCompIdx].empty()) {
|
|
|
|
this->flowsj_[oilCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(oilCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->flowsj_[waterCompIdx].empty()) {
|
|
|
|
this->flowsj_[waterCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(waterCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
if (flowsInfo.faceId == 5) {
|
|
|
|
if (!this->flowsk_[gasCompIdx].empty()) {
|
|
|
|
this->flowsk_[gasCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(gasCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->flowsk_[oilCompIdx].empty()) {
|
|
|
|
this->flowsk_[oilCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(oilCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->flowsk_[waterCompIdx].empty()) {
|
|
|
|
this->flowsk_[waterCompIdx][globalDofIdx]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(waterCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
if (flowsInfo.faceId == -2) {
|
|
|
|
if (!this->flowsn_[gasCompIdx].second.first.empty()) {
|
|
|
|
this->flowsn_[gasCompIdx].second.first[flowsInfo.nncId] = flowsInfo.nncId;
|
|
|
|
this->flowsn_[gasCompIdx].second.second[flowsInfo.nncId]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(gasCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->flowsn_[oilCompIdx].second.first.empty()) {
|
|
|
|
this->flowsn_[oilCompIdx].second.first[flowsInfo.nncId] = flowsInfo.nncId;
|
|
|
|
this->flowsn_[oilCompIdx].second.second[flowsInfo.nncId]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(oilCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->flowsn_[waterCompIdx].second.first.empty()) {
|
|
|
|
this->flowsn_[waterCompIdx].second.first[flowsInfo.nncId] = flowsInfo.nncId;
|
|
|
|
this->flowsn_[waterCompIdx].second.second[flowsInfo.nncId]
|
2023-02-06 17:41:41 -06:00
|
|
|
= flowsInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(waterCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
// flores
|
|
|
|
if (!problem.model().linearizer().getFloresInfo().empty()) {
|
|
|
|
const auto& floresInf = problem.model().linearizer().getFloresInfo();
|
|
|
|
auto floresInfos =floresInf[globalDofIdx];
|
|
|
|
for (auto& floresInfo : floresInfos) {
|
|
|
|
if (floresInfo.faceId == 1) {
|
2023-02-06 17:41:41 -06:00
|
|
|
if (!this->floresi_[gasCompIdx].empty()) {
|
|
|
|
this->floresi_[gasCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(gasCompIdx)];
|
|
|
|
}
|
|
|
|
if (!this->floresi_[oilCompIdx].empty()) {
|
|
|
|
this->floresi_[oilCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(oilCompIdx)];
|
|
|
|
}
|
|
|
|
if (!this->floresi_[waterCompIdx].empty()) {
|
|
|
|
this->floresi_[waterCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(waterCompIdx)];
|
|
|
|
}
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (floresInfo.faceId == 3) {
|
2023-02-06 17:41:41 -06:00
|
|
|
if (!this->floresj_[gasCompIdx].empty()) {
|
|
|
|
this->floresj_[gasCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(gasCompIdx)];
|
|
|
|
}
|
|
|
|
if (!this->floresj_[oilCompIdx].empty()) {
|
|
|
|
this->floresj_[oilCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(oilCompIdx)];
|
|
|
|
}
|
|
|
|
if (!this->floresj_[waterCompIdx].empty()) {
|
|
|
|
this->floresj_[waterCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(waterCompIdx)];
|
|
|
|
}
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (floresInfo.faceId == 5) {
|
2023-02-06 17:41:41 -06:00
|
|
|
if (!this->floresk_[gasCompIdx].empty()) {
|
|
|
|
this->floresk_[gasCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(gasCompIdx)];
|
|
|
|
}
|
|
|
|
if (!this->floresk_[oilCompIdx].empty()) {
|
|
|
|
this->floresk_[oilCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(oilCompIdx)];
|
|
|
|
}
|
|
|
|
if (!this->floresk_[waterCompIdx].empty()) {
|
|
|
|
this->floresk_[waterCompIdx][globalDofIdx]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(waterCompIdx)];
|
|
|
|
}
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (floresInfo.faceId == -2) {
|
|
|
|
if (!this->floresn_[gasCompIdx].second.first.empty()) {
|
|
|
|
this->floresn_[gasCompIdx].second.first[floresInfo.nncId] = floresInfo.nncId;
|
2023-02-06 17:41:41 -06:00
|
|
|
this->floresn_[gasCompIdx].second.second[floresInfo.nncId]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(gasCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->floresn_[oilCompIdx].second.first.empty()) {
|
|
|
|
this->floresn_[oilCompIdx].second.first[floresInfo.nncId] = floresInfo.nncId;
|
2023-02-06 17:41:41 -06:00
|
|
|
this->floresn_[oilCompIdx].second.second[floresInfo.nncId]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(oilCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
if (!this->floresn_[waterCompIdx].second.first.empty()) {
|
|
|
|
this->floresn_[waterCompIdx].second.first[floresInfo.nncId] = floresInfo.nncId;
|
2023-02-06 17:41:41 -06:00
|
|
|
this->floresn_[waterCompIdx].second.second[floresInfo.nncId]
|
|
|
|
= floresInfo.flow[conti0EqIdx + Indices::canonicalToActiveComponentIndex(waterCompIdx)];
|
2022-11-17 08:01:14 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2018-01-12 08:32:43 -06:00
|
|
|
}
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2023-03-23 07:44:46 -05:00
|
|
|
void processElementBlockData(const ElementContext& elemCtx)
|
|
|
|
{
|
|
|
|
OPM_TIMEBLOCK_LOCAL(processElementBlockData);
|
|
|
|
if (!std::is_same<Discretization, EcfvDiscretization<TypeTag>>::value)
|
|
|
|
return;
|
|
|
|
|
|
|
|
const auto& problem = elemCtx.simulator().problem();
|
|
|
|
for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
|
|
|
|
// Adding block data
|
|
|
|
unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
|
|
|
const auto cartesianIdx = elemCtx.simulator().vanguard().cartesianIndex(globalDofIdx);
|
|
|
|
const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
|
|
|
|
const auto& fs = intQuants.fluidState();
|
|
|
|
for (auto& val : this->blockData_) {
|
|
|
|
const auto& key = val.first;
|
|
|
|
assert(key.second > 0);
|
|
|
|
unsigned int cartesianIdxBlock = key.second - 1;
|
|
|
|
if (cartesianIdx == cartesianIdxBlock) {
|
|
|
|
if ((key.first == "BWSAT") || (key.first == "BSWAT"))
|
|
|
|
val.second = getValue(fs.saturation(waterPhaseIdx));
|
|
|
|
else if ((key.first == "BGSAT") || (key.first == "BSGAS"))
|
|
|
|
val.second = getValue(fs.saturation(gasPhaseIdx));
|
|
|
|
else if ((key.first == "BOSAT") || (key.first == "BSOIL"))
|
|
|
|
val.second = getValue(fs.saturation(oilPhaseIdx));
|
|
|
|
else if (key.first == "BNSAT")
|
|
|
|
val.second = intQuants.solventSaturation().value();
|
|
|
|
else if ((key.first == "BPR") || (key.first == "BPRESSUR")) {
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx))
|
|
|
|
val.second = getValue(fs.pressure(oilPhaseIdx));
|
|
|
|
else if (FluidSystem::phaseIsActive(gasPhaseIdx))
|
|
|
|
val.second = getValue(fs.pressure(gasPhaseIdx));
|
|
|
|
else if (FluidSystem::phaseIsActive(waterPhaseIdx))
|
|
|
|
val.second = getValue(fs.pressure(waterPhaseIdx));
|
|
|
|
} else if ((key.first == "BTCNFHEA") || (key.first == "BTEMP")) {
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx))
|
|
|
|
val.second = getValue(fs.temperature(oilPhaseIdx));
|
|
|
|
else if (FluidSystem::phaseIsActive(gasPhaseIdx))
|
|
|
|
val.second = getValue(fs.temperature(gasPhaseIdx));
|
|
|
|
else if (FluidSystem::phaseIsActive(waterPhaseIdx))
|
|
|
|
val.second = getValue(fs.temperature(waterPhaseIdx));
|
|
|
|
} else if (key.first == "BWKR" || key.first == "BKRW")
|
|
|
|
val.second = getValue(intQuants.relativePermeability(waterPhaseIdx));
|
|
|
|
else if (key.first == "BGKR" || key.first == "BKRG")
|
|
|
|
val.second = getValue(intQuants.relativePermeability(gasPhaseIdx));
|
|
|
|
else if (key.first == "BOKR" || key.first == "BKRO")
|
|
|
|
val.second = getValue(intQuants.relativePermeability(oilPhaseIdx));
|
|
|
|
else if (key.first == "BKROG") {
|
|
|
|
const auto& materialParams = problem.materialLawParams(elemCtx, dofIdx, /* timeIdx = */ 0);
|
|
|
|
const auto krog
|
|
|
|
= MaterialLaw::template relpermOilInOilGasSystem<Evaluation>(materialParams, fs);
|
|
|
|
val.second = getValue(krog);
|
|
|
|
} else if (key.first == "BKROW") {
|
|
|
|
const auto& materialParams = problem.materialLawParams(elemCtx, dofIdx, /* timeIdx = */ 0);
|
|
|
|
const auto krow
|
|
|
|
= MaterialLaw::template relpermOilInOilWaterSystem<Evaluation>(materialParams, fs);
|
|
|
|
val.second = getValue(krow);
|
|
|
|
} else if (key.first == "BWPC")
|
|
|
|
val.second = getValue(fs.pressure(oilPhaseIdx)) - getValue(fs.pressure(waterPhaseIdx));
|
|
|
|
else if (key.first == "BGPC")
|
|
|
|
val.second = getValue(fs.pressure(gasPhaseIdx)) - getValue(fs.pressure(oilPhaseIdx));
|
|
|
|
else if (key.first == "BWPR")
|
|
|
|
val.second = getValue(fs.pressure(waterPhaseIdx));
|
|
|
|
else if (key.first == "BGPR")
|
|
|
|
val.second = getValue(fs.pressure(gasPhaseIdx));
|
|
|
|
else if (key.first == "BVWAT" || key.first == "BWVIS")
|
|
|
|
val.second = getValue(fs.viscosity(waterPhaseIdx));
|
|
|
|
else if (key.first == "BVGAS" || key.first == "BGVIS")
|
|
|
|
val.second = getValue(fs.viscosity(gasPhaseIdx));
|
|
|
|
else if (key.first == "BVOIL" || key.first == "BOVIS")
|
|
|
|
val.second = getValue(fs.viscosity(oilPhaseIdx));
|
|
|
|
else if ((key.first == "BRPV") || (key.first == "BOPV") || (key.first == "BWPV")
|
|
|
|
|| (key.first == "BGPV")) {
|
|
|
|
if (key.first == "BRPV") {
|
|
|
|
val.second = 1.0;
|
|
|
|
} else if (key.first == "BOPV") {
|
|
|
|
val.second = getValue(fs.saturation(oilPhaseIdx));
|
|
|
|
} else if (key.first == "BWPV") {
|
|
|
|
val.second = getValue(fs.saturation(waterPhaseIdx));
|
|
|
|
} else {
|
|
|
|
val.second = getValue(fs.saturation(gasPhaseIdx));
|
|
|
|
}
|
|
|
|
|
|
|
|
// Include active pore-volume.
|
|
|
|
val.second *= elemCtx.simulator().model().dofTotalVolume(globalDofIdx)
|
|
|
|
* getValue(intQuants.porosity());
|
|
|
|
} else if (key.first == "BRS")
|
|
|
|
val.second = getValue(fs.Rs());
|
|
|
|
else if (key.first == "BRV")
|
|
|
|
val.second = getValue(fs.Rv());
|
|
|
|
else if ((key.first == "BOIP") || (key.first == "BOIPL") || (key.first == "BOIPG")
|
|
|
|
|| (key.first == "BGIP") || (key.first == "BGIPL") || (key.first == "BGIPG")
|
|
|
|
|| (key.first == "BWIP")) {
|
|
|
|
if ((key.first == "BOIP") || (key.first == "BOIPL")) {
|
|
|
|
val.second = getValue(fs.invB(oilPhaseIdx)) * getValue(fs.saturation(oilPhaseIdx));
|
|
|
|
|
|
|
|
if (key.first == "BOIP") {
|
|
|
|
val.second += getValue(fs.Rv()) * getValue(fs.invB(gasPhaseIdx))
|
|
|
|
* getValue(fs.saturation(gasPhaseIdx));
|
|
|
|
}
|
|
|
|
} else if (key.first == "BOIPG") {
|
|
|
|
val.second = getValue(fs.Rv()) * getValue(fs.invB(gasPhaseIdx))
|
|
|
|
* getValue(fs.saturation(gasPhaseIdx));
|
|
|
|
} else if ((key.first == "BGIP") || (key.first == "BGIPG")) {
|
|
|
|
val.second = getValue(fs.invB(gasPhaseIdx)) * getValue(fs.saturation(gasPhaseIdx));
|
|
|
|
|
|
|
|
if (key.first == "BGIP") {
|
|
|
|
val.second += getValue(fs.Rs()) * getValue(fs.invB(oilPhaseIdx))
|
|
|
|
* getValue(fs.saturation(oilPhaseIdx));
|
|
|
|
}
|
|
|
|
} else if (key.first == "BGIPL") {
|
|
|
|
val.second = getValue(fs.Rs()) * getValue(fs.invB(oilPhaseIdx))
|
|
|
|
* getValue(fs.saturation(oilPhaseIdx));
|
|
|
|
} else { // BWIP
|
|
|
|
val.second = getValue(fs.invB(waterPhaseIdx)) * getValue(fs.saturation(waterPhaseIdx));
|
|
|
|
}
|
|
|
|
|
|
|
|
// Include active pore-volume.
|
|
|
|
val.second *= elemCtx.simulator().model().dofTotalVolume(globalDofIdx)
|
|
|
|
* getValue(intQuants.porosity());
|
|
|
|
} else {
|
|
|
|
std::string logstring = "Keyword '";
|
|
|
|
logstring.append(key.first);
|
|
|
|
logstring.append("' is unhandled for output to file.");
|
|
|
|
OpmLog::warning("Unhandled output keyword", logstring);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2023-03-24 08:56:23 -05:00
|
|
|
|
2022-02-14 05:30:32 -06:00
|
|
|
/*!
|
|
|
|
* \brief Capture connection fluxes, particularly to account for inter-region flows.
|
|
|
|
*
|
|
|
|
* \tparam ActiveIndex Callable type, typically a lambda, that enables
|
|
|
|
* retrieving the active index, on the local MPI rank, of a
|
|
|
|
* particular cell/element. Must support a function call operator of
|
|
|
|
* the form
|
|
|
|
\code
|
|
|
|
int operator()(const Element& elem) const
|
|
|
|
\endcode
|
|
|
|
*
|
|
|
|
* \tparam CartesianIndex Callable type, typically a lambda, that
|
|
|
|
* enables retrieving the globally unique Cartesian index of a
|
|
|
|
* particular cell/element given its active index on the local MPI
|
|
|
|
* rank. Must support a function call operator of the form
|
|
|
|
\code
|
|
|
|
int operator()(const int activeIndex) const
|
|
|
|
\endcode
|
|
|
|
*
|
|
|
|
* \param[in] elemCtx Primary lookup structure for per-cell/element
|
|
|
|
* dynamic information.
|
|
|
|
*
|
|
|
|
* \param[in] activeIndex Mapping from cell/elements to linear indices
|
|
|
|
* on local MPI rank.
|
|
|
|
*
|
|
|
|
* \param[in] cartesianIndex Mapping from active index on local MPI rank
|
|
|
|
* to globally unique Cartesian cell/element index.
|
|
|
|
*/
|
|
|
|
template <class ActiveIndex, class CartesianIndex>
|
2023-01-16 01:53:57 -06:00
|
|
|
void processFluxes(const ElementContext& elemCtx, ActiveIndex&& activeIndex, CartesianIndex&& cartesianIndex)
|
2022-02-14 05:30:32 -06:00
|
|
|
{
|
2023-02-15 04:05:45 -06:00
|
|
|
OPM_TIMEBLOCK_LOCAL(processFluxes);
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto identifyCell = [&activeIndex, &cartesianIndex](const Element& elem) -> EclInterRegFlowMap::Cell {
|
2022-02-14 05:30:32 -06:00
|
|
|
const auto cellIndex = activeIndex(elem);
|
|
|
|
|
|
|
|
return {
|
2023-01-16 01:53:57 -06:00
|
|
|
static_cast<int>(cellIndex), cartesianIndex(cellIndex), elem.partitionType() == Dune::InteriorEntity};
|
2022-02-14 05:30:32 -06:00
|
|
|
};
|
|
|
|
|
|
|
|
const auto timeIdx = 0u;
|
|
|
|
const auto& stencil = elemCtx.stencil(timeIdx);
|
|
|
|
const auto numInteriorFaces = elemCtx.numInteriorFaces(timeIdx);
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
for (auto scvfIdx = 0 * numInteriorFaces; scvfIdx < numInteriorFaces; ++scvfIdx) {
|
2022-02-14 05:30:32 -06:00
|
|
|
const auto& face = stencil.interiorFace(scvfIdx);
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto left = identifyCell(stencil.element(face.interiorIndex()));
|
2022-02-14 05:30:32 -06:00
|
|
|
const auto right = identifyCell(stencil.element(face.exteriorIndex()));
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto rates = this->getComponentSurfaceRates(elemCtx, face.area(), scvfIdx, timeIdx);
|
2022-02-14 05:30:32 -06:00
|
|
|
|
|
|
|
this->interRegionFlows_.addConnection(left, right, rates);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Prepare for capturing connection fluxes, particularly to
|
|
|
|
* account for inter-region flows.
|
|
|
|
*/
|
|
|
|
void initializeFluxData()
|
|
|
|
{
|
|
|
|
// Inter-region flow rates. Note: ".clear()" prepares to accumulate
|
|
|
|
// contributions per bulk connection between FIP regions.
|
|
|
|
this->interRegionFlows_.clear();
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Finalize capturing connection fluxes.
|
|
|
|
*/
|
|
|
|
void finalizeFluxData()
|
|
|
|
{
|
|
|
|
this->interRegionFlows_.compress();
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Get read-only access to collection of inter-region flows.
|
|
|
|
*/
|
|
|
|
const EclInterRegFlowMap& getInterRegFlows() const
|
|
|
|
{
|
|
|
|
return this->interRegionFlows_;
|
|
|
|
}
|
|
|
|
|
2017-12-07 05:38:48 -06:00
|
|
|
template <class FluidState>
|
2018-04-30 10:17:22 -05:00
|
|
|
void assignToFluidState(FluidState& fs, unsigned elemIdx) const
|
2017-12-07 05:38:48 -06:00
|
|
|
{
|
2023-01-16 01:53:57 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->saturation_[phaseIdx].empty())
|
2018-01-09 07:01:30 -06:00
|
|
|
continue;
|
2015-01-25 10:27:08 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
fs.setSaturation(phaseIdx, this->saturation_[phaseIdx][elemIdx]);
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
2018-01-09 07:01:30 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
if (!this->fluidPressure_.empty()) {
|
2017-12-07 05:38:48 -06:00
|
|
|
// this assumes that capillary pressures only depend on the phase saturations
|
|
|
|
// and possibly on temperature. (this is always the case for ECL problems.)
|
2022-08-05 03:35:52 -05:00
|
|
|
std::array<Scalar, numPhases> pc = {0};
|
2017-12-07 05:38:48 -06:00
|
|
|
const MaterialLawParams& matParams = simulator_.problem().materialLawParams(elemIdx);
|
|
|
|
MaterialLaw::capillaryPressures(pc, matParams, fs);
|
2023-01-16 01:53:57 -06:00
|
|
|
Valgrind::CheckDefined(this->fluidPressure_[elemIdx]);
|
2021-05-05 04:22:44 -05:00
|
|
|
Valgrind::CheckDefined(pc);
|
2017-12-07 05:38:48 -06:00
|
|
|
assert(FluidSystem::phaseIsActive(oilPhaseIdx));
|
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
|
|
|
if (!FluidSystem::phaseIsActive(phaseIdx))
|
|
|
|
continue;
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
fs.setPressure(phaseIdx, this->fluidPressure_[elemIdx] + (pc[phaseIdx] - pc[oilPhaseIdx]));
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
2015-01-25 10:27:08 -06:00
|
|
|
}
|
2017-12-07 05:38:48 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->temperature_.empty())
|
|
|
|
fs.setTemperature(this->temperature_[elemIdx]);
|
|
|
|
if (!this->rs_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
fs.setRs(this->rs_[elemIdx]);
|
2022-11-24 07:32:41 -06:00
|
|
|
if (!this->rsw_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
fs.setRsw(this->rsw_[elemIdx]);
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->rv_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
fs.setRv(this->rv_[elemIdx]);
|
2022-03-01 10:54:22 -06:00
|
|
|
if (!this->rvw_.empty())
|
2023-01-16 01:53:57 -06:00
|
|
|
fs.setRvw(this->rvw_[elemIdx]);
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
|
|
|
|
2019-02-01 10:33:30 -06:00
|
|
|
void initHysteresisParams(Simulator& simulator, unsigned elemIdx) const
|
|
|
|
{
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->soMax_.empty())
|
|
|
|
simulator.problem().setMaxOilSaturation(elemIdx, this->soMax_[elemIdx]);
|
2017-12-07 05:38:48 -06:00
|
|
|
|
|
|
|
if (simulator.problem().materialLawManager()->enableHysteresis()) {
|
|
|
|
auto matLawManager = simulator.problem().materialLawManager();
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->pcSwMdcOw_.empty() && !this->krnSwMdcOw_.empty()) {
|
2017-12-07 05:38:48 -06:00
|
|
|
matLawManager->setOilWaterHysteresisParams(
|
2023-01-16 01:53:57 -06:00
|
|
|
this->pcSwMdcOw_[elemIdx], this->krnSwMdcOw_[elemIdx], elemIdx);
|
2018-01-09 07:01:30 -06:00
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->pcSwMdcGo_.empty() && !this->krnSwMdcGo_.empty()) {
|
2017-12-07 05:38:48 -06:00
|
|
|
matLawManager->setGasOilHysteresisParams(
|
2023-01-16 01:53:57 -06:00
|
|
|
this->pcSwMdcGo_[elemIdx], this->krnSwMdcGo_[elemIdx], elemIdx);
|
2018-01-09 07:01:30 -06:00
|
|
|
}
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2020-01-13 08:46:50 -06:00
|
|
|
if (simulator_.vanguard().eclState().fieldProps().has_double("SWATINIT")) {
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto& oilWaterScaledEpsInfoDrainage
|
|
|
|
= simulator.problem().materialLawManager()->oilWaterScaledEpsInfoDrainage(elemIdx);
|
2021-12-14 06:40:47 -06:00
|
|
|
const_cast<EclEpsScalingPointsInfo<Scalar>&>(oilWaterScaledEpsInfoDrainage).maxPcow = this->ppcw_[elemIdx];
|
2018-09-07 03:19:00 -05:00
|
|
|
}
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
|
|
|
|
2023-03-24 08:56:23 -05:00
|
|
|
void updateFluidInPlace(const ElementContext& elemCtx)
|
|
|
|
{
|
2023-03-21 16:29:56 -05:00
|
|
|
for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
|
|
|
|
updateFluidInPlace_(elemCtx, dofIdx);
|
|
|
|
}
|
|
|
|
}
|
2023-03-24 08:56:23 -05:00
|
|
|
void updateFluidInPlace(unsigned globalDofIdx,const IntensiveQuantities& intQuants, double totVolume)
|
|
|
|
{
|
2023-03-23 07:44:46 -05:00
|
|
|
this->updateFluidInPlace_(globalDofIdx, intQuants, totVolume);
|
|
|
|
}
|
2014-11-28 05:58:48 -06:00
|
|
|
private:
|
2021-05-17 16:18:36 -05:00
|
|
|
bool isDefunctParallelWell(std::string wname) const override
|
2020-09-22 07:12:15 -05:00
|
|
|
{
|
2023-01-16 01:53:57 -06:00
|
|
|
if (simulator_.gridView().comm().size() == 1)
|
2020-09-22 07:12:15 -05:00
|
|
|
return false;
|
|
|
|
const auto& parallelWells = simulator_.vanguard().parallelWells();
|
2023-01-16 01:53:57 -06:00
|
|
|
std::pair<std::string, bool> value {wname, true};
|
|
|
|
auto candidate = std::lower_bound(parallelWells.begin(), parallelWells.end(), value);
|
2020-09-22 07:12:15 -05:00
|
|
|
return candidate == parallelWells.end() || *candidate != value;
|
|
|
|
}
|
|
|
|
|
2023-03-24 08:56:23 -05:00
|
|
|
void updateFluidInPlace_(const ElementContext& elemCtx, unsigned dofIdx)
|
|
|
|
{
|
2023-03-23 07:44:46 -05:00
|
|
|
const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
|
|
|
|
unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
|
|
|
const auto totVolume = elemCtx.simulator().model().dofTotalVolume(globalDofIdx);
|
|
|
|
this->updateFluidInPlace_(globalDofIdx, intQuants, totVolume);
|
|
|
|
}
|
2023-03-24 08:56:23 -05:00
|
|
|
|
2023-03-23 07:44:46 -05:00
|
|
|
void updateFluidInPlace_(unsigned globalDofIdx,const IntensiveQuantities& intQuants, double totVolume)
|
2018-01-29 01:43:42 -06:00
|
|
|
{
|
2023-02-15 04:05:45 -06:00
|
|
|
OPM_TIMEBLOCK_LOCAL(updateFluidInPlace);
|
2023-03-23 07:44:46 -05:00
|
|
|
//const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
|
2018-01-29 01:43:42 -06:00
|
|
|
const auto& fs = intQuants.fluidState();
|
2023-03-23 07:44:46 -05:00
|
|
|
//unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
2018-01-29 01:43:42 -06:00
|
|
|
|
|
|
|
// Fluid in Place calculations
|
|
|
|
|
|
|
|
// 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.
|
2023-03-23 07:44:46 -05:00
|
|
|
//const auto totVolume = elemCtx.simulator().model().dofTotalVolume(globalDofIdx);
|
2021-06-21 09:15:10 -05:00
|
|
|
const double pv = totVolume * intQuants.porosity().value();
|
2018-01-29 01:43:42 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->pressureTimesHydrocarbonVolume_.empty() && !this->pressureTimesPoreVolume_.empty()) {
|
2023-01-16 01:53:57 -06:00
|
|
|
assert(this->hydrocarbonPoreVolume_.size() == this->pressureTimesHydrocarbonVolume_.size());
|
2021-05-17 16:18:36 -05:00
|
|
|
assert(this->fip_[Inplace::Phase::PoreVolume].size() == this->pressureTimesPoreVolume_.size());
|
2018-01-29 01:43:42 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
this->fip_[Inplace::Phase::PoreVolume][globalDofIdx] = totVolume * intQuants.referencePorosity();
|
2021-06-21 09:15:10 -05:00
|
|
|
|
|
|
|
this->dynamicPoreVolume_[globalDofIdx] = pv;
|
2018-01-29 01:43:42 -06:00
|
|
|
|
|
|
|
Scalar hydrocarbon = 0.0;
|
|
|
|
|
2023-03-21 03:08:10 -05:00
|
|
|
if (!this->eclState_.runspec().co2Storage()) {
|
|
|
|
// Common case. Hydrocarbon volume is fraction occupied by actual hydrocarbons.
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx))
|
|
|
|
hydrocarbon += getValue(fs.saturation(oilPhaseIdx));
|
|
|
|
if (FluidSystem::phaseIsActive(gasPhaseIdx))
|
|
|
|
hydrocarbon += getValue(fs.saturation(gasPhaseIdx));
|
|
|
|
} else {
|
|
|
|
// CO2 storage: Hydrocarbon volume is full pore-volume.
|
|
|
|
hydrocarbon = 1.0;
|
|
|
|
}
|
2021-05-17 16:18:36 -05:00
|
|
|
this->hydrocarbonPoreVolume_[globalDofIdx] = pv * hydrocarbon;
|
2018-01-29 01:43:42 -06:00
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
|
2021-05-17 16:18:36 -05:00
|
|
|
this->pressureTimesPoreVolume_[globalDofIdx] = getValue(fs.pressure(oilPhaseIdx)) * pv;
|
2023-01-16 01:53:57 -06:00
|
|
|
this->pressureTimesHydrocarbonVolume_[globalDofIdx]
|
|
|
|
= this->pressureTimesPoreVolume_[globalDofIdx] * hydrocarbon;
|
2021-09-09 06:19:06 -05:00
|
|
|
} else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
|
|
|
|
this->pressureTimesPoreVolume_[globalDofIdx] = getValue(fs.pressure(gasPhaseIdx)) * pv;
|
2023-01-16 01:53:57 -06:00
|
|
|
this->pressureTimesHydrocarbonVolume_[globalDofIdx]
|
|
|
|
= this->pressureTimesPoreVolume_[globalDofIdx] * hydrocarbon;
|
2021-09-09 06:19:06 -05:00
|
|
|
} else if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
|
2021-10-06 12:32:35 -05:00
|
|
|
this->pressureTimesPoreVolume_[globalDofIdx] = getValue(fs.pressure(waterPhaseIdx)) * pv;
|
2021-09-09 06:19:06 -05:00
|
|
|
}
|
2018-01-29 01:43:42 -06:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (this->computeFip_) {
|
2018-01-29 01:43:42 -06:00
|
|
|
Scalar fip[FluidSystem::numPhases];
|
2021-08-25 11:24:05 -05:00
|
|
|
Scalar fipr[FluidSystem::numPhases]; // at reservoir condition
|
2023-02-22 05:24:00 -06:00
|
|
|
|
2018-04-30 10:17:22 -05:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
|
|
|
|
fip[phaseIdx] = 0.0;
|
|
|
|
|
|
|
|
if (!FluidSystem::phaseIsActive(phaseIdx))
|
2018-01-29 01:43:42 -06:00
|
|
|
continue;
|
|
|
|
|
2021-05-05 04:22:44 -05:00
|
|
|
const double b = getValue(fs.invB(phaseIdx));
|
2021-10-06 12:32:35 -05:00
|
|
|
const double s = getValue(fs.saturation(phaseIdx));
|
2021-08-25 11:24:05 -05:00
|
|
|
fipr[phaseIdx] = s * pv;
|
|
|
|
fip[phaseIdx] = b * fipr[phaseIdx];
|
2018-01-29 01:43:42 -06:00
|
|
|
}
|
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx) && !this->fip_[Inplace::Phase::OIL].empty())
|
|
|
|
this->fip_[Inplace::Phase::OIL][globalDofIdx] = fip[oilPhaseIdx];
|
|
|
|
if (FluidSystem::phaseIsActive(gasPhaseIdx) && !this->fip_[Inplace::Phase::GAS].empty())
|
|
|
|
this->fip_[Inplace::Phase::GAS][globalDofIdx] = fip[gasPhaseIdx];
|
|
|
|
if (FluidSystem::phaseIsActive(waterPhaseIdx) && !this->fip_[Inplace::Phase::WATER].empty())
|
|
|
|
this->fip_[Inplace::Phase::WATER][globalDofIdx] = fip[waterPhaseIdx];
|
2018-01-29 01:43:42 -06:00
|
|
|
|
2021-08-25 11:24:05 -05:00
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx) && !this->fip_[Inplace::Phase::OilResVolume].empty())
|
|
|
|
this->fip_[Inplace::Phase::OilResVolume][globalDofIdx] = fipr[oilPhaseIdx];
|
|
|
|
if (FluidSystem::phaseIsActive(gasPhaseIdx) && !this->fip_[Inplace::Phase::GasResVolume].empty())
|
|
|
|
this->fip_[Inplace::Phase::GasResVolume][globalDofIdx] = fipr[gasPhaseIdx];
|
|
|
|
if (FluidSystem::phaseIsActive(waterPhaseIdx) && !this->fip_[Inplace::Phase::WaterResVolume].empty())
|
|
|
|
this->fip_[Inplace::Phase::WaterResVolume][globalDofIdx] = fipr[waterPhaseIdx];
|
2021-10-06 12:32:35 -05:00
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(waterPhaseIdx) && !this->fip_[Inplace::Phase::SALT].empty())
|
2021-08-26 11:33:49 -05:00
|
|
|
this->fip_[Inplace::Phase::SALT][globalDofIdx] = fipr[waterPhaseIdx] * fs.saltConcentration().value();
|
2021-08-25 11:24:05 -05:00
|
|
|
|
2018-01-29 05:23:23 -06:00
|
|
|
// Store the pure oil and gas Fip
|
2021-05-17 16:18:36 -05:00
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx) && !this->fip_[Inplace::Phase::OilInLiquidPhase].empty())
|
|
|
|
this->fip_[Inplace::Phase::OilInLiquidPhase][globalDofIdx] = fip[oilPhaseIdx];
|
2018-01-29 01:43:42 -06:00
|
|
|
|
2021-05-17 16:18:36 -05:00
|
|
|
if (FluidSystem::phaseIsActive(gasPhaseIdx) && !this->fip_[Inplace::Phase::GasInGasPhase].empty())
|
|
|
|
this->fip_[Inplace::Phase::GasInGasPhase][globalDofIdx] = fip[gasPhaseIdx];
|
2018-01-29 01:43:42 -06:00
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx) && FluidSystem::phaseIsActive(gasPhaseIdx)) {
|
|
|
|
// Gas dissolved in oil and vaporized oil
|
2021-05-05 04:22:44 -05:00
|
|
|
Scalar gasInPlaceLiquid = getValue(fs.Rs()) * fip[oilPhaseIdx];
|
|
|
|
Scalar oilInPlaceGas = getValue(fs.Rv()) * fip[gasPhaseIdx];
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->fip_[Inplace::Phase::GasInLiquidPhase].empty())
|
|
|
|
this->fip_[Inplace::Phase::GasInLiquidPhase][globalDofIdx] = gasInPlaceLiquid;
|
|
|
|
if (!this->fip_[Inplace::Phase::OilInGasPhase].empty())
|
|
|
|
this->fip_[Inplace::Phase::OilInGasPhase][globalDofIdx] = oilInPlaceGas;
|
2018-01-29 05:23:23 -06:00
|
|
|
|
|
|
|
// Add dissolved gas and vaporized oil to total Fip
|
2021-05-17 16:18:36 -05:00
|
|
|
if (!this->fip_[Inplace::Phase::OIL].empty())
|
|
|
|
this->fip_[Inplace::Phase::OIL][globalDofIdx] += oilInPlaceGas;
|
|
|
|
if (!this->fip_[Inplace::Phase::GAS].empty())
|
|
|
|
this->fip_[Inplace::Phase::GAS][globalDofIdx] += gasInPlaceLiquid;
|
2018-01-29 01:43:42 -06:00
|
|
|
}
|
2023-02-22 05:24:00 -06:00
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(waterPhaseIdx) && FluidSystem::phaseIsActive(gasPhaseIdx)) {
|
|
|
|
// Gas dissolved in water and vaporized water
|
|
|
|
Scalar gasInPlaceWater = getValue(fs.Rsw()) * fip[waterPhaseIdx];
|
|
|
|
Scalar waterInPlaceGas = getValue(fs.Rvw()) * fip[gasPhaseIdx];
|
|
|
|
|
|
|
|
if (!this->fip_[Inplace::Phase::WaterInGasPhase].empty())
|
|
|
|
this->fip_[Inplace::Phase::WaterInGasPhase][globalDofIdx] = waterInPlaceGas;
|
|
|
|
|
|
|
|
if (!this->fip_[Inplace::Phase::WaterInWaterPhase].empty())
|
|
|
|
this->fip_[Inplace::Phase::WaterInWaterPhase][globalDofIdx] = fip[waterPhaseIdx];
|
|
|
|
|
2023-03-24 09:44:16 -05:00
|
|
|
// For water+gas cases the gas in water is added to the GIPL value
|
|
|
|
if (!this->fip_[Inplace::Phase::GasInLiquidPhase].empty() && !FluidSystem::phaseIsActive(oilPhaseIdx))
|
|
|
|
this->fip_[Inplace::Phase::GasInLiquidPhase][globalDofIdx] = gasInPlaceWater;
|
|
|
|
|
2023-02-22 05:24:00 -06:00
|
|
|
// Add dissolved gas and vaporized water to total Fip
|
|
|
|
if (!this->fip_[Inplace::Phase::WATER].empty())
|
|
|
|
this->fip_[Inplace::Phase::WATER][globalDofIdx] += waterInPlaceGas;
|
|
|
|
if (!this->fip_[Inplace::Phase::GAS].empty())
|
|
|
|
this->fip_[Inplace::Phase::GAS][globalDofIdx] += gasInPlaceWater;
|
|
|
|
}
|
2023-05-10 07:05:00 -05:00
|
|
|
|
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(gasPhaseIdx) &&
|
|
|
|
(!this->fip_[Inplace::Phase::CO2InGasPhaseInMob].empty() || !this->fip_[Inplace::Phase::CO2InGasPhaseMob].empty())) {
|
|
|
|
|
|
|
|
const auto& scaledDrainageInfo
|
|
|
|
= simulator_.problem().materialLawManager()->oilWaterScaledEpsInfoDrainage(globalDofIdx);
|
|
|
|
Scalar sgcr = scaledDrainageInfo.Sgcr;
|
|
|
|
if(simulator_.problem().materialLawManager()->enableHysteresis()) {
|
|
|
|
const MaterialLawParams& matParams = simulator_.problem().materialLawParams(globalDofIdx);
|
|
|
|
sgcr = MaterialLaw::trappedGasSaturation(matParams);
|
|
|
|
}
|
|
|
|
|
2023-02-22 05:24:00 -06:00
|
|
|
const double sg = getValue(fs.saturation(gasPhaseIdx));
|
2023-03-16 02:52:15 -05:00
|
|
|
const double rhog = getValue(fs.density(gasPhaseIdx));
|
|
|
|
const double xgW = FluidSystem::phaseIsActive(waterPhaseIdx)?
|
|
|
|
FluidSystem::convertRvwToXgW(getValue(fs.Rvw()), fs.pvtRegionIndex()):
|
|
|
|
FluidSystem::convertRvToXgO(getValue(fs.Rv()), fs.pvtRegionIndex());
|
|
|
|
Scalar mM = FluidSystem::molarMass(gasCompIdx, fs.pvtRegionIndex());
|
2023-05-10 07:05:00 -05:00
|
|
|
if (!this->fip_[Inplace::Phase::CO2InGasPhaseInMob].empty()) {
|
|
|
|
Scalar imMobileGas = (1 - xgW) * pv * rhog / mM * std::min(sgcr , sg);
|
|
|
|
this->fip_[Inplace::Phase::CO2InGasPhaseInMob][globalDofIdx] = imMobileGas;
|
|
|
|
}
|
|
|
|
if (!this->fip_[Inplace::Phase::CO2InGasPhaseMob].empty()) {
|
|
|
|
Scalar mobileGas = (1 - xgW) * pv * rhog / mM * std::max(0.0, sg - sgcr);
|
|
|
|
this->fip_[Inplace::Phase::CO2InGasPhaseMob][globalDofIdx] = mobileGas;
|
|
|
|
}
|
2023-03-16 02:52:15 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(waterPhaseIdx) && !this->fip_[Inplace::Phase::CO2InWaterPhase].empty()) {
|
|
|
|
const double rhow = getValue(fs.density(waterPhaseIdx));
|
|
|
|
const double sw = getValue(fs.saturation(waterPhaseIdx));
|
|
|
|
const double xwG = FluidSystem::convertRswToXwG(getValue(fs.Rsw()), fs.pvtRegionIndex());
|
|
|
|
Scalar mM = FluidSystem::molarMass(gasCompIdx, fs.pvtRegionIndex());
|
|
|
|
Scalar co2inWater = xwG * pv * rhow * sw / mM;
|
|
|
|
this->fip_[Inplace::Phase::CO2InWaterPhase][globalDofIdx] = co2inWater;
|
|
|
|
}
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx) && !this->fip_[Inplace::Phase::CO2InWaterPhase].empty()) {
|
|
|
|
const double rhoo = getValue(fs.density(oilPhaseIdx));
|
|
|
|
const double so = getValue(fs.saturation(oilPhaseIdx));
|
|
|
|
const double xoG = FluidSystem::convertRsToXoG(getValue(fs.Rs()), fs.pvtRegionIndex());
|
|
|
|
Scalar mM = FluidSystem::molarMass(gasCompIdx, fs.pvtRegionIndex());
|
|
|
|
Scalar co2inWater = xoG * pv * rhoo * so / mM;
|
|
|
|
this->fip_[Inplace::Phase::CO2InWaterPhase][globalDofIdx] = co2inWater;
|
2023-02-22 05:24:00 -06:00
|
|
|
}
|
2018-01-29 01:43:42 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2020-09-07 03:22:42 -05:00
|
|
|
void createLocalRegion_(std::vector<int>& region)
|
2018-01-12 08:32:43 -06:00
|
|
|
{
|
2020-02-25 07:27:52 -06:00
|
|
|
size_t elemIdx = 0;
|
2022-10-12 07:27:20 -05:00
|
|
|
for (const auto& elem : elements(simulator_.gridView())) {
|
2020-02-25 07:27:52 -06:00
|
|
|
if (elem.partitionType() != Dune::InteriorEntity)
|
2020-09-07 03:22:42 -05:00
|
|
|
region[elemIdx] = 0;
|
2022-10-12 07:27:20 -05:00
|
|
|
++elemIdx;
|
2018-01-12 08:32:43 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2022-02-14 05:30:32 -06:00
|
|
|
/*!
|
|
|
|
* \brief Compute surface level component flow rates across a single
|
|
|
|
* intersection.
|
|
|
|
*
|
|
|
|
* \param[in] elemCtx Primary lookup structure for per-cell/element
|
|
|
|
* dynamic information.
|
|
|
|
*
|
|
|
|
* \param[in] scvfIdx Linear index of current interior bulk connection.
|
|
|
|
*
|
|
|
|
* \param[in] timeIdx Historical time-point at which to evaluate dynamic
|
|
|
|
* quantities (e.g., reciprocal FVF or dissolved gas concentration).
|
|
|
|
* Zero for the current time.
|
|
|
|
*
|
|
|
|
* \return Surface level component flow rates.
|
|
|
|
*/
|
2023-01-16 01:53:57 -06:00
|
|
|
data::InterRegFlowMap::FlowRates getComponentSurfaceRates(const ElementContext& elemCtx,
|
|
|
|
const Scalar faceArea,
|
|
|
|
const std::size_t scvfIdx,
|
|
|
|
const std::size_t timeIdx) const
|
2022-02-14 05:30:32 -06:00
|
|
|
{
|
|
|
|
using Component = data::InterRegFlowMap::Component;
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
auto rates = data::InterRegFlowMap::FlowRates {};
|
2022-02-14 05:30:32 -06:00
|
|
|
|
|
|
|
const auto& extQuant = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
|
|
|
|
|
|
|
|
const auto alpha = getValue(extQuant.extrusionFactor()) * faceArea;
|
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto& up = elemCtx.intensiveQuantities(extQuant.upstreamIndex(oilPhaseIdx), timeIdx);
|
2022-02-14 05:30:32 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
using FluidState = std::remove_cv_t<std::remove_reference_t<decltype(up.fluidState())>>;
|
2022-02-14 05:30:32 -06:00
|
|
|
|
|
|
|
const auto pvtReg = up.pvtRegionIndex();
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto bO = getValue(getInvB_<FluidSystem, FluidState, Scalar>(up.fluidState(), oilPhaseIdx, pvtReg));
|
2022-02-14 05:30:32 -06:00
|
|
|
|
|
|
|
const auto qO = alpha * bO * getValue(extQuant.volumeFlux(oilPhaseIdx));
|
|
|
|
|
|
|
|
rates[Component::Oil] += qO;
|
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto Rs = getValue(BlackOil::getRs_<FluidSystem, FluidState, Scalar>(up.fluidState(), pvtReg));
|
2022-02-14 05:30:32 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
rates[Component::Gas] += qO * Rs;
|
2022-02-14 05:30:32 -06:00
|
|
|
rates[Component::Disgas] += qO * Rs;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto& up = elemCtx.intensiveQuantities(extQuant.upstreamIndex(gasPhaseIdx), timeIdx);
|
2022-02-14 05:30:32 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
using FluidState = std::remove_cv_t<std::remove_reference_t<decltype(up.fluidState())>>;
|
2022-02-14 05:30:32 -06:00
|
|
|
|
|
|
|
const auto pvtReg = up.pvtRegionIndex();
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto bG = getValue(getInvB_<FluidSystem, FluidState, Scalar>(up.fluidState(), gasPhaseIdx, pvtReg));
|
2022-02-14 05:30:32 -06:00
|
|
|
|
|
|
|
const auto qG = alpha * bG * getValue(extQuant.volumeFlux(gasPhaseIdx));
|
|
|
|
|
|
|
|
rates[Component::Gas] += qG;
|
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto Rv = getValue(BlackOil::getRv_<FluidSystem, FluidState, Scalar>(up.fluidState(), pvtReg));
|
2022-02-14 05:30:32 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
rates[Component::Oil] += qG * Rv;
|
2022-02-14 05:30:32 -06:00
|
|
|
rates[Component::Vapoil] += qG * Rv;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto& up = elemCtx.intensiveQuantities(extQuant.upstreamIndex(waterPhaseIdx), timeIdx);
|
2022-02-14 05:30:32 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
using FluidState = std::remove_cv_t<std::remove_reference_t<decltype(up.fluidState())>>;
|
2022-02-14 05:30:32 -06:00
|
|
|
|
|
|
|
const auto pvtReg = up.pvtRegionIndex();
|
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
const auto bW = getValue(getInvB_<FluidSystem, FluidState, Scalar>(up.fluidState(), waterPhaseIdx, pvtReg));
|
2022-02-14 05:30:32 -06:00
|
|
|
|
2023-01-16 01:53:57 -06:00
|
|
|
rates[Component::Water] += alpha * bW * getValue(extQuant.volumeFlux(waterPhaseIdx));
|
2022-02-14 05:30:32 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
return rates;
|
|
|
|
}
|
|
|
|
|
2017-12-07 05:38:48 -06:00
|
|
|
const Simulator& simulator_;
|
2023-03-24 08:56:23 -05:00
|
|
|
};
|
2021-05-17 16:18:36 -05:00
|
|
|
|
2019-09-05 10:04:39 -05:00
|
|
|
} // namespace Opm
|
2014-11-28 05:58:48 -06:00
|
|
|
|
|
|
|
#endif
|