25#ifndef DUMUX_CC_MPFA_SCV_GRADIENTS_HH
26#define DUMUX_CC_MPFA_SCV_GRADIENTS_HH
32#include <dune/common/fvector.hh>
47 template<
class Gr
idGeometry,
class Scalar>
48 using ResultPair = std::pair< std::vector<typename GridGeometry::SubControlVolume::GlobalPosition>,
49 std::vector<Dune::FieldVector<Scalar, GridGeometry::GridView::dimension>> >;
61 template<
class Gr
idGeometry,
class Gr
idVariables,
class SolutionVector>
62 static ResultPair<GridGeometry, typename GridVariables::Scalar>
64 const GridVariables& gridVariables,
65 const SolutionVector& x,
66 unsigned int phaseIdx)
68 auto gradToVelocity = [phaseIdx] (
const auto& grad,
const auto& volVars)
70 auto vel =
mv(volVars.permeability(), grad);
71 vel *= volVars.mobility(phaseIdx);
86 template<
class Gr
idGeometry,
class Gr
idVariables,
class SolutionVector>
87 static ResultPair<GridGeometry, typename GridVariables::Scalar>
89 const GridVariables& gridVariables,
90 const SolutionVector& x,
91 unsigned int phaseIdx)
93 auto f = [] (
const auto& grad,
const auto& volVars) {
return grad; };
113 template<
class Gr
idGeometry,
class Gr
idVariables,
class SolutionVector,
class F>
114 static ResultPair<GridGeometry, typename GridVariables::Scalar>
116 const GridVariables& gridVariables,
117 const SolutionVector& x,
118 unsigned int phaseIdx,
121 using ElemVolVars =
typename GridVariables::GridVolumeVariables::LocalView;
122 using FVElementGeometry =
typename GridGeometry::LocalView;
124 auto handleFunction = [&] (
auto& result,
127 const FVElementGeometry& fvGeometry,
128 const ElemVolVars& elemVolVars)
130 handle.advectionHandle().setPhaseIndex(phaseIdx);
132 for (
unsigned int i = 0; i < iv.numScvs(); ++i)
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) );
141 return computeGradients_(gridGeometry, gridVariables, x, handleFunction);
149 template<
class Gr
idGeometry,
class Gr
idVariables,
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)
156 using GridView =
typename GridGeometry::GridView;
157 static constexpr int dim = GridView::dimension;
160 std::size_t numScvs = 0;
161 const auto& gridView = gridGeometry.gridView();
162 for (
const auto& element : elements(gridView))
163 numScvs += element.subEntities(dim);
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);
170 for (
const auto& element : elements(gridView))
172 bool allFinished =
true;
173 for (
int i = 0; i < element.subEntities(dim); ++i)
174 if (!vertexHandled[gridGeometry.vertexMapper().subIndex(element, i, dim)])
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);
189 for (
const auto& scvf : scvfs(fvGeometry))
191 if (vertexHandled[scvf.vertexIndex()])
194 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
196 const auto& iv = elemFluxVarsCache.secondaryInteractionVolume(scvf);
197 const auto& handle = elemFluxVarsCache.secondaryDataHandle(scvf);
198 handleFunction(result, handle, iv, fvGeometry, elemVolVars);
202 const auto& iv = elemFluxVarsCache.primaryInteractionVolume(scvf);
203 const auto& handle = elemFluxVarsCache.primaryDataHandle(scvf);
204 handleFunction(result, handle, iv, fvGeometry, elemVolVars);
207 vertexHandled[scvf.vertexIndex()] =
true;
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:798
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