3.6-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#include <memory>
31
32#include <dune/localfunctions/lagrange/pqkfactory.hh>
33
40
41namespace Dumux {
42
43// some implementation details
44namespace Detail {
45
47template<class ElementSolution>
48bool allStatesEqual(const ElementSolution& elemSol, std::true_type hasState)
49{
50 // determine if all states are the same at all vertices
51 auto firstState = elemSol[0].state();
52 for (std::size_t i = 1; i < elemSol.size(); ++i)
53 if (elemSol[i].state() != firstState)
54 return false;
55
56 return true;
57}
58
60template<class ElementSolution>
61bool allStatesEqual(const ElementSolution& elemSol, std::false_type hasState)
62{
63 return true;
64}
65
67template<class Geometry, class ElementSolution>
68auto minDistVertexSol(const Geometry& geometry, const typename Geometry::GlobalCoordinate& globalPos,
69 const ElementSolution& elemSol)
70{
71 // calculate the distances from the evaluation point to the vertices
72 std::vector<typename Geometry::ctype> distances(geometry.corners());
73 for (int i = 0; i < geometry.corners(); ++i)
74 distances[i] = (geometry.corner(i) - globalPos).two_norm2();
75
76 // determine the minimum distance and corresponding vertex
77 auto minDistanceIt = std::min_element(distances.begin(), distances.end());
78 return elemSol[std::distance(distances.begin(), minDistanceIt)];
79}
80
94template<class Element, class GridGeometry, class CVFEElemSol>
95typename CVFEElemSol::PrimaryVariables
96evalCVFESolution(const Element& element,
97 const typename Element::Geometry& geometry,
98 const GridGeometry& gridGeometry,
99 const CVFEElemSol& elemSol,
100 const typename Element::Geometry::GlobalCoordinate& globalPos,
101 bool ignoreState = false)
102{
103 // determine if all states are the same at all vertices
104 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
105 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
106
107 if (allStatesEqual)
108 {
109 using Scalar = typename CVFEElemSol::PrimaryVariables::value_type;
110
111 // interpolate the solution
112 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
113
114 // evaluate the shape functions at the scv center
115 const auto localPos = geometry.local(globalPos);
116 std::vector< Dune::FieldVector<Scalar, 1> > shapeValues;
117 localBasis.evaluateFunction(localPos, shapeValues);
118
119 typename CVFEElemSol::PrimaryVariables result(0.0);
120 for (int i = 0; i < shapeValues.size(); ++i)
121 {
122 auto value = elemSol[i];
123 value *= shapeValues[i][0];
124 result += value;
125 }
126
127 // set an arbitrary state if the model requires a state (models constexpr if)
128 if constexpr (HasState{})
129 if (!ignoreState)
130 result.setState(elemSol[0].state());
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
149} // end namespace Detail
150
164template<class Element, class FVElementGeometry, class PrimaryVariables>
165PrimaryVariables evalSolution(const Element& element,
166 const typename Element::Geometry& geometry,
167 const typename FVElementGeometry::GridGeometry& gridGeometry,
169 const typename Element::Geometry::GlobalCoordinate& globalPos,
170 bool ignoreState = false)
171{
172 return Detail::evalCVFESolution(element, geometry, gridGeometry, elemSol, globalPos, ignoreState);
173}
174
189template<class Element, class FVElementGeometry, class PrimaryVariables>
190PrimaryVariables evalSolution(const Element& element,
191 const typename Element::Geometry& geometry,
193 const typename Element::Geometry::GlobalCoordinate& globalPos,
194 bool ignoreState = false)
195{
196 // determine if all states are the same at all vertices
197 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
198 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
199
200 if (allStatesEqual)
201 {
202 using Scalar = typename PrimaryVariables::value_type;
203 using CoordScalar = typename Element::Geometry::GlobalCoordinate::value_type;
204 static constexpr int dim = Element::Geometry::mydimension;
205
206 // The box scheme always uses linear Ansatz functions
207 using FEFactory = Dune::PQkLocalFiniteElementFactory<CoordScalar, Scalar, dim, 1>;
208 using FiniteElement = typename FEFactory::FiniteElementType;
209 std::unique_ptr<FiniteElement> fe(FEFactory::create(geometry.type()));
210 const auto& localBasis = fe->localBasis();
211
212 // evaluate the shape functions at the scv center
213 const auto localPos = geometry.local(globalPos);
214 using ShapeValue = typename FiniteElement::Traits::LocalBasisType::Traits::RangeType;
215 std::vector< ShapeValue > shapeValues;
216 localBasis.evaluateFunction(localPos, shapeValues);
217
218 PrimaryVariables result(0.0);
219 for (int i = 0; i < geometry.corners(); ++i)
220 {
221 auto value = elemSol[i];
222 value *= shapeValues[i][0];
223 result += value;
224 }
225
226 // set an arbitrary state if the model requires a state (models constexpr if)
227 if constexpr (HasState{})
228 if (!ignoreState)
229 result.setState(elemSol[0].state());
230
231 return result;
232 }
233 else
234 {
235 static bool warnedAboutUsingMinDist = false;
236 if (!warnedAboutUsingMinDist)
237 {
238 std::cout << "Warning: Using nearest-neighbor interpolation in evalSolution"
239 << "\nbecause not all states are equal and ignoreState is false!"
240 << std::endl;
241 warnedAboutUsingMinDist = true;
242 }
243
244 return Detail::minDistVertexSol(geometry, globalPos, elemSol);
245 }
246}
247
260template<class Element, class FVElementGeometry, class PrimaryVariables>
261PrimaryVariables evalSolution(const Element& element,
262 const typename Element::Geometry& geometry,
263 const typename FVElementGeometry::GridGeometry& gridGeometry,
265 const typename Element::Geometry::GlobalCoordinate& globalPos,
266 bool ignoreState = false)
267{
268 return elemSol[0];
269}
270
284template< class Element, class FVElementGeometry, class PrimaryVariables>
285PrimaryVariables evalSolution(const Element& element,
286 const typename Element::Geometry& geometry,
288 const typename Element::Geometry::GlobalCoordinate& globalPos,
289 bool ignoreState = false)
290{
291 return elemSol[0];
292}
293
307template<class Element, class FVElementGeometry, class PrimaryVariables>
308PrimaryVariables evalSolution(const Element& element,
309 const typename Element::Geometry& geometry,
310 const typename FVElementGeometry::GridGeometry& gridGeometry,
312 const typename Element::Geometry::GlobalCoordinate& globalPos,
313 bool ignoreState = false)
314{
315 return Detail::evalCVFESolution(element, geometry, gridGeometry, elemSol, globalPos, ignoreState);
316}
317
331template<class Element, class FVElementGeometry, class PrimaryVariables>
332PrimaryVariables evalSolution(const Element& element,
333 const typename Element::Geometry& geometry,
334 const typename FVElementGeometry::GridGeometry& gridGeometry,
336 const typename Element::Geometry::GlobalCoordinate& globalPos,
337 bool ignoreState = false)
338{
339 return Detail::evalCVFESolution(element, geometry, gridGeometry, elemSol, globalPos, ignoreState);
340}
341
342} // namespace Dumux
343
344#endif
A helper function for class member function introspection.
Type traits to be used with matrix types.
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:294
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:165
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:96
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
auto minDistVertexSol(const Geometry &geometry, const typename Geometry::GlobalCoordinate &globalPos, const ElementSolution &elemSol)
return the solution at the closest dof
Definition: evalsolution.hh:68
bool allStatesEqual(const ElementSolution &elemSol, std::true_type hasState)
returns true if all states in an element solution are the same
Definition: evalsolution.hh:48
bool allStatesEqual(const ElementSolution &elemSol, std::false_type hasState)
overload if the solution is stateless
Definition: evalsolution.hh:61
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:33
The element solution vector.
Definition: box/elementsolution.hh:39
The element solution vector.
Definition: cellcentered/elementsolution.hh:40
The global face variables class for face-centered diamond models.
Definition: facecentered/diamond/elementsolution.hh:42
The element solution vector.
Definition: pq1bubble/elementsolution.hh:39
The local element solution class for the box method.
The local element solution class for cell-centered methods.
The local element solution class for the pq1bubble method.