3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/cellcentered/mpfa/localassembler.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 2 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_LOCAL_ASSEMBLER_HH
25#define DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_LOCAL_ASSEMBLER_HH
26
27#include <dune/common/exceptions.hh>
28#include <dune/common/float_cmp.hh>
29
30#include <dumux/common/math.hh>
36
37namespace Dumux {
38
49template< class P, class EG, class EV >
51: public InteractionVolumeAssemblerBase< P, EG, EV >
52{
55
56public:
58 using ParentType::ParentType;
59
68 template< class DataHandle, class IV, class TensorFunc >
69 void assembleMatrices(DataHandle& handle, IV& iv, const TensorFunc& getT)
70 {
71 assembleLocalMatrices_(handle.A(), handle.AB(), handle.CA(), handle.T(), handle.omegas(), iv, getT);
72
73 // maybe solve the local system
74 if (iv.numUnknowns() > 0)
75 Helper::solveLocalSystem(this->fvGeometry(), handle, iv);
76 }
77
86 template< class DataHandle, class IV, class GetU >
87 void assembleU(DataHandle& handle, const IV& iv, const GetU& getU)
88 {
89 auto& u = handle.uj();
90 Helper::resizeVector(u, iv.numKnowns());
91
92 // put the cell unknowns first ...
93 typename IV::Traits::IndexSet::LocalIndexType i = 0;
94 for (; i < iv.numScvs(); i++)
95 u[i] = getU( this->elemVolVars()[iv.localScv(i).gridScvIndex()] );
96
97 // ... then put facet unknowns ...
98 for (const auto& data : iv.interiorBoundaryData())
99 {
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) );
103 }
104
105 // ... then put the Dirichlet unknowns
106 for (const auto& data : iv.dirichletData())
107 u[i++] = getU( this->elemVolVars()[data.volVarIndex()] );
108 }
109
118 template< class DataHandle, class IV, class GetRho >
119 void assembleGravity(DataHandle& handle, const IV& iv, const GetRho& getRho)
120 {
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;
125
126 // resize the gravity vectors
127 auto& g = handle.g();
128 auto& deltaG = handle.deltaG();
129 auto& outsideG = handle.gOutside();
130 Helper::resizeVector(g, iv.numFaces());
131 Helper::resizeVector(deltaG, iv.numUnknowns());
132 if (isSurfaceGrid)
133 Helper::resizeVector(outsideG, iv.numFaces());
134
139 using Scalar = typename IV::Traits::MatVecTraits::TMatrix::value_type;
140 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
141
142 // xi factor for coupling conditions
143 static const Scalar xi = getParamFromGroup<Scalar>(this->problem().paramGroup(), "FacetCoupling.Xi", 1.0);
144
145 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
146 {
147 // gravitational acceleration on this face
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;
153
154 // get permeability tensor in "positive" sub volume
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(),
160 gravity);
161
162 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
163 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
164 std::vector<Scalar>,
165 Dune::ReservedVector<Scalar, 1> >;
166 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
167 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
168 Scalar rho;
169
170 if (isSurfaceGrid && !curIsInteriorBoundary)
171 Helper::resizeVector(outsideG[faceIdx], numOutsideFaces);
172
173 if (!curLocalScvf.isDirichlet())
174 {
175 const auto localDofIdx = curLocalScvf.localDofIndex();
176 const Scalar curOneMinusXi = curIsInteriorBoundary ? -(1.0 - xi) : 1.0;
177
178 rho = getRho(posVolVars);
179 deltaG[localDofIdx] = 0.0;
180
181 if (!curGlobalScvf.boundary())
182 {
183 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
184 {
185 // obtain outside tensor
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);
191
192 alpha_outside[idxInOutside] = negVolVars.extrusionFactor()*vtmv(flipScvf.unitOuterNormal(),
193 negVolVars.permeability(),
194 gravity);
195 if (isSurfaceGrid)
196 alpha_outside[idxInOutside] *= -1.0;
197
198 if (!curIsInteriorBoundary)
199 rho += getRho(negVolVars);
200
201 deltaG[localDofIdx] += curOneMinusXi*alpha_outside[idxInOutside];
202 }
203 }
204
205 if (curIsInteriorBoundary)
206 {
207 const auto posElement = this->fvGeometry().gridGeometry().element(posGlobalScv.elementIndex());
208 const auto& facetVolVars = this->problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
209 rho += getRho(facetVolVars);
210 rho /= 2.0;
211 const auto alphaFacet = posVolVars.extrusionFactor()*vtmv(curGlobalScvf.unitOuterNormal(),
212 facetVolVars.permeability(),
213 gravity);
214 deltaG[localDofIdx] += alphaFacet;
215 }
216 else
217 rho /= numOutsideFaces + 1;
218
219 deltaG[localDofIdx] -= curXiFactor*alpha_inside;
220 deltaG[localDofIdx] *= rho*curGlobalScvf.area();
221 }
222 // use average density between facet and cell density on interior Dirichlet boundaries
223 else if (curIsInteriorBoundary)
224 {
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));
228 }
229 // use density resulting from Dirichlet BCs
230 else
231 rho = getRho(this->elemVolVars()[curGlobalScvf.outsideScvIdx()]);
232
233 // add "inside" & "outside" alphas to gravity containers
234 g[faceIdx] = alpha_inside*rho*curGlobalScvf.area();
235
236 if (isSurfaceGrid && !curIsInteriorBoundary)
237 {
238 unsigned int i = 0;
239 for (const auto& alpha : alpha_outside)
240 outsideG[faceIdx][i++] = alpha*rho*curGlobalScvf.area();
241 }
242 }
243
244 // add iv-wide contributions to gravity vectors
245 handle.CA().umv(deltaG, g);
246 if (isSurfaceGrid)
247 {
248 using FaceVector = typename IV::Traits::MatVecTraits::FaceVector;
249 FaceVector AG;
250 Helper::resizeVector(AG, iv.numUnknowns());
251 handle.A().mv(deltaG, AG);
252
253 // compute gravitational accelerations
254 for (const auto& localFaceData : iv.localFaceData())
255 {
256 // continue only for "outside" faces
257 if (!localFaceData.isOutsideFace())
258 continue;
259
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];
265
266 // add contributions from all local directions
267 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
268 {
269 // the scvf corresponding to this local direction in the scv
270 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
271 if (!curLocalScvf.isDirichlet())
272 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
273 }
274 }
275 }
276 }
277
278private:
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)
289 {
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;
294
295 // xi factor for coupling conditions
296 static const Scalar xi = getParamFromGroup<Scalar>(this->problem().paramGroup(), "FacetCoupling.Xi", 1.0);
297
298 // On surface grids only xi = 1.0 can be used, as the coupling condition
299 // for xi != 1.0 does not generalize for surface grids where the normal
300 // vectors of the inside/outside elements have different orientations.
301 if (dim < dimWorld)
302 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
303 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
304
305 // resize omegas
306 Helper::resizeVector(wijk, iv.numFaces());
307
308 // if only Dirichlet faces are present in the iv,
309 // the matrices A, B & C are undefined and D = T
310 if (iv.numUnknowns() == 0)
311 {
312 // resize & reset D matrix
313 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
314
315 // Loop over all the faces, in this case these are all dirichlet boundaries
316 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
317 {
318 const auto& curLocalScvf = iv.localScvf(faceIdx);
319 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
320 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
321
322 // get tensor in "positive" sub volume
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);
328
329 // the omega factors of the "positive" sub volume
330 Helper::resizeVector(wijk[faceIdx], /*no outside scvs present*/1);
331 wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
332
333 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
334 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
335 {
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];
340 }
341 }
342 }
343 else
344 {
345 // resize & reset matrices
346 Helper::resizeMatrix(A, iv.numUnknowns(), iv.numUnknowns()); A = 0.0;
347 Helper::resizeMatrix(B, iv.numUnknowns(), iv.numKnowns()); B = 0.0;
348 Helper::resizeMatrix(C, iv.numFaces(), iv.numUnknowns()); C = 0.0;
349 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
350
351 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
352 {
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;
359
360 // get tensor in "positive" sub volume
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);
367
368 // the omega factors of the "positive" sub volume
369 Helper::resizeVector(wijk[faceIdx], neighborScvIndices.size());
370 wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
371
372 // go over the coordinate directions in the positive sub volume
373 for (unsigned int localDir = 0; localDir < dim; localDir++)
374 {
375 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
376 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
377
378 // if we are not on a Dirichlet face, add entries associated with unknown face pressures
379 // i.e. in matrix C and maybe A (if current face is not a Dirichlet face)
380 if (!otherLocalScvf.isDirichlet())
381 {
382 C[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
383 if (!curIsDirichlet)
384 A[curLocalDofIdx][otherLocalDofIdx] -= curXiFactor*wijk[faceIdx][0][localDir];
385 }
386 // the current face is a Dirichlet face and creates entries in D & maybe B
387 else
388 {
389 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
390 if (!curIsDirichlet)
391 B[curLocalDofIdx][otherLocalDofIdx] += curXiFactor*wijk[faceIdx][0][localDir];
392 }
393
394 // add entries related to pressures at the scv centers (dofs)
395 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
396 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
397
398 if (!curIsDirichlet)
399 B[curLocalDofIdx][posScvLocalDofIdx] -= curXiFactor*wijk[faceIdx][0][localDir];
400 }
401
402 // handle the entries related to the coupled facet element
403 if (curIsInteriorBoundary && !curIsDirichlet)
404 {
405 const auto facetVolVars = this->problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
406 const auto facetTensor = this->problem().couplingManager().getLowDimTensor(posElement, curGlobalScvf, getT);
407
408 // On surface grids we use the square root of the extrusion factor as approximation of the aperture
409 using std::sqrt;
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());
413
414 A[curLocalDofIdx][curLocalDofIdx] -= wFacet;
415 B[curLocalDofIdx][curLocalScvf.coupledFacetLocalDofIndex()] -= wFacet;
416 }
417
418 // If we are on an interior face (which isn't coupled via Dirichlet Bs), add values from negative sub volume
419 if (!curGlobalScvf.boundary() && !curIsDirichlet)
420 {
421 // the minus ensures the right sign of the coupling condition in the matrices
422 const Scalar curOneMinusXi = curIsInteriorBoundary ? -(1.0 - xi) : 1.0;
423
424 // loop over all the outside neighbors of this face and add entries
425 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
426 {
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);
433
434 // On surface grids, use outside face for "negative" transmissibility calculation
435 const auto& scvf = dim < dimWorld ? this->fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside)
436 : curGlobalScvf;
437 wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility(negLocalScv, scvf, negTensor, negVolVars.extrusionFactor());
438
439 // flip sign on surface grids (since we used the "outside" normal)
440 if (dim < dimWorld)
441 wijk[faceIdx][idxOnScvf] *= -1.0;
442
443 // go over the coordinate directions in the negative sub volume
444 for (int localDir = 0; localDir < dim; localDir++)
445 {
446 const auto otherLocalScvfIdx = negLocalScv.localScvfIndex(localDir);
447 const auto& otherLocalScvf = iv.localScvf(otherLocalScvfIdx);
448 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
449
450 if (!otherLocalScvf.isDirichlet())
451 A[curLocalDofIdx][otherLocalDofIdx] += curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
452 else
453 B[curLocalDofIdx][otherLocalDofIdx] -= curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
454
455 // add entries to matrix B
456 B[curLocalDofIdx][negLocalScv.localDofIndex()] += curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
457 }
458 }
459 }
460 }
461 }
462 }
463};
464
465} // end namespace
466
467#endif
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...