24#ifndef DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_LOCAL_ASSEMBLER_HH
25#define DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_LOCAL_ASSEMBLER_HH
29#include <dune/common/exceptions.hh>
30#include <dune/common/float_cmp.hh>
52template<
class P,
class EG,
class EV >
61 using Scalar =
typename IV::Traits::MatVecTraits::FaceVector::value_type;
65 using ParentType::ParentType;
80 template<
class DataHandle,
class IV,
class TensorFunc >
81 void assembleMatrices(DataHandle& handle, IV& iv,
const TensorFunc& getT, Scalar<IV> wijZeroThresh = 0.0)
83 assembleLocalMatrices_(handle.A(), handle.AB(), handle.CA(), handle.T(), handle.omegas(), iv, getT, wijZeroThresh);
86 if (iv.numUnknowns() > 0)
98 template<
class DataHandle,
class IV,
class GetU >
99 void assembleU(DataHandle& handle,
const IV& iv,
const GetU& getU)
101 auto& u = handle.uj();
105 typename IV::Traits::IndexSet::LocalIndexType i = 0;
106 for (; i < iv.numScvs(); i++)
107 u[i] = getU( this->
elemVolVars()[iv.localScv(i).gridScvIndex()] );
110 for (
const auto& data : iv.interiorBoundaryData())
112 const auto& scvf = this->
fvGeometry().scvf(data.scvfIndex());
113 const auto element = this->
fvGeometry().gridGeometry().element(scvf.insideScvIdx());
114 u[i++] = getU( this->
problem().couplingManager().getLowDimVolVars(element, scvf) );
118 for (
const auto& data : iv.dirichletData())
119 u[i++] = getU( this->
elemVolVars()[data.volVarIndex()] );
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;
155 static const Scalar xi = getParamFromGroup<Scalar>(this->
problem().paramGroup(),
"FacetCoupling.Xi", 1.0);
157 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
160 const auto& curLocalScvf = iv.localScvf(faceIdx);
161 const auto& curGlobalScvf = this->
fvGeometry().scvf(curLocalScvf.gridScvfIndex());
162 const auto& gravity = this->
problem().spatialParams().gravity(curGlobalScvf.ipGlobal());
163 const auto curIsInteriorBoundary = curLocalScvf.isOnInteriorBoundary();
164 const Scalar curXiFactor = curIsInteriorBoundary ? (curGlobalScvf.boundary() ? 1.0 : xi) : 1.0;
167 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
168 const auto& posGlobalScv = this->
fvGeometry().scv(iv.localScv(neighborScvIndices[0]).gridScvIndex());
169 const auto& posVolVars = this->
elemVolVars()[posGlobalScv];
170 const auto alpha_inside = posVolVars.extrusionFactor()*
vtmv(curGlobalScvf.unitOuterNormal(),
171 posVolVars.permeability(),
174 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
175 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
177 Dune::ReservedVector<Scalar, 1> >;
178 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
179 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
182 if (isSurfaceGrid && !curIsInteriorBoundary)
185 if (!curLocalScvf.isDirichlet())
187 const auto localDofIdx = curLocalScvf.localDofIndex();
188 const Scalar curOneMinusXi = curIsInteriorBoundary ? -(1.0 - xi) : 1.0;
190 rho = getRho(posVolVars);
191 deltaG[localDofIdx] = 0.0;
193 if (!curGlobalScvf.boundary())
195 for (
unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
198 const auto negLocalScvIdx = neighborScvIndices[idxInOutside+1];
199 const auto& negGlobalScv = this->
fvGeometry().scv(iv.localScv(negLocalScvIdx).gridScvIndex());
200 const auto& negVolVars = this->
elemVolVars()[negGlobalScv];
201 const auto& flipScvf = !isSurfaceGrid ? curGlobalScvf
202 : this->
fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside);
204 alpha_outside[idxInOutside] = negVolVars.extrusionFactor()*
vtmv(flipScvf.unitOuterNormal(),
205 negVolVars.permeability(),
208 alpha_outside[idxInOutside] *= -1.0;
210 if (!curIsInteriorBoundary)
211 rho += getRho(negVolVars);
213 deltaG[localDofIdx] += curOneMinusXi*alpha_outside[idxInOutside];
217 if (curIsInteriorBoundary)
219 const auto posElement = this->
fvGeometry().gridGeometry().element(posGlobalScv.elementIndex());
220 const auto& facetVolVars = this->
problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
221 rho += getRho(facetVolVars);
223 const auto alphaFacet = posVolVars.extrusionFactor()*
vtmv(curGlobalScvf.unitOuterNormal(),
224 facetVolVars.permeability(),
226 deltaG[localDofIdx] += alphaFacet;
229 rho /= numOutsideFaces + 1;
231 deltaG[localDofIdx] -= curXiFactor*alpha_inside;
232 deltaG[localDofIdx] *= rho*Extrusion::area(curGlobalScvf);
235 else if (curIsInteriorBoundary)
237 const auto posElement = this->
fvGeometry().gridGeometry().element(posGlobalScv.elementIndex());
238 const auto& facetVolVars = this->
problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
239 rho = 0.5*(getRho(facetVolVars) + getRho(posVolVars));
243 rho = getRho(this->
elemVolVars()[curGlobalScvf.outsideScvIdx()]);
246 g[faceIdx] = alpha_inside*rho*Extrusion::area(curGlobalScvf);
248 if (isSurfaceGrid && !curIsInteriorBoundary)
251 for (
const auto& alpha : alpha_outside)
252 outsideG[faceIdx][i++] = alpha*rho*Extrusion::area(curGlobalScvf);
257 handle.CA().umv(deltaG, g);
260 using FaceVector =
typename IV::Traits::MatVecTraits::FaceVector;
263 handle.A().mv(deltaG, AG);
266 for (
const auto& localFaceData : iv.localFaceData())
269 if (!localFaceData.isOutsideFace())
272 const auto localScvIdx = localFaceData.ivLocalInsideScvIndex();
273 const auto localScvfIdx = localFaceData.ivLocalScvfIndex();
274 const auto idxInOutside = localFaceData.scvfLocalOutsideScvfIndex();
275 const auto& posLocalScv = iv.localScv(localScvIdx);
276 const auto& wijk = handle.omegas()[localScvfIdx][idxInOutside+1];
279 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
282 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
283 if (!curLocalScvf.isDirichlet())
284 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
294 template<
class IV,
class TensorFunc >
295 void assembleLocalMatrices_(
typename IV::Traits::MatVecTraits::AMatrix& A,
296 typename IV::Traits::MatVecTraits::BMatrix& B,
297 typename IV::Traits::MatVecTraits::CMatrix& C,
298 typename IV::Traits::MatVecTraits::DMatrix& D,
299 typename IV::Traits::MatVecTraits::OmegaStorage& wijk,
300 IV& iv,
const TensorFunc& getT,
301 Scalar<IV> wijZeroThresh)
303 using Scalar =
typename IV::Traits::MatVecTraits::TMatrix::value_type;
304 using LocalIndexType =
typename IV::Traits::IndexSet::LocalIndexType;
305 static constexpr int dim = IV::Traits::GridView::dimension;
306 static constexpr int dimWorld = IV::Traits::GridView::dimensionworld;
309 static const Scalar xi = getParamFromGroup<Scalar>(this->
problem().paramGroup(),
"FacetCoupling.Xi", 1.0);
315 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
316 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
318 std::vector< std::pair<LocalIndexType, LocalIndexType> > faceMarkers;
323 if (iv.numUnknowns() == 0)
329 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
331 const auto& curLocalScvf = iv.localScvf(faceIdx);
332 const auto& curGlobalScvf = this->
fvGeometry().scvf(curLocalScvf.gridScvfIndex());
333 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
336 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
337 const auto& posGlobalScv = this->
fvGeometry().scv(posLocalScv.gridScvIndex());
338 const auto& posVolVars = this->
elemVolVars()[posGlobalScv];
339 const auto tensor = getT(posVolVars);
343 wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
345 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
346 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
348 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
349 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
350 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
351 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
363 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
365 const auto& curLocalScvf = iv.localScvf(faceIdx);
366 const auto& curGlobalScvf = this->
fvGeometry().scvf(curLocalScvf.gridScvfIndex());
367 const auto curIsDirichlet = curLocalScvf.isDirichlet();
368 const auto curLocalDofIdx = curLocalScvf.localDofIndex();
369 const auto curIsInteriorBoundary = curLocalScvf.isOnInteriorBoundary();
370 const Scalar curXiFactor = curIsInteriorBoundary ? (curGlobalScvf.boundary() ? 1.0 : xi) : 1.0;
373 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
374 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
375 const auto& posGlobalScv = this->
fvGeometry().scv(posLocalScv.gridScvIndex());
376 const auto& posVolVars = this->
elemVolVars()[posGlobalScv];
377 const auto& posElement = iv.element(neighborScvIndices[0]);
378 const auto tensor = getT(posVolVars);
382 wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
385 bool isZeroWij =
false;
388 for (
unsigned int localDir = 0; localDir < dim; localDir++)
390 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
391 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
393 if (otherLocalDofIdx == curLocalDofIdx)
395 if (abs(wijk[faceIdx][0][localDir]) <= wijZeroThresh)
400 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
407 if (!otherLocalScvf.isDirichlet())
409 C[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
411 A[curLocalDofIdx][otherLocalDofIdx] -= curXiFactor*wijk[faceIdx][0][localDir];
416 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
418 B[curLocalDofIdx][otherLocalDofIdx] += curXiFactor*wijk[faceIdx][0][localDir];
422 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
423 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
426 B[curLocalDofIdx][posScvLocalDofIdx] -= curXiFactor*wijk[faceIdx][0][localDir];
430 if (curIsInteriorBoundary && !curIsDirichlet)
432 const auto facetVolVars = this->
problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
433 const auto facetTensor = getT(facetVolVars);
437 const auto wFacet = 2.0*Extrusion::area(curGlobalScvf)*posVolVars.extrusionFactor()
438 *
vtmv(curGlobalScvf.unitOuterNormal(), facetTensor, curGlobalScvf.unitOuterNormal())
439 / (dim < dimWorld ? sqrt(facetVolVars.extrusionFactor()) : facetVolVars.extrusionFactor());
441 A[curLocalDofIdx][curLocalDofIdx] -= wFacet;
442 B[curLocalDofIdx][curLocalScvf.coupledFacetLocalDofIndex()] -= wFacet;
445 if (!isZeroWij && abs(wFacet) <= wijZeroThresh)
448 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
453 if (!curGlobalScvf.boundary() && !curIsDirichlet)
456 const Scalar curOneMinusXi = curIsInteriorBoundary ? -(1.0 - xi) : 1.0;
459 for (
unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
461 const auto idxOnScvf = idxInOutside+1;
462 const auto& negLocalScv = iv.localScv( neighborScvIndices[idxOnScvf] );
463 const auto& negGlobalScv = this->
fvGeometry().scv(negLocalScv.gridScvIndex());
464 const auto& negVolVars = this->
elemVolVars()[negGlobalScv];
465 const auto negTensor = getT(negVolVars);
468 const auto& scvf = dim < dimWorld ? this->
fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside)
470 wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility<EG>(negLocalScv, scvf, negTensor, negVolVars.extrusionFactor());
474 wijk[faceIdx][idxOnScvf] *= -1.0;
477 for (
int localDir = 0; localDir < dim; localDir++)
479 const auto otherLocalScvfIdx = negLocalScv.localScvfIndex(localDir);
480 const auto& otherLocalScvf = iv.localScvf(otherLocalScvfIdx);
481 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
484 if (otherLocalDofIdx == curLocalDofIdx && !isZeroWij)
485 if (abs(wijk[faceIdx][idxOnScvf][localDir]) <= wijZeroThresh)
486 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
488 if (!otherLocalScvf.isDirichlet())
489 A[curLocalDofIdx][otherLocalDofIdx] += curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
491 B[curLocalDofIdx][otherLocalDofIdx] -= curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
494 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.
A class that contains helper functions as well as functionality which is common to different mpfa sch...
The available mpfa schemes in Dumux.
Defines the general interface of classes used for the assembly of the local systems of equations invo...
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:849
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Defines the general interface of the local assembler classes for the assembly of the interaction volu...
Definition: localassemblerbase.hh:57
const FVElementGeometry & fvGeometry() const
Definition: localassemblerbase.hh:87
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:88
const Problem & problem() const
Definition: localassemblerbase.hh:86
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:55
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:99
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 i...
Definition: multidomain/facet/cellcentered/mpfa/localassembler.hh:81
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:131
This file contains free functions to evaluate the transmissibilities associated with flux evaluations...