3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
evalgradients.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_GRADIENTS_HH
25#define DUMUX_DISCRETIZATION_EVAL_GRADIENTS_HH
26
27#include <dune/localfunctions/lagrange/pqkfactory.hh>
28#include <dune/common/exceptions.hh>
29
36
37#include "evalsolution.hh"
38
39namespace Dumux {
40
41// some implementation details
42namespace Detail{
58template<class Element, class GridGeometry, class CVFEElemSol>
59auto evalCVFEGradients(const Element& element,
60 const typename Element::Geometry& geometry,
61 const GridGeometry& gridGeometry,
62 const CVFEElemSol& elemSol,
63 const typename Element::Geometry::GlobalCoordinate& globalPos,
64 bool ignoreState = false)
65{
66 // determine if all states are the same at all vertices
67 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
68 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
69
71 {
72 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
73
74 // evaluate gradients using the local finite element basis
75 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
76
77 // evaluate the shape function gradients at the scv center
78 using ShapeJacobian = typename std::decay_t< decltype(localBasis) >::Traits::JacobianType;
79 const auto localPos = geometry.local(globalPos);
80 std::vector< ShapeJacobian > shapeJacobian;
81 localBasis.evaluateJacobian(localPos, shapeJacobian);
82
83 // the inverse transposed of the jacobian matrix
84 const auto jacInvT = geometry.jacobianInverseTransposed(localPos);
85
86 // interpolate the gradients
87 Dune::FieldVector<GlobalPosition, CVFEElemSol::PrimaryVariables::dimension> result( GlobalPosition(0.0) );
88 for (int i = 0; i < shapeJacobian.size(); ++i)
89 {
90 // the global shape function gradient
91 GlobalPosition gradN;
92 jacInvT.mv(shapeJacobian[i][0], gradN);
93
94 // add gradient to global privar gradients
95 for (unsigned int pvIdx = 0; pvIdx < CVFEElemSol::PrimaryVariables::dimension; ++pvIdx)
96 {
97 GlobalPosition tmp(gradN);
98 tmp *= elemSol[i][pvIdx];
99 result[pvIdx] += tmp;
100 }
101 }
102
103 return result;
104 }
105 else
106 {
107 DUNE_THROW(Dune::NotImplemented, "Element dofs have different phase states. Enforce calculation by setting ignoreState to true.");
108 }
109}
110
111} // end namespace Detail
112
124template<class Element, class FVElementGeometry, class PrimaryVariables>
125auto evalGradients(const Element& element,
126 const typename Element::Geometry& geometry,
127 const typename FVElementGeometry::GridGeometry& gridGeometry,
129 const typename Element::Geometry::GlobalCoordinate& globalPos,
130 bool ignoreState = false)
131{
132 return Detail::evalCVFEGradients(element, geometry, gridGeometry, elemSol, globalPos, ignoreState);
133}
134
146template<class Element, class FVElementGeometry, class PrimaryVariables>
147auto evalGradients(const Element& element,
148 const typename Element::Geometry& geometry,
149 const typename FVElementGeometry::GridGeometry& gridGeometry,
151 const typename Element::Geometry::GlobalCoordinate& globalPos,
152 bool ignoreState = false)
153{
154 return Detail::evalCVFEGradients(element, geometry, gridGeometry, elemSol, globalPos, ignoreState);
155}
156
168template<class Element, class FVElementGeometry, class PrimaryVariables>
169auto evalGradients(const Element& element,
170 const typename Element::Geometry& geometry,
171 const typename FVElementGeometry::GridGeometry& gridGeometry,
173 const typename Element::Geometry::GlobalCoordinate& globalPos,
174 bool ignoreState = false)
175{
176 return Detail::evalCVFEGradients(element, geometry, gridGeometry, elemSol, globalPos, ignoreState);
177}
178
198template<class Element, class FVElementGeometry, class PrimaryVariables>
199Dune::FieldVector<typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension>
200evalGradients(const Element& element,
201 const typename Element::Geometry& geometry,
202 const typename FVElementGeometry::GridGeometry& gridGeometry,
204 const typename Element::Geometry::GlobalCoordinate& globalPos)
205{ DUNE_THROW(Dune::NotImplemented, "General gradient evaluation for cell-centered methods"); }
206
207} // namespace Dumux
208
209#endif
A helper function for class member function introspection.
Type traits to be used with matrix types.
free functions for the evaluation of primary variables inside elements.
auto evalGradients(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)
Evaluates the gradient of a given box element solution to a given global position.
Definition: evalgradients.hh:125
auto evalCVFEGradients(const Element &element, const typename Element::Geometry &geometry, const GridGeometry &gridGeometry, const CVFEElemSol &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Evaluates the gradient of a control-volume finite element solution to a given global position.
Definition: evalgradients.hh:59
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
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
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.