3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
25#define DUMUX_DISCRETIZATION_EVAL_SOLUTION_HH
26
27#include <iterator>
28#include <algorithm>
29#include <type_traits>
30
31#include <dune/localfunctions/lagrange/pqkfactory.hh>
32
37
38namespace Dumux {
39
40// some implementation details
41namespace Detail {
42
44template<class ElementSolution>
45bool allStatesEqual(const ElementSolution& elemSol, std::true_type hasState)
46{
47 // determine if all states are the same at all vertices
48 auto firstState = elemSol[0].state();
49 for (std::size_t i = 1; i < elemSol.size(); ++i)
50 if (elemSol[i].state() != firstState)
51 return false;
52
53 return true;
54}
55
57template<class ElementSolution>
58bool allStatesEqual(const ElementSolution& elemSol, std::false_type hasState)
59{
60 return true;
61}
62
64template<class Geometry, class ElementSolution>
65auto minDistVertexSol(const Geometry& geometry, const typename Geometry::GlobalCoordinate& globalPos,
66 const ElementSolution& elemSol)
67{
68 // calculate the distances from the evaluation point to the vertices
69 std::vector<typename Geometry::ctype> distances(geometry.corners());
70 for (int i = 0; i < geometry.corners(); ++i)
71 distances[i] = (geometry.corner(i) - globalPos).two_norm2();
72
73 // determine the minimum distance and corresponding vertex
74 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
75 return elemSol[std::distance(distances.begin(), minDistanceIt)];
76}
77
78} // end namespace Detail
79
93template<class Element, class FVElementGeometry, class PrimaryVariables>
94PrimaryVariables evalSolution(const Element& element,
95 const typename Element::Geometry& geometry,
96 const typename FVElementGeometry::GridGeometry& gridGeometry,
98 const typename Element::Geometry::GlobalCoordinate& globalPos,
99 bool ignoreState = false)
100{
101 // determine if all states are the same at all vertices
102 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
103 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
104
105 if (allStatesEqual)
106 {
107 using Scalar = typename PrimaryVariables::value_type;
108
109 // interpolate the solution
110 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
111
112 // evaluate the shape functions at the scv center
113 const auto localPos = geometry.local(globalPos);
114 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
115 localBasis.evaluateFunction(localPos, shapeValues);
116
117 PrimaryVariables result(0.0);
118 for (int i = 0; i < geometry.corners(); ++i)
119 {
120 auto value = elemSol[i];
121 value *= shapeValues[i][0];
122 result += value;
123 }
124
125 // set an arbitrary state if the model requires a state (models constexpr if)
126 if constexpr (HasState{})
127 if (!ignoreState)
128 result.setState(elemSol[0].state());
129
130 return result;
131 }
132 else
133 {
134 static bool warnedAboutUsingMinDist = false;
135 if (!warnedAboutUsingMinDist)
136 {
137 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
138 << "\nbecause not all states are equal and ignoreState is false!"
139 << std::endl;
140 warnedAboutUsingMinDist = true;
141 }
142
143 return Detail::minDistVertexSol(geometry, globalPos, elemSol);
144 }
145}
146
161template<class Element, class FVElementGeometry, class PrimaryVariables>
162PrimaryVariables evalSolution(const Element& element,
163 const typename Element::Geometry& geometry,
165 const typename Element::Geometry::GlobalCoordinate& globalPos,
166 bool ignoreState = false)
167{
168 // determine if all states are the same at all vertices
169 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
170 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
171
172 if (allStatesEqual)
173 {
174 using Scalar = typename PrimaryVariables::value_type;
175 using CoordScalar = typename Element::Geometry::GlobalCoordinate::value_type;
176 static constexpr int dim = Element::Geometry::mydimension;
177
179 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
180 using ShapeValue = typename FeCache::FiniteElementType::Traits::LocalBasisType::Traits::RangeType;
181
182 // obtain local finite element basis
183 FeCache feCache;
184 const auto& localBasis = feCache.get(geometry.type()).localBasis();
185
186 // evaluate the shape functions at the scv center
187 const auto localPos = geometry.local(globalPos);
188 std::vector< ShapeValue > shapeValues;
189 localBasis.evaluateFunction(localPos, shapeValues);
190
191 PrimaryVariables result(0.0);
192 for (int i = 0; i < geometry.corners(); ++i)
193 {
194 auto value = elemSol[i];
195 value *= shapeValues[i][0];
196 result += value;
197 }
198
199 // set an arbitrary state if the model requires a state (models constexpr if)
200 if constexpr (HasState{})
201 if (!ignoreState)
202 result.setState(elemSol[0].state());
203
204 return result;
205 }
206 else
207 {
208 static bool warnedAboutUsingMinDist = false;
209 if (!warnedAboutUsingMinDist)
210 {
211 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
212 << "\nbecause not all states are equal and ignoreState is false!"
213 << std::endl;
214 warnedAboutUsingMinDist = true;
215 }
216
217 return Detail::minDistVertexSol(geometry, globalPos, elemSol);
218 }
219}
220
233template<class Element, class FVElementGeometry, class PrimaryVariables>
234PrimaryVariables evalSolution(const Element& element,
235 const typename Element::Geometry& geometry,
236 const typename FVElementGeometry::GridGeometry& gridGeometry,
238 const typename Element::Geometry::GlobalCoordinate& globalPos,
239 bool ignoreState = false)
240{
241 return elemSol[0];
242}
243
257template< class Element, class FVElementGeometry, class PrimaryVariables>
258PrimaryVariables evalSolution(const Element& element,
259 const typename Element::Geometry& geometry,
261 const typename Element::Geometry::GlobalCoordinate& globalPos,
262 bool ignoreState = false)
263{
264 return elemSol[0];
265}
266
267} // namespace Dumux
268
269#endif
Type traits to be used with matrix types.
A helper function for class member function introspection.
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: geometry/distance.hh:138
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:94
Definition: adapt.hh:29
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
auto minDistVertexSol(const Geometry &geometry, const typename Geometry::GlobalCoordinate &globalPos, const ElementSolution &elemSol)
return the solution at the closest dof
Definition: evalsolution.hh:65
bool allStatesEqual(const ElementSolution &elemSol, std::true_type hasState)
returns true if all states in an element solution are the same
Definition: evalsolution.hh:45
bool allStatesEqual(const ElementSolution &elemSol, std::false_type hasState)
overload if the solution is stateless
Definition: evalsolution.hh:58
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:32
The element solution vector.
Definition: box/elementsolution.hh:39
The element solution vector.
Definition: cellcentered/elementsolution.hh:40
The local element solution class for the box method.
The local element solution class for cell-centered methods.