26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_LOCAL_ASSEMBLER_BASE_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_LOCAL_ASSEMBLER_BASE_HH
33#include <dune/common/exceptions.hh>
34#include <dune/common/reservedvector.hh>
53template<
class P,
class EG,
class EV >
58 using FVElementGeometry = EG;
59 using ElementVolumeVariables = EV;
82 const Problem&
problem()
const {
return *problemPtr_; }
83 const FVElementGeometry&
fvGeometry()
const {
return *fvGeometryPtr_; }
84 const ElementVolumeVariables&
elemVolVars()
const {
return *elemVolVarsPtr_; }
99 template<
class DataHandle,
class IV,
class TensorFunc >
102 DUNE_THROW(Dune::NotImplemented,
"Implementation does not provide a assembleMatrices() function");
116 template<
class DataHandle,
class IV,
class GetU >
117 void assembleU(DataHandle& handle,
const IV& iv,
const GetU& getU)
119 DUNE_THROW(Dune::NotImplemented,
"Implementation does not provide a assemble() function for the cell/Dirichlet unknowns");
130 template<
class DataHandle,
class IV,
class GetRho >
133 using GridView =
typename IV::Traits::GridView;
134 static constexpr int dim = GridView::dimension;
135 static constexpr int dimWorld = GridView::dimensionworld;
136 static constexpr bool isSurfaceGrid = dim < dimWorld;
139 auto& g = handle.g();
140 auto& deltaG = handle.deltaG();
141 auto& outsideG = handle.gOutside();
151 using Scalar =
typename IV::Traits::MatVecTraits::TMatrix::value_type;
152 using LocalIndexType =
typename IV::Traits::IndexSet::LocalIndexType;
154 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
157 const auto& curLocalScvf = iv.localScvf(faceIdx);
158 const auto& curGlobalScvf =
fvGeometry().scvf(curLocalScvf.gridScvfIndex());
159 const auto& gravity =
problem().spatialParams().gravity(curGlobalScvf.ipGlobal());
162 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
163 const auto& posGlobalScv =
fvGeometry().scv(iv.localScv(neighborScvIndices[0]).gridScvIndex());
164 const auto& posVolVars =
elemVolVars()[posGlobalScv];
165 const auto alpha_inside = posVolVars.extrusionFactor()*
vtmv(curGlobalScvf.unitOuterNormal(),
166 posVolVars.permeability(),
169 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
170 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
172 Dune::ReservedVector<Scalar, 1> >;
173 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
174 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
180 if (!curLocalScvf.isDirichlet())
182 const auto localDofIdx = curLocalScvf.localDofIndex();
184 rho = getRho(posVolVars);
185 deltaG[localDofIdx] = 0.0;
187 if (!curGlobalScvf.boundary())
189 for (
unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
192 const auto negLocalScvIdx = neighborScvIndices[idxInOutside+1];
193 const auto& negGlobalScv =
fvGeometry().scv(iv.localScv(negLocalScvIdx).gridScvIndex());
194 const auto& negVolVars =
elemVolVars()[negGlobalScv];
195 const auto& flipScvf = !isSurfaceGrid ? curGlobalScvf
196 :
fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside);
198 alpha_outside[idxInOutside] = negVolVars.extrusionFactor()*
vtmv(flipScvf.unitOuterNormal(),
199 negVolVars.permeability(),
202 alpha_outside[idxInOutside] *= -1.0;
204 rho += getRho(negVolVars);
205 deltaG[localDofIdx] += alpha_outside[idxInOutside];
209 rho /= numOutsideFaces + 1;
210 deltaG[localDofIdx] -= alpha_inside;
211 deltaG[localDofIdx] *= rho*curGlobalScvf.area();
215 rho = getRho(
elemVolVars()[curGlobalScvf.outsideScvIdx()]);
218 g[faceIdx] = alpha_inside*rho*curGlobalScvf.area();
223 for (
const auto& alpha : alpha_outside)
224 outsideG[faceIdx][i++] = alpha*rho*curGlobalScvf.area();
229 handle.CA().umv(deltaG, g);
232 using FaceVector =
typename IV::Traits::MatVecTraits::FaceVector;
235 handle.A().mv(deltaG, AG);
238 for (
const auto& localFaceData : iv.localFaceData())
241 if (!localFaceData.isOutsideFace())
244 const auto localScvIdx = localFaceData.ivLocalInsideScvIndex();
245 const auto localScvfIdx = localFaceData.ivLocalScvfIndex();
246 const auto idxInOutside = localFaceData.scvfLocalOutsideScvfIndex();
247 const auto& posLocalScv = iv.localScv(localScvIdx);
248 const auto& wijk = handle.omegas()[localScvfIdx][idxInOutside+1];
251 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
254 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
255 if (!curLocalScvf.isDirichlet())
256 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
264 const Problem* problemPtr_;
265 const FVElementGeometry* fvGeometryPtr_;
266 const ElementVolumeVariables* elemVolVarsPtr_;
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...
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:840
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Defines the general interface of the local assembler classes for the assembly of the interaction volu...
Definition: localassemblerbase.hh:56
void assembleMatrices(DataHandle &handle, IV &iv, const TensorFunc &getT)
Assembles the matrices involved in the flux expressions and the local system of equations within an m...
Definition: localassemblerbase.hh:100
const FVElementGeometry & fvGeometry() const
Definition: localassemblerbase.hh:83
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:84
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:72
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:117
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:131
const Problem & problem() const
Definition: localassemblerbase.hh:82
A class that contains helper functions as well as functionality which is common to different mpfa sch...
Definition: localassemblerhelper.hh:46
static void resizeVector(Vector &v, size_type size)
resizes a vector to the given size (specialization for dynamic matrix type)
Definition: localassemblerhelper.hh:238