24#ifndef DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_LOCAL_ASSEMBLER_HH
25#define DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_LOCAL_ASSEMBLER_HH
27#include <dune/common/exceptions.hh>
28#include <dune/common/float_cmp.hh>
49template<
class P,
class EG,
class EV >
58 using ParentType::ParentType;
68 template<
class DataHandle,
class IV,
class TensorFunc >
71 assembleLocalMatrices_(handle.A(), handle.AB(), handle.CA(), handle.T(), handle.omegas(), iv, getT);
74 if (iv.numUnknowns() > 0)
86 template<
class DataHandle,
class IV,
class GetU >
87 void assembleU(DataHandle& handle,
const IV& iv,
const GetU& getU)
89 auto& u = handle.uj();
93 typename IV::Traits::IndexSet::LocalIndexType i = 0;
94 for (; i < iv.numScvs(); i++)
95 u[i] = getU( this->
elemVolVars()[iv.localScv(i).gridScvIndex()] );
98 for (
const auto& data : iv.interiorBoundaryData())
100 const auto& scvf = this->
fvGeometry().scvf(data.scvfIndex());
101 const auto element = this->
fvGeometry().gridGeometry().element(scvf.insideScvIdx());
102 u[i++] = getU( this->
problem().couplingManager().getLowDimVolVars(element, scvf) );
106 for (
const auto& data : iv.dirichletData())
107 u[i++] = getU( this->
elemVolVars()[data.volVarIndex()] );
118 template<
class DataHandle,
class IV,
class GetRho >
121 using GridView =
typename IV::Traits::GridView;
122 static constexpr int dim = GridView::dimension;
123 static constexpr int dimWorld = GridView::dimensionworld;
124 static constexpr bool isSurfaceGrid = dim < dimWorld;
127 auto& g = handle.g();
128 auto& deltaG = handle.deltaG();
129 auto& outsideG = handle.gOutside();
139 using Scalar =
typename IV::Traits::MatVecTraits::TMatrix::value_type;
140 using LocalIndexType =
typename IV::Traits::IndexSet::LocalIndexType;
143 static const Scalar xi = getParamFromGroup<Scalar>(this->
problem().paramGroup(),
"FacetCoupling.Xi", 1.0);
145 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
148 const auto& curLocalScvf = iv.localScvf(faceIdx);
149 const auto& curGlobalScvf = this->
fvGeometry().scvf(curLocalScvf.gridScvfIndex());
150 const auto& gravity = this->
problem().spatialParams().gravity(curGlobalScvf.ipGlobal());
151 const auto curIsInteriorBoundary = curLocalScvf.isOnInteriorBoundary();
152 const Scalar curXiFactor = curIsInteriorBoundary ? (curGlobalScvf.boundary() ? 1.0 : xi) : 1.0;
155 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
156 const auto& posGlobalScv = this->
fvGeometry().scv(iv.localScv(neighborScvIndices[0]).gridScvIndex());
157 const auto& posVolVars = this->
elemVolVars()[posGlobalScv];
158 const auto alpha_inside = posVolVars.extrusionFactor()*
vtmv(curGlobalScvf.unitOuterNormal(),
159 posVolVars.permeability(),
162 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
163 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
165 Dune::ReservedVector<Scalar, 1> >;
166 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
167 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
170 if (isSurfaceGrid && !curIsInteriorBoundary)
173 if (!curLocalScvf.isDirichlet())
175 const auto localDofIdx = curLocalScvf.localDofIndex();
176 const Scalar curOneMinusXi = curIsInteriorBoundary ? -(1.0 - xi) : 1.0;
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 = this->
fvGeometry().scv(iv.localScv(negLocalScvIdx).gridScvIndex());
188 const auto& negVolVars = this->
elemVolVars()[negGlobalScv];
189 const auto& flipScvf = !isSurfaceGrid ? curGlobalScvf
190 : this->
fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside);
192 alpha_outside[idxInOutside] = negVolVars.extrusionFactor()*
vtmv(flipScvf.unitOuterNormal(),
193 negVolVars.permeability(),
196 alpha_outside[idxInOutside] *= -1.0;
198 if (!curIsInteriorBoundary)
199 rho += getRho(negVolVars);
201 deltaG[localDofIdx] += curOneMinusXi*alpha_outside[idxInOutside];
205 if (curIsInteriorBoundary)
207 const auto posElement = this->
fvGeometry().gridGeometry().element(posGlobalScv.elementIndex());
208 const auto& facetVolVars = this->
problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
209 rho += getRho(facetVolVars);
211 const auto alphaFacet = posVolVars.extrusionFactor()*
vtmv(curGlobalScvf.unitOuterNormal(),
212 facetVolVars.permeability(),
214 deltaG[localDofIdx] += alphaFacet;
217 rho /= numOutsideFaces + 1;
219 deltaG[localDofIdx] -= curXiFactor*alpha_inside;
220 deltaG[localDofIdx] *= rho*curGlobalScvf.area();
223 else if (curIsInteriorBoundary)
225 const auto posElement = this->
fvGeometry().gridGeometry().element(posGlobalScv.elementIndex());
226 const auto& facetVolVars = this->
problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
227 rho = 0.5*(getRho(facetVolVars) + getRho(posVolVars));
231 rho = getRho(this->
elemVolVars()[curGlobalScvf.outsideScvIdx()]);
234 g[faceIdx] = alpha_inside*rho*curGlobalScvf.area();
236 if (isSurfaceGrid && !curIsInteriorBoundary)
239 for (
const auto& alpha : alpha_outside)
240 outsideG[faceIdx][i++] = alpha*rho*curGlobalScvf.area();
245 handle.CA().umv(deltaG, g);
248 using FaceVector =
typename IV::Traits::MatVecTraits::FaceVector;
251 handle.A().mv(deltaG, AG);
254 for (
const auto& localFaceData : iv.localFaceData())
257 if (!localFaceData.isOutsideFace())
260 const auto localScvIdx = localFaceData.ivLocalInsideScvIndex();
261 const auto localScvfIdx = localFaceData.ivLocalScvfIndex();
262 const auto idxInOutside = localFaceData.scvfLocalOutsideScvfIndex();
263 const auto& posLocalScv = iv.localScv(localScvIdx);
264 const auto& wijk = handle.omegas()[localScvfIdx][idxInOutside+1];
267 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
270 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
271 if (!curLocalScvf.isDirichlet())
272 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
282 template<
class IV,
class TensorFunc >
283 void assembleLocalMatrices_(
typename IV::Traits::MatVecTraits::AMatrix& A,
284 typename IV::Traits::MatVecTraits::BMatrix& B,
285 typename IV::Traits::MatVecTraits::CMatrix& C,
286 typename IV::Traits::MatVecTraits::DMatrix& D,
287 typename IV::Traits::MatVecTraits::OmegaStorage& wijk,
288 IV& iv,
const TensorFunc& getT)
290 using Scalar =
typename IV::Traits::MatVecTraits::TMatrix::value_type;
291 using LocalIndexType =
typename IV::Traits::IndexSet::LocalIndexType;
292 static constexpr int dim = IV::Traits::GridView::dimension;
293 static constexpr int dimWorld = IV::Traits::GridView::dimensionworld;
296 static const Scalar xi = getParamFromGroup<Scalar>(this->
problem().paramGroup(),
"FacetCoupling.Xi", 1.0);
302 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
303 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
310 if (iv.numUnknowns() == 0)
316 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
318 const auto& curLocalScvf = iv.localScvf(faceIdx);
319 const auto& curGlobalScvf = this->
fvGeometry().scvf(curLocalScvf.gridScvfIndex());
320 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
323 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
324 const auto& posGlobalScv = this->
fvGeometry().scv(posLocalScv.gridScvIndex());
325 const auto& posVolVars = this->
elemVolVars()[posGlobalScv];
326 const auto& posElement = iv.element(neighborScvIndices[0]);
327 const auto tensor = getT(this->
problem(), posElement, posVolVars, this->
fvGeometry(), posGlobalScv);
333 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
334 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
336 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
337 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
338 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
339 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
351 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
353 const auto& curLocalScvf = iv.localScvf(faceIdx);
354 const auto& curGlobalScvf = this->
fvGeometry().scvf(curLocalScvf.gridScvfIndex());
355 const auto curIsDirichlet = curLocalScvf.isDirichlet();
356 const auto curLocalDofIdx = curLocalScvf.localDofIndex();
357 const auto curIsInteriorBoundary = curLocalScvf.isOnInteriorBoundary();
358 const Scalar curXiFactor = curIsInteriorBoundary ? (curGlobalScvf.boundary() ? 1.0 : xi) : 1.0;
361 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
362 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
363 const auto& posGlobalScv = this->
fvGeometry().scv(posLocalScv.gridScvIndex());
364 const auto& posVolVars = this->
elemVolVars()[posGlobalScv];
365 const auto& posElement = iv.element(neighborScvIndices[0]);
366 const auto tensor = getT(this->
problem(), posElement, posVolVars, this->
fvGeometry(), posGlobalScv);
373 for (
unsigned int localDir = 0; localDir < dim; localDir++)
375 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
376 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
380 if (!otherLocalScvf.isDirichlet())
382 C[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
384 A[curLocalDofIdx][otherLocalDofIdx] -= curXiFactor*wijk[faceIdx][0][localDir];
389 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
391 B[curLocalDofIdx][otherLocalDofIdx] += curXiFactor*wijk[faceIdx][0][localDir];
395 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
396 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
399 B[curLocalDofIdx][posScvLocalDofIdx] -= curXiFactor*wijk[faceIdx][0][localDir];
403 if (curIsInteriorBoundary && !curIsDirichlet)
405 const auto facetVolVars = this->
problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
406 const auto facetTensor = this->
problem().couplingManager().getLowDimTensor(posElement, curGlobalScvf, getT);
410 const auto wFacet = 2.0*curGlobalScvf.area()*posVolVars.extrusionFactor()
411 *
vtmv(curGlobalScvf.unitOuterNormal(), facetTensor, curGlobalScvf.unitOuterNormal())
412 / (dim < dimWorld ? sqrt(facetVolVars.extrusionFactor()) : facetVolVars.extrusionFactor());
414 A[curLocalDofIdx][curLocalDofIdx] -= wFacet;
415 B[curLocalDofIdx][curLocalScvf.coupledFacetLocalDofIndex()] -= wFacet;
419 if (!curGlobalScvf.boundary() && !curIsDirichlet)
422 const Scalar curOneMinusXi = curIsInteriorBoundary ? -(1.0 - xi) : 1.0;
425 for (
unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
427 const auto idxOnScvf = idxInOutside+1;
428 const auto& negLocalScv = iv.localScv( neighborScvIndices[idxOnScvf] );
429 const auto& negGlobalScv = this->
fvGeometry().scv(negLocalScv.gridScvIndex());
430 const auto& negVolVars = this->
elemVolVars()[negGlobalScv];
431 const auto& negElement = iv.element( neighborScvIndices[idxOnScvf] );
432 const auto negTensor = getT(this->
problem(), negElement, negVolVars, this->
fvGeometry(), negGlobalScv);
435 const auto& scvf = dim < dimWorld ? this->
fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside)
441 wijk[faceIdx][idxOnScvf] *= -1.0;
444 for (
int localDir = 0; localDir < dim; localDir++)
446 const auto otherLocalScvfIdx = negLocalScv.localScvfIndex(localDir);
447 const auto& otherLocalScvf = iv.localScvf(otherLocalScvfIdx);
448 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
450 if (!otherLocalScvf.isDirichlet())
451 A[curLocalDofIdx][otherLocalDofIdx] += curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
453 B[curLocalDofIdx][otherLocalDofIdx] -= curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
456 B[curLocalDofIdx][negLocalScv.localDofIndex()] += curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Defines the general interface of classes used for the assembly of the local systems of equations invo...
A class that contains helper functions as well as functionality which is common to different mpfa sch...
The available mpfa schemes in Dumux.
Dune::FieldVector< typename Tensor::field_type, Scv::myDimension > computeMpfaTransmissibility(const Scv &scv, const Scvf &scvf, const Tensor &T, typename Scv::ctype extrusionFactor)
Free function to evaluate the Mpfa transmissibility associated with the flux (in the form of flux = T...
Definition: mpfa/computetransmissibility.hh:51
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
const FVElementGeometry & fvGeometry() const
Definition: localassemblerbase.hh:83
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:84
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 solveLocalSystem(const FVElementGeometry &fvGeometry, DataHandle &handle, IV &iv)
Solves a previously assembled iv-local system of equations and stores the resulting transmissibilitie...
Definition: localassemblerhelper.hh:82
static void resizeMatrix(Matrix &M, size_type rows, size_type cols)
resizes a matrix to the given sizes (specialization for dynamic matrix type)
Definition: localassemblerhelper.hh:222
static void resizeVector(Vector &v, size_type size)
resizes a vector to the given size (specialization for dynamic matrix type)
Definition: localassemblerhelper.hh:238
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: multidomain/facet/cellcentered/mpfa/localassembler.hh:52
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: multidomain/facet/cellcentered/mpfa/localassembler.hh:87
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 i...
Definition: multidomain/facet/cellcentered/mpfa/localassembler.hh:69
void assembleGravity(DataHandle &handle, const IV &iv, const GetRho &getRho)
Assembles the gravitational flux contributions on the scvfs within an interaction volume.
Definition: multidomain/facet/cellcentered/mpfa/localassembler.hh:119
This file contains free functions to evaluate the transmissibilities associated with flux evaluations...