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)];
81template<
class Element,
class Gr
idGeometry,
class CVFEElemSol>
82typename CVFEElemSol::PrimaryVariables
84 const typename Element::Geometry& geometry,
85 const GridGeometry& gridGeometry,
86 const CVFEElemSol& elemSol,
87 const typename Element::Geometry::GlobalCoordinate& globalPos,
88 bool ignoreState =
false)
96 using Scalar =
typename CVFEElemSol::PrimaryVariables::value_type;
99 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
102 const auto localPos = geometry.local(globalPos);
103 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
104 localBasis.evaluateFunction(localPos, shapeValues);
106 typename CVFEElemSol::PrimaryVariables result(0.0);
107 for (
int i = 0; i < shapeValues.size(); ++i)
109 auto value = elemSol[i];
110 value *= shapeValues[i][0];
115 if constexpr (HasState{})
117 result.setState(elemSol[0].state());
123 static bool warnedAboutUsingMinDist =
false;
124 if (!warnedAboutUsingMinDist)
126 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
127 <<
"\nbecause not all states are equal and ignoreState is false!"
129 warnedAboutUsingMinDist =
true;
151template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
153 const typename Element::Geometry& geometry,
154 const typename FVElementGeometry::GridGeometry& gridGeometry,
156 const typename Element::Geometry::GlobalCoordinate& globalPos,
157 bool ignoreState =
false)
176template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
178 const typename Element::Geometry& geometry,
180 const typename Element::Geometry::GlobalCoordinate& globalPos,
181 bool ignoreState =
false)
183 static_assert(FVElementGeometry::GridGeometry::discMethod ==
DiscretizationMethods::box,
"Only implemented for the Box method");
191 using Scalar =
typename PrimaryVariables::value_type;
192 using CoordScalar =
typename Element::Geometry::GlobalCoordinate::value_type;
193 static constexpr int dim = Element::Geometry::mydimension;
196 using FEFactory = Dune::PQkLocalFiniteElementFactory<CoordScalar, Scalar, dim, 1>;
197 using FiniteElement =
typename FEFactory::FiniteElementType;
198 std::unique_ptr<FiniteElement> fe(FEFactory::create(geometry.type()));
199 const auto& localBasis = fe->localBasis();
202 const auto localPos = geometry.local(globalPos);
203 using ShapeValue =
typename FiniteElement::Traits::LocalBasisType::Traits::RangeType;
204 std::vector< ShapeValue > shapeValues;
205 localBasis.evaluateFunction(localPos, shapeValues);
207 PrimaryVariables result(0.0);
208 for (
int i = 0; i < geometry.corners(); ++i)
210 auto value = elemSol[i];
211 value *= shapeValues[i][0];
216 if constexpr (HasState{})
218 result.setState(elemSol[0].state());
224 static bool warnedAboutUsingMinDist =
false;
225 if (!warnedAboutUsingMinDist)
227 std::cout <<
"Warning: Using nearest-neighbor interpolation in evalSolution"
228 <<
"\nbecause not all states are equal and ignoreState is false!"
230 warnedAboutUsingMinDist =
true;
249template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
251 const typename Element::Geometry& geometry,
252 const typename FVElementGeometry::GridGeometry& gridGeometry,
254 const typename Element::Geometry::GlobalCoordinate& globalPos,
255 bool ignoreState =
false)
273template<
class Element,
class FVElementGeometry,
class PrimaryVariables>
275 const typename Element::Geometry& geometry,
277 const typename Element::Geometry::GlobalCoordinate& globalPos,
278 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:28
The local element solution class for control-volume finite element methods.
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:83
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:152
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 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
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