2015-06-18 06:43:59 -05:00
|
|
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
|
|
|
|
// vi: set et ts=4 sw=4 sts=4:
|
2014-04-15 11:30:06 -05:00
|
|
|
/*
|
2023-09-09 16:02:09 -05:00
|
|
|
Copyright 2023 INRIA
|
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
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.
|
|
|
|
|
|
|
|
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.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
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-04-15 11:30:06 -05:00
|
|
|
*/
|
|
|
|
/*!
|
|
|
|
* \file
|
|
|
|
*
|
2024-02-05 12:22:15 -06:00
|
|
|
* \copydoc Opm::FlowProblem
|
2014-04-15 11:30:06 -05:00
|
|
|
*/
|
2024-02-05 12:22:15 -06:00
|
|
|
#ifndef OPM_FLOW_PROBLEM_HPP
|
|
|
|
#define OPM_FLOW_PROBLEM_HPP
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <dune/common/version.hh>
|
|
|
|
#include <dune/common/fvector.hh>
|
|
|
|
#include <dune/common/fmatrix.hh>
|
2020-11-10 01:54:44 -06:00
|
|
|
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <opm/common/utility/TimeService.hpp>
|
2021-01-26 07:59:22 -06:00
|
|
|
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
|
|
|
|
#include <opm/input/eclipse/Schedule/Schedule.hpp>
|
2024-10-15 11:09:25 -05:00
|
|
|
#include <opm/input/eclipse/Units/Units.hpp>
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <opm/material/common/ConditionalStorage.hpp>
|
|
|
|
#include <opm/material/common/Valgrind.hpp>
|
2022-09-26 10:13:03 -05:00
|
|
|
#include <opm/material/densead/Evaluation.hpp>
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
|
|
|
|
#include <opm/material/thermal/EclThermalLawManager.hpp>
|
2019-03-12 09:51:41 -05:00
|
|
|
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <opm/models/common/directionalmobility.hh>
|
|
|
|
#include <opm/models/utils/pffgridvector.hh>
|
|
|
|
#include <opm/models/discretization/ecfv/ecfvdiscretization.hh>
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2017-12-07 05:38:48 -06:00
|
|
|
#include <opm/output/eclipse/EclipseIO.hpp>
|
|
|
|
|
2024-02-02 03:46:44 -06:00
|
|
|
#include <opm/simulators/flow/CpGridVanguard.hpp>
|
2024-02-02 03:46:44 -06:00
|
|
|
#include <opm/simulators/flow/DummyGradientCalculator.hpp>
|
2024-02-05 12:22:15 -06:00
|
|
|
#include <opm/simulators/flow/EclWriter.hpp>
|
2024-02-02 03:46:44 -06:00
|
|
|
#include <opm/simulators/flow/EquilInitializer.hpp>
|
2024-02-02 03:46:44 -06:00
|
|
|
#include <opm/simulators/flow/FlowGenericProblem.hpp>
|
2024-09-11 07:25:34 -05:00
|
|
|
// TODO: maybe we can name it FlowProblemProperties.hpp
|
|
|
|
#include <opm/simulators/flow/FlowBaseProblemProperties.hpp>
|
2024-10-03 08:34:35 -05:00
|
|
|
#include <opm/simulators/flow/FlowUtils.hpp>
|
2024-02-05 12:22:15 -06:00
|
|
|
#include <opm/simulators/flow/TracerModel.hpp>
|
2024-02-02 03:46:44 -06:00
|
|
|
#include <opm/simulators/flow/Transmissibility.hpp>
|
2024-01-31 07:14:50 -06:00
|
|
|
#include <opm/simulators/timestepping/AdaptiveTimeStepping.hpp>
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <opm/simulators/timestepping/SimulatorReport.hpp>
|
2024-06-25 03:20:41 -05:00
|
|
|
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
|
|
|
|
#include <opm/simulators/utils/ParallelSerialization.hpp>
|
2024-06-25 03:20:41 -05:00
|
|
|
#include <opm/simulators/utils/satfunc/RelpermDiagnostics.hpp>
|
2023-08-02 05:42:33 -05:00
|
|
|
|
|
|
|
#include <opm/utility/CopyablePtr.hpp>
|
|
|
|
|
2017-12-20 05:22:53 -06:00
|
|
|
#include <algorithm>
|
2024-10-15 11:09:25 -05:00
|
|
|
#include <cstddef>
|
2022-08-12 07:33:15 -05:00
|
|
|
#include <functional>
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <set>
|
2024-10-15 11:09:25 -05:00
|
|
|
#include <stdexcept>
|
2023-08-02 05:42:33 -05:00
|
|
|
#include <string>
|
|
|
|
#include <vector>
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2019-09-05 10:04:39 -05:00
|
|
|
namespace Opm {
|
2014-04-15 11:30:06 -05:00
|
|
|
|
|
|
|
/*!
|
2024-02-05 12:22:15 -06:00
|
|
|
* \ingroup BlackOilSimulator
|
2014-04-15 11:30:06 -05:00
|
|
|
*
|
2014-11-28 05:58:48 -06:00
|
|
|
* \brief This problem simulates an input file given in the data format used by the
|
|
|
|
* commercial ECLiPSE simulator.
|
2014-04-15 11:30:06 -05:00
|
|
|
*/
|
|
|
|
template <class TypeTag>
|
2024-02-05 12:22:15 -06:00
|
|
|
class FlowProblem : public GetPropType<TypeTag, Properties::BaseProblem>
|
|
|
|
, public FlowGenericProblem<GetPropType<TypeTag, Properties::GridView>,
|
2024-02-22 08:17:09 -06:00
|
|
|
GetPropType<TypeTag, Properties::FluidSystem>>
|
2014-04-15 11:30:06 -05:00
|
|
|
{
|
2024-09-11 07:25:34 -05:00
|
|
|
protected:
|
2024-02-02 03:46:44 -06:00
|
|
|
using BaseType = FlowGenericProblem<GetPropType<TypeTag, Properties::GridView>,
|
2024-02-22 08:17:09 -06:00
|
|
|
GetPropType<TypeTag, Properties::FluidSystem>>;
|
2020-08-26 03:49:52 -05:00
|
|
|
using ParentType = GetPropType<TypeTag, Properties::BaseProblem>;
|
|
|
|
using Implementation = GetPropType<TypeTag, Properties::Problem>;
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2020-08-26 03:49:52 -05:00
|
|
|
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
|
|
|
|
using GridView = GetPropType<TypeTag, Properties::GridView>;
|
|
|
|
using Stencil = GetPropType<TypeTag, Properties::Stencil>;
|
|
|
|
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
|
|
|
|
using GlobalEqVector = GetPropType<TypeTag, Properties::GlobalEqVector>;
|
|
|
|
using EqVector = GetPropType<TypeTag, Properties::EqVector>;
|
2021-05-05 05:13:25 -05:00
|
|
|
using Vanguard = GetPropType<TypeTag, Properties::Vanguard>;
|
2014-04-15 11:30:06 -05:00
|
|
|
|
|
|
|
// Grid and world dimension
|
|
|
|
enum { dim = GridView::dimension };
|
|
|
|
enum { dimWorld = GridView::dimensionworld };
|
|
|
|
|
|
|
|
// copy some indices for convenience
|
2020-08-27 02:13:30 -05:00
|
|
|
enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
|
2014-04-15 11:30:06 -05:00
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
|
|
|
enum { numComponents = FluidSystem::numComponents };
|
2024-09-11 07:25:34 -05:00
|
|
|
|
|
|
|
enum { enableConvectiveMixing = getPropValue<TypeTag, Properties::EnableConvectiveMixing>() };
|
|
|
|
enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
|
|
|
|
enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
|
|
|
|
enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
|
|
|
|
enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
|
2020-08-27 02:13:30 -05:00
|
|
|
enum { enableExperiments = getPropValue<TypeTag, Properties::EnableExperiments>() };
|
2024-09-11 07:25:34 -05:00
|
|
|
enum { enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>() };
|
|
|
|
enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
|
|
|
|
enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
|
2020-08-27 02:13:30 -05:00
|
|
|
enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
|
|
|
|
enum { enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>() };
|
2024-09-11 07:25:34 -05:00
|
|
|
enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
|
|
|
|
enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
|
2020-08-27 02:13:30 -05:00
|
|
|
enum { enableTemperature = getPropValue<TypeTag, Properties::EnableTemperature>() };
|
|
|
|
enum { enableThermalFluxBoundaries = getPropValue<TypeTag, Properties::EnableThermalFluxBoundaries>() };
|
2024-09-11 07:25:34 -05:00
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
|
|
|
|
enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
|
|
|
|
enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
|
2024-09-11 07:25:34 -05:00
|
|
|
|
|
|
|
// TODO: later, gasCompIdx, oilCompIdx and waterCompIdx should go to the FlowProblemBlackoil in the future
|
|
|
|
// we do not want them in the compositional setting
|
2014-04-15 11:30:06 -05:00
|
|
|
enum { gasCompIdx = FluidSystem::gasCompIdx };
|
|
|
|
enum { oilCompIdx = FluidSystem::oilCompIdx };
|
|
|
|
enum { waterCompIdx = FluidSystem::waterCompIdx };
|
|
|
|
|
2020-08-26 03:49:52 -05:00
|
|
|
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
|
|
|
|
using RateVector = GetPropType<TypeTag, Properties::RateVector>;
|
|
|
|
using Simulator = GetPropType<TypeTag, Properties::Simulator>;
|
|
|
|
using Element = typename GridView::template Codim<0>::Entity;
|
|
|
|
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
|
2020-08-27 03:30:29 -05:00
|
|
|
using EclMaterialLawManager = typename GetProp<TypeTag, Properties::MaterialLaw>::EclMaterialLawManager;
|
|
|
|
using EclThermalLawManager = typename GetProp<TypeTag, Properties::SolidEnergyLaw>::EclThermalLawManager;
|
2020-08-26 03:49:52 -05:00
|
|
|
using MaterialLawParams = typename EclMaterialLawManager::MaterialLawParams;
|
|
|
|
using SolidEnergyLawParams = typename EclThermalLawManager::SolidEnergyLawParams;
|
|
|
|
using ThermalConductionLawParams = typename EclThermalLawManager::ThermalConductionLawParams;
|
|
|
|
using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
|
|
|
|
using DofMapper = GetPropType<TypeTag, Properties::DofMapper>;
|
|
|
|
using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
|
|
|
|
using Indices = GetPropType<TypeTag, Properties::Indices>;
|
|
|
|
using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
|
2024-02-06 01:18:32 -06:00
|
|
|
using WellModel = GetPropType<TypeTag, Properties::WellModel>;
|
2024-02-06 01:18:32 -06:00
|
|
|
using AquiferModel = GetPropType<TypeTag, Properties::AquiferModel>;
|
2017-06-01 00:51:44 -05:00
|
|
|
|
2021-05-26 07:24:16 -05:00
|
|
|
using Toolbox = MathToolbox<Evaluation>;
|
|
|
|
using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2017-02-10 06:43:17 -06:00
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
using TracerModel = GetPropType<TypeTag, Properties::TracerModel>;
|
2024-02-05 12:22:15 -06:00
|
|
|
using DirectionalMobilityPtr = Utility::CopyablePtr<DirectionalMobility<TypeTag, Evaluation>>;
|
2018-11-14 06:10:01 -06:00
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
public:
|
2024-02-22 08:17:09 -06:00
|
|
|
using BaseType::briefDescription;
|
|
|
|
using BaseType::helpPreamble;
|
|
|
|
using BaseType::shouldWriteOutput;
|
|
|
|
using BaseType::shouldWriteRestartFile;
|
|
|
|
using BaseType::rockCompressibility;
|
|
|
|
using BaseType::rockReferencePressure;
|
|
|
|
using BaseType::porosity;
|
2021-05-25 08:49:14 -05:00
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::registerParameters
|
|
|
|
*/
|
|
|
|
static void registerParameters()
|
|
|
|
{
|
|
|
|
ParentType::registerParameters();
|
2023-09-09 16:02:09 -05:00
|
|
|
|
2024-08-15 02:26:02 -05:00
|
|
|
registerFlowProblemParameters<Scalar>();
|
2014-04-15 11:30:06 -05:00
|
|
|
}
|
|
|
|
|
2018-08-14 03:26:36 -05:00
|
|
|
|
2018-06-13 09:44:50 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::handlePositionalParameter
|
|
|
|
*/
|
2024-09-04 06:13:55 -05:00
|
|
|
static int handlePositionalParameter(std::function<void(const std::string&,
|
|
|
|
const std::string&)> addKey,
|
|
|
|
std::set<std::string>& seenParams,
|
2018-08-10 10:32:46 -05:00
|
|
|
std::string& errorMsg,
|
2021-05-26 07:18:44 -05:00
|
|
|
int,
|
2018-06-13 09:44:50 -05:00
|
|
|
const char** argv,
|
|
|
|
int paramIdx,
|
2021-05-26 07:18:44 -05:00
|
|
|
int)
|
2018-06-13 09:44:50 -05:00
|
|
|
{
|
2024-10-03 08:34:35 -05:00
|
|
|
return detail::eclPositionalParameter(addKey,
|
|
|
|
seenParams,
|
|
|
|
errorMsg,
|
|
|
|
argv,
|
|
|
|
paramIdx);
|
2018-06-13 09:44:50 -05:00
|
|
|
}
|
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc Doxygen::defaultProblemConstructor
|
|
|
|
*/
|
2024-09-11 07:25:34 -05:00
|
|
|
explicit FlowProblem(Simulator& simulator)
|
2014-04-15 11:30:06 -05:00
|
|
|
: ParentType(simulator)
|
2024-02-22 08:17:09 -06:00
|
|
|
, BaseType(simulator.vanguard().eclState(),
|
|
|
|
simulator.vanguard().schedule(),
|
|
|
|
simulator.vanguard().gridView())
|
2021-05-05 05:13:25 -05:00
|
|
|
, transmissibilities_(simulator.vanguard().eclState(),
|
|
|
|
simulator.vanguard().gridView(),
|
|
|
|
simulator.vanguard().cartesianIndexMapper(),
|
|
|
|
simulator.vanguard().grid(),
|
2021-05-15 06:35:21 -05:00
|
|
|
simulator.vanguard().cellCentroids(),
|
|
|
|
enableEnergy,
|
2023-10-25 12:46:55 -05:00
|
|
|
enableDiffusion,
|
|
|
|
enableDispersion)
|
2018-07-09 05:13:57 -05:00
|
|
|
, wellModel_(simulator)
|
2018-10-23 10:47:38 -05:00
|
|
|
, aquiferModel_(simulator)
|
2016-11-17 12:41:20 -06:00
|
|
|
, pffDofData_(simulator.gridView(), this->elementMapper())
|
2018-11-14 06:10:01 -06:00
|
|
|
, tracerModel_(simulator)
|
2024-09-11 07:25:34 -05:00
|
|
|
{
|
2024-10-15 11:09:25 -05:00
|
|
|
this->enableDriftCompensation_ = Parameters::Get<Parameters::EnableDriftCompensation>();
|
|
|
|
this->enableVtkOutput_ = Parameters::Get<Parameters::EnableVtkOutput>();
|
2024-07-06 03:22:47 -05:00
|
|
|
this->enableTuning_ = Parameters::Get<Parameters::EnableTuning>();
|
2024-10-15 11:09:25 -05:00
|
|
|
|
2024-07-06 03:22:47 -05:00
|
|
|
this->initialTimeStepSize_ = Parameters::Get<Parameters::InitialTimeStepSize<Scalar>>();
|
2024-10-15 11:09:25 -05:00
|
|
|
this->maxTimeStepAfterWellEvent_ = unit::convert::from
|
|
|
|
(Parameters::Get<Parameters::TimeStepAfterEventInDays<Scalar>>(), unit::day);
|
2023-07-27 00:02:37 -05:00
|
|
|
|
2024-10-15 11:09:25 -05:00
|
|
|
// The value N for this parameter is defined in the following order of precedence:
|
|
|
|
//
|
2023-06-27 06:40:24 -05:00
|
|
|
// 1. Command line value (--num-pressure-points-equil=N)
|
2024-10-15 11:09:25 -05:00
|
|
|
//
|
|
|
|
// 2. EQLDIMS item 2. Default value from
|
|
|
|
// opm-common/opm/input/eclipse/share/keywords/000_Eclipse100/E/EQLDIMS
|
2021-01-26 07:59:22 -06:00
|
|
|
|
2024-10-15 11:09:25 -05:00
|
|
|
this->numPressurePointsEquil_ = Parameters::IsSet<Parameters::NumPressurePointsEquil>()
|
|
|
|
? Parameters::Get<Parameters::NumPressurePointsEquil>()
|
|
|
|
: simulator.vanguard().eclState().getTableManager().getEqldims().getNumDepthNodesP();
|
2024-04-26 09:46:27 -05:00
|
|
|
|
2024-10-15 11:09:25 -05:00
|
|
|
this->explicitRockCompaction_ = Parameters::Get<Parameters::ExplicitRockCompaction>();
|
2024-04-26 09:46:27 -05:00
|
|
|
|
2024-10-15 11:27:57 -05:00
|
|
|
if (! Parameters::Get<Parameters::CheckSatfuncConsistency>()) {
|
|
|
|
// User did not enable the "new" saturation function consistency
|
|
|
|
// check module. Run the original checker instead. This is a
|
|
|
|
// temporary measure.
|
|
|
|
RelpermDiagnostics relpermDiagnostics{};
|
|
|
|
relpermDiagnostics.diagnosis(simulator.vanguard().eclState(),
|
2024-10-01 09:22:46 -05:00
|
|
|
simulator.vanguard().levelCartesianIndexMapper());
|
2024-10-15 11:27:57 -05:00
|
|
|
}
|
2014-11-28 05:58:48 -06:00
|
|
|
}
|
2014-07-25 08:31:01 -05:00
|
|
|
|
2024-09-17 00:51:56 -05:00
|
|
|
virtual ~FlowProblem() = default;
|
|
|
|
|
2016-10-30 12:42:06 -05:00
|
|
|
void prefetch(const Element& elem) const
|
2024-10-15 11:09:25 -05:00
|
|
|
{ this->pffDofData_.prefetch(elem); }
|
2016-10-30 12:42:06 -05:00
|
|
|
|
2015-01-19 09:40:26 -06:00
|
|
|
/*!
|
2018-10-23 10:47:38 -05:00
|
|
|
* \brief This method restores the complete state of the problem and its sub-objects
|
2015-01-19 09:40:26 -06:00
|
|
|
* from disk.
|
|
|
|
*
|
2018-10-23 10:47:38 -05:00
|
|
|
* The serialization format used by this method is ad-hoc. It is the inverse of the
|
|
|
|
* serialize() method.
|
2015-01-19 09:40:26 -06:00
|
|
|
*
|
|
|
|
* \tparam Restarter The deserializer type
|
|
|
|
*
|
|
|
|
* \param res The deserializer object
|
|
|
|
*/
|
|
|
|
template <class Restarter>
|
2016-11-07 08:14:07 -06:00
|
|
|
void deserialize(Restarter& res)
|
2015-03-03 07:48:11 -06:00
|
|
|
{
|
|
|
|
// reload the current episode/report step from the deck
|
2024-09-11 07:25:34 -05:00
|
|
|
this->beginEpisode();
|
2015-03-03 07:48:11 -06:00
|
|
|
|
|
|
|
// deserialize the wells
|
2018-07-09 05:13:57 -05:00
|
|
|
wellModel_.deserialize(res);
|
2018-10-23 10:47:38 -05:00
|
|
|
|
2024-02-06 01:55:42 -06:00
|
|
|
// deserialize the aquifer
|
|
|
|
aquiferModel_.deserialize(res);
|
2015-03-03 07:48:11 -06:00
|
|
|
}
|
2015-01-19 09:40:26 -06:00
|
|
|
|
|
|
|
/*!
|
2018-10-23 10:47:38 -05:00
|
|
|
* \brief This method writes the complete state of the problem and its subobjects to
|
|
|
|
* disk.
|
|
|
|
*
|
|
|
|
* The file format used here is ad-hoc.
|
2015-01-19 09:40:26 -06:00
|
|
|
*/
|
|
|
|
template <class Restarter>
|
2016-11-07 08:14:07 -06:00
|
|
|
void serialize(Restarter& res)
|
2018-10-23 10:47:38 -05:00
|
|
|
{
|
|
|
|
wellModel_.serialize(res);
|
2019-06-03 04:11:41 -05:00
|
|
|
|
2024-02-06 01:55:42 -06:00
|
|
|
aquiferModel_.serialize(res);
|
2018-10-23 10:47:38 -05:00
|
|
|
}
|
2015-01-19 09:40:26 -06:00
|
|
|
|
2021-05-25 08:49:14 -05:00
|
|
|
int episodeIndex() const
|
|
|
|
{
|
|
|
|
return std::max(this->simulator().episodeIndex(), 0);
|
|
|
|
}
|
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
/*!
|
2014-07-09 05:01:07 -05:00
|
|
|
* \brief Called by the simulator before an episode begins.
|
2014-04-15 11:30:06 -05:00
|
|
|
*/
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void beginEpisode()
|
2015-03-03 07:48:11 -06:00
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(beginEpisode);
|
2015-03-03 07:48:11 -06:00
|
|
|
// Proceed to the next report step
|
2019-04-03 10:26:57 -05:00
|
|
|
auto& simulator = this->simulator();
|
2019-10-07 07:51:23 -05:00
|
|
|
int episodeIdx = simulator.episodeIndex();
|
2019-04-03 10:26:57 -05:00
|
|
|
auto& eclState = simulator.vanguard().eclState();
|
|
|
|
const auto& schedule = simulator.vanguard().schedule();
|
2021-01-10 14:46:45 -06:00
|
|
|
const auto& events = schedule[episodeIdx].events();
|
2016-10-17 07:02:09 -05:00
|
|
|
|
2021-05-05 04:22:44 -05:00
|
|
|
if (episodeIdx >= 0 && events.hasEvent(ScheduleEvents::GEO_MODIFIER)) {
|
2016-10-17 07:02:09 -05:00
|
|
|
// bring the contents of the keywords to the current state of the SCHEDULE
|
2019-04-03 10:26:57 -05:00
|
|
|
// section.
|
2016-10-17 07:02:09 -05:00
|
|
|
//
|
|
|
|
// TODO (?): make grid topology changes possible (depending on what exactly
|
|
|
|
// has changed, the grid may need be re-created which has some serious
|
|
|
|
// implications on e.g., the solution of the simulation.)
|
2021-01-10 14:46:45 -06:00
|
|
|
const auto& miniDeck = schedule[episodeIdx].geo_keywords();
|
2021-11-08 05:02:51 -06:00
|
|
|
const auto& cc = simulator.vanguard().grid().comm();
|
2021-11-07 12:58:20 -06:00
|
|
|
eclState.apply_schedule_keywords( miniDeck );
|
2021-11-08 05:02:51 -06:00
|
|
|
eclBroadcast(cc, eclState.getTransMult() );
|
2016-10-17 07:02:09 -05:00
|
|
|
|
2021-12-01 07:00:21 -06:00
|
|
|
// Re-ordering in case of ALUGrid
|
2023-08-01 06:41:51 -05:00
|
|
|
std::function<unsigned int(unsigned int)> equilGridToGrid = [&simulator](unsigned int i) {
|
2022-11-14 08:10:03 -06:00
|
|
|
return simulator.vanguard().gridEquilIdxToGridIdx(i);
|
|
|
|
};
|
2021-12-01 07:00:21 -06:00
|
|
|
|
2016-10-17 07:02:09 -05:00
|
|
|
// re-compute all quantities which may possibly be affected.
|
2024-06-06 08:15:56 -05:00
|
|
|
using TransUpdateQuantities = typename Vanguard::TransmissibilityType::TransUpdateQuantities;
|
|
|
|
transmissibilities_.update(true, TransUpdateQuantities::All, equilGridToGrid);
|
2021-05-25 08:49:14 -05:00
|
|
|
this->referencePorosity_[1] = this->referencePorosity_[0];
|
2019-03-12 09:51:41 -05:00
|
|
|
updateReferencePorosity_();
|
2016-12-08 04:11:56 -06:00
|
|
|
updatePffDofData_();
|
2022-10-27 02:54:40 -05:00
|
|
|
this->model().linearizer().updateDiscretizationParameters();
|
2016-10-17 07:02:09 -05:00
|
|
|
}
|
2015-03-03 07:48:11 -06:00
|
|
|
|
2024-08-21 02:43:19 -05:00
|
|
|
bool tuningEvent = this->beginEpisode_(enableExperiments, this->episodeIndex());
|
2019-03-11 04:18:27 -05:00
|
|
|
|
2019-02-15 06:39:58 -06:00
|
|
|
// set up the wells for the next episode.
|
2019-04-03 10:26:57 -05:00
|
|
|
wellModel_.beginEpisode();
|
2017-12-19 11:00:14 -06:00
|
|
|
|
2019-04-03 10:26:57 -05:00
|
|
|
// set up the aquifers for the next episode.
|
2024-02-06 01:55:42 -06:00
|
|
|
aquiferModel_.beginEpisode();
|
2018-10-23 10:47:38 -05:00
|
|
|
|
2019-04-03 10:26:57 -05:00
|
|
|
// set the size of the initial time step of the episode
|
|
|
|
Scalar dt = limitNextTimeStepSize_(simulator.episodeLength());
|
2023-08-23 09:06:46 -05:00
|
|
|
// negative value of initialTimeStepSize_ indicates no active limit from TSINIT or NEXTSTEP
|
|
|
|
if ( (episodeIdx == 0 || tuningEvent) && this->initialTimeStepSize_ > 0)
|
2019-04-03 10:26:57 -05:00
|
|
|
// allow the size of the initial time step to be set via an external parameter
|
|
|
|
// if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
|
2021-05-25 08:49:14 -05:00
|
|
|
dt = std::min(dt, this->initialTimeStepSize_);
|
2019-04-03 10:26:57 -05:00
|
|
|
simulator.setTimeStepSize(dt);
|
2015-03-03 07:48:11 -06:00
|
|
|
}
|
2014-07-09 05:01:07 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Called by the simulator before each time integration.
|
|
|
|
*/
|
2014-05-07 08:07:39 -05:00
|
|
|
void beginTimeStep()
|
2015-09-10 08:35:17 -05:00
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(beginTimeStep);
|
2024-08-16 03:40:49 -05:00
|
|
|
const int episodeIdx = this->episodeIndex();
|
|
|
|
const int timeStepSize = this->simulator().timeStepSize();
|
2019-10-07 07:51:23 -05:00
|
|
|
|
2021-05-25 08:49:14 -05:00
|
|
|
this->beginTimeStep_(enableExperiments,
|
|
|
|
episodeIdx,
|
|
|
|
this->simulator().timeStepIndex(),
|
|
|
|
this->simulator().startTime(),
|
|
|
|
this->simulator().time(),
|
2024-08-16 03:40:49 -05:00
|
|
|
timeStepSize,
|
2021-05-25 08:49:14 -05:00
|
|
|
this->simulator().endTime());
|
2019-03-12 09:51:41 -05:00
|
|
|
|
2019-10-08 08:49:48 -05:00
|
|
|
// update maximum water saturation and minimum pressure
|
|
|
|
// used when ROCKCOMP is activated
|
2024-08-28 09:28:51 -05:00
|
|
|
// Do not update max RS first step after a restart
|
2024-09-11 07:25:34 -05:00
|
|
|
this->updateExplicitQuantities_(episodeIdx, timeStepSize, first_step_ && (episodeIdx > 0));
|
2024-08-28 09:28:51 -05:00
|
|
|
first_step_ = false;
|
2020-09-28 02:42:10 -05:00
|
|
|
|
2023-04-21 08:17:58 -05:00
|
|
|
if (nonTrivialBoundaryConditions()) {
|
|
|
|
this->model().linearizer().updateBoundaryConditionData();
|
|
|
|
}
|
|
|
|
|
2019-02-15 06:39:58 -06:00
|
|
|
wellModel_.beginTimeStep();
|
2024-02-06 01:55:42 -06:00
|
|
|
aquiferModel_.beginTimeStep();
|
2018-11-14 06:10:01 -06:00
|
|
|
tracerModel_.beginTimeStep();
|
|
|
|
|
2015-09-10 08:35:17 -05:00
|
|
|
}
|
2015-07-28 10:24:44 -05:00
|
|
|
|
2014-05-07 08:07:39 -05:00
|
|
|
/*!
|
|
|
|
* \brief Called by the simulator before each Newton-Raphson iteration.
|
|
|
|
*/
|
|
|
|
void beginIteration()
|
2016-06-07 10:43:45 -05:00
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(beginIteration);
|
2019-02-15 06:39:58 -06:00
|
|
|
wellModel_.beginIteration();
|
2024-02-06 01:55:42 -06:00
|
|
|
aquiferModel_.beginIteration();
|
2016-06-07 10:43:45 -05:00
|
|
|
}
|
2014-05-07 08:07:39 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Called by the simulator after each Newton-Raphson iteration.
|
|
|
|
*/
|
|
|
|
void endIteration()
|
2016-06-07 10:43:45 -05:00
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(endIteration);
|
2019-02-15 06:39:58 -06:00
|
|
|
wellModel_.endIteration();
|
2024-02-06 01:55:42 -06:00
|
|
|
aquiferModel_.endIteration();
|
2016-06-07 10:43:45 -05:00
|
|
|
}
|
2014-05-07 08:07:39 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Called by the simulator after each time integration.
|
|
|
|
*/
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void endTimeStep()
|
2014-07-22 05:41:56 -05:00
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(endTimeStep);
|
2024-05-03 03:29:40 -05:00
|
|
|
|
2014-05-07 08:07:39 -05:00
|
|
|
#ifndef NDEBUG
|
2024-08-14 07:05:43 -05:00
|
|
|
if constexpr (getPropValue<TypeTag, Properties::EnableDebuggingChecks>()) {
|
2024-05-03 03:29:40 -05:00
|
|
|
// in debug mode, we don't care about performance, so we check
|
|
|
|
// if the model does the right thing (i.e., the mass change
|
|
|
|
// inside the whole reservoir must be equivalent to the fluxes
|
|
|
|
// over the grid's boundaries plus the source rates specified by
|
|
|
|
// the problem).
|
|
|
|
const int rank = this->simulator().gridView().comm().rank();
|
|
|
|
if (rank == 0) {
|
2019-06-03 04:11:41 -05:00
|
|
|
std::cout << "checking conservativeness of solution\n";
|
2024-05-03 03:29:40 -05:00
|
|
|
}
|
|
|
|
|
2016-06-07 11:21:14 -05:00
|
|
|
this->model().checkConservativeness(/*tolerance=*/-1, /*verbose=*/true);
|
2024-05-03 03:29:40 -05:00
|
|
|
if (rank == 0) {
|
2019-06-03 04:11:41 -05:00
|
|
|
std::cout << "solution is sufficiently conservative\n";
|
2024-05-03 03:29:40 -05:00
|
|
|
}
|
2016-06-07 11:21:14 -05:00
|
|
|
}
|
2014-07-22 05:41:56 -05:00
|
|
|
#endif // NDEBUG
|
2016-02-11 09:25:06 -06:00
|
|
|
|
2019-12-06 06:12:36 -06:00
|
|
|
auto& simulator = this->simulator();
|
2024-06-14 03:35:41 -05:00
|
|
|
simulator.setTimeStepIndex(simulator.timeStepIndex()+1);
|
2018-11-14 06:10:01 -06:00
|
|
|
|
2024-05-03 03:29:40 -05:00
|
|
|
this->wellModel_.endTimeStep();
|
|
|
|
this->aquiferModel_.endTimeStep();
|
|
|
|
this->tracerModel_.endTimeStep();
|
2023-12-05 06:08:23 -06:00
|
|
|
|
|
|
|
// Compute flux for output
|
|
|
|
this->model().linearizer().updateFlowsInfo();
|
|
|
|
|
2024-05-03 03:29:40 -05:00
|
|
|
if (this->enableDriftCompensation_) {
|
|
|
|
OPM_TIMEBLOCK(driftCompansation);
|
|
|
|
|
2019-03-05 04:32:11 -06:00
|
|
|
const auto& residual = this->model().linearizer().residual();
|
2024-05-03 03:29:40 -05:00
|
|
|
|
2019-03-05 04:32:11 -06:00
|
|
|
for (unsigned globalDofIdx = 0; globalDofIdx < residual.size(); globalDofIdx ++) {
|
2024-09-25 03:09:39 -05:00
|
|
|
int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(globalDofIdx);
|
|
|
|
this->drift_[sfcdofIdx] = residual[sfcdofIdx] * simulator.timeStepSize();
|
2024-05-03 03:29:40 -05:00
|
|
|
|
|
|
|
if constexpr (getPropValue<TypeTag, Properties::UseVolumetricResidual>()) {
|
2024-09-25 03:09:39 -05:00
|
|
|
this->drift_[sfcdofIdx] *= this->model().dofTotalVolume(sfcdofIdx);
|
2024-05-03 03:29:40 -05:00
|
|
|
}
|
2019-03-05 04:32:11 -06:00
|
|
|
}
|
|
|
|
}
|
2014-07-22 05:41:56 -05:00
|
|
|
}
|
2014-07-09 05:01:07 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Called by the simulator after the end of an episode.
|
|
|
|
*/
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void endEpisode()
|
2014-04-15 11:30:06 -05:00
|
|
|
{
|
2024-05-23 07:53:03 -05:00
|
|
|
const int episodeIdx = this->episodeIndex();
|
|
|
|
|
|
|
|
this->wellModel_.endEpisode();
|
|
|
|
this->aquiferModel_.endEpisode();
|
|
|
|
|
|
|
|
const auto& schedule = this->simulator().vanguard().schedule();
|
|
|
|
|
|
|
|
// End simulation when completed.
|
|
|
|
if (episodeIdx + 1 >= static_cast<int>(schedule.size()) - 1) {
|
|
|
|
this->simulator().setFinished(true);
|
2015-04-01 07:33:54 -05:00
|
|
|
return;
|
|
|
|
}
|
2019-04-03 10:26:57 -05:00
|
|
|
|
2024-05-23 07:53:03 -05:00
|
|
|
// Otherwise, start next episode (report step).
|
|
|
|
this->simulator().startNextEpisode(schedule.stepLength(episodeIdx + 1));
|
2014-04-15 11:30:06 -05:00
|
|
|
}
|
|
|
|
|
2014-11-28 05:58:48 -06:00
|
|
|
/*!
|
|
|
|
* \brief Write the requested quantities of the current solution into the output
|
|
|
|
* files.
|
|
|
|
*/
|
2024-10-01 16:31:24 -05:00
|
|
|
virtual void writeOutput(bool verbose)
|
2014-11-28 05:58:48 -06:00
|
|
|
{
|
2023-02-15 04:05:45 -06:00
|
|
|
OPM_TIMEBLOCK(problemWriteOutput);
|
2018-07-12 06:47:00 -05:00
|
|
|
// use the generic code to prepare the output fields and to
|
|
|
|
// write the desired VTK files.
|
2024-07-06 03:22:47 -05:00
|
|
|
if (Parameters::Get<Parameters::EnableWriteAllSolutions>() ||
|
2024-04-05 05:53:20 -05:00
|
|
|
this->simulator().episodeWillBeOver()) {
|
2023-10-09 02:04:22 -05:00
|
|
|
ParentType::writeOutput(verbose);
|
|
|
|
}
|
2017-12-07 05:38:48 -06:00
|
|
|
}
|
2014-11-28 05:58:48 -06:00
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseMultiPhaseProblem::intrinsicPermeability
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
const DimMatrix& intrinsicPermeability(const Context& context,
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned spaceIdx,
|
|
|
|
unsigned timeIdx) const
|
2014-04-15 11:30:06 -05:00
|
|
|
{
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
2016-12-05 12:49:58 -06:00
|
|
|
return transmissibilities_.permeability(globalSpaceIdx);
|
2014-04-15 11:30:06 -05:00
|
|
|
}
|
|
|
|
|
2014-07-02 10:50:35 -05:00
|
|
|
/*!
|
2014-12-27 08:19:15 -06:00
|
|
|
* \brief This method returns the intrinsic permeability tensor
|
|
|
|
* given a global element index.
|
|
|
|
*
|
|
|
|
* Its main (only?) usage is the ECL transmissibility calculation code...
|
2014-07-02 10:50:35 -05:00
|
|
|
*/
|
2016-11-07 08:14:07 -06:00
|
|
|
const DimMatrix& intrinsicPermeability(unsigned globalElemIdx) const
|
2016-12-05 12:49:58 -06:00
|
|
|
{ return transmissibilities_.permeability(globalElemIdx); }
|
2014-07-02 10:50:35 -05:00
|
|
|
|
2014-12-27 08:19:15 -06:00
|
|
|
/*!
|
2018-01-30 04:46:23 -06:00
|
|
|
* \copydoc EclTransmissiblity::transmissibility
|
2014-12-27 08:19:15 -06:00
|
|
|
*/
|
2016-10-30 12:41:58 -05:00
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
Scalar transmissibility(const Context& context,
|
2021-05-26 07:19:23 -05:00
|
|
|
[[maybe_unused]] unsigned fromDofLocalIdx,
|
2016-10-30 12:41:58 -05:00
|
|
|
unsigned toDofLocalIdx) const
|
|
|
|
{
|
2016-10-30 12:42:02 -05:00
|
|
|
assert(fromDofLocalIdx == 0);
|
|
|
|
return pffDofData_.get(context.element(), toDofLocalIdx).transmissibility;
|
2016-10-30 12:41:58 -05:00
|
|
|
}
|
2014-07-02 10:50:35 -05:00
|
|
|
|
2022-07-05 02:48:38 -05:00
|
|
|
/*!
|
|
|
|
* \brief Direct access to the transmissibility between two elements.
|
|
|
|
*/
|
|
|
|
Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
|
|
|
|
{
|
|
|
|
return transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
|
|
|
|
}
|
|
|
|
|
2021-01-14 08:33:28 -06:00
|
|
|
/*!
|
|
|
|
* \copydoc EclTransmissiblity::diffusivity
|
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
Scalar diffusivity(const Context& context,
|
2021-05-26 07:19:23 -05:00
|
|
|
[[maybe_unused]] unsigned fromDofLocalIdx,
|
2021-01-14 08:33:28 -06:00
|
|
|
unsigned toDofLocalIdx) const
|
|
|
|
{
|
|
|
|
assert(fromDofLocalIdx == 0);
|
|
|
|
return *pffDofData_.get(context.element(), toDofLocalIdx).diffusivity;
|
|
|
|
}
|
|
|
|
|
2023-06-30 05:54:18 -05:00
|
|
|
/*!
|
|
|
|
* give the transmissibility for a face i.e. pair. should be symmetric?
|
|
|
|
*/
|
|
|
|
Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
|
|
|
|
return transmissibilities_.diffusivity(globalCellIn, globalCellOut);
|
|
|
|
}
|
|
|
|
|
2023-10-25 12:46:55 -05:00
|
|
|
/*!
|
|
|
|
* give the dispersivity for a face i.e. pair.
|
|
|
|
*/
|
|
|
|
Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
|
|
|
|
return transmissibilities_.dispersivity(globalCellIn, globalCellOut);
|
|
|
|
}
|
|
|
|
|
2023-06-30 05:54:18 -05:00
|
|
|
/*!
|
|
|
|
* \brief Direct access to a boundary transmissibility.
|
|
|
|
*/
|
|
|
|
Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx,
|
|
|
|
const unsigned boundaryFaceIdx) const
|
|
|
|
{
|
|
|
|
return transmissibilities_.thermalTransmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2018-01-30 04:46:23 -06:00
|
|
|
/*!
|
|
|
|
* \copydoc EclTransmissiblity::transmissibilityBoundary
|
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
Scalar transmissibilityBoundary(const Context& elemCtx,
|
|
|
|
unsigned boundaryFaceIdx) const
|
|
|
|
{
|
|
|
|
unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
|
|
|
|
return transmissibilities_.transmissibilityBoundary(elemIdx, boundaryFaceIdx);
|
|
|
|
}
|
|
|
|
|
2022-09-14 03:33:47 -05:00
|
|
|
/*!
|
|
|
|
* \brief Direct access to a boundary transmissibility.
|
|
|
|
*/
|
|
|
|
Scalar transmissibilityBoundary(const unsigned globalSpaceIdx,
|
|
|
|
const unsigned boundaryFaceIdx) const
|
|
|
|
{
|
|
|
|
return transmissibilities_.transmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
|
|
|
|
}
|
|
|
|
|
2023-06-30 05:54:18 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \copydoc EclTransmissiblity::thermalHalfTransmissibility
|
|
|
|
*/
|
|
|
|
Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn,
|
|
|
|
const unsigned globalSpaceIdxOut) const
|
|
|
|
{
|
|
|
|
return transmissibilities_.thermalHalfTrans(globalSpaceIdxIn,globalSpaceIdxOut);
|
|
|
|
}
|
|
|
|
|
2018-01-30 04:46:23 -06:00
|
|
|
/*!
|
|
|
|
* \copydoc EclTransmissiblity::thermalHalfTransmissibility
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2021-03-26 09:27:14 -05:00
|
|
|
Scalar thermalHalfTransmissibilityIn(const Context& context,
|
|
|
|
unsigned faceIdx,
|
|
|
|
unsigned timeIdx) const
|
2018-01-30 04:46:23 -06:00
|
|
|
{
|
2018-03-26 03:53:43 -05:00
|
|
|
const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
|
|
|
|
unsigned toDofLocalIdx = face.exteriorIndex();
|
2021-03-26 09:27:14 -05:00
|
|
|
return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransIn;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \copydoc EclTransmissiblity::thermalHalfTransmissibility
|
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
Scalar thermalHalfTransmissibilityOut(const Context& context,
|
|
|
|
unsigned faceIdx,
|
|
|
|
unsigned timeIdx) const
|
|
|
|
{
|
|
|
|
const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
|
|
|
|
unsigned toDofLocalIdx = face.exteriorIndex();
|
|
|
|
return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransOut;
|
2018-01-30 04:46:23 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \copydoc EclTransmissiblity::thermalHalfTransmissibility
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2018-01-30 04:46:23 -06:00
|
|
|
Scalar thermalHalfTransmissibilityBoundary(const Context& elemCtx,
|
2018-01-30 04:46:23 -06:00
|
|
|
unsigned boundaryFaceIdx) const
|
|
|
|
{
|
2018-01-30 04:46:23 -06:00
|
|
|
unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
|
2018-01-30 04:46:23 -06:00
|
|
|
return transmissibilities_.thermalHalfTransBoundary(elemIdx, boundaryFaceIdx);
|
|
|
|
}
|
|
|
|
|
2016-12-16 11:47:52 -06:00
|
|
|
/*!
|
|
|
|
* \brief Return a reference to the object that handles the "raw" transmissibilities.
|
|
|
|
*/
|
2021-05-05 05:13:25 -05:00
|
|
|
const typename Vanguard::TransmissibilityType& eclTransmissibilities() const
|
2016-12-16 11:47:52 -06:00
|
|
|
{ return transmissibilities_; }
|
|
|
|
|
2018-04-28 06:31:27 -05:00
|
|
|
|
2024-02-05 12:22:15 -06:00
|
|
|
const TracerModel& tracerModel() const
|
2018-11-14 06:10:01 -06:00
|
|
|
{ return tracerModel_; }
|
|
|
|
|
2024-02-05 12:22:15 -06:00
|
|
|
TracerModel& tracerModel()
|
2021-11-23 05:57:20 -06:00
|
|
|
{ return tracerModel_; }
|
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseMultiPhaseProblem::porosity
|
2019-03-12 09:51:41 -05:00
|
|
|
*
|
2024-02-05 12:22:15 -06:00
|
|
|
* For the FlowProblem, this method is identical to referencePorosity(). The intensive
|
2019-03-12 09:51:41 -05:00
|
|
|
* quantities object may apply various multipliers (e.g. ones which model rock
|
|
|
|
* compressibility and water induced rock compaction) to it which depend on the
|
|
|
|
* current physical conditions.
|
2014-04-15 11:30:06 -05:00
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
Scalar porosity(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2014-04-15 11:30:06 -05:00
|
|
|
{
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
2022-07-05 02:48:38 -05:00
|
|
|
return this->porosity(globalSpaceIdx, timeIdx);
|
|
|
|
}
|
|
|
|
|
2016-10-25 08:54:14 -05:00
|
|
|
/*!
|
|
|
|
* \brief Returns the depth of an degree of freedom [m]
|
|
|
|
*
|
|
|
|
* For ECL problems this is defined as the average of the depth of an element and is
|
|
|
|
* thus slightly different from the depth of an element's centroid.
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
Scalar dofCenterDepth(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2016-10-25 08:54:14 -05:00
|
|
|
{
|
|
|
|
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
2022-07-05 02:48:38 -05:00
|
|
|
return this->dofCenterDepth(globalSpaceIdx);
|
2021-12-01 07:00:21 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
2022-07-05 02:48:38 -05:00
|
|
|
* \brief Direct indexed acces to the depth of an degree of freedom [m]
|
2021-12-01 07:00:21 -06:00
|
|
|
*
|
|
|
|
* For ECL problems this is defined as the average of the depth of an element and is
|
|
|
|
* thus slightly different from the depth of an element's centroid.
|
|
|
|
*/
|
|
|
|
Scalar dofCenterDepth(unsigned globalSpaceIdx) const
|
|
|
|
{
|
2020-12-08 09:09:01 -06:00
|
|
|
return this->simulator().vanguard().cellCenterDepth(globalSpaceIdx);
|
2016-10-25 08:54:14 -05:00
|
|
|
}
|
|
|
|
|
2014-12-04 12:22:56 -06:00
|
|
|
/*!
|
|
|
|
* \copydoc BlackoilProblem::rockCompressibility
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
Scalar rockCompressibility(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2014-12-04 12:22:56 -06:00
|
|
|
{
|
2022-08-15 03:16:11 -05:00
|
|
|
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
|
|
|
return this->rockCompressibility(globalSpaceIdx);
|
2022-07-05 02:48:38 -05:00
|
|
|
}
|
|
|
|
|
2014-12-04 12:22:56 -06:00
|
|
|
/*!
|
|
|
|
* \copydoc BlackoilProblem::rockReferencePressure
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
Scalar rockReferencePressure(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2014-12-04 12:22:56 -06:00
|
|
|
{
|
2022-08-15 03:16:11 -05:00
|
|
|
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
|
|
|
return this->rockReferencePressure(globalSpaceIdx);
|
2022-07-05 02:48:38 -05:00
|
|
|
}
|
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseMultiPhaseProblem::materialLawParams
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
const MaterialLawParams& materialLawParams(const Context& context,
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned spaceIdx, unsigned timeIdx) const
|
2014-04-15 11:30:06 -05:00
|
|
|
{
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
2022-07-05 02:48:38 -05:00
|
|
|
return this->materialLawParams(globalSpaceIdx);
|
2014-04-15 11:30:06 -05:00
|
|
|
}
|
|
|
|
|
2015-11-18 04:54:35 -06:00
|
|
|
const MaterialLawParams& materialLawParams(unsigned globalDofIdx) const
|
2022-07-05 02:48:38 -05:00
|
|
|
{
|
|
|
|
return materialLawManager_->materialLawParams(globalDofIdx);
|
|
|
|
}
|
2015-07-29 06:43:17 -05:00
|
|
|
|
2023-01-08 05:33:16 -06:00
|
|
|
const MaterialLawParams& materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
|
|
|
|
{
|
|
|
|
return materialLawManager_->materialLawParams(globalDofIdx, facedir);
|
|
|
|
}
|
|
|
|
|
2018-01-30 04:46:23 -06:00
|
|
|
/*!
|
2018-01-30 04:46:23 -06:00
|
|
|
* \brief Return the parameters for the energy storage law of the rock
|
2018-01-30 04:46:23 -06:00
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
const SolidEnergyLawParams&
|
2021-05-26 07:18:44 -05:00
|
|
|
solidEnergyLawParams(const Context& context,
|
|
|
|
unsigned spaceIdx,
|
|
|
|
unsigned timeIdx) const
|
2018-01-30 04:46:23 -06:00
|
|
|
{
|
|
|
|
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
|
|
|
return thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseMultiPhaseProblem::thermalConductionParams
|
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
const ThermalConductionLawParams &
|
|
|
|
thermalConductionLawParams(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
|
|
|
{
|
|
|
|
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
|
|
|
return thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
|
|
|
|
}
|
|
|
|
|
2016-11-30 07:24:15 -06:00
|
|
|
/*!
|
|
|
|
* \brief Returns the ECL material law manager
|
|
|
|
*
|
|
|
|
* Note that this method is *not* part of the generic eWoms problem API because it
|
|
|
|
* would force all problens use the ECL material laws.
|
|
|
|
*/
|
2016-12-10 03:53:28 -06:00
|
|
|
std::shared_ptr<const EclMaterialLawManager> materialLawManager() const
|
2016-11-22 07:07:14 -06:00
|
|
|
{ return materialLawManager_; }
|
2023-03-24 08:56:23 -05:00
|
|
|
|
2022-09-26 10:13:03 -05:00
|
|
|
template <class FluidState>
|
|
|
|
void updateRelperms(
|
|
|
|
std::array<Evaluation,numPhases> &mobility,
|
|
|
|
DirectionalMobilityPtr &dirMob,
|
|
|
|
FluidState &fluidState,
|
|
|
|
unsigned globalSpaceIdx) const
|
2022-09-22 16:54:29 -05:00
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK_LOCAL(updateRelperms);
|
2023-01-08 05:33:16 -06:00
|
|
|
{
|
|
|
|
// calculate relative permeabilities. note that we store the result into the
|
|
|
|
// mobility_ class attribute. the division by the phase viscosity happens later.
|
|
|
|
const auto& materialParams = materialLawParams(globalSpaceIdx);
|
|
|
|
MaterialLaw::relativePermeabilities(mobility, materialParams, fluidState);
|
|
|
|
Valgrind::CheckDefined(mobility);
|
|
|
|
}
|
|
|
|
if (materialLawManager_->hasDirectionalRelperms()
|
|
|
|
|| materialLawManager_->hasDirectionalImbnum())
|
|
|
|
{
|
2022-09-26 10:13:03 -05:00
|
|
|
using Dir = FaceDir::DirEnum;
|
|
|
|
constexpr int ndim = 3;
|
|
|
|
dirMob = std::make_unique<DirectionalMobility<TypeTag, Evaluation>>();
|
|
|
|
Dir facedirs[ndim] = {Dir::XPlus, Dir::YPlus, Dir::ZPlus};
|
|
|
|
for (int i = 0; i<ndim; i++) {
|
2023-01-08 05:33:16 -06:00
|
|
|
const auto& materialParams = materialLawParams(globalSpaceIdx, facedirs[i]);
|
2022-09-26 10:13:03 -05:00
|
|
|
auto& mob_array = dirMob->getArray(i);
|
2023-01-08 05:33:16 -06:00
|
|
|
MaterialLaw::relativePermeabilities(mob_array, materialParams, fluidState);
|
2022-09-26 10:13:03 -05:00
|
|
|
}
|
|
|
|
}
|
2022-09-22 16:54:29 -05:00
|
|
|
}
|
|
|
|
|
2016-11-30 07:24:15 -06:00
|
|
|
/*!
|
|
|
|
* \copydoc materialLawManager()
|
|
|
|
*/
|
2016-11-22 07:07:14 -06:00
|
|
|
std::shared_ptr<EclMaterialLawManager> materialLawManager()
|
|
|
|
{ return materialLawManager_; }
|
|
|
|
|
2024-02-22 08:17:09 -06:00
|
|
|
using BaseType::pvtRegionIndex;
|
2014-08-05 09:52:52 -05:00
|
|
|
/*!
|
|
|
|
* \brief Returns the index of the relevant region for thermodynmic properties
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
unsigned pvtRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2015-08-26 09:04:10 -05:00
|
|
|
{ return pvtRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
|
|
|
|
|
2024-02-22 08:17:09 -06:00
|
|
|
using BaseType::satnumRegionIndex;
|
2017-05-05 03:46:53 -05:00
|
|
|
/*!
|
|
|
|
* \brief Returns the index of the relevant region for thermodynmic properties
|
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
unsigned satnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2021-05-25 08:49:14 -05:00
|
|
|
{ return this->satnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
|
2017-05-05 03:46:53 -05:00
|
|
|
|
2024-02-22 08:17:09 -06:00
|
|
|
using BaseType::miscnumRegionIndex;
|
2017-05-15 05:21:31 -05:00
|
|
|
/*!
|
|
|
|
* \brief Returns the index of the relevant region for thermodynmic properties
|
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
unsigned miscnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2021-05-25 08:49:14 -05:00
|
|
|
{ return this->miscnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
|
2017-05-15 05:21:31 -05:00
|
|
|
|
2024-02-22 08:17:09 -06:00
|
|
|
using BaseType::plmixnumRegionIndex;
|
2017-06-21 03:27:25 -05:00
|
|
|
/*!
|
|
|
|
* \brief Returns the index of the relevant region for thermodynmic properties
|
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
unsigned plmixnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2021-05-25 08:49:14 -05:00
|
|
|
{ return this->plmixnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
|
2017-06-21 03:27:25 -05:00
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
// TODO: polymer related might need to go to the blackoil side
|
2024-02-22 08:17:09 -06:00
|
|
|
using BaseType::maxPolymerAdsorption;
|
2017-06-01 00:51:44 -05:00
|
|
|
/*!
|
|
|
|
* \brief Returns the max polymer adsorption value
|
|
|
|
*/
|
|
|
|
template <class Context>
|
|
|
|
Scalar maxPolymerAdsorption(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2021-05-25 08:49:14 -05:00
|
|
|
{ return this->maxPolymerAdsorption(context.globalSpaceIndex(spaceIdx, timeIdx)); }
|
2017-06-01 00:51:44 -05:00
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::name
|
|
|
|
*/
|
2014-04-25 10:22:28 -05:00
|
|
|
std::string name() const
|
2018-02-01 09:26:58 -06:00
|
|
|
{ return this->simulator().vanguard().caseName(); }
|
2014-04-15 11:30:06 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseMultiPhaseProblem::temperature
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
2015-02-04 04:19:05 -06:00
|
|
|
{
|
2019-03-21 06:29:42 -05:00
|
|
|
// use the initial temperature of the DOF if temperature is not a primary
|
|
|
|
// variable
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
2024-09-11 07:25:34 -05:00
|
|
|
return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
|
2015-02-04 04:19:05 -06:00
|
|
|
}
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2023-06-30 05:54:18 -05:00
|
|
|
|
2023-09-04 07:35:20 -05:00
|
|
|
Scalar temperature(unsigned globalDofIdx, unsigned /*timeIdx*/) const
|
2023-06-30 05:54:18 -05:00
|
|
|
{
|
|
|
|
// use the initial temperature of the DOF if temperature is not a primary
|
|
|
|
// variable
|
2024-09-11 07:25:34 -05:00
|
|
|
return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
|
2023-06-30 05:54:18 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
const SolidEnergyLawParams&
|
|
|
|
solidEnergyLawParams(unsigned globalSpaceIdx,
|
|
|
|
unsigned /*timeIdx*/) const
|
|
|
|
{
|
|
|
|
return this->thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
|
|
|
|
}
|
|
|
|
const ThermalConductionLawParams &
|
|
|
|
thermalConductionLawParams(unsigned globalSpaceIdx,
|
|
|
|
unsigned /*timeIdx*/)const
|
|
|
|
{
|
|
|
|
return this->thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
|
|
|
|
}
|
|
|
|
|
2021-05-25 08:49:14 -05:00
|
|
|
/*!
|
|
|
|
* \brief Returns an element's historic maximum oil phase saturation that was
|
|
|
|
* observed during the simulation.
|
|
|
|
*
|
|
|
|
* In this context, "historic" means the the time before the current timestep began.
|
|
|
|
*
|
|
|
|
* This is a bit of a hack from the conceptional point of view, but it is required to
|
|
|
|
* match the results of the 'flow' and ECLIPSE 100 simulators.
|
|
|
|
*/
|
|
|
|
Scalar maxOilSaturation(unsigned globalDofIdx) const
|
|
|
|
{
|
|
|
|
if (!this->vapparsActive(this->episodeIndex()))
|
|
|
|
return 0.0;
|
|
|
|
|
|
|
|
return this->maxOilSaturation_[globalDofIdx];
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Sets an element's maximum oil phase saturation observed during the
|
|
|
|
* simulation.
|
|
|
|
*
|
|
|
|
* In this context, "historic" means the the time before the current timestep began.
|
|
|
|
*
|
|
|
|
* This a hack on top of the maxOilSaturation() hack but it is currently required to
|
|
|
|
* do restart externally. i.e. from the flow code.
|
|
|
|
*/
|
|
|
|
void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
|
|
|
|
{
|
|
|
|
if (!this->vapparsActive(this->episodeIndex()))
|
|
|
|
return;
|
|
|
|
|
|
|
|
this->maxOilSaturation_[globalDofIdx] = value;
|
|
|
|
}
|
|
|
|
|
2016-02-01 10:25:23 -06:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::initialSolutionApplied()
|
|
|
|
*/
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void initialSolutionApplied()
|
2015-07-28 10:24:44 -05:00
|
|
|
{
|
2023-05-22 10:40:27 -05:00
|
|
|
// Calculate all intensive quantities.
|
2023-03-23 16:37:44 -05:00
|
|
|
this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/0);
|
2023-05-22 10:40:27 -05:00
|
|
|
|
|
|
|
// We also need the intensive quantities for timeIdx == 1
|
|
|
|
// corresponding to the start of the current timestep, if we
|
|
|
|
// do not use the storage cache, or if we cannot recycle the
|
|
|
|
// first iteration storage.
|
|
|
|
if (!this->model().enableStorageCache() || !this->recycleFirstIterationStorage()) {
|
|
|
|
this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/1);
|
|
|
|
}
|
|
|
|
|
2019-02-15 06:39:58 -06:00
|
|
|
// initialize the wells. Note that this needs to be done after initializing the
|
|
|
|
// intrinsic permeabilities and the after applying the initial solution because
|
|
|
|
// the well model uses these...
|
2019-04-03 10:26:57 -05:00
|
|
|
wellModel_.init();
|
2016-02-01 10:25:23 -06:00
|
|
|
|
2024-02-06 01:55:42 -06:00
|
|
|
aquiferModel_.initialSolutionApplied();
|
2023-11-09 05:38:25 -06:00
|
|
|
|
2024-06-25 05:55:02 -05:00
|
|
|
const bool invalidateFromHyst = updateHysteresis_();
|
|
|
|
if (invalidateFromHyst) {
|
|
|
|
OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
|
|
|
|
this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
|
|
|
|
}
|
2015-07-28 10:24:44 -05:00
|
|
|
}
|
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
/*!
|
|
|
|
* \copydoc FvBaseProblem::source
|
|
|
|
*
|
|
|
|
* For this problem, the source term of all components is 0 everywhere.
|
|
|
|
*/
|
|
|
|
template <class Context>
|
2016-11-07 08:14:07 -06:00
|
|
|
void source(RateVector& rate,
|
|
|
|
const Context& context,
|
2015-11-18 04:54:35 -06:00
|
|
|
unsigned spaceIdx,
|
|
|
|
unsigned timeIdx) const
|
2022-07-05 02:48:38 -05:00
|
|
|
{
|
|
|
|
const unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
|
|
|
source(rate, globalDofIdx, timeIdx);
|
|
|
|
}
|
|
|
|
|
|
|
|
void source(RateVector& rate,
|
|
|
|
unsigned globalDofIdx,
|
|
|
|
unsigned timeIdx) const
|
2014-04-15 11:30:06 -05:00
|
|
|
{
|
2023-03-02 02:20:07 -06:00
|
|
|
OPM_TIMEBLOCK_LOCAL(eclProblemSource);
|
2015-12-31 06:20:37 -06:00
|
|
|
rate = 0.0;
|
2015-05-21 09:18:45 -05:00
|
|
|
|
2022-10-11 09:10:00 -05:00
|
|
|
// Add well contribution to source here.
|
2022-08-10 03:01:54 -05:00
|
|
|
wellModel_.computeTotalRatesForDof(rate, globalDofIdx);
|
2014-05-07 08:07:39 -05:00
|
|
|
|
2019-02-15 06:39:58 -06:00
|
|
|
// convert the source term from the total mass rate of the
|
|
|
|
// cell to the one per unit of volume as used by the model.
|
|
|
|
for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
|
|
|
|
rate[eqIdx] /= this->model().dofTotalVolume(globalDofIdx);
|
2018-01-30 04:46:23 -06:00
|
|
|
|
2021-05-05 04:22:44 -05:00
|
|
|
Valgrind::CheckDefined(rate[eqIdx]);
|
|
|
|
assert(isfinite(rate[eqIdx]));
|
2016-06-07 10:43:45 -05:00
|
|
|
}
|
2018-10-23 10:47:38 -05:00
|
|
|
|
2022-10-11 09:10:00 -05:00
|
|
|
// Add non-well sources.
|
|
|
|
addToSourceDense(rate, globalDofIdx, timeIdx);
|
|
|
|
}
|
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void addToSourceDense(RateVector& rate,
|
|
|
|
unsigned globalDofIdx,
|
|
|
|
unsigned timeIdx) const = 0;
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2015-11-13 06:24:39 -06:00
|
|
|
/*!
|
|
|
|
* \brief Returns a reference to the ECL well manager used by the problem.
|
|
|
|
*
|
|
|
|
* This can be used for inspecting wells outside of the problem.
|
|
|
|
*/
|
2024-02-06 01:18:32 -06:00
|
|
|
const WellModel& wellModel() const
|
2018-07-09 05:13:57 -05:00
|
|
|
{ return wellModel_; }
|
|
|
|
|
2024-02-06 01:18:32 -06:00
|
|
|
WellModel& wellModel()
|
2018-07-09 05:13:57 -05:00
|
|
|
{ return wellModel_; }
|
2015-11-13 06:24:39 -06:00
|
|
|
|
2024-02-06 01:18:32 -06:00
|
|
|
const AquiferModel& aquiferModel() const
|
2021-05-12 10:12:03 -05:00
|
|
|
{ return aquiferModel_; }
|
|
|
|
|
2024-02-06 01:18:32 -06:00
|
|
|
AquiferModel& mutableAquiferModel()
|
2019-11-28 11:26:22 -06:00
|
|
|
{ return aquiferModel_; }
|
|
|
|
|
2019-03-19 06:31:28 -05:00
|
|
|
bool nonTrivialBoundaryConditions() const
|
|
|
|
{ return nonTrivialBoundaryConditions_; }
|
2019-02-06 04:18:51 -06:00
|
|
|
|
2019-02-15 07:05:45 -06:00
|
|
|
/*!
|
|
|
|
* \brief Propose the size of the next time step to the simulator.
|
|
|
|
*
|
|
|
|
* This method is only called if the Newton solver does converge, the simulator
|
|
|
|
* automatically cuts the time step in half without consultating this method again.
|
|
|
|
*/
|
|
|
|
Scalar nextTimeStepSize() const
|
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(nexTimeStepSize);
|
2019-02-15 09:58:50 -06:00
|
|
|
// allow external code to do the timestepping
|
|
|
|
if (this->nextTimeStepSize_ > 0.0)
|
|
|
|
return this->nextTimeStepSize_;
|
|
|
|
|
2019-04-03 10:26:57 -05:00
|
|
|
const auto& simulator = this->simulator();
|
|
|
|
int episodeIdx = simulator.episodeIndex();
|
2019-02-15 07:05:45 -06:00
|
|
|
|
|
|
|
// for the initial episode, we use a fixed time step size
|
|
|
|
if (episodeIdx < 0)
|
2021-05-25 08:49:14 -05:00
|
|
|
return this->initialTimeStepSize_;
|
2019-02-15 07:05:45 -06:00
|
|
|
|
|
|
|
// ask the newton method for a suggestion. This suggestion will be based on how
|
|
|
|
// well the previous time step converged. After that, apply the runtime time
|
|
|
|
// stepping constraints.
|
|
|
|
const auto& newtonMethod = this->model().newtonMethod();
|
|
|
|
return limitNextTimeStepSize_(newtonMethod.suggestTimeStepSize(simulator.timeStepSize()));
|
|
|
|
}
|
|
|
|
|
2019-03-12 09:51:41 -05:00
|
|
|
/*!
|
|
|
|
* \brief Calculate the porosity multiplier due to water induced rock compaction.
|
|
|
|
*
|
|
|
|
* TODO: The API of this is a bit ad-hoc, it would be better to use context objects.
|
|
|
|
*/
|
|
|
|
template <class LhsEval>
|
|
|
|
LhsEval rockCompPoroMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
|
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK_LOCAL(rockCompPoroMultiplier);
|
2021-05-25 08:49:14 -05:00
|
|
|
if (this->rockCompPoroMult_.empty() && this->rockCompPoroMultWc_.empty())
|
2019-03-12 09:51:41 -05:00
|
|
|
return 1.0;
|
|
|
|
|
|
|
|
unsigned tableIdx = 0;
|
2021-05-25 08:49:14 -05:00
|
|
|
if (!this->rockTableIdx_.empty())
|
|
|
|
tableIdx = this->rockTableIdx_[elementIdx];
|
2019-03-12 09:51:41 -05:00
|
|
|
|
|
|
|
const auto& fs = intQuants.fluidState();
|
2024-01-11 03:24:25 -06:00
|
|
|
LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
|
|
|
|
if (!this->minRefPressure_.empty())
|
2019-03-12 09:51:41 -05:00
|
|
|
// The pore space change is irreversible
|
2024-01-11 03:24:25 -06:00
|
|
|
effectivePressure =
|
2023-12-21 01:40:45 -06:00
|
|
|
min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
|
2024-01-11 03:24:25 -06:00
|
|
|
this->minRefPressure_[elementIdx]);
|
2019-03-12 09:51:41 -05:00
|
|
|
|
2021-05-25 08:49:14 -05:00
|
|
|
if (!this->overburdenPressure_.empty())
|
2024-01-11 03:24:25 -06:00
|
|
|
effectivePressure -= this->overburdenPressure_[elementIdx];
|
2019-03-12 09:51:41 -05:00
|
|
|
|
2020-09-25 08:10:39 -05:00
|
|
|
|
2021-05-25 08:49:14 -05:00
|
|
|
if (!this->rockCompPoroMult_.empty()) {
|
2024-01-11 03:24:25 -06:00
|
|
|
return this->rockCompPoroMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
|
2020-09-25 08:10:39 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
// water compaction
|
2021-05-25 08:49:14 -05:00
|
|
|
assert(!this->rockCompPoroMultWc_.empty());
|
|
|
|
LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
|
2024-09-11 07:25:34 -05:00
|
|
|
LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
|
2020-09-25 08:10:39 -05:00
|
|
|
|
2024-01-11 03:24:25 -06:00
|
|
|
return this->rockCompPoroMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
|
2019-03-12 09:51:41 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Calculate the transmissibility multiplier due to water induced rock compaction.
|
|
|
|
*
|
|
|
|
* TODO: The API of this is a bit ad-hoc, it would be better to use context objects.
|
|
|
|
*/
|
|
|
|
template <class LhsEval>
|
|
|
|
LhsEval rockCompTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
|
|
|
|
{
|
2024-04-26 09:46:27 -05:00
|
|
|
const bool implicit = !this->explicitRockCompaction_;
|
2024-01-04 06:24:36 -06:00
|
|
|
return implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx)
|
|
|
|
: this->simulator().problem().getRockCompTransMultVal(elementIdx);
|
2019-03-12 09:51:41 -05:00
|
|
|
}
|
|
|
|
|
2023-10-09 07:08:18 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Return the well transmissibility multiplier due to rock changues.
|
|
|
|
*/
|
|
|
|
template <class LhsEval>
|
|
|
|
LhsEval wellTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
|
|
|
|
{
|
|
|
|
OPM_TIMEBLOCK_LOCAL(wellTransMultiplier);
|
|
|
|
|
2024-04-26 09:46:27 -05:00
|
|
|
const bool implicit = !this->explicitRockCompaction_;
|
2024-01-04 06:24:36 -06:00
|
|
|
double trans_mult = implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<double>(intQuants, elementIdx)
|
|
|
|
: this->simulator().problem().getRockCompTransMultVal(elementIdx);
|
2023-10-09 07:08:18 -05:00
|
|
|
trans_mult *= this->simulator().problem().template permFactTransMultiplier<double>(intQuants);
|
|
|
|
|
|
|
|
return trans_mult;
|
|
|
|
}
|
|
|
|
|
2023-04-21 08:17:58 -05:00
|
|
|
std::pair<BCType, RateVector> boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
|
2022-09-16 09:17:37 -05:00
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK_LOCAL(boundaryCondition);
|
2022-09-22 02:31:54 -05:00
|
|
|
if (!nonTrivialBoundaryConditions_) {
|
2023-04-21 08:17:58 -05:00
|
|
|
return { BCType::NONE, RateVector(0.0) };
|
2022-09-16 09:17:37 -05:00
|
|
|
}
|
2022-10-17 04:31:46 -05:00
|
|
|
FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
|
2023-04-21 08:17:58 -05:00
|
|
|
const auto& schedule = this->simulator().vanguard().schedule();
|
|
|
|
if (bcindex_(dir)[globalSpaceIdx] == 0) {
|
|
|
|
return { BCType::NONE, RateVector(0.0) };
|
|
|
|
}
|
2023-09-05 07:57:54 -05:00
|
|
|
if (schedule[this->episodeIndex()].bcprop.size() == 0) {
|
|
|
|
return { BCType::NONE, RateVector(0.0) };
|
|
|
|
}
|
2023-04-21 08:17:58 -05:00
|
|
|
const auto& bc = schedule[this->episodeIndex()].bcprop[bcindex_(dir)[globalSpaceIdx]];
|
|
|
|
if (bc.bctype!=BCType::RATE) {
|
|
|
|
return { bc.bctype, RateVector(0.0) };
|
|
|
|
}
|
|
|
|
|
|
|
|
RateVector rate = 0.0;
|
|
|
|
switch (bc.component) {
|
|
|
|
case BCComponent::OIL:
|
|
|
|
rate[Indices::canonicalToActiveComponentIndex(oilCompIdx)] = bc.rate;
|
|
|
|
break;
|
|
|
|
case BCComponent::GAS:
|
|
|
|
rate[Indices::canonicalToActiveComponentIndex(gasCompIdx)] = bc.rate;
|
|
|
|
break;
|
|
|
|
case BCComponent::WATER:
|
|
|
|
rate[Indices::canonicalToActiveComponentIndex(waterCompIdx)] = bc.rate;
|
|
|
|
break;
|
|
|
|
case BCComponent::SOLVENT:
|
2024-09-11 07:25:34 -05:00
|
|
|
this->handleSolventBC(bc, rate);
|
2023-04-21 08:17:58 -05:00
|
|
|
break;
|
|
|
|
case BCComponent::POLYMER:
|
2024-09-11 07:25:34 -05:00
|
|
|
this->handlePolymerBC(bc, rate);
|
2023-04-21 08:17:58 -05:00
|
|
|
break;
|
|
|
|
case BCComponent::NONE:
|
|
|
|
throw std::logic_error("you need to specify the component when RATE type is set in BC");
|
|
|
|
break;
|
|
|
|
}
|
|
|
|
//TODO add support for enthalpy rate
|
|
|
|
return {bc.bctype, rate};
|
2022-09-16 09:17:37 -05:00
|
|
|
}
|
|
|
|
|
2024-06-07 01:58:27 -05:00
|
|
|
|
2023-02-02 04:52:08 -06:00
|
|
|
template<class Serializer>
|
|
|
|
void serializeOp(Serializer& serializer)
|
|
|
|
{
|
|
|
|
serializer(static_cast<BaseType&>(*this));
|
|
|
|
serializer(drift_);
|
|
|
|
serializer(wellModel_);
|
|
|
|
serializer(aquiferModel_);
|
|
|
|
serializer(tracerModel_);
|
2023-02-24 02:58:53 -06:00
|
|
|
serializer(*materialLawManager_);
|
2023-03-24 08:56:23 -05:00
|
|
|
}
|
2024-08-16 03:40:49 -05:00
|
|
|
|
2023-03-23 04:12:02 -05:00
|
|
|
private:
|
|
|
|
Implementation& asImp_()
|
|
|
|
{ return *static_cast<Implementation *>(this); }
|
2024-01-04 06:24:36 -06:00
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
const Implementation& asImp_() const
|
|
|
|
{ return *static_cast<const Implementation *>(this); }
|
2023-02-02 04:52:08 -06:00
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
protected:
|
2022-10-05 04:25:32 -05:00
|
|
|
template<class UpdateFunc>
|
|
|
|
void updateProperty_(const std::string& failureMsg,
|
|
|
|
UpdateFunc func)
|
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(updateProperty);
|
2023-03-23 04:12:02 -05:00
|
|
|
const auto& model = this->simulator().model();
|
|
|
|
const auto& primaryVars = model.solution(/*timeIdx*/0);
|
2022-10-05 04:25:32 -05:00
|
|
|
const auto& vanguard = this->simulator().vanguard();
|
2023-08-15 02:32:10 -05:00
|
|
|
std::size_t numGridDof = primaryVars.size();
|
2022-10-05 04:25:32 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
2023-03-23 04:12:02 -05:00
|
|
|
#ifdef _OPENMP
|
|
|
|
#pragma omp parallel for
|
|
|
|
#endif
|
|
|
|
for (unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
|
|
|
|
const auto& iq = *model.cachedIntensiveQuantities(dofIdx, /*timeIdx=*/ 0);
|
|
|
|
func(dofIdx, iq);
|
2022-10-05 04:25:32 -05:00
|
|
|
}
|
|
|
|
OPM_END_PARALLEL_TRY_CATCH(failureMsg, vanguard.grid().comm());
|
|
|
|
}
|
|
|
|
|
2017-12-19 11:00:14 -06:00
|
|
|
bool updateMaxOilSaturation_()
|
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(updateMaxOilSaturation);
|
2021-05-25 08:49:14 -05:00
|
|
|
int episodeIdx = this->episodeIndex();
|
2019-04-03 10:26:57 -05:00
|
|
|
|
2017-12-19 11:00:14 -06:00
|
|
|
// we use VAPPARS
|
2021-05-25 08:49:14 -05:00
|
|
|
if (this->vapparsActive(episodeIdx)) {
|
2024-02-05 12:22:15 -06:00
|
|
|
this->updateProperty_("FlowProblem::updateMaxOilSaturation_() failed:",
|
2023-03-23 04:12:02 -05:00
|
|
|
[this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
|
|
|
|
{
|
|
|
|
this->updateMaxOilSaturation_(compressedDofIdx,iq);
|
|
|
|
});
|
2017-12-19 11:00:14 -06:00
|
|
|
return true;
|
2017-12-19 05:42:10 -06:00
|
|
|
}
|
2017-12-19 11:00:14 -06:00
|
|
|
|
|
|
|
return false;
|
2017-12-19 05:42:10 -06:00
|
|
|
}
|
|
|
|
|
2023-03-24 08:56:23 -05:00
|
|
|
bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
|
2023-03-23 04:12:02 -05:00
|
|
|
{
|
|
|
|
OPM_TIMEBLOCK_LOCAL(updateMaxOilSaturation);
|
|
|
|
const auto& fs = iq.fluidState();
|
2023-12-21 01:40:45 -06:00
|
|
|
const Scalar So = decay<Scalar>(fs.saturation(refPressurePhaseIdx_()));
|
2023-03-23 04:12:02 -05:00
|
|
|
auto& mos = this->maxOilSaturation_;
|
2023-03-23 07:53:14 -05:00
|
|
|
if(mos[compressedDofIdx] < So){
|
2023-03-23 04:12:02 -05:00
|
|
|
mos[compressedDofIdx] = So;
|
|
|
|
return true;
|
|
|
|
}else{
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2019-03-12 09:51:41 -05:00
|
|
|
bool updateMaxWaterSaturation_()
|
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(updateMaxWaterSaturation);
|
2019-03-12 09:51:41 -05:00
|
|
|
// water compaction is activated in ROCKCOMP
|
2021-05-25 08:49:14 -05:00
|
|
|
if (this->maxWaterSaturation_.empty())
|
2019-03-12 09:51:41 -05:00
|
|
|
return false;
|
|
|
|
|
2021-05-25 08:49:14 -05:00
|
|
|
this->maxWaterSaturation_[/*timeIdx=*/1] = this->maxWaterSaturation_[/*timeIdx=*/0];
|
2024-02-05 12:22:15 -06:00
|
|
|
this->updateProperty_("FlowProblem::updateMaxWaterSaturation_() failed:",
|
2022-10-05 04:26:44 -05:00
|
|
|
[this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
|
|
|
|
{
|
2023-03-23 04:12:02 -05:00
|
|
|
this->updateMaxWaterSaturation_(compressedDofIdx,iq);
|
2022-10-05 04:26:44 -05:00
|
|
|
});
|
2019-03-12 09:51:41 -05:00
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2023-03-23 04:12:02 -05:00
|
|
|
|
|
|
|
bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
|
|
|
|
{
|
|
|
|
OPM_TIMEBLOCK_LOCAL(updateMaxWaterSaturation);
|
|
|
|
const auto& fs = iq.fluidState();
|
|
|
|
const Scalar Sw = decay<Scalar>(fs.saturation(waterPhaseIdx));
|
|
|
|
auto& mow = this->maxWaterSaturation_;
|
|
|
|
if(mow[compressedDofIdx]< Sw){
|
|
|
|
mow[compressedDofIdx] = Sw;
|
|
|
|
return true;
|
|
|
|
}else{
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2019-03-12 09:51:41 -05:00
|
|
|
bool updateMinPressure_()
|
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(updateMinPressure);
|
2019-03-12 09:51:41 -05:00
|
|
|
// IRREVERS option is used in ROCKCOMP
|
2024-01-11 03:24:25 -06:00
|
|
|
if (this->minRefPressure_.empty())
|
2019-03-12 09:51:41 -05:00
|
|
|
return false;
|
|
|
|
|
2024-02-05 12:22:15 -06:00
|
|
|
this->updateProperty_("FlowProblem::updateMinPressure_() failed:",
|
2022-10-05 04:26:44 -05:00
|
|
|
[this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
|
|
|
|
{
|
2023-03-23 04:12:02 -05:00
|
|
|
this->updateMinPressure_(compressedDofIdx,iq);
|
2022-10-05 04:26:44 -05:00
|
|
|
});
|
2019-03-12 09:51:41 -05:00
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2023-03-23 04:12:02 -05:00
|
|
|
bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities& iq){
|
|
|
|
OPM_TIMEBLOCK_LOCAL(updateMinPressure);
|
|
|
|
const auto& fs = iq.fluidState();
|
2024-01-11 03:24:25 -06:00
|
|
|
const Scalar min_pressure = getValue(fs.pressure(refPressurePhaseIdx_()));
|
|
|
|
auto& min_pressures = this->minRefPressure_;
|
|
|
|
if(min_pressures[compressedDofIdx]> min_pressure){
|
|
|
|
min_pressures[compressedDofIdx] = min_pressure;
|
2023-03-23 04:12:02 -05:00
|
|
|
return true;
|
|
|
|
}else{
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2023-12-11 10:54:31 -06:00
|
|
|
// \brief Function to assign field properties of type double, on the leaf grid view.
|
|
|
|
//
|
|
|
|
// For CpGrid with local grid refinement, the field property of a cell on the leaf
|
|
|
|
// is inherited from its parent or equivalent (when has no parent) cell on level zero.
|
2024-01-26 02:40:05 -06:00
|
|
|
std::function<std::vector<double>(const FieldPropsManager&, const std::string&)>
|
2023-12-12 11:48:58 -06:00
|
|
|
fieldPropDoubleOnLeafAssigner_()
|
2023-12-11 10:54:31 -06:00
|
|
|
{
|
|
|
|
const auto& lookup = this->lookUpData_;
|
2024-01-26 02:40:05 -06:00
|
|
|
return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString)
|
2023-12-11 10:54:31 -06:00
|
|
|
{
|
2024-01-26 02:40:05 -06:00
|
|
|
return lookup.assignFieldPropsDoubleOnLeaf(fieldPropManager, propString);
|
2023-12-11 10:54:31 -06:00
|
|
|
};
|
|
|
|
}
|
|
|
|
|
|
|
|
// \brief Function to assign field properties of type int, unsigned int, ..., on the leaf grid view.
|
|
|
|
//
|
|
|
|
// For CpGrid with local grid refinement, the field property of a cell on the leaf
|
|
|
|
// is inherited from its parent or equivalent (when has no parent) cell on level zero.
|
|
|
|
template<typename IntType>
|
2024-01-26 02:40:05 -06:00
|
|
|
std::function<std::vector<IntType>(const FieldPropsManager&, const std::string&, bool)>
|
2023-12-12 11:48:58 -06:00
|
|
|
fieldPropIntTypeOnLeafAssigner_()
|
2023-12-11 10:54:31 -06:00
|
|
|
{
|
|
|
|
const auto& lookup = this->lookUpData_;
|
2024-01-26 02:40:05 -06:00
|
|
|
return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString, bool needsTranslation)
|
2023-12-11 10:54:31 -06:00
|
|
|
{
|
2024-01-26 02:40:05 -06:00
|
|
|
return lookup.template assignFieldPropsIntOnLeaf<IntType>(fieldPropManager, propString, needsTranslation);
|
2023-12-11 10:54:31 -06:00
|
|
|
};
|
|
|
|
}
|
|
|
|
|
2014-06-04 11:05:12 -05:00
|
|
|
void readMaterialParameters_()
|
2014-04-15 11:30:06 -05:00
|
|
|
{
|
2023-02-15 02:41:37 -06:00
|
|
|
OPM_TIMEBLOCK(readMaterialParameters);
|
2019-04-03 10:26:57 -05:00
|
|
|
const auto& simulator = this->simulator();
|
|
|
|
const auto& vanguard = simulator.vanguard();
|
2018-02-01 09:26:58 -06:00
|
|
|
const auto& eclState = vanguard.eclState();
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2017-05-05 03:46:53 -05:00
|
|
|
// the PVT and saturation region numbers
|
2023-06-13 08:44:17 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
2021-05-25 08:49:14 -05:00
|
|
|
this->updatePvtnum_();
|
|
|
|
this->updateSatnum_();
|
2016-10-30 05:22:38 -05:00
|
|
|
|
2017-06-21 03:27:25 -05:00
|
|
|
// the MISC region numbers (solvent model)
|
2021-05-25 08:49:14 -05:00
|
|
|
this->updateMiscnum_();
|
2017-06-21 03:27:25 -05:00
|
|
|
// the PLMIX region numbers (polymer model)
|
2021-05-25 08:49:14 -05:00
|
|
|
this->updatePlmixnum_();
|
2017-06-21 03:27:25 -05:00
|
|
|
|
2023-06-13 08:44:17 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH("Invalid region numbers: ", vanguard.gridView().comm());
|
2014-07-02 10:50:35 -05:00
|
|
|
////////////////////////////////
|
2016-09-06 08:16:03 -05:00
|
|
|
// porosity
|
2019-03-12 09:51:41 -05:00
|
|
|
updateReferencePorosity_();
|
2021-05-25 08:49:14 -05:00
|
|
|
this->referencePorosity_[1] = this->referencePorosity_[0];
|
2016-10-17 07:02:09 -05:00
|
|
|
////////////////////////////////
|
|
|
|
|
|
|
|
////////////////////////////////
|
|
|
|
// fluid-matrix interactions (saturation functions; relperm/capillary pressure)
|
|
|
|
materialLawManager_ = std::make_shared<EclMaterialLawManager>();
|
2020-03-04 02:59:59 -06:00
|
|
|
materialLawManager_->initFromState(eclState);
|
2023-12-12 11:48:58 -06:00
|
|
|
materialLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
|
2023-12-14 08:01:12 -06:00
|
|
|
this-> template fieldPropIntTypeOnLeafAssigner_<int>(),
|
|
|
|
this-> lookupIdxOnLevelZeroAssigner_());
|
2016-10-17 07:02:09 -05:00
|
|
|
////////////////////////////////
|
|
|
|
}
|
|
|
|
|
2018-01-30 04:46:23 -06:00
|
|
|
void readThermalParameters_()
|
|
|
|
{
|
2021-05-26 07:31:57 -05:00
|
|
|
if constexpr (enableEnergy)
|
|
|
|
{
|
|
|
|
const auto& simulator = this->simulator();
|
|
|
|
const auto& vanguard = simulator.vanguard();
|
|
|
|
const auto& eclState = vanguard.eclState();
|
2018-01-30 04:46:23 -06:00
|
|
|
|
2021-05-26 07:31:57 -05:00
|
|
|
// fluid-matrix interactions (saturation functions; relperm/capillary pressure)
|
|
|
|
thermalLawManager_ = std::make_shared<EclThermalLawManager>();
|
2023-12-11 10:54:31 -06:00
|
|
|
thermalLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
|
2023-12-12 11:48:58 -06:00
|
|
|
this-> fieldPropDoubleOnLeafAssigner_(),
|
|
|
|
this-> template fieldPropIntTypeOnLeafAssigner_<unsigned int>());
|
2021-05-26 07:31:57 -05:00
|
|
|
}
|
2018-01-30 04:46:23 -06:00
|
|
|
}
|
|
|
|
|
2019-03-12 09:51:41 -05:00
|
|
|
void updateReferencePorosity_()
|
2016-10-17 07:02:09 -05:00
|
|
|
{
|
2019-04-03 10:26:57 -05:00
|
|
|
const auto& simulator = this->simulator();
|
|
|
|
const auto& vanguard = simulator.vanguard();
|
2018-02-01 09:26:58 -06:00
|
|
|
const auto& eclState = vanguard.eclState();
|
2016-10-17 07:02:09 -05:00
|
|
|
|
2023-08-15 02:32:10 -05:00
|
|
|
std::size_t numDof = this->model().numGridDof();
|
2016-10-17 07:02:09 -05:00
|
|
|
|
2021-05-25 08:49:14 -05:00
|
|
|
this->referencePorosity_[/*timeIdx=*/0].resize(numDof);
|
2016-02-14 16:44:21 -06:00
|
|
|
|
2020-01-13 08:46:50 -06:00
|
|
|
const auto& fp = eclState.fieldProps();
|
2024-01-29 08:05:52 -06:00
|
|
|
const std::vector<double> porvData = this -> fieldPropDoubleOnLeafAssigner_()(fp, "PORV");
|
2023-08-15 02:32:10 -05:00
|
|
|
for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
|
2024-09-25 03:09:39 -05:00
|
|
|
int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(dofIdx);
|
2021-02-28 08:48:02 -06:00
|
|
|
Scalar poreVolume = porvData[dofIdx];
|
2016-09-06 08:16:03 -05:00
|
|
|
|
|
|
|
// we define the porosity as the accumulated pore volume divided by the
|
|
|
|
// geometric volume of the element. Note that -- in pathetic cases -- it can
|
|
|
|
// be larger than 1.0!
|
2024-09-25 03:09:39 -05:00
|
|
|
Scalar dofVolume = simulator.model().dofTotalVolume(sfcdofIdx);
|
2016-09-06 08:16:03 -05:00
|
|
|
assert(dofVolume > 0.0);
|
2024-09-25 03:09:39 -05:00
|
|
|
this->referencePorosity_[/*timeIdx=*/0][sfcdofIdx] = poreVolume/dofVolume;
|
2016-02-14 16:44:21 -06:00
|
|
|
}
|
2014-04-15 11:30:06 -05:00
|
|
|
}
|
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void readInitialCondition_()
|
2015-07-29 06:43:17 -05:00
|
|
|
{
|
2024-09-11 07:25:34 -05:00
|
|
|
// TODO: whether we should move this to FlowProblemBlackoil
|
2019-04-03 10:26:57 -05:00
|
|
|
const auto& simulator = this->simulator();
|
|
|
|
const auto& vanguard = simulator.vanguard();
|
2020-01-22 04:45:03 -06:00
|
|
|
const auto& eclState = vanguard.eclState();
|
2015-07-29 06:43:17 -05:00
|
|
|
|
2020-01-22 04:45:03 -06:00
|
|
|
if (eclState.getInitConfig().hasEquil())
|
2015-07-29 06:43:17 -05:00
|
|
|
readEquilInitialCondition_();
|
2020-01-22 04:45:03 -06:00
|
|
|
else
|
|
|
|
readExplicitInitialCondition_();
|
2017-06-01 00:51:44 -05:00
|
|
|
|
2019-03-12 09:51:41 -05:00
|
|
|
//initialize min/max values
|
2023-08-15 02:32:10 -05:00
|
|
|
std::size_t numElems = this->model().numGridDof();
|
|
|
|
for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
|
2024-09-11 07:25:34 -05:00
|
|
|
const auto& fs = asImp_().initialFluidStates()[elemIdx];
|
2024-09-12 03:55:56 -05:00
|
|
|
if (!this->maxWaterSaturation_.empty() && waterPhaseIdx > -1)
|
2021-05-25 08:49:14 -05:00
|
|
|
this->maxWaterSaturation_[elemIdx] = std::max(this->maxWaterSaturation_[elemIdx], fs.saturation(waterPhaseIdx));
|
2024-09-12 03:55:56 -05:00
|
|
|
if (!this->maxOilSaturation_.empty() && oilPhaseIdx > -1)
|
2021-05-25 08:49:14 -05:00
|
|
|
this->maxOilSaturation_[elemIdx] = std::max(this->maxOilSaturation_[elemIdx], fs.saturation(oilPhaseIdx));
|
2024-09-12 03:55:56 -05:00
|
|
|
if (!this->minRefPressure_.empty() && refPressurePhaseIdx_() > -1)
|
2024-01-11 03:24:25 -06:00
|
|
|
this->minRefPressure_[elemIdx] = std::min(this->minRefPressure_[elemIdx], fs.pressure(refPressurePhaseIdx_()));
|
2019-03-12 09:51:41 -05:00
|
|
|
}
|
2015-07-29 06:43:17 -05:00
|
|
|
}
|
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void readEquilInitialCondition_() = 0;
|
|
|
|
virtual void readExplicitInitialCondition_() = 0;
|
2017-12-19 11:00:14 -06:00
|
|
|
|
2015-07-28 10:24:44 -05:00
|
|
|
// update the hysteresis parameters of the material laws for the whole grid
|
2016-02-11 09:25:06 -06:00
|
|
|
bool updateHysteresis_()
|
2015-07-28 10:24:44 -05:00
|
|
|
{
|
2015-07-29 06:43:17 -05:00
|
|
|
if (!materialLawManager_->enableHysteresis())
|
2016-02-11 09:25:06 -06:00
|
|
|
return false;
|
2015-07-28 10:24:44 -05:00
|
|
|
|
2017-04-15 08:09:34 -05:00
|
|
|
// we need to update the hysteresis data for _all_ elements (i.e., not just the
|
|
|
|
// interior ones) to avoid desynchronization of the processes in the parallel case!
|
2024-02-05 12:22:15 -06:00
|
|
|
this->updateProperty_("FlowProblem::updateHysteresis_() failed:",
|
2022-10-05 04:26:44 -05:00
|
|
|
[this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
|
|
|
|
{
|
|
|
|
materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
|
|
|
|
});
|
2016-02-11 09:25:06 -06:00
|
|
|
return true;
|
2015-07-28 10:24:44 -05:00
|
|
|
}
|
|
|
|
|
2023-03-23 04:12:02 -05:00
|
|
|
|
|
|
|
bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
|
|
|
|
{
|
|
|
|
OPM_TIMEBLOCK_LOCAL(updateHysteresis_);
|
|
|
|
materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
|
|
|
|
//TODO change materials to give a bool
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2024-01-04 06:24:36 -06:00
|
|
|
Scalar getRockCompTransMultVal(std::size_t dofIdx) const
|
|
|
|
{
|
|
|
|
if (this->rockCompTransMultVal_.empty())
|
|
|
|
return 1.0;
|
|
|
|
|
|
|
|
return this->rockCompTransMultVal_[dofIdx];
|
|
|
|
}
|
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
protected:
|
2016-10-30 12:42:02 -05:00
|
|
|
struct PffDofData_
|
|
|
|
{
|
2021-05-05 04:22:44 -05:00
|
|
|
ConditionalStorage<enableEnergy, Scalar> thermalHalfTransIn;
|
|
|
|
ConditionalStorage<enableEnergy, Scalar> thermalHalfTransOut;
|
|
|
|
ConditionalStorage<enableDiffusion, Scalar> diffusivity;
|
2023-10-25 12:46:55 -05:00
|
|
|
ConditionalStorage<enableDispersion, Scalar> dispersivity;
|
2016-10-30 12:42:02 -05:00
|
|
|
Scalar transmissibility;
|
|
|
|
};
|
|
|
|
|
|
|
|
// update the prefetch friendly data object
|
|
|
|
void updatePffDofData_()
|
|
|
|
{
|
2016-11-07 08:14:07 -06:00
|
|
|
const auto& distFn =
|
2016-10-30 12:42:02 -05:00
|
|
|
[this](PffDofData_& dofData,
|
|
|
|
const Stencil& stencil,
|
|
|
|
unsigned localDofIdx)
|
|
|
|
-> void
|
|
|
|
{
|
|
|
|
const auto& elementMapper = this->model().elementMapper();
|
|
|
|
|
|
|
|
unsigned globalElemIdx = elementMapper.index(stencil.entity(localDofIdx));
|
|
|
|
if (localDofIdx != 0) {
|
2016-10-30 16:38:41 -05:00
|
|
|
unsigned globalCenterElemIdx = elementMapper.index(stencil.entity(/*dofIdx=*/0));
|
2016-10-30 12:42:02 -05:00
|
|
|
dofData.transmissibility = transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
|
2018-01-30 04:46:23 -06:00
|
|
|
|
2021-05-26 07:31:57 -05:00
|
|
|
if constexpr (enableEnergy) {
|
2021-03-26 09:27:14 -05:00
|
|
|
*dofData.thermalHalfTransIn = transmissibilities_.thermalHalfTrans(globalCenterElemIdx, globalElemIdx);
|
|
|
|
*dofData.thermalHalfTransOut = transmissibilities_.thermalHalfTrans(globalElemIdx, globalCenterElemIdx);
|
|
|
|
}
|
2021-05-26 07:31:57 -05:00
|
|
|
if constexpr (enableDiffusion)
|
2021-01-14 08:33:28 -06:00
|
|
|
*dofData.diffusivity = transmissibilities_.diffusivity(globalCenterElemIdx, globalElemIdx);
|
2023-10-25 12:46:55 -05:00
|
|
|
if (enableDispersion)
|
|
|
|
dofData.dispersivity = transmissibilities_.dispersivity(globalCenterElemIdx, globalElemIdx);
|
2016-10-30 12:42:02 -05:00
|
|
|
}
|
|
|
|
};
|
|
|
|
|
|
|
|
pffDofData_.update(distFn);
|
|
|
|
}
|
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart) = 0;
|
|
|
|
|
2019-02-06 04:18:51 -06:00
|
|
|
void readBoundaryConditions_()
|
|
|
|
{
|
2019-04-03 10:26:57 -05:00
|
|
|
const auto& simulator = this->simulator();
|
|
|
|
const auto& vanguard = simulator.vanguard();
|
2020-01-25 02:40:50 -06:00
|
|
|
const auto& bcconfig = vanguard.eclState().getSimulationConfig().bcconfig();
|
|
|
|
if (bcconfig.size() > 0) {
|
2019-03-19 06:31:28 -05:00
|
|
|
nonTrivialBoundaryConditions_ = true;
|
2019-05-28 04:36:05 -05:00
|
|
|
|
2023-08-15 02:32:10 -05:00
|
|
|
std::size_t numCartDof = vanguard.cartesianSize();
|
2019-02-19 08:49:50 -06:00
|
|
|
unsigned numElems = vanguard.gridView().size(/*codim=*/0);
|
2020-01-25 02:40:50 -06:00
|
|
|
std::vector<int> cartesianToCompressedElemIdx(numCartDof, -1);
|
2019-02-19 08:49:50 -06:00
|
|
|
|
|
|
|
for (unsigned elemIdx = 0; elemIdx < numElems; ++elemIdx)
|
|
|
|
cartesianToCompressedElemIdx[vanguard.cartesianIndex(elemIdx)] = elemIdx;
|
|
|
|
|
2023-04-21 08:17:58 -05:00
|
|
|
bcindex_.resize(numElems, 0);
|
2022-10-17 04:31:46 -05:00
|
|
|
auto loopAndApply = [&cartesianToCompressedElemIdx,
|
|
|
|
&vanguard](const auto& bcface,
|
|
|
|
auto apply)
|
|
|
|
{
|
|
|
|
for (int i = bcface.i1; i <= bcface.i2; ++i) {
|
|
|
|
for (int j = bcface.j1; j <= bcface.j2; ++j) {
|
|
|
|
for (int k = bcface.k1; k <= bcface.k2; ++k) {
|
|
|
|
std::array<int, 3> tmp = {i,j,k};
|
|
|
|
auto elemIdx = cartesianToCompressedElemIdx[vanguard.cartesianIndex(tmp)];
|
|
|
|
if (elemIdx >= 0)
|
|
|
|
apply(elemIdx);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
};
|
2020-01-25 02:40:50 -06:00
|
|
|
for (const auto& bcface : bcconfig) {
|
2023-04-21 08:17:58 -05:00
|
|
|
std::vector<int>& data = bcindex_(bcface.dir);
|
|
|
|
const int index = bcface.index;
|
2022-10-17 04:31:46 -05:00
|
|
|
loopAndApply(bcface,
|
2023-04-21 08:17:58 -05:00
|
|
|
[&data,index](int elemIdx)
|
|
|
|
{ data[elemIdx] = index; });
|
2019-02-19 08:49:50 -06:00
|
|
|
}
|
|
|
|
}
|
2019-02-06 04:18:51 -06:00
|
|
|
}
|
|
|
|
|
2019-02-15 07:05:45 -06:00
|
|
|
// this method applies the runtime constraints specified via the deck and/or command
|
2019-04-03 10:26:57 -05:00
|
|
|
// line parameters for the size of the next time step.
|
2019-02-15 07:05:45 -06:00
|
|
|
Scalar limitNextTimeStepSize_(Scalar dtNext) const
|
|
|
|
{
|
2021-05-26 07:31:57 -05:00
|
|
|
if constexpr (enableExperiments) {
|
|
|
|
const auto& simulator = this->simulator();
|
2023-08-30 03:51:04 -05:00
|
|
|
const auto& schedule = simulator.vanguard().schedule();
|
2021-05-26 07:31:57 -05:00
|
|
|
int episodeIdx = simulator.episodeIndex();
|
|
|
|
|
|
|
|
// first thing in the morning, limit the time step size to the maximum size
|
2024-07-06 03:22:47 -05:00
|
|
|
Scalar maxTimeStepSize = Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60;
|
2023-08-30 03:51:04 -05:00
|
|
|
int reportStepIdx = std::max(episodeIdx, 0);
|
|
|
|
if (this->enableTuning_) {
|
|
|
|
const auto& tuning = schedule[reportStepIdx].tuning();
|
|
|
|
maxTimeStepSize = tuning.TSMAXZ;
|
|
|
|
}
|
|
|
|
|
|
|
|
dtNext = std::min(dtNext, maxTimeStepSize);
|
2021-05-26 07:31:57 -05:00
|
|
|
|
|
|
|
Scalar remainingEpisodeTime =
|
|
|
|
simulator.episodeStartTime() + simulator.episodeLength()
|
|
|
|
- (simulator.startTime() + simulator.time());
|
|
|
|
assert(remainingEpisodeTime >= 0.0);
|
|
|
|
|
|
|
|
// if we would have a small amount of time left over in the current episode, make
|
|
|
|
// two equal time steps instead of a big and a small one
|
|
|
|
if (remainingEpisodeTime/2.0 < dtNext && dtNext < remainingEpisodeTime*(1.0 - 1e-5))
|
|
|
|
// note: limiting to the maximum time step size here is probably not strictly
|
|
|
|
// necessary, but it should not hurt and is more fool-proof
|
2023-08-30 03:51:04 -05:00
|
|
|
dtNext = std::min(maxTimeStepSize, remainingEpisodeTime/2.0);
|
2021-05-26 07:31:57 -05:00
|
|
|
|
|
|
|
if (simulator.episodeStarts()) {
|
2021-07-02 07:56:19 -05:00
|
|
|
// if a well event occurred, respect the limit for the maximum time step after
|
2021-05-26 07:31:57 -05:00
|
|
|
// that, too
|
|
|
|
const auto& events = simulator.vanguard().schedule()[reportStepIdx].events();
|
|
|
|
bool wellEventOccured =
|
|
|
|
events.hasEvent(ScheduleEvents::NEW_WELL)
|
|
|
|
|| events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE)
|
|
|
|
|| events.hasEvent(ScheduleEvents::INJECTION_UPDATE)
|
|
|
|
|| events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
|
2021-05-25 08:49:14 -05:00
|
|
|
if (episodeIdx >= 0 && wellEventOccured && this->maxTimeStepAfterWellEvent_ > 0)
|
|
|
|
dtNext = std::min(dtNext, this->maxTimeStepAfterWellEvent_);
|
2021-05-26 07:31:57 -05:00
|
|
|
}
|
2019-06-20 08:42:48 -05:00
|
|
|
}
|
2019-02-15 07:05:45 -06:00
|
|
|
|
|
|
|
return dtNext;
|
|
|
|
}
|
|
|
|
|
2023-12-21 01:40:45 -06:00
|
|
|
int refPressurePhaseIdx_() const {
|
|
|
|
if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
|
|
|
|
return oilPhaseIdx;
|
|
|
|
}
|
|
|
|
else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
|
|
|
|
return gasPhaseIdx;
|
|
|
|
}
|
|
|
|
else {
|
|
|
|
return waterPhaseIdx;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2024-01-04 06:24:36 -06:00
|
|
|
void updateRockCompTransMultVal_()
|
|
|
|
{
|
|
|
|
const auto& model = this->simulator().model();
|
|
|
|
std::size_t numGridDof = this->model().numGridDof();
|
|
|
|
this->rockCompTransMultVal_.resize(numGridDof, 1.0);
|
|
|
|
for (std::size_t elementIdx = 0; elementIdx < numGridDof; ++elementIdx) {
|
|
|
|
const auto& iq = *model.cachedIntensiveQuantities(elementIdx, /*timeIdx=*/ 0);
|
|
|
|
Scalar trans_mult = computeRockCompTransMultiplier_<Scalar>(iq, elementIdx);
|
|
|
|
this->rockCompTransMultVal_[elementIdx] = trans_mult;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Calculate the transmissibility multiplier due to water induced rock compaction.
|
|
|
|
*
|
|
|
|
* TODO: The API of this is a bit ad-hoc, it would be better to use context objects.
|
|
|
|
*/
|
|
|
|
template <class LhsEval>
|
|
|
|
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities& intQuants, unsigned elementIdx) const
|
|
|
|
{
|
|
|
|
OPM_TIMEBLOCK_LOCAL(computeRockCompTransMultiplier);
|
|
|
|
if (this->rockCompTransMult_.empty() && this->rockCompTransMultWc_.empty())
|
|
|
|
return 1.0;
|
|
|
|
|
|
|
|
unsigned tableIdx = 0;
|
|
|
|
if (!this->rockTableIdx_.empty())
|
|
|
|
tableIdx = this->rockTableIdx_[elementIdx];
|
|
|
|
|
|
|
|
const auto& fs = intQuants.fluidState();
|
|
|
|
LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
|
|
|
|
|
|
|
|
if (!this->minRefPressure_.empty())
|
|
|
|
// The pore space change is irreversible
|
|
|
|
effectivePressure =
|
|
|
|
min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
|
|
|
|
this->minRefPressure_[elementIdx]);
|
|
|
|
|
|
|
|
if (!this->overburdenPressure_.empty())
|
|
|
|
effectivePressure -= this->overburdenPressure_[elementIdx];
|
|
|
|
|
|
|
|
if (!this->rockCompTransMult_.empty())
|
|
|
|
return this->rockCompTransMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
|
|
|
|
|
|
|
|
// water compaction
|
|
|
|
assert(!this->rockCompTransMultWc_.empty());
|
|
|
|
LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
|
2024-09-11 07:25:34 -05:00
|
|
|
LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
|
2024-01-04 06:24:36 -06:00
|
|
|
|
|
|
|
return this->rockCompTransMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
|
|
|
|
}
|
|
|
|
|
2021-05-05 05:13:25 -05:00
|
|
|
typename Vanguard::TransmissibilityType transmissibilities_;
|
2014-07-02 10:50:35 -05:00
|
|
|
|
2015-07-29 06:43:17 -05:00
|
|
|
std::shared_ptr<EclMaterialLawManager> materialLawManager_;
|
2018-01-30 04:46:23 -06:00
|
|
|
std::shared_ptr<EclThermalLawManager> thermalLawManager_;
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2019-03-05 04:32:11 -06:00
|
|
|
bool enableDriftCompensation_;
|
|
|
|
GlobalEqVector drift_;
|
|
|
|
|
2024-02-06 01:18:32 -06:00
|
|
|
WellModel wellModel_;
|
2024-02-06 01:18:32 -06:00
|
|
|
AquiferModel aquiferModel_;
|
2019-04-03 10:26:57 -05:00
|
|
|
|
2024-07-31 06:35:56 -05:00
|
|
|
bool enableVtkOutput_;
|
2024-09-11 07:25:34 -05:00
|
|
|
|
2016-10-30 12:42:02 -05:00
|
|
|
|
2016-10-30 14:09:58 -05:00
|
|
|
PffGridVector<GridView, Stencil, PffDofData_, DofMapper> pffDofData_;
|
2018-11-14 06:10:01 -06:00
|
|
|
TracerModel tracerModel_;
|
2022-08-12 07:20:33 -05:00
|
|
|
|
2022-10-17 04:31:46 -05:00
|
|
|
template<class T>
|
|
|
|
struct BCData
|
|
|
|
{
|
|
|
|
std::array<std::vector<T>,6> data;
|
|
|
|
|
2023-08-15 02:32:10 -05:00
|
|
|
void resize(std::size_t size, T defVal)
|
2022-10-17 04:31:46 -05:00
|
|
|
{
|
|
|
|
for (auto& d : data)
|
|
|
|
d.resize(size, defVal);
|
|
|
|
}
|
|
|
|
|
|
|
|
const std::vector<T>& operator()(FaceDir::DirEnum dir) const
|
|
|
|
{
|
|
|
|
if (dir == FaceDir::DirEnum::Unknown)
|
2022-11-01 15:44:59 -05:00
|
|
|
throw std::runtime_error("Tried to access BC data for the 'Unknown' direction");
|
2022-10-17 04:31:46 -05:00
|
|
|
int idx = 0;
|
|
|
|
int div = static_cast<int>(dir);
|
|
|
|
while ((div /= 2) >= 1)
|
|
|
|
++idx;
|
|
|
|
assert(idx >= 0 && idx <= 5);
|
|
|
|
return data[idx];
|
|
|
|
}
|
|
|
|
|
|
|
|
std::vector<T>& operator()(FaceDir::DirEnum dir)
|
|
|
|
{
|
|
|
|
return const_cast<std::vector<T>&>(std::as_const(*this)(dir));
|
|
|
|
}
|
|
|
|
};
|
|
|
|
|
2024-09-11 07:25:34 -05:00
|
|
|
virtual void handleSolventBC(const BCProp::BCFace&, RateVector&) const = 0;
|
|
|
|
|
|
|
|
virtual void handlePolymerBC(const BCProp::BCFace&, RateVector&) const = 0;
|
|
|
|
|
2023-04-21 08:17:58 -05:00
|
|
|
BCData<int> bcindex_;
|
2022-10-17 04:31:46 -05:00
|
|
|
bool nonTrivialBoundaryConditions_ = false;
|
2024-04-26 09:46:27 -05:00
|
|
|
bool explicitRockCompaction_ = false;
|
2024-08-28 09:28:51 -05:00
|
|
|
bool first_step_ = true;
|
2024-08-27 08:51:56 -05:00
|
|
|
|
2014-04-15 11:30:06 -05:00
|
|
|
};
|
2018-08-06 03:43:28 -05:00
|
|
|
|
2019-09-05 10:04:39 -05:00
|
|
|
} // namespace Opm
|
2014-04-15 11:30:06 -05:00
|
|
|
|
2024-02-05 12:22:15 -06:00
|
|
|
#endif // OPM_FLOW_PROBLEM_HPP
|