version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_CC_MPFA_SCV_GRADIENTS_HH
14#define DUMUX_CC_MPFA_SCV_GRADIENTS_HH
15
16#include <vector>
17#include <utility>
18#include <type_traits>
19
20#include <dune/common/fvector.hh>
21
22#include <dumux/common/math.hh>
24
25namespace Dumux {
26
33{
35 template<class GridGeometry, class Scalar>
36 using ResultPair = std::pair< std::vector<typename GridGeometry::SubControlVolume::GlobalPosition>,
37 std::vector<Dune::FieldVector<Scalar, GridGeometry::GridView::dimension>> >;
38
39public:
40
49 template<class GridGeometry, class GridVariables, class SolutionVector>
50 static ResultPair<GridGeometry, typename GridVariables::Scalar>
51 computeVelocities(const GridGeometry& gridGeometry,
52 const GridVariables& gridVariables,
53 const SolutionVector& x,
54 unsigned int phaseIdx)
55 {
56 auto gradToVelocity = [phaseIdx] (const auto& grad, const auto& volVars)
57 {
58 auto vel = mv(volVars.permeability(), grad);
59 vel *= volVars.mobility(phaseIdx);
60 return vel;
61 };
62
63 return computePressureGradients(gridGeometry, gridVariables, x, phaseIdx, gradToVelocity);
64 }
65
74 template<class GridGeometry, class GridVariables, class SolutionVector>
75 static ResultPair<GridGeometry, typename GridVariables::Scalar>
76 computePressureGradients(const GridGeometry& gridGeometry,
77 const GridVariables& gridVariables,
78 const SolutionVector& x,
79 unsigned int phaseIdx)
80 {
81 auto f = [] (const auto& grad, const auto& volVars) { return grad; };
82 return computePressureGradients(gridGeometry, gridVariables, x, phaseIdx, f);
83 }
84
101 template<class GridGeometry, class GridVariables, class SolutionVector, class F>
102 static ResultPair<GridGeometry, typename GridVariables::Scalar>
103 computePressureGradients(const GridGeometry& gridGeometry,
104 const GridVariables& gridVariables,
105 const SolutionVector& x,
106 unsigned int phaseIdx,
107 F& f)
108 {
109 using ElemVolVars = typename GridVariables::GridVolumeVariables::LocalView;
110 using FVElementGeometry = typename GridGeometry::LocalView;
111
112 auto handleFunction = [&] (auto& result,
113 const auto& handle,
114 const auto& iv,
115 const FVElementGeometry& fvGeometry,
116 const ElemVolVars& elemVolVars)
117 {
118 handle.advectionHandle().setPhaseIndex(phaseIdx);
119 const auto grads = InteractionVolumeAssemblerHelper::assembleScvGradients(handle.advectionHandle(), iv);
120 for (unsigned int i = 0; i < iv.numScvs(); ++i)
121 {
122 const auto& volVars = elemVolVars[iv.localScv(i).gridScvIndex()];
123 const auto scvGeometry = iv.getScvGeometry(i, fvGeometry);
124 result.first.push_back(scvGeometry.center());
125 result.second.push_back( f(grads[i], volVars) );
126 }
127 };
128
129 return computeGradients_(gridGeometry, gridVariables, x, handleFunction);
130 }
131
132private:
137 template<class GridGeometry, class GridVariables, class SolutionVector, class HandleFunction>
138 static ResultPair<GridGeometry, typename GridVariables::Scalar>
139 computeGradients_(const GridGeometry& gridGeometry,
140 const GridVariables& gridVariables,
141 const SolutionVector& x,
142 const HandleFunction& handleFunction)
143 {
144 using GridView = typename GridGeometry::GridView;
145 static constexpr int dim = GridView::dimension;
146
147 // first, find out how many scvs live on this grid
148 std::size_t numScvs = 0;
149 const auto& gridView = gridGeometry.gridView();
150 for (const auto& element : elements(gridView))
151 numScvs += element.subEntities(dim);
152
153 ResultPair<GridGeometry, typename GridVariables::Scalar> result;
154 result.first.reserve(numScvs);
155 result.second.reserve(numScvs);
156 std::vector<bool> vertexHandled(gridView.size(dim), false);
157
158 auto fvGeometry = localView(gridGeometry);
159 auto elemVolVars = localView(gridVariables.curGridVolVars());
160 auto elemFluxVarsCache = localView(gridVariables.gridFluxVarsCache());
161
162 for (const auto& element : elements(gridView))
163 {
164 bool allFinished = true;
165 for (int i = 0; i < element.subEntities(dim); ++i)
166 if (!vertexHandled[gridGeometry.vertexMapper().subIndex(element, i, dim)])
167 allFinished = false;
168
169 // bind views only if there is unfinished business
170 if (allFinished)
171 continue;
172
173 // compute gradients in all scvs of all interaction volumes in this element
174 fvGeometry.bind(element);
175 elemVolVars.bind(element, fvGeometry, x);
176 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
177
178 for (const auto& scvf : scvfs(fvGeometry))
179 {
180 if (vertexHandled[scvf.vertexIndex()])
181 continue;
182
183 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
184 {
185 const auto& iv = elemFluxVarsCache.secondaryInteractionVolume(scvf);
186 const auto& handle = elemFluxVarsCache.secondaryDataHandle(scvf);
187 handleFunction(result, handle, iv, fvGeometry, elemVolVars);
188 }
189 else
190 {
191 const auto& iv = elemFluxVarsCache.primaryInteractionVolume(scvf);
192 const auto& handle = elemFluxVarsCache.primaryDataHandle(scvf);
193 handleFunction(result, handle, iv, fvGeometry, elemVolVars);
194 }
195
196 vertexHandled[scvf.vertexIndex()] = true;
197 }
198 }
199
200 return result;
201 }
202};
203
204} // end namespace Dumux
205
206#endif
Class providing functionality for the reconstruction of the gradients in the sub-control volumes invo...
Definition: scvgradients.hh:33
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:51
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:103
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:76
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:173
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:829
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
A class that contains helper functions as well as functionality which is common to different mpfa sch...
Define some often used mathematical functions.
Definition: adapt.hh:17