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>
54template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
56 const typename Element::Geometry& geometry,
57 const typename FVElementGeometry::GridGeometry& gridGeometry,
59 const typename Element::Geometry::GlobalCoordinate& globalPos,
60 bool ignoreState =
false)
68 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
71 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
74 using ShapeJacobian =
typename std::decay_t<
decltype(localBasis) >::Traits::JacobianType;
75 const auto localPos = geometry.local(globalPos);
76 std::vector< ShapeJacobian > shapeJacobian;
77 localBasis.evaluateJacobian(localPos, shapeJacobian);
80 const auto jacInvT = geometry.jacobianInverseTransposed(localPos);
83 Dune::FieldVector<GlobalPosition, PrimaryVariables::dimension> result( GlobalPosition(0.0) );
84 for (
int i = 0; i < element.subEntities(Element::Geometry::mydimension); ++i)
88 jacInvT.mv(shapeJacobian[i][0], gradN);
91 for (
unsigned int pvIdx = 0; pvIdx < PrimaryVariables::dimension; ++pvIdx)
93 GlobalPosition tmp(gradN);
94 tmp *= elemSol[i][pvIdx];
103 DUNE_THROW(Dune::NotImplemented,
"Element vertices have different phase states. Enforce calculation by setting ignoreState to true.");
126template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
127Dune::FieldVector<typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension>
129 const typename Element::Geometry& geometry,
130 const typename FVElementGeometry::GridGeometry& gridGeometry,
132 const typename Element::Geometry::GlobalCoordinate& globalPos)
133{ 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:55
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
bool allStatesEqual(const ElementSolution &elemSol, std::true_type hasState)
returns true if all states in an element solution are the same
Definition: evalsolution.hh:45
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:32
The element solution vector.
Definition: box/elementsolution.hh:39
The element solution vector.
Definition: cellcentered/elementsolution.hh:40
The local element solution class for the box method.
The local element solution class for cell-centered methods.