24#ifndef DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
25#define DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
31#include <dune/common/hybridutilities.hh>
32#include <dune/localfunctions/lagrange/pqkfactory.hh>
45template<
class ElementSolution>
49 auto firstState = elemSol[0].state();
50 for (std::size_t i = 1; i < elemSol.size(); ++i)
51 if (elemSol[i].state() != firstState)
58template<
class ElementSolution>
65template<
class Geometry,
class ElementSolution>
66auto minDistVertexSol(
const Geometry& geometry,
const typename Geometry::GlobalCoordinate& globalPos,
67 const ElementSolution& elemSol)
70 std::vector<typename Geometry::ctype> distances(geometry.corners());
71 for (
int i = 0; i < geometry.corners(); ++i)
72 distances[i] = (geometry.corner(i) - globalPos).two_norm2();
75 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
76 return elemSol[std::distance(distances.begin(), minDistanceIt)];
94template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
96 const typename Element::Geometry& geometry,
97 const typename FVElementGeometry::GridGeometry& gridGeometry,
99 const typename Element::Geometry::GlobalCoordinate& globalPos,
100 bool ignoreState =
false)
108 using Scalar =
typename PrimaryVariables::value_type;
111 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
114 const auto localPos = geometry.local(globalPos);
115 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
116 localBasis.evaluateFunction(localPos, shapeValues);
118 PrimaryVariables result(0.0);
119 for (
int i = 0; i < geometry.corners(); ++i)
121 auto value = elemSol[i];
122 value *= shapeValues[i][0];
127 Dune::Hybrid::ifElse(HasState{}, [&](
auto&& evalLazy){
129 evalLazy(result).setState(evalLazy(elemSol)[0].state());
136 static bool warnedAboutUsingMinDist =
false;
137 if (!warnedAboutUsingMinDist)
139 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
140 <<
"\nbecause not all states are equal and ignoreState is false!"
142 warnedAboutUsingMinDist =
true;
163template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
165 const typename Element::Geometry& geometry,
167 const typename Element::Geometry::GlobalCoordinate& globalPos,
168 bool ignoreState =
false)
176 using Scalar =
typename PrimaryVariables::value_type;
177 using CoordScalar =
typename Element::Geometry::GlobalCoordinate::value_type;
178 static constexpr int dim = Element::Geometry::mydimension;
181 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
182 using ShapeValue =
typename FeCache::FiniteElementType::Traits::LocalBasisType::Traits::RangeType;
186 const auto& localBasis = feCache.get(geometry.type()).localBasis();
189 const auto localPos = geometry.local(globalPos);
190 std::vector< ShapeValue > shapeValues;
191 localBasis.evaluateFunction(localPos, shapeValues);
193 PrimaryVariables result(0.0);
194 for (
int i = 0; i < geometry.corners(); ++i)
196 auto value = elemSol[i];
197 value *= shapeValues[i][0];
202 Dune::Hybrid::ifElse(HasState{}, [&](
auto&& evalLazy){
204 evalLazy(result).setState(evalLazy(elemSol)[0].state());
211 static bool warnedAboutUsingMinDist =
false;
212 if (!warnedAboutUsingMinDist)
214 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
215 <<
"\nbecause not all states are equal and ignoreState is false!"
217 warnedAboutUsingMinDist =
true;
236template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
238 const typename Element::Geometry& geometry,
239 const typename FVElementGeometry::GridGeometry& gridGeometry,
241 const typename Element::Geometry::GlobalCoordinate& globalPos,
242 bool ignoreState =
false)
260template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
262 const typename Element::Geometry& geometry,
264 const typename Element::Geometry::GlobalCoordinate& globalPos,
265 bool ignoreState =
false)
A helper function for class member function introspection.
Type traits to be used with matrix types.
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:95
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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:66
bool allStatesEqual(const ElementSolution &elemSol, std::true_type hasState)
returns true if all states in an element solution are the same
Definition: evalsolution.hh:46
bool allStatesEqual(const ElementSolution &elemSol, std::false_type hasState)
overload if the solution is stateless
Definition: evalsolution.hh:59
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.