24#ifndef DUMUX_DISCRETIZATION_EVAL_GRADIENTS_HH
25#define DUMUX_DISCRETIZATION_EVAL_GRADIENTS_HH
27#include <dune/localfunctions/lagrange/pqkfactory.hh>
28#include <dune/common/exceptions.hh>
58template<
class Element,
class Gr
idGeometry,
class CVFEElemSol>
60 const typename Element::Geometry& geometry,
61 const GridGeometry& gridGeometry,
62 const CVFEElemSol& elemSol,
63 const typename Element::Geometry::GlobalCoordinate& globalPos,
64 bool ignoreState =
false)
72 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
75 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
78 using ShapeJacobian =
typename std::decay_t<
decltype(localBasis) >::Traits::JacobianType;
79 const auto localPos = geometry.local(globalPos);
80 std::vector< ShapeJacobian > shapeJacobian;
81 localBasis.evaluateJacobian(localPos, shapeJacobian);
84 const auto jacInvT = geometry.jacobianInverseTransposed(localPos);
87 Dune::FieldVector<GlobalPosition, CVFEElemSol::PrimaryVariables::dimension> result( GlobalPosition(0.0) );
88 for (
int i = 0; i < shapeJacobian.size(); ++i)
92 jacInvT.mv(shapeJacobian[i][0], gradN);
95 for (
unsigned int pvIdx = 0; pvIdx < CVFEElemSol::PrimaryVariables::dimension; ++pvIdx)
97 GlobalPosition tmp(gradN);
98 tmp *= elemSol[i][pvIdx];
107 DUNE_THROW(Dune::NotImplemented,
"Element dofs have different phase states. Enforce calculation by setting ignoreState to true.");
124template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
126 const typename Element::Geometry& geometry,
127 const typename FVElementGeometry::GridGeometry& gridGeometry,
129 const typename Element::Geometry::GlobalCoordinate& globalPos,
130 bool ignoreState =
false)
146template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
148 const typename Element::Geometry& geometry,
149 const typename FVElementGeometry::GridGeometry& gridGeometry,
151 const typename Element::Geometry::GlobalCoordinate& globalPos,
152 bool ignoreState =
false)
168template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
170 const typename Element::Geometry& geometry,
171 const typename FVElementGeometry::GridGeometry& gridGeometry,
173 const typename Element::Geometry::GlobalCoordinate& globalPos,
174 bool ignoreState =
false)
198template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
199Dune::FieldVector<typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension>
201 const typename Element::Geometry& geometry,
202 const typename FVElementGeometry::GridGeometry& gridGeometry,
204 const typename Element::Geometry::GlobalCoordinate& globalPos)
205{ DUNE_THROW(Dune::NotImplemented,
"General gradient evaluation for cell-centered methods"); }
A helper function for class member function introspection.
Type traits to be used with matrix types.
free functions for the evaluation of primary variables inside elements.
auto evalGradients(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const BoxElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Evaluates the gradient of a given box element solution to a given global position.
Definition: evalgradients.hh:125
auto evalCVFEGradients(const Element &element, const typename Element::Geometry &geometry, const GridGeometry &gridGeometry, const CVFEElemSol &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Evaluates the gradient of a control-volume finite element solution to a given global position.
Definition: evalgradients.hh:59
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
bool allStatesEqual(const ElementSolution &elemSol, std::true_type hasState)
returns true if all states in an element solution are the same
Definition: evalsolution.hh:48
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:33
The element solution vector.
Definition: box/elementsolution.hh:39
The element solution vector.
Definition: cellcentered/elementsolution.hh:40
The global face variables class for face-centered diamond models.
Definition: facecentered/diamond/elementsolution.hh:42
The element solution vector.
Definition: pq1bubble/elementsolution.hh:39
The local element solution class for the box method.
The local element solution class for cell-centered methods.
The local element solution class for the pq1bubble method.