3.2-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 <vector>
28
29#include <dune/common/exceptions.hh>
30#include <dune/common/float_cmp.hh>
31
32#include <dumux/common/math.hh>
38
39namespace Dumux {
40
51template< class P, class EG, class EV >
53: public InteractionVolumeAssemblerBase< P, EG, EV >
54{
57
58 template< class IV >
59 using Scalar = typename IV::Traits::MatVecTraits::FaceVector::value_type;
60
61public:
63 using ParentType::ParentType;
64
78 template< class DataHandle, class IV, class TensorFunc >
79 void assembleMatrices(DataHandle& handle, IV& iv, const TensorFunc& getT, Scalar<IV> wijZeroThresh = 0.0)
80 {
81 assembleLocalMatrices_(handle.A(), handle.AB(), handle.CA(), handle.T(), handle.omegas(), iv, getT, wijZeroThresh);
82
83 // maybe solve the local system
84 if (iv.numUnknowns() > 0)
85 Helper::solveLocalSystem(this->fvGeometry(), handle, iv);
86 }
87
96 template< class DataHandle, class IV, class GetU >
97 void assembleU(DataHandle& handle, const IV& iv, const GetU& getU)
98 {
99 auto& u = handle.uj();
100 Helper::resizeVector(u, iv.numKnowns());
101
102 // put the cell unknowns first ...
103 typename IV::Traits::IndexSet::LocalIndexType i = 0;
104 for (; i < iv.numScvs(); i++)
105 u[i] = getU( this->elemVolVars()[iv.localScv(i).gridScvIndex()] );
106
107 // ... then put facet unknowns ...
108 for (const auto& data : iv.interiorBoundaryData())
109 {
110 const auto& scvf = this->fvGeometry().scvf(data.scvfIndex());
111 const auto element = this->fvGeometry().gridGeometry().element(scvf.insideScvIdx());
112 u[i++] = getU( this->problem().couplingManager().getLowDimVolVars(element, scvf) );
113 }
114
115 // ... then put the Dirichlet unknowns
116 for (const auto& data : iv.dirichletData())
117 u[i++] = getU( this->elemVolVars()[data.volVarIndex()] );
118 }
119
128 template< class DataHandle, class IV, class GetRho >
129 void assembleGravity(DataHandle& handle, const IV& iv, const GetRho& getRho)
130 {
131 using GridView = typename IV::Traits::GridView;
132 static constexpr int dim = GridView::dimension;
133 static constexpr int dimWorld = GridView::dimensionworld;
134 static constexpr bool isSurfaceGrid = dim < dimWorld;
135
136 // resize the gravity vectors
137 auto& g = handle.g();
138 auto& deltaG = handle.deltaG();
139 auto& outsideG = handle.gOutside();
140 Helper::resizeVector(g, iv.numFaces());
141 Helper::resizeVector(deltaG, iv.numUnknowns());
142 if (isSurfaceGrid)
143 Helper::resizeVector(outsideG, iv.numFaces());
144
149 using Scalar = typename IV::Traits::MatVecTraits::TMatrix::value_type;
150 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
151
152 // xi factor for coupling conditions
153 static const Scalar xi = getParamFromGroup<Scalar>(this->problem().paramGroup(), "FacetCoupling.Xi", 1.0);
154
155 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
156 {
157 // gravitational acceleration on this face
158 const auto& curLocalScvf = iv.localScvf(faceIdx);
159 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
160 const auto& gravity = this->problem().spatialParams().gravity(curGlobalScvf.ipGlobal());
161 const auto curIsInteriorBoundary = curLocalScvf.isOnInteriorBoundary();
162 const Scalar curXiFactor = curIsInteriorBoundary ? (curGlobalScvf.boundary() ? 1.0 : xi) : 1.0;
163
164 // get permeability tensor in "positive" sub volume
165 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
166 const auto& posGlobalScv = this->fvGeometry().scv(iv.localScv(neighborScvIndices[0]).gridScvIndex());
167 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
168 const auto alpha_inside = posVolVars.extrusionFactor()*vtmv(curGlobalScvf.unitOuterNormal(),
169 posVolVars.permeability(),
170 gravity);
171
172 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
173 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
174 std::vector<Scalar>,
175 Dune::ReservedVector<Scalar, 1> >;
176 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
177 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
178 Scalar rho;
179
180 if (isSurfaceGrid && !curIsInteriorBoundary)
181 Helper::resizeVector(outsideG[faceIdx], numOutsideFaces);
182
183 if (!curLocalScvf.isDirichlet())
184 {
185 const auto localDofIdx = curLocalScvf.localDofIndex();
186 const Scalar curOneMinusXi = curIsInteriorBoundary ? -(1.0 - xi) : 1.0;
187
188 rho = getRho(posVolVars);
189 deltaG[localDofIdx] = 0.0;
190
191 if (!curGlobalScvf.boundary())
192 {
193 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
194 {
195 // obtain outside tensor
196 const auto negLocalScvIdx = neighborScvIndices[idxInOutside+1];
197 const auto& negGlobalScv = this->fvGeometry().scv(iv.localScv(negLocalScvIdx).gridScvIndex());
198 const auto& negVolVars = this->elemVolVars()[negGlobalScv];
199 const auto& flipScvf = !isSurfaceGrid ? curGlobalScvf
200 : this->fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside);
201
202 alpha_outside[idxInOutside] = negVolVars.extrusionFactor()*vtmv(flipScvf.unitOuterNormal(),
203 negVolVars.permeability(),
204 gravity);
205 if (isSurfaceGrid)
206 alpha_outside[idxInOutside] *= -1.0;
207
208 if (!curIsInteriorBoundary)
209 rho += getRho(negVolVars);
210
211 deltaG[localDofIdx] += curOneMinusXi*alpha_outside[idxInOutside];
212 }
213 }
214
215 if (curIsInteriorBoundary)
216 {
217 const auto posElement = this->fvGeometry().gridGeometry().element(posGlobalScv.elementIndex());
218 const auto& facetVolVars = this->problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
219 rho += getRho(facetVolVars);
220 rho /= 2.0;
221 const auto alphaFacet = posVolVars.extrusionFactor()*vtmv(curGlobalScvf.unitOuterNormal(),
222 facetVolVars.permeability(),
223 gravity);
224 deltaG[localDofIdx] += alphaFacet;
225 }
226 else
227 rho /= numOutsideFaces + 1;
228
229 deltaG[localDofIdx] -= curXiFactor*alpha_inside;
230 deltaG[localDofIdx] *= rho*curGlobalScvf.area();
231 }
232 // use average density between facet and cell density on interior Dirichlet boundaries
233 else if (curIsInteriorBoundary)
234 {
235 const auto posElement = this->fvGeometry().gridGeometry().element(posGlobalScv.elementIndex());
236 const auto& facetVolVars = this->problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
237 rho = 0.5*(getRho(facetVolVars) + getRho(posVolVars));
238 }
239 // use density resulting from Dirichlet BCs
240 else
241 rho = getRho(this->elemVolVars()[curGlobalScvf.outsideScvIdx()]);
242
243 // add "inside" & "outside" alphas to gravity containers
244 g[faceIdx] = alpha_inside*rho*curGlobalScvf.area();
245
246 if (isSurfaceGrid && !curIsInteriorBoundary)
247 {
248 unsigned int i = 0;
249 for (const auto& alpha : alpha_outside)
250 outsideG[faceIdx][i++] = alpha*rho*curGlobalScvf.area();
251 }
252 }
253
254 // add iv-wide contributions to gravity vectors
255 handle.CA().umv(deltaG, g);
256 if (isSurfaceGrid)
257 {
258 using FaceVector = typename IV::Traits::MatVecTraits::FaceVector;
259 FaceVector AG;
260 Helper::resizeVector(AG, iv.numUnknowns());
261 handle.A().mv(deltaG, AG);
262
263 // compute gravitational accelerations
264 for (const auto& localFaceData : iv.localFaceData())
265 {
266 // continue only for "outside" faces
267 if (!localFaceData.isOutsideFace())
268 continue;
269
270 const auto localScvIdx = localFaceData.ivLocalInsideScvIndex();
271 const auto localScvfIdx = localFaceData.ivLocalScvfIndex();
272 const auto idxInOutside = localFaceData.scvfLocalOutsideScvfIndex();
273 const auto& posLocalScv = iv.localScv(localScvIdx);
274 const auto& wijk = handle.omegas()[localScvfIdx][idxInOutside+1];
275
276 // add contributions from all local directions
277 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
278 {
279 // the scvf corresponding to this local direction in the scv
280 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
281 if (!curLocalScvf.isDirichlet())
282 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
283 }
284 }
285 }
286 }
287
288private:
292 template< class IV, class TensorFunc >
293 void assembleLocalMatrices_(typename IV::Traits::MatVecTraits::AMatrix& A,
294 typename IV::Traits::MatVecTraits::BMatrix& B,
295 typename IV::Traits::MatVecTraits::CMatrix& C,
296 typename IV::Traits::MatVecTraits::DMatrix& D,
297 typename IV::Traits::MatVecTraits::OmegaStorage& wijk,
298 IV& iv, const TensorFunc& getT,
299 Scalar<IV> wijZeroThresh)
300 {
301 using Scalar = typename IV::Traits::MatVecTraits::TMatrix::value_type;
302 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
303 static constexpr int dim = IV::Traits::GridView::dimension;
304 static constexpr int dimWorld = IV::Traits::GridView::dimensionworld;
305
306 // xi factor for coupling conditions
307 static const Scalar xi = getParamFromGroup<Scalar>(this->problem().paramGroup(), "FacetCoupling.Xi", 1.0);
308
309 // On surface grids only xi = 1.0 can be used, as the coupling condition
310 // for xi != 1.0 does not generalize for surface grids where the normal
311 // vectors of the inside/outside elements have different orientations.
312 if (dim < dimWorld)
313 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
314 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
315
316 std::vector< std::pair<LocalIndexType, LocalIndexType> > faceMarkers;
317 Helper::resizeVector(wijk, iv.numFaces());
318
319 // if only Dirichlet faces are present in the iv,
320 // the matrices A, B & C are undefined and D = T
321 if (iv.numUnknowns() == 0)
322 {
323 // resize & reset D matrix
324 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
325
326 // Loop over all the faces, in this case these are all dirichlet boundaries
327 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
328 {
329 const auto& curLocalScvf = iv.localScvf(faceIdx);
330 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
331 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
332
333 // get tensor in "positive" sub volume
334 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
335 const auto& posGlobalScv = this->fvGeometry().scv(posLocalScv.gridScvIndex());
336 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
337 const auto tensor = getT(posVolVars);
338
339 // the omega factors of the "positive" sub volume
340 Helper::resizeVector(wijk[faceIdx], /*no outside scvs present*/1);
341 wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
342
343 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
344 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
345 {
346 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
347 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
348 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
349 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
350 }
351 }
352 }
353 else
354 {
355 // resize & reset matrices
356 Helper::resizeMatrix(A, iv.numUnknowns(), iv.numUnknowns()); A = 0.0;
357 Helper::resizeMatrix(B, iv.numUnknowns(), iv.numKnowns()); B = 0.0;
358 Helper::resizeMatrix(C, iv.numFaces(), iv.numUnknowns()); C = 0.0;
359 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
360
361 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
362 {
363 const auto& curLocalScvf = iv.localScvf(faceIdx);
364 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
365 const auto curIsDirichlet = curLocalScvf.isDirichlet();
366 const auto curLocalDofIdx = curLocalScvf.localDofIndex();
367 const auto curIsInteriorBoundary = curLocalScvf.isOnInteriorBoundary();
368 const Scalar curXiFactor = curIsInteriorBoundary ? (curGlobalScvf.boundary() ? 1.0 : xi) : 1.0;
369
370 // get tensor in "positive" sub volume
371 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
372 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
373 const auto& posGlobalScv = this->fvGeometry().scv(posLocalScv.gridScvIndex());
374 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
375 const auto& posElement = iv.element(neighborScvIndices[0]);
376 const auto tensor = getT(posVolVars);
377
378 // the omega factors of the "positive" sub volume
379 Helper::resizeVector(wijk[faceIdx], neighborScvIndices.size());
380 wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
381
382 using std::abs;
383 bool isZeroWij = false;
384
385 // go over the coordinate directions in the positive sub volume
386 for (unsigned int localDir = 0; localDir < dim; localDir++)
387 {
388 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
389 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
390
391 if (otherLocalDofIdx == curLocalDofIdx)
392 {
393 if (abs(wijk[faceIdx][0][localDir]) <= wijZeroThresh)
394 {
395 if (!curIsDirichlet)
396 {
397 isZeroWij = true;
398 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
399 }
400 }
401 }
402
403 // if we are not on a Dirichlet face, add entries associated with unknown face pressures
404 // i.e. in matrix C and maybe A (if current face is not a Dirichlet face)
405 if (!otherLocalScvf.isDirichlet())
406 {
407 C[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
408 if (!curIsDirichlet)
409 A[curLocalDofIdx][otherLocalDofIdx] -= curXiFactor*wijk[faceIdx][0][localDir];
410 }
411 // the current face is a Dirichlet face and creates entries in D & maybe B
412 else
413 {
414 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
415 if (!curIsDirichlet)
416 B[curLocalDofIdx][otherLocalDofIdx] += curXiFactor*wijk[faceIdx][0][localDir];
417 }
418
419 // add entries related to pressures at the scv centers (dofs)
420 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
421 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
422
423 if (!curIsDirichlet)
424 B[curLocalDofIdx][posScvLocalDofIdx] -= curXiFactor*wijk[faceIdx][0][localDir];
425 }
426
427 // handle the entries related to the coupled facet element
428 if (curIsInteriorBoundary && !curIsDirichlet)
429 {
430 const auto facetVolVars = this->problem().couplingManager().getLowDimVolVars(posElement, curGlobalScvf);
431 const auto facetTensor = getT(facetVolVars);
432
433 // On surface grids we use the square root of the extrusion factor as approximation of the aperture
434 using std::sqrt;
435 const auto wFacet = 2.0*curGlobalScvf.area()*posVolVars.extrusionFactor()
436 *vtmv(curGlobalScvf.unitOuterNormal(), facetTensor, curGlobalScvf.unitOuterNormal())
437 / (dim < dimWorld ? sqrt(facetVolVars.extrusionFactor()) : facetVolVars.extrusionFactor());
438
439 A[curLocalDofIdx][curLocalDofIdx] -= wFacet;
440 B[curLocalDofIdx][curLocalScvf.coupledFacetLocalDofIndex()] -= wFacet;
441
442 // check for zero transmissibilities (skip if inside has been zero already)
443 if (!isZeroWij && abs(wFacet) <= wijZeroThresh)
444 {
445 isZeroWij = true;
446 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
447 }
448 }
449
450 // If we are on an interior face (which isn't coupled via Dirichlet Bs), add values from negative sub volume
451 if (!curGlobalScvf.boundary() && !curIsDirichlet)
452 {
453 // the minus ensures the right sign of the coupling condition in the matrices
454 const Scalar curOneMinusXi = curIsInteriorBoundary ? -(1.0 - xi) : 1.0;
455
456 // loop over all the outside neighbors of this face and add entries
457 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
458 {
459 const auto idxOnScvf = idxInOutside+1;
460 const auto& negLocalScv = iv.localScv( neighborScvIndices[idxOnScvf] );
461 const auto& negGlobalScv = this->fvGeometry().scv(negLocalScv.gridScvIndex());
462 const auto& negVolVars = this->elemVolVars()[negGlobalScv];
463 const auto negTensor = getT(negVolVars);
464
465 // On surface grids, use outside face for "negative" transmissibility calculation
466 const auto& scvf = dim < dimWorld ? this->fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside)
467 : curGlobalScvf;
468 wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility(negLocalScv, scvf, negTensor, negVolVars.extrusionFactor());
469
470 // flip sign on surface grids (since we used the "outside" normal)
471 if (dim < dimWorld)
472 wijk[faceIdx][idxOnScvf] *= -1.0;
473
474 // go over the coordinate directions in the negative sub volume
475 for (int localDir = 0; localDir < dim; localDir++)
476 {
477 const auto otherLocalScvfIdx = negLocalScv.localScvfIndex(localDir);
478 const auto& otherLocalScvf = iv.localScvf(otherLocalScvfIdx);
479 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
480
481 // check for zero transmissibilities (skip if inside has been zero already)
482 if (otherLocalDofIdx == curLocalDofIdx && !isZeroWij)
483 if (abs(wijk[faceIdx][idxOnScvf][localDir]) <= wijZeroThresh)
484 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
485
486 if (!otherLocalScvf.isDirichlet())
487 A[curLocalDofIdx][otherLocalDofIdx] += curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
488 else
489 B[curLocalDofIdx][otherLocalDofIdx] -= curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
490
491 // add entries to matrix B
492 B[curLocalDofIdx][negLocalScv.localDofIndex()] += curOneMinusXi*wijk[faceIdx][idxOnScvf][localDir];
493 }
494 }
495 }
496 }
497 }
498 }
499};
500
501} // end namespace
502
503#endif
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...
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
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:85
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:86
const Problem & problem() const
Definition: localassemblerbase.hh:84
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:54
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:97
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:79
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:129
This file contains free functions to evaluate the transmissibilities associated with flux evaluations...