14#ifndef DUMUX_DISCRETIZATION_CC_MPFA_LOCAL_ASSEMBLER_BASE_HH
15#define DUMUX_DISCRETIZATION_CC_MPFA_LOCAL_ASSEMBLER_BASE_HH
21#include <dune/common/exceptions.hh>
22#include <dune/common/reservedvector.hh>
42template<
class P,
class EG,
class EV >
47 using FVElementGeometry = EG;
48 using ElementVolumeVariables = EV;
53 using Scalar =
typename IV::Traits::MatVecTraits::FaceVector::value_type;
74 const Problem&
problem()
const {
return *problemPtr_; }
75 const FVElementGeometry&
fvGeometry()
const {
return *fvGeometryPtr_; }
76 const ElementVolumeVariables&
elemVolVars()
const {
return *elemVolVarsPtr_; }
93 template<
class DataHandle,
class IV,
class TensorFunc >
94 void assembleMatrices(DataHandle& handle, IV& iv,
const TensorFunc& getT, Scalar<IV> wijZeroThresh = 0.0)
96 DUNE_THROW(Dune::NotImplemented,
"Implementation does not provide an assembleMatrices() function");
110 template<
class DataHandle,
class IV,
class GetU >
111 void assembleU(DataHandle& handle,
const IV& iv,
const GetU& getU)
113 DUNE_THROW(Dune::NotImplemented,
"Implementation does not provide an assemble() function for the cell/Dirichlet unknowns");
124 template<
class DataHandle,
class IV,
class GetRho >
127 using GridView =
typename IV::Traits::GridView;
128 static constexpr int dim = GridView::dimension;
129 static constexpr int dimWorld = GridView::dimensionworld;
130 static constexpr bool isSurfaceGrid = dim < dimWorld;
133 auto& g = handle.g();
134 auto& deltaG = handle.deltaG();
135 auto& outsideG = handle.gOutside();
145 using Scalar =
typename IV::Traits::MatVecTraits::TMatrix::value_type;
146 using LocalIndexType =
typename IV::Traits::IndexSet::LocalIndexType;
148 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
151 const auto& curLocalScvf = iv.localScvf(faceIdx);
152 const auto& curGlobalScvf =
fvGeometry().scvf(curLocalScvf.gridScvfIndex());
153 const auto& gravity =
problem().spatialParams().gravity(curGlobalScvf.ipGlobal());
156 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
157 const auto& posGlobalScv =
fvGeometry().scv(iv.localScv(neighborScvIndices[0]).gridScvIndex());
158 const auto& posVolVars =
elemVolVars()[posGlobalScv];
159 const auto alpha_inside = posVolVars.extrusionFactor()*
vtmv(curGlobalScvf.unitOuterNormal(),
160 posVolVars.permeability(),
163 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
164 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
166 Dune::ReservedVector<Scalar, 1> >;
167 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
168 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
174 if (!curLocalScvf.isDirichlet())
176 const auto localDofIdx = curLocalScvf.localDofIndex();
178 rho = getRho(posVolVars);
179 deltaG[localDofIdx] = 0.0;
181 if (!curGlobalScvf.boundary())
183 for (
unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
186 const auto negLocalScvIdx = neighborScvIndices[idxInOutside+1];
187 const auto& negGlobalScv =
fvGeometry().scv(iv.localScv(negLocalScvIdx).gridScvIndex());
188 const auto& negVolVars =
elemVolVars()[negGlobalScv];
189 const auto& flipScvf = !isSurfaceGrid ? curGlobalScvf
190 :
fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside);
192 alpha_outside[idxInOutside] = negVolVars.extrusionFactor()*
vtmv(flipScvf.unitOuterNormal(),
193 negVolVars.permeability(),
196 alpha_outside[idxInOutside] *= -1.0;
198 rho += getRho(negVolVars);
199 deltaG[localDofIdx] += alpha_outside[idxInOutside];
203 rho /= numOutsideFaces + 1;
204 deltaG[localDofIdx] -= alpha_inside;
205 deltaG[localDofIdx] *= rho*Extrusion::area(
fvGeometry(), curGlobalScvf);
209 rho = getRho(
elemVolVars()[curGlobalScvf.outsideScvIdx()]);
212 g[faceIdx] = alpha_inside*rho*Extrusion::area(
fvGeometry(), curGlobalScvf);
217 for (
const auto& alpha : alpha_outside)
218 outsideG[faceIdx][i++] = alpha*rho*Extrusion::area(
fvGeometry(), curGlobalScvf);
223 handle.CA().umv(deltaG, g);
226 using FaceVector =
typename IV::Traits::MatVecTraits::FaceVector;
229 handle.A().mv(deltaG, AG);
232 for (
const auto& localFaceData : iv.localFaceData())
235 if (!localFaceData.isOutsideFace())
238 const auto localScvIdx = localFaceData.ivLocalInsideScvIndex();
239 const auto localScvfIdx = localFaceData.ivLocalScvfIndex();
240 const auto idxInOutside = localFaceData.scvfLocalOutsideScvfIndex();
241 const auto& posLocalScv = iv.localScv(localScvIdx);
242 const auto& wijk = handle.omegas()[localScvfIdx][idxInOutside+1];
245 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
248 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
249 if (!curLocalScvf.isDirichlet())
250 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
258 const Problem* problemPtr_;
259 const FVElementGeometry* fvGeometryPtr_;
260 const ElementVolumeVariables* elemVolVarsPtr_;
Defines the general interface of the local assembler classes for the assembly of the interaction volu...
Definition: localassemblerbase.hh:45
const FVElementGeometry & fvGeometry() const
Definition: localassemblerbase.hh:75
void assembleMatrices(DataHandle &handle, IV &iv, const TensorFunc &getT, Scalar< IV > wijZeroThresh=0.0)
Assembles the matrices involved in the flux expressions and the local system of equations within an m...
Definition: localassemblerbase.hh:94
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:76
InteractionVolumeAssemblerBase(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars)
The constructor. Sets pointers to the objects required for a subsequent call to assemble().
Definition: localassemblerbase.hh:64
void assembleU(DataHandle &handle, const IV &iv, const GetU &getU)
Assembles the vector of primary (cell) unknowns and (maybe) Dirichlet boundary conditions within an i...
Definition: localassemblerbase.hh:111
void assembleGravity(DataHandle &handle, const IV &iv, const GetRho &getRho)
Assembles the gravitational flux contributions on the scvfs within an interaction volume.
Definition: localassemblerbase.hh:125
const Problem & problem() const
Definition: localassemblerbase.hh:74
A class that contains helper functions as well as functionality which is common to different mpfa sch...
Definition: localassemblerhelper.hh:34
static void resizeVector(Vector &v, size_type size)
resizes a vector to the given size (specialization for dynamic matrix type)
Definition: localassemblerhelper.hh:226
Helper classes to compute the integration elements.
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:851
A helper function for class member function introspection.
A class that contains helper functions as well as functionality which is common to different mpfa sch...
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166