12#ifndef DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
13#define DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
20#include <dune/localfunctions/lagrange/pqkfactory.hh>
34template<
class ElementSolution>
38 auto firstState = elemSol[0].state();
39 for (std::size_t i = 1; i < elemSol.size(); ++i)
40 if (elemSol[i].state() != firstState)
47template<
class ElementSolution>
54template<
class Geometry,
class ElementSolution>
55auto minDistVertexSol(
const Geometry& geometry,
const typename Geometry::GlobalCoordinate& globalPos,
56 const ElementSolution& elemSol)
59 std::vector<typename Geometry::ctype> distances(geometry.corners());
60 for (
int i = 0; i < geometry.corners(); ++i)
61 distances[i] = (geometry.corner(i) - globalPos).two_norm2();
64 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
65 return elemSol[
std::distance(distances.begin(), minDistanceIt)];
69template<
class Element,
class Gr
idGeometry,
class ElementSolution>
71 const GridGeometry& gridGeometry,
72 const typename Element::Geometry::LocalCoordinate& localPos,
73 const ElementSolution& elemSol)
75 using GeometryHelper = GridGeometry::Cache::GeometryHelper;
76 const auto& localCoeffs = gridGeometry.feCache().get(element.type()).localCoefficients();
78 std::vector<typename Element::Geometry::ctype> distances(localCoeffs.size());
79 for (
int idx = 0; idx < localCoeffs.size(); ++idx)
81 const auto& localDofPos = GeometryHelper::localDofPosition(element.type(), localCoeffs.localKey(idx));
82 distances[idx] = (localPos - localDofPos).two_norm2();
86 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
87 return elemSol[
std::distance(distances.begin(), minDistanceIt)];
103template<
class Element,
class Gr
idGeometry,
class CVFEElemSol>
104typename CVFEElemSol::PrimaryVariables
106 const typename Element::Geometry& geometry,
107 const GridGeometry& gridGeometry,
108 const CVFEElemSol& elemSol,
109 const typename Element::Geometry::LocalCoordinate& localPos,
110 bool ignoreState =
false)
118 using Scalar =
typename CVFEElemSol::PrimaryVariables::value_type;
121 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
124 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
125 localBasis.evaluateFunction(localPos, shapeValues);
127 typename CVFEElemSol::PrimaryVariables result(0.0);
128 for (
int i = 0; i < shapeValues.size(); ++i)
130 auto value = elemSol[i];
131 value *= shapeValues[i][0];
136 if constexpr (HasState{})
138 result.setState(elemSol[0].state());
144 static bool warnedAboutUsingMinDist =
false;
145 if (!warnedAboutUsingMinDist)
147 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
148 <<
"\nbecause not all states are equal and ignoreState is false!"
150 warnedAboutUsingMinDist =
true;
172template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
174 const typename Element::Geometry& geometry,
175 const typename FVElementGeometry::GridGeometry& gridGeometry,
177 const typename Element::Geometry::LocalCoordinate& localPos,
178 bool ignoreState =
false)
197template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
199 const typename Element::Geometry& geometry,
200 const typename FVElementGeometry::GridGeometry& gridGeometry,
202 const typename Element::Geometry::GlobalCoordinate& globalPos,
203 bool ignoreState =
false)
205 const auto& localPos = geometry.local(globalPos);
223template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
225 const typename Element::Geometry& geometry,
227 const typename Element::Geometry::GlobalCoordinate& globalPos,
228 bool ignoreState =
false)
230 static_assert(FVElementGeometry::GridGeometry::discMethod ==
DiscretizationMethods::box,
"Only implemented for the Box method");
238 using Scalar =
typename PrimaryVariables::value_type;
239 using CoordScalar =
typename Element::Geometry::GlobalCoordinate::value_type;
240 static constexpr int dim = Element::Geometry::mydimension;
243 using FEFactory = Dune::PQkLocalFiniteElementFactory<CoordScalar, Scalar, dim, 1>;
244 using FiniteElement =
typename FEFactory::FiniteElementType;
245 std::unique_ptr<FiniteElement> fe(FEFactory::create(geometry.type()));
246 const auto& localBasis = fe->localBasis();
249 const auto localPos = geometry.local(globalPos);
250 using ShapeValue =
typename FiniteElement::Traits::LocalBasisType::Traits::RangeType;
251 std::vector< ShapeValue > shapeValues;
252 localBasis.evaluateFunction(localPos, shapeValues);
254 PrimaryVariables result(0.0);
255 for (
int i = 0; i < geometry.corners(); ++i)
257 auto value = elemSol[i];
258 value *= shapeValues[i][0];
263 if constexpr (HasState{})
265 result.setState(elemSol[0].state());
271 static bool warnedAboutUsingMinDist =
false;
272 if (!warnedAboutUsingMinDist)
274 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
275 <<
"\nbecause not all states are equal and ignoreState is false!"
277 warnedAboutUsingMinDist =
true;
296template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
298 const typename Element::Geometry& geometry,
299 const typename FVElementGeometry::GridGeometry& gridGeometry,
301 const typename Element::Geometry::GlobalCoordinate& globalPos,
302 bool ignoreState =
false)
319template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
321 const typename Element::Geometry& geometry,
322 const typename FVElementGeometry::GridGeometry& gridGeometry,
324 const typename Element::Geometry::LocalCoordinate& localPos,
325 bool ignoreState =
false)
344template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
346 const typename Element::Geometry& geometry,
348 const typename Element::Geometry::GlobalCoordinate& globalPos,
349 bool ignoreState =
false)
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.
PrimaryVariables evalSolution(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)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:198
CVFEElemSol::PrimaryVariables evalCVFESolutionAtLocalPos(const Element &element, const typename Element::Geometry &geometry, const GridGeometry &gridGeometry, const CVFEElemSol &elemSol, const typename Element::Geometry::LocalCoordinate &localPos, bool ignoreState=false)
Interpolates a given control-volume finite element solution at a given global position....
Definition: evalsolution.hh:105
PrimaryVariables evalSolutionAtLocalPos(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)
Interpolates a given cvfe element solution at a given local position. Uses the finite element cache o...
Definition: evalsolution.hh:173
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:282
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.
The available discretization methods in Dumux.
auto minDistVertexSol(const Geometry &geometry, const typename Geometry::GlobalCoordinate &globalPos, const ElementSolution &elemSol)
return the solution at the closest vertex dof
Definition: evalsolution.hh:55
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
bool allStatesEqual(const ElementSolution &elemSol, std::false_type hasState)
overload if the solution is stateless
Definition: evalsolution.hh:48
auto minDistDofSol(const Element &element, const GridGeometry &gridGeometry, const typename Element::Geometry::LocalCoordinate &localPos, const ElementSolution &elemSol)
return the solution at the closest dof
Definition: evalsolution.hh:70
constexpr Box box
Definition: method.hh:147
Type traits to be used with matrix types.
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:21