3.6-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
39
40namespace Dumux {
41
54template< class P, class EG, class EV >
57{
58 using Problem = P;
59 using FVElementGeometry = EG;
60 using ElementVolumeVariables = EV;
63
64 template< class IV >
65 using Scalar = typename IV::Traits::MatVecTraits::FaceVector::value_type;
66
67 public:
77 const FVElementGeometry& fvGeometry,
78 const ElementVolumeVariables& elemVolVars)
79 {
80 problemPtr_ = &problem;
81 fvGeometryPtr_ = &fvGeometry;
82 elemVolVarsPtr_ = &elemVolVars;
83 }
84
85 // return functions to the local views & problem
86 const Problem& problem() const { return *problemPtr_; }
87 const FVElementGeometry& fvGeometry() const { return *fvGeometryPtr_; }
88 const ElementVolumeVariables& elemVolVars() const { return *elemVolVarsPtr_; }
89
105 template< class DataHandle, class IV, class TensorFunc >
106 void assembleMatrices(DataHandle& handle, IV& iv, const TensorFunc& getT, Scalar<IV> wijZeroThresh = 0.0)
107 {
108 DUNE_THROW(Dune::NotImplemented, "Implementation does not provide an assembleMatrices() function");
109 }
110
122 template< class DataHandle, class IV, class GetU >
123 void assembleU(DataHandle& handle, const IV& iv, const GetU& getU)
124 {
125 DUNE_THROW(Dune::NotImplemented, "Implementation does not provide an assemble() function for the cell/Dirichlet unknowns");
126 }
127
136 template< class DataHandle, class IV, class GetRho >
137 void assembleGravity(DataHandle& handle, const IV& iv, const GetRho& getRho)
138 {
139 using GridView = typename IV::Traits::GridView;
140 static constexpr int dim = GridView::dimension;
141 static constexpr int dimWorld = GridView::dimensionworld;
142 static constexpr bool isSurfaceGrid = dim < dimWorld;
143
144 // resize the gravity vectors
145 auto& g = handle.g();
146 auto& deltaG = handle.deltaG();
147 auto& outsideG = handle.gOutside();
148 Helper::resizeVector(g, iv.numFaces());
149 Helper::resizeVector(deltaG, iv.numUnknowns());
150 if (isSurfaceGrid)
151 Helper::resizeVector(outsideG, iv.numFaces());
152
157 using Scalar = typename IV::Traits::MatVecTraits::TMatrix::value_type;
158 using LocalIndexType = typename IV::Traits::IndexSet::LocalIndexType;
159
160 for (LocalIndexType faceIdx = 0; faceIdx < iv.numFaces(); ++faceIdx)
161 {
162 // gravitational acceleration on this face
163 const auto& curLocalScvf = iv.localScvf(faceIdx);
164 const auto& curGlobalScvf = fvGeometry().scvf(curLocalScvf.gridScvfIndex());
165 const auto& gravity = problem().spatialParams().gravity(curGlobalScvf.ipGlobal());
166
167 // get permeability tensor in "positive" sub volume
168 const auto& neighborScvIndices = curLocalScvf.neighboringLocalScvIndices();
169 const auto& posGlobalScv = fvGeometry().scv(iv.localScv(neighborScvIndices[0]).gridScvIndex());
170 const auto& posVolVars = elemVolVars()[posGlobalScv];
171 const auto alpha_inside = posVolVars.extrusionFactor()*vtmv(curGlobalScvf.unitOuterNormal(),
172 posVolVars.permeability(),
173 gravity);
174
175 const auto numOutsideFaces = !curGlobalScvf.boundary() ? curGlobalScvf.numOutsideScvs() : 0;
176 using OutsideAlphaStorage = std::conditional_t< isSurfaceGrid,
177 std::vector<Scalar>,
178 Dune::ReservedVector<Scalar, 1> >;
179 OutsideAlphaStorage alpha_outside; alpha_outside.resize(numOutsideFaces);
180 std::fill(alpha_outside.begin(), alpha_outside.end(), 0.0);
181 Scalar rho;
182
183 if (isSurfaceGrid)
184 Helper::resizeVector(outsideG[faceIdx], numOutsideFaces);
185
186 if (!curLocalScvf.isDirichlet())
187 {
188 const auto localDofIdx = curLocalScvf.localDofIndex();
189
190 rho = getRho(posVolVars);
191 deltaG[localDofIdx] = 0.0;
192
193 if (!curGlobalScvf.boundary())
194 {
195 for (unsigned int idxInOutside = 0; idxInOutside < curGlobalScvf.numOutsideScvs(); ++idxInOutside)
196 {
197 // obtain outside tensor
198 const auto negLocalScvIdx = neighborScvIndices[idxInOutside+1];
199 const auto& negGlobalScv = fvGeometry().scv(iv.localScv(negLocalScvIdx).gridScvIndex());
200 const auto& negVolVars = elemVolVars()[negGlobalScv];
201 const auto& flipScvf = !isSurfaceGrid ? curGlobalScvf
202 : fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside);
203
204 alpha_outside[idxInOutside] = negVolVars.extrusionFactor()*vtmv(flipScvf.unitOuterNormal(),
205 negVolVars.permeability(),
206 gravity);
207 if (isSurfaceGrid)
208 alpha_outside[idxInOutside] *= -1.0;
209
210 rho += getRho(negVolVars);
211 deltaG[localDofIdx] += alpha_outside[idxInOutside];
212 }
213 }
214
215 rho /= numOutsideFaces + 1;
216 deltaG[localDofIdx] -= alpha_inside;
217 deltaG[localDofIdx] *= rho*Extrusion::area(fvGeometry(), curGlobalScvf);
218 }
219 // use density resulting from Dirichlet BCs
220 else
221 rho = getRho(elemVolVars()[curGlobalScvf.outsideScvIdx()]);
222
223 // add "inside" & "outside" alphas to gravity containers
224 g[faceIdx] = alpha_inside*rho*Extrusion::area(fvGeometry(), curGlobalScvf);
225
226 if (isSurfaceGrid)
227 {
228 unsigned int i = 0;
229 for (const auto& alpha : alpha_outside)
230 outsideG[faceIdx][i++] = alpha*rho*Extrusion::area(fvGeometry(), curGlobalScvf);
231 }
232 }
233
234 // add iv-wide contributions to gravity vectors
235 handle.CA().umv(deltaG, g);
236 if (isSurfaceGrid)
237 {
238 using FaceVector = typename IV::Traits::MatVecTraits::FaceVector;
239 FaceVector AG;
240 Helper::resizeVector(AG, iv.numUnknowns());
241 handle.A().mv(deltaG, AG);
242
243 // compute gravitational accelerations
244 for (const auto& localFaceData : iv.localFaceData())
245 {
246 // continue only for "outside" faces
247 if (!localFaceData.isOutsideFace())
248 continue;
249
250 const auto localScvIdx = localFaceData.ivLocalInsideScvIndex();
251 const auto localScvfIdx = localFaceData.ivLocalScvfIndex();
252 const auto idxInOutside = localFaceData.scvfLocalOutsideScvfIndex();
253 const auto& posLocalScv = iv.localScv(localScvIdx);
254 const auto& wijk = handle.omegas()[localScvfIdx][idxInOutside+1];
255
256 // add contributions from all local directions
257 for (LocalIndexType localDir = 0; localDir < dim; localDir++)
258 {
259 // the scvf corresponding to this local direction in the scv
260 const auto& curLocalScvf = iv.localScvf(posLocalScv.localScvfIndex(localDir));
261 if (!curLocalScvf.isDirichlet())
262 outsideG[localScvfIdx][idxInOutside] -= wijk[localDir]*AG[curLocalScvf.localDofIndex()];
263 }
264 }
265 }
266 }
267
268private:
269 // pointers to the data required for assembly
270 const Problem* problemPtr_;
271 const FVElementGeometry* fvGeometryPtr_;
272 const ElementVolumeVariables* elemVolVarsPtr_;
273};
274
275} // end namespace Dumux
276
277#endif
A helper function for class member function introspection.
Helper classes to compute the integration elements.
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:863
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
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
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:106
const ElementVolumeVariables & elemVolVars() const
Definition: localassemblerbase.hh:88
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:76
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:123
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:137
const Problem & problem() const
Definition: localassemblerbase.hh:86
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