24#ifndef DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
25#define DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
32#include <dune/localfunctions/lagrange/pqkfactory.hh>
47template<
class ElementSolution>
51 auto firstState = elemSol[0].state();
52 for (std::size_t i = 1; i < elemSol.size(); ++i)
53 if (elemSol[i].state() != firstState)
60template<
class ElementSolution>
67template<
class Geometry,
class ElementSolution>
68auto minDistVertexSol(
const Geometry& geometry,
const typename Geometry::GlobalCoordinate& globalPos,
69 const ElementSolution& elemSol)
72 std::vector<typename Geometry::ctype> distances(geometry.corners());
73 for (
int i = 0; i < geometry.corners(); ++i)
74 distances[i] = (geometry.corner(i) - globalPos).two_norm2();
77 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
78 return elemSol[
std::distance(distances.begin(), minDistanceIt)];
94template<
class Element,
class Gr
idGeometry,
class CVFEElemSol>
95typename CVFEElemSol::PrimaryVariables
97 const typename Element::Geometry& geometry,
98 const GridGeometry& gridGeometry,
99 const CVFEElemSol& elemSol,
100 const typename Element::Geometry::GlobalCoordinate& globalPos,
101 bool ignoreState =
false)
109 using Scalar =
typename CVFEElemSol::PrimaryVariables::value_type;
112 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
115 const auto localPos = geometry.local(globalPos);
116 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
117 localBasis.evaluateFunction(localPos, shapeValues);
119 typename CVFEElemSol::PrimaryVariables result(0.0);
120 for (
int i = 0; i < shapeValues.size(); ++i)
122 auto value = elemSol[i];
123 value *= shapeValues[i][0];
128 if constexpr (HasState{})
130 result.setState(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;
164template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
166 const typename Element::Geometry& geometry,
167 const typename FVElementGeometry::GridGeometry& gridGeometry,
169 const typename Element::Geometry::GlobalCoordinate& globalPos,
170 bool ignoreState =
false)
189template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
191 const typename Element::Geometry& geometry,
193 const typename Element::Geometry::GlobalCoordinate& globalPos,
194 bool ignoreState =
false)
202 using Scalar =
typename PrimaryVariables::value_type;
203 using CoordScalar =
typename Element::Geometry::GlobalCoordinate::value_type;
204 static constexpr int dim = Element::Geometry::mydimension;
207 using FEFactory = Dune::PQkLocalFiniteElementFactory<CoordScalar, Scalar, dim, 1>;
208 using FiniteElement =
typename FEFactory::FiniteElementType;
209 std::unique_ptr<FiniteElement> fe(FEFactory::create(geometry.type()));
210 const auto& localBasis = fe->localBasis();
213 const auto localPos = geometry.local(globalPos);
214 using ShapeValue =
typename FiniteElement::Traits::LocalBasisType::Traits::RangeType;
215 std::vector< ShapeValue > shapeValues;
216 localBasis.evaluateFunction(localPos, shapeValues);
218 PrimaryVariables result(0.0);
219 for (
int i = 0; i < geometry.corners(); ++i)
221 auto value = elemSol[i];
222 value *= shapeValues[i][0];
227 if constexpr (HasState{})
229 result.setState(elemSol[0].state());
235 static bool warnedAboutUsingMinDist =
false;
236 if (!warnedAboutUsingMinDist)
238 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
239 <<
"\nbecause not all states are equal and ignoreState is false!"
241 warnedAboutUsingMinDist =
true;
260template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
262 const typename Element::Geometry& geometry,
263 const typename FVElementGeometry::GridGeometry& gridGeometry,
265 const typename Element::Geometry::GlobalCoordinate& globalPos,
266 bool ignoreState =
false)
284template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
286 const typename Element::Geometry& geometry,
288 const typename Element::Geometry::GlobalCoordinate& globalPos,
289 bool ignoreState =
false)
307template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
309 const typename Element::Geometry& geometry,
310 const typename FVElementGeometry::GridGeometry& gridGeometry,
312 const typename Element::Geometry::GlobalCoordinate& globalPos,
313 bool ignoreState =
false)
331template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
333 const typename Element::Geometry& geometry,
334 const typename FVElementGeometry::GridGeometry& gridGeometry,
336 const typename Element::Geometry::GlobalCoordinate& globalPos,
337 bool ignoreState =
false)
A helper function for class member function introspection.
Type traits to be used with matrix types.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
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:165
CVFEElemSol::PrimaryVariables evalCVFESolution(const Element &element, const typename Element::Geometry &geometry, const GridGeometry &gridGeometry, const CVFEElemSol &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given control-volume finite element solution at a given global position....
Definition: evalsolution.hh:96
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
auto minDistVertexSol(const Geometry &geometry, const typename Geometry::GlobalCoordinate &globalPos, const ElementSolution &elemSol)
return the solution at the closest dof
Definition: evalsolution.hh:68
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
bool allStatesEqual(const ElementSolution &elemSol, std::false_type hasState)
overload if the solution is stateless
Definition: evalsolution.hh:61
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.