3.2-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;
61
62 template< class IV >
63 using Scalar = typename IV::Traits::MatVecTraits::FaceVector::value_type;
64
65 public:
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars)
77 {
78 problemPtr_ = &problem;
79 fvGeometryPtr_ = &fvGeometry;
80 elemVolVarsPtr_ = &elemVolVars;
81 }
82
83 // return functions to the local views & problem
84 const Problem& problem() const { return *problemPtr_; }
85 const FVElementGeometry& fvGeometry() const { return *fvGeometryPtr_; }
86 const ElementVolumeVariables& elemVolVars() const { return *elemVolVarsPtr_; }
87
103 template< class DataHandle, class IV, class TensorFunc >
104 void assembleMatrices(DataHandle& handle, IV& iv, const TensorFunc& getT, Scalar<IV> wijZeroThresh = 0.0)
105 {
106 DUNE_THROW(Dune::NotImplemented, "Implementation does not provide an assembleMatrices() function");
107 }
108
120 template< class DataHandle, class IV, class GetU >
121 void assembleU(DataHandle& handle, const IV& iv, const GetU& getU)
122 {
123 DUNE_THROW(Dune::NotImplemented, "Implementation does not provide an assemble() function for the cell/Dirichlet unknowns");
124 }
125
134 template< class DataHandle, class IV, class GetRho >
135 void assembleGravity(DataHandle& handle, const IV& iv, const GetRho& getRho)
136 {
137 using GridView = typename IV::Traits::GridView;
138 static constexpr int dim = GridView::dimension;
139 static constexpr int dimWorld = GridView::dimensionworld;
140 static constexpr bool isSurfaceGrid = dim < dimWorld;
141
142 // resize the gravity vectors
143 auto& g = handle.g();
144 auto& deltaG = handle.deltaG();
145 auto& outsideG = handle.gOutside();
146 Helper::resizeVector(g, iv.numFaces());
147 Helper::resizeVector(deltaG, iv.numUnknowns());
148 if (isSurfaceGrid)
149 Helper::resizeVector(outsideG, iv.numFaces());
150
155 using Scalar = typename IV::Traits::MatVecTraits::TMatrix::value_type;
156 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
157
158 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
159 {
160 // gravitational acceleration on this face
161 const auto& curLocalScvf = iv.localScvf(faceIdx);
162 const auto& curGlobalScvf = fvGeometry().scvf(curLocalScvf.gridScvfIndex());
163 const auto& gravity = problem().spatialParams().gravity(curGlobalScvf.ipGlobal());
164
165 // get permeability tensor in "positive" sub volume
166 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
167 const auto& posGlobalScv = fvGeometry().scv(iv.localScv(neighborScvIndices[0]).gridScvIndex());
168 const auto& posVolVars = elemVolVars()[posGlobalScv];
169 const auto alpha_inside = posVolVars.extrusionFactor()*vtmv(curGlobalScvf.unitOuterNormal(),
170 posVolVars.permeability(),
171 gravity);
172
173 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
174 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
175 std::vector<Scalar>,
176 Dune::ReservedVector<Scalar, 1> >;
177 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
178 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
179 Scalar rho;
180
181 if (isSurfaceGrid)
182 Helper::resizeVector(outsideG[faceIdx], numOutsideFaces);
183
184 if (!curLocalScvf.isDirichlet())
185 {
186 const auto localDofIdx = curLocalScvf.localDofIndex();
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 = fvGeometry().scv(iv.localScv(negLocalScvIdx).gridScvIndex());
198 const auto& negVolVars = elemVolVars()[negGlobalScv];
199 const auto& flipScvf = !isSurfaceGrid ? curGlobalScvf
200 : 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 rho += getRho(negVolVars);
209 deltaG[localDofIdx] += alpha_outside[idxInOutside];
210 }
211 }
212
213 rho /= numOutsideFaces + 1;
214 deltaG[localDofIdx] -= alpha_inside;
215 deltaG[localDofIdx] *= rho*curGlobalScvf.area();
216 }
217 // use density resulting from Dirichlet BCs
218 else
219 rho = getRho(elemVolVars()[curGlobalScvf.outsideScvIdx()]);
220
221 // add "inside" & "outside" alphas to gravity containers
222 g[faceIdx] = alpha_inside*rho*curGlobalScvf.area();
223
224 if (isSurfaceGrid)
225 {
226 unsigned int i = 0;
227 for (const auto& alpha : alpha_outside)
228 outsideG[faceIdx][i++] = alpha*rho*curGlobalScvf.area();
229 }
230 }
231
232 // add iv-wide contributions to gravity vectors
233 handle.CA().umv(deltaG, g);
234 if (isSurfaceGrid)
235 {
236 using FaceVector = typename IV::Traits::MatVecTraits::FaceVector;
237 FaceVector AG;
238 Helper::resizeVector(AG, iv.numUnknowns());
239 handle.A().mv(deltaG, AG);
240
241 // compute gravitational accelerations
242 for (const auto& localFaceData : iv.localFaceData())
243 {
244 // continue only for "outside" faces
245 if (!localFaceData.isOutsideFace())
246 continue;
247
248 const auto localScvIdx = localFaceData.ivLocalInsideScvIndex();
249 const auto localScvfIdx = localFaceData.ivLocalScvfIndex();
250 const auto idxInOutside = localFaceData.scvfLocalOutsideScvfIndex();
251 const auto& posLocalScv = iv.localScv(localScvIdx);
252 const auto& wijk = handle.omegas()[localScvfIdx][idxInOutside+1];
253
254 // add contributions from all local directions
255 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
256 {
257 // the scvf corresponding to this local direction in the scv
258 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
259 if (!curLocalScvf.isDirichlet())
260 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
261 }
262 }
263 }
264 }
265
266private:
267 // pointers to the data required for assembly
268 const Problem* problemPtr_;
269 const FVElementGeometry* fvGeometryPtr_;
270 const ElementVolumeVariables* elemVolVarsPtr_;
271};
272
273} // end namespace Dumux
274
275#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
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
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 m...
Definition: localassemblerbase.hh:104
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:86
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: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: localassemblerbase.hh:121
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:135
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 resizeVector(Vector &v, size_type size)
resizes a vector to the given size (specialization for dynamic matrix type)
Definition: localassemblerhelper.hh:238