version 3.11-dev
evalsolution.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
13#define DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
14
15#include <iterator>
16#include <algorithm>
17#include <type_traits>
18#include <memory>
19
20#include <dune/localfunctions/lagrange/pqkfactory.hh>
21
27
28namespace Dumux {
29
30// some implementation details
31namespace Detail {
32
34template<class ElementSolution>
35bool allStatesEqual(const ElementSolution& elemSol, std::true_type hasState)
36{
37 // determine if all states are the same at all vertices
38 auto firstState = elemSol[0].state();
39 for (std::size_t i = 1; i < elemSol.size(); ++i)
40 if (elemSol[i].state() != firstState)
41 return false;
42
43 return true;
44}
45
47template<class ElementSolution>
48bool allStatesEqual(const ElementSolution& elemSol, std::false_type hasState)
49{
50 return true;
51}
52
54template<class Geometry, class ElementSolution>
55auto minDistVertexSol(const Geometry& geometry, const typename Geometry::GlobalCoordinate& globalPos,
56 const ElementSolution& elemSol)
57{
58 // calculate the distances from the evaluation point to the vertices
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();
62
63 // determine the minimum distance and corresponding vertex
64 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
65 return elemSol[std::distance(distances.begin(), minDistanceIt)];
66}
67
69template<class Element, class GridGeometry, class ElementSolution>
70auto minDistDofSol(const Element& element,
71 const GridGeometry& gridGeometry,
72 const typename Element::Geometry::LocalCoordinate& localPos,
73 const ElementSolution& elemSol)
74{
75 using GeometryHelper = GridGeometry::Cache::GeometryHelper;
76 const auto& localCoeffs = gridGeometry.feCache().get(element.type()).localCoefficients();
77 // calculate the distances from the evaluation point to the local positions of dofs
78 std::vector<typename Element::Geometry::ctype> distances(localCoeffs.size());
79 for (int idx = 0; idx < localCoeffs.size(); ++idx)
80 {
81 const auto& localDofPos = GeometryHelper::localDofPosition(element.type(), localCoeffs.localKey(idx));
82 distances[idx] = (localPos - localDofPos).two_norm2();
83 }
84
85 // determine the minimum distance and corresponding vertex
86 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
87 return elemSol[std::distance(distances.begin(), minDistanceIt)];
88}
89
103template<class Element, class GridGeometry, class CVFEElemSol>
104typename CVFEElemSol::PrimaryVariables
105evalCVFESolutionAtLocalPos(const Element& element,
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)
111{
112 // determine if all states are the same at all vertices
113 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
114 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
115
116 if (allStatesEqual)
117 {
118 using Scalar = typename CVFEElemSol::PrimaryVariables::value_type;
119
120 // interpolate the solution
121 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
122
123 // evaluate the shape functions at the scv center
124 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
125 localBasis.evaluateFunction(localPos, shapeValues);
126
127 typename CVFEElemSol::PrimaryVariables result(0.0);
128 for (int i = 0; i < shapeValues.size(); ++i)
129 {
130 auto value = elemSol[i];
131 value *= shapeValues[i][0];
132 result += value;
133 }
134
135 // set an arbitrary state if the model requires a state (models constexpr if)
136 if constexpr (HasState{})
137 if (!ignoreState)
138 result.setState(elemSol[0].state());
139
140 return result;
141 }
142 else
143 {
144 static bool warnedAboutUsingMinDist = false;
145 if (!warnedAboutUsingMinDist)
146 {
147 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
148 << "\nbecause not all states are equal and ignoreState is false!"
149 << std::endl;
150 warnedAboutUsingMinDist = true;
151 }
152
153 return Detail::minDistDofSol(element, gridGeometry, localPos, elemSol);
154 }
155}
156
157} // end namespace Detail
158
172template<class Element, class FVElementGeometry, class PrimaryVariables>
173PrimaryVariables evalSolutionAtLocalPos(const Element& element,
174 const typename Element::Geometry& geometry,
175 const typename FVElementGeometry::GridGeometry& gridGeometry,
177 const typename Element::Geometry::LocalCoordinate& localPos,
178 bool ignoreState = false)
179{
180 return Detail::evalCVFESolutionAtLocalPos(element, geometry, gridGeometry, elemSol, localPos, ignoreState);
181}
182
183
197template<class Element, class FVElementGeometry, class PrimaryVariables>
198PrimaryVariables evalSolution(const Element& element,
199 const typename Element::Geometry& geometry,
200 const typename FVElementGeometry::GridGeometry& gridGeometry,
202 const typename Element::Geometry::GlobalCoordinate& globalPos,
203 bool ignoreState = false)
204{
205 const auto& localPos = geometry.local(globalPos);
206 return evalSolutionAtLocalPos(element, geometry, gridGeometry, elemSol, localPos, ignoreState);
207}
208
223template<class Element, class FVElementGeometry, class PrimaryVariables>
224PrimaryVariables evalSolution(const Element& element,
225 const typename Element::Geometry& geometry,
227 const typename Element::Geometry::GlobalCoordinate& globalPos,
228 bool ignoreState = false)
229{
230 static_assert(FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::box, "Only implemented for the Box method");
231
232 // determine if all states are the same at all vertices
233 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
234 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
235
236 if (allStatesEqual)
237 {
238 using Scalar = typename PrimaryVariables::value_type;
239 using CoordScalar = typename Element::Geometry::GlobalCoordinate::value_type;
240 static constexpr int dim = Element::Geometry::mydimension;
241
242 // The box scheme always uses linear Ansatz functions
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();
247
248 // evaluate the shape functions at the scv center
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);
253
254 PrimaryVariables result(0.0);
255 for (int i = 0; i < geometry.corners(); ++i)
256 {
257 auto value = elemSol[i];
258 value *= shapeValues[i][0];
259 result += value;
260 }
261
262 // set an arbitrary state if the model requires a state (models constexpr if)
263 if constexpr (HasState{})
264 if (!ignoreState)
265 result.setState(elemSol[0].state());
266
267 return result;
268 }
269 else
270 {
271 static bool warnedAboutUsingMinDist = false;
272 if (!warnedAboutUsingMinDist)
273 {
274 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
275 << "\nbecause not all states are equal and ignoreState is false!"
276 << std::endl;
277 warnedAboutUsingMinDist = true;
278 }
279
280 return Detail::minDistVertexSol(geometry, globalPos, elemSol);
281 }
282}
283
296template<class Element, class FVElementGeometry, class PrimaryVariables>
297PrimaryVariables evalSolution(const Element& element,
298 const typename Element::Geometry& geometry,
299 const typename FVElementGeometry::GridGeometry& gridGeometry,
301 const typename Element::Geometry::GlobalCoordinate& globalPos,
302 bool ignoreState = false)
303{
304 return elemSol[0];
305}
306
319template<class Element, class FVElementGeometry, class PrimaryVariables>
320PrimaryVariables evalSolutionAtLocalPos(const Element& element,
321 const typename Element::Geometry& geometry,
322 const typename FVElementGeometry::GridGeometry& gridGeometry,
324 const typename Element::Geometry::LocalCoordinate& localPos,
325 bool ignoreState = false)
326{
327 return elemSol[0];
328}
329
330
344template< class Element, class FVElementGeometry, class PrimaryVariables>
345PrimaryVariables evalSolution(const Element& element,
346 const typename Element::Geometry& geometry,
348 const typename Element::Geometry::GlobalCoordinate& globalPos,
349 bool ignoreState = false)
350{
351 return elemSol[0];
352}
353
354} // namespace Dumux
355
356#endif
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
Definition: adapt.hh:17
Type traits to be used with matrix types.
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:21