12#ifndef DUMUX_DISCRETIZATION_EVAL_GRADIENTS_HH
13#define DUMUX_DISCRETIZATION_EVAL_GRADIENTS_HH
15#include <dune/common/exceptions.hh>
43template<
class Element,
class Gr
idGeometry,
class CVFEElemSol>
45 const typename Element::Geometry& geometry,
46 const GridGeometry& gridGeometry,
47 const CVFEElemSol& elemSol,
48 const typename Element::Geometry::LocalCoordinate& localPos,
49 bool ignoreState =
false)
57 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
60 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
63 using ShapeJacobian =
typename std::decay_t<
decltype(localBasis) >::Traits::JacobianType;
64 std::vector< ShapeJacobian > shapeJacobian;
65 localBasis.evaluateJacobian(localPos, shapeJacobian);
68 const auto jacInvT = geometry.jacobianInverseTransposed(localPos);
71 Dune::FieldVector<GlobalPosition, CVFEElemSol::PrimaryVariables::dimension> result( GlobalPosition(0.0) );
72 for (
int i = 0; i < shapeJacobian.size(); ++i)
76 jacInvT.mv(shapeJacobian[i][0], gradN);
79 for (
unsigned int pvIdx = 0; pvIdx < CVFEElemSol::PrimaryVariables::dimension; ++pvIdx)
81 GlobalPosition tmp(gradN);
82 tmp *= elemSol[i][pvIdx];
91 DUNE_THROW(Dune::NotImplemented,
"Element dofs have different phase states. Enforce calculation by setting ignoreState to true.");
108template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
110 const typename Element::Geometry& geometry,
111 const typename FVElementGeometry::GridGeometry& gridGeometry,
113 const typename Element::Geometry::LocalCoordinate& localPos,
114 bool ignoreState =
false)
130template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
132 const typename Element::Geometry& geometry,
133 const typename FVElementGeometry::GridGeometry& gridGeometry,
135 const typename Element::Geometry::GlobalCoordinate& globalPos,
136 bool ignoreState =
false)
138 const auto& localPos = geometry.local(globalPos);
161template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
162Dune::FieldVector<typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension>
164 const typename Element::Geometry& geometry,
165 const typename FVElementGeometry::GridGeometry& gridGeometry,
167 const typename Element::Geometry::GlobalCoordinate& globalPos)
168{ DUNE_THROW(Dune::NotImplemented,
"General gradient evaluation for cell-centered methods"); }
The local element solution class for cell-centered methods.
The element solution vector.
Definition: cellcentered/elementsolution.hh:28
The element solution vector.
Definition: cvfe/elementsolution.hh:30
The local element solution class for control-volume finite element methods.
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 CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Evaluates the gradient of a given CVFE element solution to a given global position.
Definition: evalgradients.hh:131
auto evalCVFEGradientsAtLocalPos(const Element &element, const typename Element::Geometry &geometry, const GridGeometry &gridGeometry, const CVFEElemSol &elemSol, const typename Element::Geometry::LocalCoordinate &localPos, bool ignoreState=false)
Evaluates the gradient of a control-volume finite element solution to a given global position.
Definition: evalgradients.hh:44
auto evalGradientsAtLocalPos(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::LocalCoordinate &localPos, bool ignoreState=false)
Evaluates the gradient of a given CVFE element solution to a given global position.
Definition: evalgradients.hh:109
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:81
A helper function for class member function introspection.
bool allStatesEqual(const ElementSolution &elemSol, std::true_type hasState)
returns true if all states in an element solution are the same
Definition: evalsolution.hh:35
Type traits to be used with matrix types.
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:21