3.1-git
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/common/hybridutilities.hh>
32#include <dune/localfunctions/lagrange/pqkfactory.hh>
33
38
39namespace Dumux {
40
41// some implementation details
42namespace Detail {
43
45template<class ElementSolution>
46bool allStatesEqual(const ElementSolution& elemSol, std::true_type hasState)
47{
48 // determine if all states are the same at all vertices
49 auto firstState = elemSol[0].state();
50 for (std::size_t i = 1; i < elemSol.size(); ++i)
51 if (elemSol[i].state() != firstState)
52 return false;
53
54 return true;
55}
56
58template<class ElementSolution>
59bool allStatesEqual(const ElementSolution& elemSol, std::false_type hasState)
60{
61 return true;
62}
63
65template<class Geometry, class ElementSolution>
66auto minDistVertexSol(const Geometry& geometry, const typename Geometry::GlobalCoordinate& globalPos,
67 const ElementSolution& elemSol)
68{
69 // calculate the distances from the evaluation point to the vertices
70 std::vector<typename Geometry::ctype> distances(geometry.corners());
71 for (int i = 0; i < geometry.corners(); ++i)
72 distances[i] = (geometry.corner(i) - globalPos).two_norm2();
73
74 // determine the minimum distance and corresponding vertex
75 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
76 return elemSol[std::distance(distances.begin(), minDistanceIt)];
77}
78
79} // end namespace Detail
80
94template<class Element, class FVElementGeometry, class PrimaryVariables>
95PrimaryVariables evalSolution(const Element& element,
96 const typename Element::Geometry& geometry,
97 const typename FVElementGeometry::GridGeometry& gridGeometry,
99 const typename Element::Geometry::GlobalCoordinate& globalPos,
100 bool ignoreState = false)
101{
102 // determine if all states are the same at all vertices
103 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
104 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
105
106 if (allStatesEqual)
107 {
108 using Scalar = typename PrimaryVariables::value_type;
109
110 // interpolate the solution
111 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
112
113 // evaluate the shape functions at the scv center
114 const auto localPos = geometry.local(globalPos);
115 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
116 localBasis.evaluateFunction(localPos, shapeValues);
117
118 PrimaryVariables result(0.0);
119 for (int i = 0; i < geometry.corners(); ++i)
120 {
121 auto value = elemSol[i];
122 value *= shapeValues[i][0];
123 result += value;
124 }
125
126 // set an arbitrary state if the model requires a state (models constexpr if)
127 Dune::Hybrid::ifElse(HasState{}, [&](auto&& evalLazy){
128 if (!ignoreState)
129 evalLazy(result).setState(evalLazy(elemSol)[0].state());
130 });
131
132 return result;
133 }
134 else
135 {
136 static bool warnedAboutUsingMinDist = false;
137 if (!warnedAboutUsingMinDist)
138 {
139 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
140 << "\nbecause not all states are equal and ignoreState is false!"
141 << std::endl;
142 warnedAboutUsingMinDist = true;
143 }
144
145 return Detail::minDistVertexSol(geometry, globalPos, elemSol);
146 }
147}
148
163template<class Element, class FVElementGeometry, class PrimaryVariables>
164PrimaryVariables evalSolution(const Element& element,
165 const typename Element::Geometry& geometry,
167 const typename Element::Geometry::GlobalCoordinate& globalPos,
168 bool ignoreState = false)
169{
170 // determine if all states are the same at all vertices
171 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
172 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
173
174 if (allStatesEqual)
175 {
176 using Scalar = typename PrimaryVariables::value_type;
177 using CoordScalar = typename Element::Geometry::GlobalCoordinate::value_type;
178 static constexpr int dim = Element::Geometry::mydimension;
179
181 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
182 using ShapeValue = typename FeCache::FiniteElementType::Traits::LocalBasisType::Traits::RangeType;
183
184 // obtain local finite element basis
185 FeCache feCache;
186 const auto& localBasis = feCache.get(geometry.type()).localBasis();
187
188 // evaluate the shape functions at the scv center
189 const auto localPos = geometry.local(globalPos);
190 std::vector< ShapeValue > shapeValues;
191 localBasis.evaluateFunction(localPos, shapeValues);
192
193 PrimaryVariables result(0.0);
194 for (int i = 0; i < geometry.corners(); ++i)
195 {
196 auto value = elemSol[i];
197 value *= shapeValues[i][0];
198 result += value;
199 }
200
201 // set an arbitrary state if the model requires a state (models constexpr if)
202 Dune::Hybrid::ifElse(HasState{}, [&](auto&& evalLazy){
203 if (!ignoreState)
204 evalLazy(result).setState(evalLazy(elemSol)[0].state());
205 });
206
207 return result;
208 }
209 else
210 {
211 static bool warnedAboutUsingMinDist = false;
212 if (!warnedAboutUsingMinDist)
213 {
214 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
215 << "\nbecause not all states are equal and ignoreState is false!"
216 << std::endl;
217 warnedAboutUsingMinDist = true;
218 }
219
220 return Detail::minDistVertexSol(geometry, globalPos, elemSol);
221 }
222}
223
236template<class Element, class FVElementGeometry, class PrimaryVariables>
237PrimaryVariables evalSolution(const Element& element,
238 const typename Element::Geometry& geometry,
239 const typename FVElementGeometry::GridGeometry& gridGeometry,
241 const typename Element::Geometry::GlobalCoordinate& globalPos,
242 bool ignoreState = false)
243{
244 return elemSol[0];
245}
246
260template< class Element, class FVElementGeometry, class PrimaryVariables>
261PrimaryVariables evalSolution(const Element& element,
262 const typename Element::Geometry& geometry,
264 const typename Element::Geometry::GlobalCoordinate& globalPos,
265 bool ignoreState = false)
266{
267 return elemSol[0];
268}
269
270} // namespace Dumux
271
272#endif
A helper function for class member function introspection.
Type traits to be used with matrix types.
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:95
make the local view function available whenever we use the grid geometry
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:66
bool allStatesEqual(const ElementSolution &elemSol, std::true_type hasState)
returns true if all states in an element solution are the same
Definition: evalsolution.hh:46
bool allStatesEqual(const ElementSolution &elemSol, std::false_type hasState)
overload if the solution is stateless
Definition: evalsolution.hh:59
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.