24#ifndef DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
25#define DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
31#include <dune/localfunctions/lagrange/pqkfactory.hh>
44template<
class ElementSolution>
48 auto firstState = elemSol[0].state();
49 for (std::size_t i = 1; i < elemSol.size(); ++i)
50 if (elemSol[i].state() != firstState)
57template<
class ElementSolution>
64template<
class Geometry,
class ElementSolution>
65auto minDistVertexSol(
const Geometry& geometry,
const typename Geometry::GlobalCoordinate& globalPos,
66 const ElementSolution& elemSol)
69 std::vector<typename Geometry::ctype> distances(geometry.corners());
70 for (
int i = 0; i < geometry.corners(); ++i)
71 distances[i] = (geometry.corner(i) - globalPos).two_norm2();
74 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
75 return elemSol[
std::distance(distances.begin(), minDistanceIt)];
93template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
95 const typename Element::Geometry& geometry,
96 const typename FVElementGeometry::GridGeometry& gridGeometry,
98 const typename Element::Geometry::GlobalCoordinate& globalPos,
99 bool ignoreState =
false)
107 using Scalar =
typename PrimaryVariables::value_type;
110 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
113 const auto localPos = geometry.local(globalPos);
114 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
115 localBasis.evaluateFunction(localPos, shapeValues);
117 PrimaryVariables result(0.0);
118 for (
int i = 0; i < geometry.corners(); ++i)
120 auto value = elemSol[i];
121 value *= shapeValues[i][0];
126 if constexpr (HasState{})
128 result.setState(elemSol[0].state());
134 static bool warnedAboutUsingMinDist =
false;
135 if (!warnedAboutUsingMinDist)
137 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
138 <<
"\nbecause not all states are equal and ignoreState is false!"
140 warnedAboutUsingMinDist =
true;
161template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
163 const typename Element::Geometry& geometry,
165 const typename Element::Geometry::GlobalCoordinate& globalPos,
166 bool ignoreState =
false)
174 using Scalar =
typename PrimaryVariables::value_type;
175 using CoordScalar =
typename Element::Geometry::GlobalCoordinate::value_type;
176 static constexpr int dim = Element::Geometry::mydimension;
179 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
180 using ShapeValue =
typename FeCache::FiniteElementType::Traits::LocalBasisType::Traits::RangeType;
184 const auto& localBasis = feCache.get(geometry.type()).localBasis();
187 const auto localPos = geometry.local(globalPos);
188 std::vector< ShapeValue > shapeValues;
189 localBasis.evaluateFunction(localPos, shapeValues);
191 PrimaryVariables result(0.0);
192 for (
int i = 0; i < geometry.corners(); ++i)
194 auto value = elemSol[i];
195 value *= shapeValues[i][0];
200 if constexpr (HasState{})
202 result.setState(elemSol[0].state());
208 static bool warnedAboutUsingMinDist =
false;
209 if (!warnedAboutUsingMinDist)
211 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
212 <<
"\nbecause not all states are equal and ignoreState is false!"
214 warnedAboutUsingMinDist =
true;
233template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
235 const typename Element::Geometry& geometry,
236 const typename FVElementGeometry::GridGeometry& gridGeometry,
238 const typename Element::Geometry::GlobalCoordinate& globalPos,
239 bool ignoreState =
false)
257template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
259 const typename Element::Geometry& geometry,
261 const typename Element::Geometry::GlobalCoordinate& globalPos,
262 bool ignoreState =
false)
Type traits to be used with matrix types.
A helper function for class member function introspection.
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
PrimaryVariables evalSolution(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)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:94
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
auto minDistVertexSol(const Geometry &geometry, const typename Geometry::GlobalCoordinate &globalPos, const ElementSolution &elemSol)
return the solution at the closest dof
Definition: evalsolution.hh:65
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
bool allStatesEqual(const ElementSolution &elemSol, std::false_type hasState)
overload if the solution is stateless
Definition: evalsolution.hh:58
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.