3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
localassemblerbase.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 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_LOCAL_ASSEMBLER_BASE_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_LOCAL_ASSEMBLER_BASE_HH
28
29#include <algorithm>
30#include <vector>
31#include <type_traits>
32
33#include <dune/common/exceptions.hh>
34#include <dune/common/reservedvector.hh>
35
38
39namespace Dumux {
40
53template< class P, class EG, class EV >
56{
57 using Problem = P;
58 using FVElementGeometry = EG;
59 using ElementVolumeVariables = EV;
60
62
63 public:
73 const FVElementGeometry& fvGeometry,
74 const ElementVolumeVariables& elemVolVars)
75 {
76 problemPtr_ = &problem;
77 fvGeometryPtr_ = &fvGeometry;
78 elemVolVarsPtr_ = &elemVolVars;
79 }
80
81 // return functions to the local views & problem
82 const Problem& problem() const { return *problemPtr_; }
83 const FVElementGeometry& fvGeometry() const { return *fvGeometryPtr_; }
84 const ElementVolumeVariables& elemVolVars() const { return *elemVolVarsPtr_; }
85
99 template< class DataHandle, class IV, class TensorFunc >
100 void assembleMatrices(DataHandle& handle, IV& iv, const TensorFunc& getT)
101 {
102 DUNE_THROW(Dune::NotImplemented, "Implementation does not provide a assembleMatrices() function");
103 }
104
116 template< class DataHandle, class IV, class GetU >
117 void assembleU(DataHandle& handle, const IV& iv, const GetU& getU)
118 {
119 DUNE_THROW(Dune::NotImplemented, "Implementation does not provide a assemble() function for the cell/Dirichlet unknowns");
120 }
121
130 template< class DataHandle, class IV, class GetRho >
131 void assembleGravity(DataHandle& handle, const IV& iv, const GetRho& getRho)
132 {
133 using GridView = typename IV::Traits::GridView;
134 static constexpr int dim = GridView::dimension;
135 static constexpr int dimWorld = GridView::dimensionworld;
136 static constexpr bool isSurfaceGrid = dim < dimWorld;
137
138 // resize the gravity vectors
139 auto& g = handle.g();
140 auto& deltaG = handle.deltaG();
141 auto& outsideG = handle.gOutside();
142 Helper::resizeVector(g, iv.numFaces());
143 Helper::resizeVector(deltaG, iv.numUnknowns());
144 if (isSurfaceGrid)
145 Helper::resizeVector(outsideG, iv.numFaces());
146
151 using Scalar = typename IV::Traits::MatVecTraits::TMatrix::value_type;
152 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
153
154 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
155 {
156 // gravitational acceleration on this face
157 const auto& curLocalScvf = iv.localScvf(faceIdx);
158 const auto& curGlobalScvf = fvGeometry().scvf(curLocalScvf.gridScvfIndex());
159 const auto& gravity = problem().spatialParams().gravity(curGlobalScvf.ipGlobal());
160
161 // get permeability tensor in "positive" sub volume
162 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
163 const auto& posGlobalScv = fvGeometry().scv(iv.localScv(neighborScvIndices[0]).gridScvIndex());
164 const auto& posVolVars = elemVolVars()[posGlobalScv];
165 const auto alpha_inside = posVolVars.extrusionFactor()*vtmv(curGlobalScvf.unitOuterNormal(),
166 posVolVars.permeability(),
167 gravity);
168
169 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
170 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
171 std::vector<Scalar>,
172 Dune::ReservedVector<Scalar, 1> >;
173 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
174 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
175 Scalar rho;
176
177 if (isSurfaceGrid)
178 Helper::resizeVector(outsideG[faceIdx], numOutsideFaces);
179
180 if (!curLocalScvf.isDirichlet())
181 {
182 const auto localDofIdx = curLocalScvf.localDofIndex();
183
184 rho = getRho(posVolVars);
185 deltaG[localDofIdx] = 0.0;
186
187 if (!curGlobalScvf.boundary())
188 {
189 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
190 {
191 // obtain outside tensor
192 const auto negLocalScvIdx = neighborScvIndices[idxInOutside+1];
193 const auto& negGlobalScv = fvGeometry().scv(iv.localScv(negLocalScvIdx).gridScvIndex());
194 const auto& negVolVars = elemVolVars()[negGlobalScv];
195 const auto& flipScvf = !isSurfaceGrid ? curGlobalScvf
196 : fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside);
197
198 alpha_outside[idxInOutside] = negVolVars.extrusionFactor()*vtmv(flipScvf.unitOuterNormal(),
199 negVolVars.permeability(),
200 gravity);
201 if (isSurfaceGrid)
202 alpha_outside[idxInOutside] *= -1.0;
203
204 rho += getRho(negVolVars);
205 deltaG[localDofIdx] += alpha_outside[idxInOutside];
206 }
207 }
208
209 rho /= numOutsideFaces + 1;
210 deltaG[localDofIdx] -= alpha_inside;
211 deltaG[localDofIdx] *= rho*curGlobalScvf.area();
212 }
213 // use density resulting from Dirichlet BCs
214 else
215 rho = getRho(elemVolVars()[curGlobalScvf.outsideScvIdx()]);
216
217 // add "inside" & "outside" alphas to gravity containers
218 g[faceIdx] = alpha_inside*rho*curGlobalScvf.area();
219
220 if (isSurfaceGrid)
221 {
222 unsigned int i = 0;
223 for (const auto& alpha : alpha_outside)
224 outsideG[faceIdx][i++] = alpha*rho*curGlobalScvf.area();
225 }
226 }
227
228 // add iv-wide contributions to gravity vectors
229 handle.CA().umv(deltaG, g);
230 if (isSurfaceGrid)
231 {
232 using FaceVector = typename IV::Traits::MatVecTraits::FaceVector;
233 FaceVector AG;
234 Helper::resizeVector(AG, iv.numUnknowns());
235 handle.A().mv(deltaG, AG);
236
237 // compute gravitational accelerations
238 for (const auto& localFaceData : iv.localFaceData())
239 {
240 // continue only for "outside" faces
241 if (!localFaceData.isOutsideFace())
242 continue;
243
244 const auto localScvIdx = localFaceData.ivLocalInsideScvIndex();
245 const auto localScvfIdx = localFaceData.ivLocalScvfIndex();
246 const auto idxInOutside = localFaceData.scvfLocalOutsideScvfIndex();
247 const auto& posLocalScv = iv.localScv(localScvIdx);
248 const auto& wijk = handle.omegas()[localScvfIdx][idxInOutside+1];
249
250 // add contributions from all local directions
251 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
252 {
253 // the scvf corresponding to this local direction in the scv
254 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
255 if (!curLocalScvf.isDirichlet())
256 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
257 }
258 }
259 }
260 }
261
262private:
263 // pointers to the data required for assembly
264 const Problem* problemPtr_;
265 const FVElementGeometry* fvGeometryPtr_;
266 const ElementVolumeVariables* elemVolVarsPtr_;
267};
268
269} // end namespace Dumux
270
271#endif
A helper function for class member function introspection.
A class that contains helper functions as well as functionality which is common to different mpfa sch...
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
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 m...
Definition: localassemblerbase.hh:100
const FVElementGeometry & fvGeometry() const
Definition: localassemblerbase.hh:83
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:84
InteractionVolumeAssemblerBase(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars)
The constructor. Sets pointers to the objects required for a subsequent call to assemble().
Definition: localassemblerbase.hh:72
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: localassemblerbase.hh:117
void assembleGravity(DataHandle &handle, const IV &iv, const GetRho &getRho)
Assembles the gravitational flux contributions on the scvfs within an interaction volume.
Definition: localassemblerbase.hh:131
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 resizeVector(Vector &v, size_type size)
resizes a vector to the given size (specialization for dynamic matrix type)
Definition: localassemblerhelper.hh:238