28#ifndef EWOMS_MULTI_PHASE_BASE_PROBLEM_HH
29#define EWOMS_MULTI_PHASE_BASE_PROBLEM_HH
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
34#include <dune/grid/common/partitionset.hh>
36#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
37#include <opm/material/common/Means.hpp>
38#include <opm/material/densead/Evaluation.hpp>
47#include <opm/utility/CopyablePtr.hpp>
61template<
class TypeTag>
63 :
public FvBaseProblem<TypeTag>
64 ,
public GetPropType<TypeTag, Properties::FluxModule>::FluxBaseProblem
67 using ParentType = FvBaseProblem<TypeTag>;
78 using DirectionalMobilityPtr = Utility::CopyablePtr<DirectionalMobility<TypeTag>>;
80 enum { dimWorld = GridView::dimensionworld };
83 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
84 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
100 ParentType::registerParameters();
103 (
"Use the gravity correction for the pressure gradients.");
115 template <
class Context>
117 const Context& context,
118 unsigned intersectionIdx,
119 unsigned timeIdx)
const
121 const auto& scvf = context.stencil(timeIdx).interiorFace(intersectionIdx);
123 const DimMatrix& K1 = asImp_().intrinsicPermeability(context, scvf.interiorIndex(), timeIdx);
124 const DimMatrix& K2 = asImp_().intrinsicPermeability(context, scvf.exteriorIndex(), timeIdx);
128 for (
unsigned i = 0; i < dimWorld; ++i) {
129 for (
unsigned j = 0; j < dimWorld; ++j) {
130 result[i][j] = harmonicMean(K1[i][j], K2[i][j]);
148 template <
class Context>
153 throw std::logic_error(
"Not implemented: Problem::intrinsicPermeability()");
165 template <
class Context>
170 throw std::logic_error(
"Not implemented: Problem::porosity()");
182 template <
class Context>
183 const SolidEnergyLawParams&
188 throw std::logic_error(
"Not implemented: Problem::solidEnergyParams()");
200 template <
class Context>
201 const ThermalConductionLawParams&
206 throw std::logic_error(
"Not implemented: Problem::thermalConductionParams()");
217 template <
class Context>
222 throw std::logic_error(
"Not implemented: Problem::tortuosity()");
233 template <
class Context>
238 throw std::logic_error(
"Not implemented: Problem::dispersivity()");
254 template <
class Context>
255 const MaterialLawParams &
260 static MaterialLawParams dummy;
264 template <
class Flu
idState>
265 void updateRelperms([[maybe_unused]] std::array<Evaluation,numPhases>& mobility,
266 [[maybe_unused]] DirectionalMobilityPtr& dirMob,
267 [[maybe_unused]] FluidState& fluidState,
268 [[maybe_unused]]
unsigned globalSpaceIdx)
const
279 template <
class Context>
283 {
return asImp_().temperature(); }
293 {
throw std::logic_error(
"Not implemented:temperature() method not implemented by the actual problem"); }
304 template <
class Context>
308 {
return asImp_().gravity(); }
329 using Toolbox = MathToolbox<Evaluation>;
331 unsigned numMarked = 0;
332 ElementContext elemCtx( this->
simulator() );
334 auto& grid = this->
simulator().vanguard().grid();
335 for (
const auto& element : elements(
gridView, Dune::Partitions::interior)) {
336 elemCtx.updateAll(element);
339 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
340 Scalar minSat = 1e100 ;
341 Scalar maxSat = -1e100;
342 size_t nDofs = elemCtx.numDof(0);
343 for (
unsigned dofIdx = 0; dofIdx < nDofs; ++dofIdx) {
344 const auto& intQuant = elemCtx.intensiveQuantities( dofIdx, 0);
345 minSat = std::min(minSat,
346 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
347 maxSat = std::max(maxSat,
348 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
351 const Scalar indicator =
352 (maxSat - minSat) / (std::max<Scalar>(0.01, maxSat + minSat) / 2);
353 if (indicator > 0.2 && element.level() < 2) {
354 grid.mark(1, element);
357 else if (indicator < 0.025) {
358 grid.mark(-1, element);
362 grid.mark(0, element);
368 return this->
simulator().vanguard().grid().comm().sum(numMarked);
387 for (
unsigned i = 0; i < DimMatrix::rows; ++i) {
397 Implementation& asImp_()
398 {
return *
static_cast<Implementation *
>(
this); }
401 const Implementation& asImp_()
const
402 {
return *
static_cast<const Implementation *
>(
this); }
407 if (Parameters::Get<Parameters::EnableGravity>()) {
408 gravity_[dimWorld-1] = -9.81;
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition fvbaseproblem.hh:694
const GridView & gridView() const
The GridView which used by the problem.
Definition fvbaseproblem.hh:662
void intersectionIntrinsicPermeability(DimMatrix &result, const Context &context, unsigned intersectionIdx, unsigned timeIdx) const
Returns the intrinsic permeability of an intersection.
Definition multiphasebaseproblem.hh:116
const ThermalConductionLawParams & thermalConductionParams(const Context &, unsigned, unsigned) const
Returns the parameter object for the thermal conductivity law in a sub-control volume.
Definition multiphasebaseproblem.hh:202
const MaterialLawParams & materialLawParams(const Context &, unsigned, unsigned) const
Returns the material law parameters within a control volume.
Definition multiphasebaseproblem.hh:256
DimMatrix toDimMatrix_(Scalar val) const
Converts a Scalar value to an isotropic Tensor.
Definition multiphasebaseproblem.hh:384
Scalar temperature() const
Returns the temperature for an isothermal problem.
Definition multiphasebaseproblem.hh:292
Scalar tortuosity(const Context &, unsigned, unsigned) const
Define the tortuosity.
Definition multiphasebaseproblem.hh:218
Scalar temperature(const Context &, unsigned, unsigned) const
Returns the temperature within a control volume.
Definition multiphasebaseproblem.hh:280
const DimVector & gravity() const
Returns the acceleration due to gravity .
Definition multiphasebaseproblem.hh:319
const DimVector & gravity(const Context &, unsigned, unsigned) const
Returns the acceleration due to gravity .
Definition multiphasebaseproblem.hh:305
Scalar porosity(const Context &, unsigned, unsigned) const
Returns the porosity [] of the porous medium for a given control volume.
Definition multiphasebaseproblem.hh:166
MultiPhaseBaseProblem(Simulator &simulator)
Definition multiphasebaseproblem.hh:91
const SolidEnergyLawParams & solidEnergyParams(const Context &, unsigned, unsigned) const
Returns the parameter object for the energy storage law of the solid in a sub-control volume.
Definition multiphasebaseproblem.hh:184
Scalar dispersivity(const Context &, unsigned, unsigned) const
Define the dispersivity.
Definition multiphasebaseproblem.hh:234
const DimMatrix & intrinsicPermeability(const Context &, unsigned, unsigned) const
Returns the intrinsic permeability tensor at a given position.
Definition multiphasebaseproblem.hh:149
unsigned markForGridAdaptation()
Mark grid cells for refinement or coarsening.
Definition multiphasebaseproblem.hh:327
static void registerParameters()
Register all run-time parameters for the problem and the model.
Definition multiphasebaseproblem.hh:98
This file contains definitions related to directional mobilities.
Base class for all problems which use a finite volume spatial discretization.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common parameters for the porous medium multi-phase models.
Defines the common properties required by the porous medium multi-phase models.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:233
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
void Register(const char *usageString)
Register a run-time parameter.
Definition parametersystem.hpp:292