version 3.10-dev
discretization/cellcentered/mpfa/omethod/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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_ASSEMBLER_HH
14#define DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_ASSEMBLER_HH
15
16#include <algorithm>
17
21
22namespace Dumux {
23
33template< class P, class EG, class EV >
35: public InteractionVolumeAssemblerBase< P, EG, EV >
36{
39
40 template< class IV >
41 using Scalar = typename IV::Traits::MatVecTraits::FaceVector::value_type;
42
43public:
45 using ParentType::ParentType;
46
61 template< class DataHandle, class IV, class TensorFunc >
62 void assembleMatrices(DataHandle& handle, IV& iv, const TensorFunc& getT, Scalar<IV> wijZeroThresh = 0.0)
63 {
64 const auto zeroRows = assembleLocalMatrices_(handle.A(), handle.AB(), handle.CA(), handle.T(), handle.omegas(), iv, getT, wijZeroThresh);
65
66 // maybe solve the local system
67 if (iv.numUnknowns() > 0)
68 {
69 // for now, B is carried in what will be AB after solve
70 auto& B = handle.AB();
71 auto& A = handle.A();
72
73 // If the tensor is zero in some cells, we might have zero rows in the matrix.
74 // In this case we set the diagonal entry of A to 1.0 and the row of B to zero.
75 // This will ultimatively lead to zero transmissibilities for this face.
76 for (const auto& zeroRowIndices : zeroRows)
77 {
78 const auto zeroRowDofIdx = zeroRowIndices.first;
79 for (auto& row : A)
80 row[zeroRowDofIdx] = 0.0;
81 A[zeroRowDofIdx] = 0.0;
82 A[zeroRowDofIdx][zeroRowDofIdx] = 1.0;
83 B[zeroRowDofIdx] = 0.0;
84 }
85
86 Helper::solveLocalSystem(this->fvGeometry(), handle, iv);
87
88 // make sure the inverse of A now carries zeroes in the zero rows, as
89 // well as CA and T. AB will have the zeroes already because we set the rows
90 // of B to zero above.
91 for (const auto& zeroRowIndices : zeroRows)
92 {
93 const auto faceIdx = zeroRowIndices.second;
94 A[zeroRowIndices.first] = 0.0;
95 handle.CA()[faceIdx] = 0.0;
96 handle.T()[faceIdx] = 0.0;
97
98 // reset outside transmissibilities on surface grids
99 static constexpr int dim = IV::Traits::GridView::dimension;
100 static constexpr int dimWorld = IV::Traits::GridView::dimensionworld;
101 if constexpr (dim < dimWorld)
102 std::for_each( handle.tijOutside()[faceIdx].begin(),
103 handle.tijOutside()[faceIdx].end(),
104 [] (auto& outsideTij) { outsideTij = 0.0; } );
105 }
106 }
107 }
108
117 template< class DataHandle, class IV, class GetU >
118 void assembleU(DataHandle& handle, const IV& iv, const GetU& getU)
119 {
120 auto& u = handle.uj();
121 Helper::resizeVector(u, iv.numKnowns());
122
123 // put the cell unknowns first, then Dirichlet values
124 typename IV::Traits::IndexSet::LocalIndexType i = 0;
125 for (; i < iv.numScvs(); i++)
126 u[i] = getU( this->elemVolVars()[iv.localScv(i).gridScvIndex()] );
127 for (const auto& data : iv.dirichletData())
128 u[i++] = getU( this->elemVolVars()[data.volVarIndex()] );
129 }
130
131private:
162 template< class IV, class TensorFunc >
163 auto assembleLocalMatrices_(typename IV::Traits::MatVecTraits::AMatrix& A,
164 typename IV::Traits::MatVecTraits::BMatrix& B,
165 typename IV::Traits::MatVecTraits::CMatrix& C,
166 typename IV::Traits::MatVecTraits::DMatrix& D,
167 typename IV::Traits::MatVecTraits::OmegaStorage& wijk,
168 IV& iv, const TensorFunc& getT,
169 Scalar<IV> wijZeroThresh)
170 {
171 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
172 static constexpr int dim = IV::Traits::GridView::dimension;
173 static constexpr int dimWorld = IV::Traits::GridView::dimensionworld;
174
175 std::vector< std::pair<LocalIndexType, LocalIndexType> > faceMarkers;
176 Helper::resizeVector(wijk, iv.numFaces());
177
178 // if only Dirichlet faces are present in the iv,
179 // the matrices A, B & C are undefined and D = T
180 if (iv.numUnknowns() == 0)
181 {
182 // resize & reset D matrix
183 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
184
185 // Loop over all the faces, in this case these are all dirichlet boundaries
186 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
187 {
188 const auto& curLocalScvf = iv.localScvf(faceIdx);
189 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
190 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
191
192 // get tensor in "positive" sub volume
193 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
194 const auto& posGlobalScv = this->fvGeometry().scv(posLocalScv.gridScvIndex());
195 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
196 const auto tensor = getT(posVolVars);
197
198 // the omega factors of the "positive" sub volume
199 Helper::resizeVector(wijk[faceIdx], /*no outside scvs present*/1);
200 wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(this->fvGeometry(), posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
201
202 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
203 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
204 {
205 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
206 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
207 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
208 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
209 }
210 }
211 }
212 else
213 {
214 // resize & reset matrices
215 Helper::resizeMatrix(A, iv.numUnknowns(), iv.numUnknowns()); A = 0.0;
216 Helper::resizeMatrix(B, iv.numUnknowns(), iv.numKnowns()); B = 0.0;
217 Helper::resizeMatrix(C, iv.numFaces(), iv.numUnknowns()); C = 0.0;
218 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
219
220 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
221 {
222 const auto& curLocalScvf = iv.localScvf(faceIdx);
223 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
224 const auto curIsDirichlet = curLocalScvf.isDirichlet();
225 const auto curLocalDofIdx = curLocalScvf.localDofIndex();
226
227 // get tensor in "positive" sub volume
228 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
229 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
230 const auto& posGlobalScv = this->fvGeometry().scv(posLocalScv.gridScvIndex());
231 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
232 const auto tensor = getT(posVolVars);
233
234 // the omega factors of the "positive" sub volume
235 Helper::resizeVector(wijk[faceIdx], neighborScvIndices.size());
236 wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(this->fvGeometry(), posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
237
238 using std::abs;
239 bool insideZeroWij = false;
240
241 // go over the coordinate directions in the positive sub volume
242 for (unsigned int localDir = 0; localDir < dim; localDir++)
243 {
244 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
245 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
246
247 if (otherLocalDofIdx == curLocalDofIdx)
248 {
249 if (abs(wijk[faceIdx][0][localDir]) <= wijZeroThresh)
250 {
251 if (!curIsDirichlet)
252 {
253 insideZeroWij = true;
254 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
255 }
256 }
257 }
258
259 // if we are not on a Dirichlet face, add entries associated with unknown face pressures
260 // i.e. in matrix C and maybe A (if current face is not a Dirichlet face)
261 if (!otherLocalScvf.isDirichlet())
262 {
263 C[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
264 if (!curIsDirichlet)
265 A[curLocalDofIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
266 }
267 // the current face is a Dirichlet face and creates entries in D & maybe B
268 else
269 {
270 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
271 if (!curIsDirichlet)
272 B[curLocalDofIdx][otherLocalDofIdx] += wijk[faceIdx][0][localDir];
273 }
274
275 // add entries related to pressures at the scv centers (dofs)
276 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
277 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
278
279 if (!curIsDirichlet)
280 B[curLocalDofIdx][posScvLocalDofIdx] -= wijk[faceIdx][0][localDir];
281 }
282
283 // If we are on an interior face, add values from negative sub volume
284 if (!curGlobalScvf.boundary())
285 {
286 // loop over all the outside neighbors of this face and add entries
287 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
288 {
289 const auto idxOnScvf = idxInOutside+1;
290 const auto& negLocalScv = iv.localScv( neighborScvIndices[idxOnScvf] );
291 const auto& negGlobalScv = this->fvGeometry().scv(negLocalScv.gridScvIndex());
292 const auto& negVolVars = this->elemVolVars()[negGlobalScv];
293 const auto negTensor = getT(negVolVars);
294
295 // On surface grids, use outside face for "negative" transmissibility calculation
296 const auto& scvf = dim < dimWorld ? this->fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside)
297 : curGlobalScvf;
298 wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility<EG>(this->fvGeometry(), negLocalScv, scvf, negTensor, negVolVars.extrusionFactor());
299
300 // flip sign on surface grids (since we used the "outside" normal)
301 if (dim < dimWorld)
302 wijk[faceIdx][idxOnScvf] *= -1.0;
303
304 // go over the coordinate directions in the negative sub volume
305 for (int localDir = 0; localDir < dim; localDir++)
306 {
307 const auto otherLocalScvfIdx = negLocalScv.localScvfIndex(localDir);
308 const auto& otherLocalScvf = iv.localScvf(otherLocalScvfIdx);
309 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
310
311 // check for zero transmissibilities (skip if inside has been zero already)
312 if (otherLocalDofIdx == curLocalDofIdx && !insideZeroWij)
313 if (abs(wijk[faceIdx][idxOnScvf][localDir]) <= wijZeroThresh)
314 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
315
316 if (!otherLocalScvf.isDirichlet())
317 A[curLocalDofIdx][otherLocalDofIdx] += wijk[faceIdx][idxOnScvf][localDir];
318 else
319 B[curLocalDofIdx][otherLocalDofIdx] -= wijk[faceIdx][idxOnScvf][localDir];
320
321 // add entries to matrix B
322 B[curLocalDofIdx][negLocalScv.localDofIndex()] += wijk[faceIdx][idxOnScvf][localDir];
323 }
324 }
325 }
326 }
327 }
328
329 return faceMarkers;
330 }
331};
332
333} // end namespace
334
335#endif
Defines the general interface of the local assembler classes for the assembly of the interaction volu...
Definition: localassemblerbase.hh:45
const FVElementGeometry & fvGeometry() const
Definition: localassemblerbase.hh:75
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:76
A class that contains helper functions as well as functionality which is common to different mpfa sch...
Definition: localassemblerhelper.hh:34
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:70
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:210
static void resizeVector(Vector &v, size_type size)
resizes a vector to the given size (specialization for dynamic matrix type)
Definition: localassemblerhelper.hh:226
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: discretization/cellcentered/mpfa/omethod/localassembler.hh:36
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: discretization/cellcentered/mpfa/omethod/localassembler.hh:62
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: discretization/cellcentered/mpfa/omethod/localassembler.hh:118
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...
This file contains free functions to evaluate the transmissibilities associated with flux evaluations...
Definition: adapt.hh:17