3.1-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
31
32namespace Dumux {
33
43template< class P, class EG, class EV >
45: public InteractionVolumeAssemblerBase< P, EG, EV >
46{
49
50public:
52 using ParentType::ParentType;
53
63 template< class DataHandle, class IV, class TensorFunc >
64 void assembleMatrices(DataHandle& handle, IV& iv, const TensorFunc& getT)
65 {
66 assembleLocalMatrices_(handle.A(), handle.AB(), handle.CA(), handle.T(), handle.omegas(), iv, getT);
67
68 // maybe solve the local system
69 if (iv.numUnknowns() > 0)
70 Helper::solveLocalSystem(this->fvGeometry(), handle, iv);
71 }
72
81 template< class DataHandle, class IV, class GetU >
82 void assembleU(DataHandle& handle, const IV& iv, const GetU& getU)
83 {
84 auto& u = handle.uj();
85 Helper::resizeVector(u, iv.numKnowns());
86
87 // put the cell unknowns first, then Dirichlet values
88 typename IV::Traits::IndexSet::LocalIndexType i = 0;
89 for (; i < iv.numScvs(); i++)
90 u[i] = getU( this->elemVolVars()[iv.localScv(i).gridScvIndex()] );
91 for (const auto& data : iv.dirichletData())
92 u[i++] = getU( this->elemVolVars()[data.volVarIndex()] );
93 }
94
95private:
118 template< class IV, class TensorFunc >
119 void assembleLocalMatrices_(typename IV::Traits::MatVecTraits::AMatrix& A,
120 typename IV::Traits::MatVecTraits::BMatrix& B,
121 typename IV::Traits::MatVecTraits::CMatrix& C,
122 typename IV::Traits::MatVecTraits::DMatrix& D,
123 typename IV::Traits::MatVecTraits::OmegaStorage& wijk,
124 IV& iv, const TensorFunc& getT)
125 {
126 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
127 static constexpr int dim = IV::Traits::GridView::dimension;
128 static constexpr int dimWorld = IV::Traits::GridView::dimensionworld;
129
130 // resize omegas
131 Helper::resizeVector(wijk, iv.numFaces());
132
133 // if only Dirichlet faces are present in the iv,
134 // the matrices A, B & C are undefined and D = T
135 if (iv.numUnknowns() == 0)
136 {
137 // resize & reset D matrix
138 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
139
140 // Loop over all the faces, in this case these are all dirichlet boundaries
141 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
142 {
143 const auto& curLocalScvf = iv.localScvf(faceIdx);
144 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
145 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
146
147 // get tensor in "positive" sub volume
148 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
149 const auto& posGlobalScv = this->fvGeometry().scv(posLocalScv.gridScvIndex());
150 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
151 const auto& posElement = iv.element(neighborScvIndices[0]);
152 const auto tensor = getT(this->problem(), posElement, posVolVars, this->fvGeometry(), posGlobalScv);
153
154 // the omega factors of the "positive" sub volume
155 Helper::resizeVector(wijk[faceIdx], /*no outside scvs present*/1);
156 wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
157
158 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
159 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
160 {
161 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
162 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
163 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
164 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
165 }
166 }
167 }
168 else
169 {
170 // resize & reset matrices
171 Helper::resizeMatrix(A, iv.numUnknowns(), iv.numUnknowns()); A = 0.0;
172 Helper::resizeMatrix(B, iv.numUnknowns(), iv.numKnowns()); B = 0.0;
173 Helper::resizeMatrix(C, iv.numFaces(), iv.numUnknowns()); C = 0.0;
174 Helper::resizeMatrix(D, iv.numFaces(), iv.numKnowns()); D = 0.0;
175
176 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
177 {
178 const auto& curLocalScvf = iv.localScvf(faceIdx);
179 const auto& curGlobalScvf = this->fvGeometry().scvf(curLocalScvf.gridScvfIndex());
180 const auto curIsDirichlet = curLocalScvf.isDirichlet();
181 const auto curLocalDofIdx = curLocalScvf.localDofIndex();
182
183 // get tensor in "positive" sub volume
184 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
185 const auto& posLocalScv = iv.localScv(neighborScvIndices[0]);
186 const auto& posGlobalScv = this->fvGeometry().scv(posLocalScv.gridScvIndex());
187 const auto& posVolVars = this->elemVolVars()[posGlobalScv];
188 const auto& posElement = iv.element(neighborScvIndices[0]);
189 const auto tensor = getT(this->problem(), posElement, posVolVars, this->fvGeometry(), posGlobalScv);
190
191 // the omega factors of the "positive" sub volume
192 Helper::resizeVector(wijk[faceIdx], neighborScvIndices.size());
193 wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor());
194
195 // go over the coordinate directions in the positive sub volume
196 for (unsigned int localDir = 0; localDir < dim; localDir++)
197 {
198 const auto& otherLocalScvf = iv.localScvf( posLocalScv.localScvfIndex(localDir) );
199 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
200
201 // if we are not on a Dirichlet face, add entries associated with unknown face pressures
202 // i.e. in matrix C and maybe A (if current face is not a Dirichlet face)
203 if (!otherLocalScvf.isDirichlet())
204 {
205 C[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
206 if (!curIsDirichlet)
207 A[curLocalDofIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
208 }
209 // the current face is a Dirichlet face and creates entries in D & maybe B
210 else
211 {
212 D[faceIdx][otherLocalDofIdx] -= wijk[faceIdx][0][localDir];
213 if (!curIsDirichlet)
214 B[curLocalDofIdx][otherLocalDofIdx] += wijk[faceIdx][0][localDir];
215 }
216
217 // add entries related to pressures at the scv centers (dofs)
218 const auto posScvLocalDofIdx = posLocalScv.localDofIndex();
219 D[faceIdx][posScvLocalDofIdx] += wijk[faceIdx][0][localDir];
220
221 if (!curIsDirichlet)
222 B[curLocalDofIdx][posScvLocalDofIdx] -= wijk[faceIdx][0][localDir];
223 }
224
225 // If we are on an interior face, add values from negative sub volume
226 if (!curGlobalScvf.boundary())
227 {
228 // loop over all the outside neighbors of this face and add entries
229 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
230 {
231 const auto idxOnScvf = idxInOutside+1;
232 const auto& negLocalScv = iv.localScv( neighborScvIndices[idxOnScvf] );
233 const auto& negGlobalScv = this->fvGeometry().scv(negLocalScv.gridScvIndex());
234 const auto& negVolVars = this->elemVolVars()[negGlobalScv];
235 const auto& negElement = iv.element( neighborScvIndices[idxOnScvf] );
236 const auto negTensor = getT(this->problem(), negElement, negVolVars, this->fvGeometry(), negGlobalScv);
237
238 // On surface grids, use outside face for "negative" transmissibility calculation
239 const auto& scvf = dim < dimWorld ? this->fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside)
240 : curGlobalScvf;
241 wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility(negLocalScv, scvf, negTensor, negVolVars.extrusionFactor());
242
243 // flip sign on surface grids (since we used the "outside" normal)
244 if (dim < dimWorld)
245 wijk[faceIdx][idxOnScvf] *= -1.0;
246
247 // go over the coordinate directions in the negative sub volume
248 for (int localDir = 0; localDir < dim; localDir++)
249 {
250 const auto otherLocalScvfIdx = negLocalScv.localScvfIndex(localDir);
251 const auto& otherLocalScvf = iv.localScvf(otherLocalScvfIdx);
252 const auto otherLocalDofIdx = otherLocalScvf.localDofIndex();
253
254 if (!otherLocalScvf.isDirichlet())
255 A[curLocalDofIdx][otherLocalDofIdx] += wijk[faceIdx][idxOnScvf][localDir];
256 else
257 B[curLocalDofIdx][otherLocalDofIdx] -= wijk[faceIdx][idxOnScvf][localDir];
258
259 // add entries to matrix B
260 B[curLocalDofIdx][negLocalScv.localDofIndex()] += wijk[faceIdx][idxOnScvf][localDir];
261 }
262 }
263 }
264 }
265 }
266 }
267};
268
269} // end namespace
270
271#endif
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...
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
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: discretization/cellcentered/mpfa/omethod/localassembler.hh:46
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: discretization/cellcentered/mpfa/omethod/localassembler.hh:64
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:82
This file contains free functions to evaluate the transmissibilities associated with flux evaluations...