13#ifndef DUMUX_2P_BOX_MATERIAL_INTERFACES_HH
14#define DUMUX_2P_BOX_MATERIAL_INTERFACES_HH
16#include <dune/common/exceptions.hh>
33template<
class Gr
idGeometry,
class PcKrSw>
36 using SubControlVolume =
typename GridGeometry::SubControlVolume;
39 template<
class SpatialParams,
class SolutionVector>
41 const SpatialParams& spatialParams,
42 const SolutionVector& x)
44 update(gridGeometry, spatialParams, x);
54 template<
class SpatialParams,
class SolutionVector>
55 void update(
const GridGeometry& gridGeometry,
56 const SpatialParams& spatialParams,
57 const SolutionVector& x)
61 DUNE_THROW(Dune::InvalidStateException,
"Determination of the interface material parameters with "
62 "this class only makes sense when using the box method!");
64 isOnMaterialInterface_.resize(gridGeometry.numDofs(),
false);
65 pcSwAtDof_.resize(gridGeometry.numDofs(),
nullptr);
66 auto fvGeometry =
localView(gridGeometry);
67 for (
const auto& element : elements(gridGeometry.gridView()))
70 fvGeometry.bind(element);
72 for (
const auto& scv : scvs(fvGeometry))
74 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
75 const auto& pcKrSw = fluidMatrixInteraction.pcSwCurve();
78 static_assert(std::is_lvalue_reference<
typename std::decay_t<
decltype(fluidMatrixInteraction)>::PcKrSwType>::value,
79 "In order to use the box-interface solver please provide access "
80 "to the material law parameters via returning (const) references");
83 if (!pcSwAtDof_[scv.dofIndex()])
84 pcSwAtDof_[scv.dofIndex()] = &pcKrSw;
87 else if (pcKrSw.endPointPc() < pcSwAtDof_[scv.dofIndex()]->endPointPc())
89 pcSwAtDof_[scv.dofIndex()] = &pcKrSw;
90 isOnMaterialInterface_[scv.dofIndex()] =
true;
94 else if ( !(pcKrSw == *(pcSwAtDof_[scv.dofIndex()])) )
95 isOnMaterialInterface_[scv.dofIndex()] =
true;
102 {
return isOnMaterialInterface_[scv.dofIndex()]; }
106 {
return *(pcSwAtDof_[scv.dofIndex()]); }
109 std::vector<bool> isOnMaterialInterface_;
110 std::vector<const PcKrSw*> pcSwAtDof_;
Class that determines the material with the lowest capillary pressure (under fully water-saturated co...
Definition: boxmaterialinterfaces.hh:35
void update(const GridGeometry &gridGeometry, const SpatialParams &spatialParams, const SolutionVector &x)
Updates the scv -> dofparameter map.
Definition: boxmaterialinterfaces.hh:55
bool isOnMaterialInterface(const SubControlVolume &scv) const
Returns if this scv is connected to a material interface.
Definition: boxmaterialinterfaces.hh:101
BoxMaterialInterfaces(const GridGeometry &gridGeometry, const SpatialParams &spatialParams, const SolutionVector &x)
Definition: boxmaterialinterfaces.hh:40
const PcKrSw & pcSwAtDof(const SubControlVolume &scv) const
Returns the material parameters associated with the dof.
Definition: boxmaterialinterfaces.hh:105
The local element solution class for control-volume finite element methods.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
PcKrSw(T &&) -> PcKrSw< T >
Deduction guide for the PcKrSw class. Makes sure that PcKrSw stores a copy of T if the constructor is...
The available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147