version 3.8
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-FileCopyrightInfo: 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
81template<class Element, class GridGeometry, class CVFEElemSol>
82typename CVFEElemSol::PrimaryVariables
83evalCVFESolution(const Element& element,
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)
89{
90 // determine if all states are the same at all vertices
91 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
92 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
93
95 {
96 using Scalar = typename CVFEElemSol::PrimaryVariables::value_type;
97
98 // interpolate the solution
99 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
100
101 // evaluate the shape functions at the scv center
102 const auto localPos = geometry.local(globalPos);
103 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
104 localBasis.evaluateFunction(localPos, shapeValues);
105
106 typename CVFEElemSol::PrimaryVariables result(0.0);
107 for (int i = 0; i < shapeValues.size(); ++i)
108 {
109 auto value = elemSol[i];
110 value *= shapeValues[i][0];
111 result += value;
112 }
113
114 // set an arbitrary state if the model requires a state (models constexpr if)
115 if constexpr (HasState{})
116 if (!ignoreState)
117 result.setState(elemSol[0].state());
118
119 return result;
120 }
121 else
122 {
123 static bool warnedAboutUsingMinDist = false;
124 if (!warnedAboutUsingMinDist)
125 {
126 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
127 << "\nbecause not all states are equal and ignoreState is false!"
128 << std::endl;
129 warnedAboutUsingMinDist = true;
130 }
131
132 return Detail::minDistVertexSol(geometry, globalPos, elemSol);
133 }
134}
135
136} // end namespace Detail
137
151template<class Element, class FVElementGeometry, class PrimaryVariables>
152PrimaryVariables evalSolution(const Element& element,
153 const typename Element::Geometry& geometry,
154 const typename FVElementGeometry::GridGeometry& gridGeometry,
156 const typename Element::Geometry::GlobalCoordinate& globalPos,
157 bool ignoreState = false)
158{
159 return Detail::evalCVFESolution(element, geometry, gridGeometry, elemSol, globalPos, ignoreState);
160}
161
176template<class Element, class FVElementGeometry, class PrimaryVariables>
177PrimaryVariables evalSolution(const Element& element,
178 const typename Element::Geometry& geometry,
180 const typename Element::Geometry::GlobalCoordinate& globalPos,
181 bool ignoreState = false)
182{
183 static_assert(FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::box, "Only implemented for the Box method");
184
185 // determine if all states are the same at all vertices
186 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
187 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
188
189 if (allStatesEqual)
190 {
191 using Scalar = typename PrimaryVariables::value_type;
192 using CoordScalar = typename Element::Geometry::GlobalCoordinate::value_type;
193 static constexpr int dim = Element::Geometry::mydimension;
194
195 // The box scheme always uses linear Ansatz functions
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();
200
201 // evaluate the shape functions at the scv center
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);
206
207 PrimaryVariables result(0.0);
208 for (int i = 0; i < geometry.corners(); ++i)
209 {
210 auto value = elemSol[i];
211 value *= shapeValues[i][0];
212 result += value;
213 }
214
215 // set an arbitrary state if the model requires a state (models constexpr if)
216 if constexpr (HasState{})
217 if (!ignoreState)
218 result.setState(elemSol[0].state());
219
220 return result;
221 }
222 else
223 {
224 static bool warnedAboutUsingMinDist = false;
225 if (!warnedAboutUsingMinDist)
226 {
227 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
228 << "\nbecause not all states are equal and ignoreState is false!"
229 << std::endl;
230 warnedAboutUsingMinDist = true;
231 }
232
233 return Detail::minDistVertexSol(geometry, globalPos, elemSol);
234 }
235}
236
249template<class Element, class FVElementGeometry, class PrimaryVariables>
250PrimaryVariables evalSolution(const Element& element,
251 const typename Element::Geometry& geometry,
252 const typename FVElementGeometry::GridGeometry& gridGeometry,
254 const typename Element::Geometry::GlobalCoordinate& globalPos,
255 bool ignoreState = false)
256{
257 return elemSol[0];
258}
259
273template< class Element, class FVElementGeometry, class PrimaryVariables>
274PrimaryVariables evalSolution(const Element& element,
275 const typename Element::Geometry& geometry,
277 const typename Element::Geometry::GlobalCoordinate& globalPos,
278 bool ignoreState = false)
279{
280 return elemSol[0];
281}
282
283} // namespace Dumux
284
285#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: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
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