3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 3 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 *****************************************************************************/
25#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_ASSEMBLER_HH
26#define DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_ASSEMBLER_HH
27
28#include <algorithm>
29
33
34namespace Dumux {
35
45template< class P, class EG, class EV >
47: public InteractionVolumeAssemblerBase< P, EG, EV >
48{
51
52 template< class IV >
53 using Scalar = typename IV::Traits::MatVecTraits::FaceVector::value_type;
54
55public:
57 using ParentType::ParentType;
58
73 template< class DataHandle, class IV, class TensorFunc >
74 void assembleMatrices(DataHandle& handle, IV& iv, const TensorFunc& getT, Scalar<IV> wijZeroThresh = 0.0)
75 {
76 const auto zeroRows = assembleLocalMatrices_(handle.A(), handle.AB(), handle.CA(), handle.T(), handle.omegas(), iv, getT, wijZeroThresh);
77
78 // maybe solve the local system
79 if (iv.numUnknowns() > 0)
80 {
81 // for now, B is carried in what will be AB after solve
82 auto& B = handle.AB();
83 auto& A = handle.A();
84
85 // If the tensor is zero in some cells, we might have zero rows in the matrix.
86 // In this case we set the diagonal entry of A to 1.0 and the row of B to zero.
87 // This will ultimatively lead to zero transmissibilities for this face.
88 for (const auto& zeroRowIndices : zeroRows)
89 {
90 const auto zeroRowDofIdx = zeroRowIndices.first;
91 for (auto& row : A)
92 row[zeroRowDofIdx] = 0.0;
93 A[zeroRowDofIdx] = 0.0;
94 A[zeroRowDofIdx][zeroRowDofIdx] = 1.0;
95 B[zeroRowDofIdx] = 0.0;
96 }
97
98 Helper::solveLocalSystem(this->fvGeometry(), handle, iv);
99
100 // make sure the inverse of A now carries zeroes in the zero rows, as
101 // well as CA and T. AB will have the zeroes already because we set the rows
102 // of B to zero above.
103 for (const auto& zeroRowIndices : zeroRows)
104 {
105 const auto faceIdx = zeroRowIndices.second;
106 A[zeroRowIndices.first] = 0.0;
107 handle.CA()[faceIdx] = 0.0;
108 handle.T()[faceIdx] = 0.0;
109
110 // reset outside transmissibilities on surface grids
111 static constexpr int dim = IV::Traits::GridView::dimension;
112 static constexpr int dimWorld = IV::Traits::GridView::dimensionworld;
113 if constexpr (dim < dimWorld)
114 std::for_each( handle.tijOutside()[faceIdx].begin(),
115 handle.tijOutside()[faceIdx].end(),
116 [] (auto& outsideTij) { outsideTij = 0.0; } );
117 }
118 }
119 }
120
129 template< class DataHandle, class IV, class GetU >
130 void assembleU(DataHandle& handle, const IV& iv, const GetU& getU)
131 {
132 auto& u = handle.uj();
133 Helper::resizeVector(u, iv.numKnowns());
134
135 // put the cell unknowns first, then Dirichlet values
136 typename IV::Traits::IndexSet::LocalIndexType i = 0;
137 for (; i < iv.numScvs(); i++)
138 u[i] = getU( this->elemVolVars()[iv.localScv(i).gridScvIndex()] );
139 for (const auto& data : iv.dirichletData())
140 u[i++] = getU( this->elemVolVars()[data.volVarIndex()] );
141 }
142
143private:
174 template< class IV, class TensorFunc >
175 auto assembleLocalMatrices_(typename IV::Traits::MatVecTraits::AMatrix& A,
176 typename IV::Traits::MatVecTraits::BMatrix& B,
177 typename IV::Traits::MatVecTraits::CMatrix& C,
178 typename IV::Traits::MatVecTraits::DMatrix& D,
179 typename IV::Traits::MatVecTraits::OmegaStorage& wijk,
180 IV& iv, const TensorFunc& getT,
181 Scalar<IV> wijZeroThresh)
182 {
183 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
184 static constexpr int dim = IV::Traits::GridView::dimension;
185 static constexpr int dimWorld = IV::Traits::GridView::dimensionworld;
186
187 std::vector< std::pair<LocalIndexType, LocalIndexType> > faceMarkers;
188 Helper::resizeVector(wijk, iv.numFaces());
189
190 // if only Dirichlet faces are present in the iv,
191 // the matrices A, B & C are undefined and D = T
192 if (iv.numUnknowns() == 0)
193 {
194 // resize & reset D matrix
195 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
196
197 // Loop over all the faces, in this case these are all dirichlet boundaries
198 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
199 {
200 const auto& curLocalScvf = iv.localScvf(faceIdx);
201 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
202 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
203
204 // get tensor in "positive" sub volume
205 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
206 const auto& posGlobalScv = this->fvGeometry().scv(posLocalScv.gridScvIndex());
207 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
208 const auto tensor = getT(posVolVars);
209
210 // the omega factors of the "positive" sub volume
211 Helper::resizeVector(wijk[faceIdx], /*no outside scvs present*/1);
212 wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(this->fvGeometry(), posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
213
214 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
215 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
216 {
217 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
218 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
219 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
220 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
221 }
222 }
223 }
224 else
225 {
226 // resize & reset matrices
227 Helper::resizeMatrix(A, iv.numUnknowns(), iv.numUnknowns()); A = 0.0;
228 Helper::resizeMatrix(B, iv.numUnknowns(), iv.numKnowns()); B = 0.0;
229 Helper::resizeMatrix(C, iv.numFaces(), iv.numUnknowns()); C = 0.0;
230 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
231
232 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
233 {
234 const auto& curLocalScvf = iv.localScvf(faceIdx);
235 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
236 const auto curIsDirichlet = curLocalScvf.isDirichlet();
237 const auto curLocalDofIdx = curLocalScvf.localDofIndex();
238
239 // get tensor in "positive" sub volume
240 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
241 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
242 const auto& posGlobalScv = this->fvGeometry().scv(posLocalScv.gridScvIndex());
243 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
244 const auto tensor = getT(posVolVars);
245
246 // the omega factors of the "positive" sub volume
247 Helper::resizeVector(wijk[faceIdx], neighborScvIndices.size());
248 wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(this->fvGeometry(), posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
249
250 using std::abs;
251 bool insideZeroWij = false;
252
253 // go over the coordinate directions in the positive sub volume
254 for (unsigned int localDir = 0; localDir < dim; localDir++)
255 {
256 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
257 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
258
259 if (otherLocalDofIdx == curLocalDofIdx)
260 {
261 if (abs(wijk[faceIdx][0][localDir]) <= wijZeroThresh)
262 {
263 if (!curIsDirichlet)
264 {
265 insideZeroWij = true;
266 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
267 }
268 }
269 }
270
271 // if we are not on a Dirichlet face, add entries associated with unknown face pressures
272 // i.e. in matrix C and maybe A (if current face is not a Dirichlet face)
273 if (!otherLocalScvf.isDirichlet())
274 {
275 C[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
276 if (!curIsDirichlet)
277 A[curLocalDofIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
278 }
279 // the current face is a Dirichlet face and creates entries in D & maybe B
280 else
281 {
282 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
283 if (!curIsDirichlet)
284 B[curLocalDofIdx][otherLocalDofIdx] += wijk[faceIdx][0][localDir];
285 }
286
287 // add entries related to pressures at the scv centers (dofs)
288 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
289 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
290
291 if (!curIsDirichlet)
292 B[curLocalDofIdx][posScvLocalDofIdx] -= wijk[faceIdx][0][localDir];
293 }
294
295 // If we are on an interior face, add values from negative sub volume
296 if (!curGlobalScvf.boundary())
297 {
298 // loop over all the outside neighbors of this face and add entries
299 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
300 {
301 const auto idxOnScvf = idxInOutside+1;
302 const auto& negLocalScv = iv.localScv( neighborScvIndices[idxOnScvf] );
303 const auto& negGlobalScv = this->fvGeometry().scv(negLocalScv.gridScvIndex());
304 const auto& negVolVars = this->elemVolVars()[negGlobalScv];
305 const auto negTensor = getT(negVolVars);
306
307 // On surface grids, use outside face for "negative" transmissibility calculation
308 const auto& scvf = dim < dimWorld ? this->fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside)
309 : curGlobalScvf;
310 wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility<EG>(this->fvGeometry(), negLocalScv, scvf, negTensor, negVolVars.extrusionFactor());
311
312 // flip sign on surface grids (since we used the "outside" normal)
313 if (dim < dimWorld)
314 wijk[faceIdx][idxOnScvf] *= -1.0;
315
316 // go over the coordinate directions in the negative sub volume
317 for (int localDir = 0; localDir < dim; localDir++)
318 {
319 const auto otherLocalScvfIdx = negLocalScv.localScvfIndex(localDir);
320 const auto& otherLocalScvf = iv.localScvf(otherLocalScvfIdx);
321 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
322
323 // check for zero transmissibilities (skip if inside has been zero already)
324 if (otherLocalDofIdx == curLocalDofIdx && !insideZeroWij)
325 if (abs(wijk[faceIdx][idxOnScvf][localDir]) <= wijZeroThresh)
326 faceMarkers.emplace_back( std::make_pair(curLocalDofIdx, faceIdx) );
327
328 if (!otherLocalScvf.isDirichlet())
329 A[curLocalDofIdx][otherLocalDofIdx] += wijk[faceIdx][idxOnScvf][localDir];
330 else
331 B[curLocalDofIdx][otherLocalDofIdx] -= wijk[faceIdx][idxOnScvf][localDir];
332
333 // add entries to matrix B
334 B[curLocalDofIdx][negLocalScv.localDofIndex()] += wijk[faceIdx][idxOnScvf][localDir];
335 }
336 }
337 }
338 }
339 }
340
341 return faceMarkers;
342 }
343};
344
345} // end namespace
346
347#endif
A class that contains helper functions as well as functionality which is common to different mpfa sch...
Defines the general interface of classes used for the assembly of the local systems of equations invo...
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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
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: discretization/cellcentered/mpfa/omethod/localassembler.hh:48
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:74
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:130
This file contains free functions to evaluate the transmissibilities associated with flux evaluations...