version 3.10-dev
boxmaterialinterfaces.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_2P_BOX_MATERIAL_INTERFACES_HH
14#define DUMUX_2P_BOX_MATERIAL_INTERFACES_HH
15
16#include <dune/common/exceptions.hh>
17
20
21namespace Dumux {
22
33template<class GridGeometry, class PcKrSw>
35{
36 using SubControlVolume = typename GridGeometry::SubControlVolume;
37
38public:
39 template<class SpatialParams, class SolutionVector>
40 BoxMaterialInterfaces(const GridGeometry& gridGeometry,
41 const SpatialParams& spatialParams,
42 const SolutionVector& x)
43 {
44 update(gridGeometry, spatialParams, x);
45 }
46
54 template<class SpatialParams, class SolutionVector>
55 void update(const GridGeometry& gridGeometry,
56 const SpatialParams& spatialParams,
57 const SolutionVector& x)
58 {
59 // make sure this is only called for geometries of the box method!
60 if (GridGeometry::discMethod != DiscretizationMethods::box)
61 DUNE_THROW(Dune::InvalidStateException, "Determination of the interface material parameters with "
62 "this class only makes sense when using the box method!");
63
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()))
68 {
69 const auto elemSol = elementSolution(element, x, gridGeometry);
70 fvGeometry.bind(element);
71
72 for (const auto& scv : scvs(fvGeometry))
73 {
74 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
75 const auto& pcKrSw = fluidMatrixInteraction.pcSwCurve();
76
77 // assert current preconditions on requiring the spatial params to store the pckrSw curve
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");
81
82 // if no parameters had been set, set them now
83 if (!pcSwAtDof_[scv.dofIndex()])
84 pcSwAtDof_[scv.dofIndex()] = &pcKrSw;
85
86 // otherwise only use the current ones if endPointPc (e.g. Brooks-Corey entry pressure) is lower
87 else if (pcKrSw.endPointPc() < pcSwAtDof_[scv.dofIndex()]->endPointPc())
88 {
89 pcSwAtDof_[scv.dofIndex()] = &pcKrSw;
90 isOnMaterialInterface_[scv.dofIndex()] = true;
91 }
92
93 // keep track of material interfaces in any case
94 else if ( !(pcKrSw == *(pcSwAtDof_[scv.dofIndex()])) )
95 isOnMaterialInterface_[scv.dofIndex()] = true;
96 }
97 }
98 }
99
101 bool isOnMaterialInterface(const SubControlVolume& scv) const
102 { return isOnMaterialInterface_[scv.dofIndex()]; }
103
105 const PcKrSw& pcSwAtDof(const SubControlVolume& scv) const
106 { return *(pcSwAtDof_[scv.dofIndex()]); }
107
108private:
109 std::vector<bool> isOnMaterialInterface_;
110 std::vector<const PcKrSw*> pcSwAtDof_;
111};
112
113} // end namespace Dumux
114
115#endif
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
Definition: adapt.hh:17