3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
scvgradients.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_CC_MPFA_SCV_GRADIENTS_HH
26#define DUMUX_CC_MPFA_SCV_GRADIENTS_HH
27
28#include <vector>
29#include <utility>
30#include <type_traits>
31
32#include <dune/common/fvector.hh>
33
34#include <dumux/common/math.hh>
36
37namespace Dumux {
38
45{
47 template<class GridGeometry, class Scalar>
48 using ResultPair = std::pair< std::vector<typename GridGeometry::SubControlVolume::GlobalPosition>,
49 std::vector<Dune::FieldVector<Scalar, GridGeometry::GridView::dimension>> >;
50
51public:
52
61 template<class GridGeometry, class GridVariables, class SolutionVector>
62 static ResultPair<GridGeometry, typename GridVariables::Scalar>
63 computeVelocities(const GridGeometry& gridGeometry,
64 const GridVariables& gridVariables,
65 const SolutionVector& x,
66 unsigned int phaseIdx)
67 {
68 auto gradToVelocity = [phaseIdx] (const auto& grad, const auto& volVars)
69 {
70 auto vel = mv(volVars.permeability(), grad);
71 vel *= volVars.mobility(phaseIdx);
72 return vel;
73 };
74
75 return computePressureGradients(gridGeometry, gridVariables, x, phaseIdx, gradToVelocity);
76 }
77
86 template<class GridGeometry, class GridVariables, class SolutionVector>
87 static ResultPair<GridGeometry, typename GridVariables::Scalar>
88 computePressureGradients(const GridGeometry& gridGeometry,
89 const GridVariables& gridVariables,
90 const SolutionVector& x,
91 unsigned int phaseIdx)
92 {
93 auto f = [] (const auto& grad, const auto& volVars) { return grad; };
94 return computePressureGradients(gridGeometry, gridVariables, x, phaseIdx, f);
95 }
96
113 template<class GridGeometry, class GridVariables, class SolutionVector, class F>
114 static ResultPair<GridGeometry, typename GridVariables::Scalar>
115 computePressureGradients(const GridGeometry& gridGeometry,
116 const GridVariables& gridVariables,
117 const SolutionVector& x,
118 unsigned int phaseIdx,
119 F& f)
120 {
121 using ElemVolVars = typename GridVariables::GridVolumeVariables::LocalView;
122 using FVElementGeometry = typename GridGeometry::LocalView;
123
124 auto handleFunction = [&] (auto& result,
125 const auto& handle,
126 const auto& iv,
127 const FVElementGeometry& fvGeometry,
128 const ElemVolVars& elemVolVars)
129 {
130 handle.advectionHandle().setPhaseIndex(phaseIdx);
131 const auto grads = InteractionVolumeAssemblerHelper::assembleScvGradients(handle.advectionHandle(), iv);
132 for (unsigned int i = 0; i < iv.numScvs(); ++i)
133 {
134 const auto& volVars = elemVolVars[iv.localScv(i).gridScvIndex()];
135 const auto scvGeometry = iv.getScvGeometry(i, fvGeometry);
136 result.first.push_back(scvGeometry.center());
137 result.second.push_back( f(grads[i], volVars) );
138 }
139 };
140
141 return computeGradients_(gridGeometry, gridVariables, x, handleFunction);
142 }
143
144private:
149 template<class GridGeometry, class GridVariables, class SolutionVector, class HandleFunction>
150 static ResultPair<GridGeometry, typename GridVariables::Scalar>
151 computeGradients_(const GridGeometry& gridGeometry,
152 const GridVariables& gridVariables,
153 const SolutionVector& x,
154 const HandleFunction& handleFunction)
155 {
156 using GridView = typename GridGeometry::GridView;
157 static constexpr int dim = GridView::dimension;
158
159 // first, find out how many scvs live on this grid
160 std::size_t numScvs = 0;
161 const auto& gridView = gridGeometry.gridView();
162 for (const auto& element : elements(gridView))
163 numScvs += element.subEntities(dim);
164
165 ResultPair<GridGeometry, typename GridVariables::Scalar> result;
166 result.first.reserve(numScvs);
167 result.second.reserve(numScvs);
168 std::vector<bool> vertexHandled(gridView.size(dim), false);
169
170 for (const auto& element : elements(gridView))
171 {
172 bool allFinished = true;
173 for (int i = 0; i < element.subEntities(dim); ++i)
174 if (!vertexHandled[gridGeometry.vertexMapper().subIndex(element, i, dim)])
175 allFinished = false;
176
177 // bind views only if there is unfinished buisness
178 if (allFinished)
179 continue;
180
181 // compute gradients in all scvs of all interaction volumes in this element
182 auto fvGeometry = localView(gridGeometry);
183 auto elemVolVars = localView(gridVariables.curGridVolVars());
184 auto elemFluxVarsCache = localView(gridVariables.gridFluxVarsCache());
185 fvGeometry.bind(element);
186 elemVolVars.bind(element, fvGeometry, x);
187 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
188
189 for (const auto& scvf : scvfs(fvGeometry))
190 {
191 if (vertexHandled[scvf.vertexIndex()])
192 continue;
193
194 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
195 {
196 const auto& iv = elemFluxVarsCache.secondaryInteractionVolume(scvf);
197 const auto& handle = elemFluxVarsCache.secondaryDataHandle(scvf);
198 handleFunction(result, handle, iv, fvGeometry, elemVolVars);
199 }
200 else
201 {
202 const auto& iv = elemFluxVarsCache.primaryInteractionVolume(scvf);
203 const auto& handle = elemFluxVarsCache.primaryDataHandle(scvf);
204 handleFunction(result, handle, iv, fvGeometry, elemVolVars);
205 }
206
207 vertexHandled[scvf.vertexIndex()] = true;
208 }
209 }
210
211 return result;
212 }
213};
214
215} // end namespace Dumux
216
217#endif
Define some often used mathematical functions.
A class that contains helper functions as well as functionality which is common to different mpfa sch...
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Dune::DenseVector< V >::derived_type mv(const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V > &v)
Returns the result of the projection of a vector v with a Matrix M.
Definition: math.hh:789
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
static std::vector< typename IV::Traits::LocalScvType::GlobalCoordinate > assembleScvGradients(const DataHandle &handle, const IV &iv)
Assembles the solution gradients in the sub-control volumes within an interaction volume.
Definition: localassemblerhelper.hh:185
Class roviding functionality for the reconstruction of the gradients in the sub-control volumes invol...
Definition: scvgradients.hh:45
static ResultPair< GridGeometry, typename GridVariables::Scalar > computeVelocities(const GridGeometry &gridGeometry, const GridVariables &gridVariables, const SolutionVector &x, unsigned int phaseIdx)
Computes the phase velocities in the scvs of the grid.
Definition: scvgradients.hh:63
static ResultPair< GridGeometry, typename GridVariables::Scalar > computePressureGradients(const GridGeometry &gridGeometry, const GridVariables &gridVariables, const SolutionVector &x, unsigned int phaseIdx, F &f)
Computes the pressure gradients in the scvs of the grid.
Definition: scvgradients.hh:115
static ResultPair< GridGeometry, typename GridVariables::Scalar > computePressureGradients(const GridGeometry &gridGeometry, const GridVariables &gridVariables, const SolutionVector &x, unsigned int phaseIdx)
Computes the pressure gradients in the scvs of the grid.
Definition: scvgradients.hh:88