3.4
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
34
35#include "evalsolution.hh"
36
37namespace Dumux {
38
54template<class Element, class FVElementGeometry, class PrimaryVariables>
55auto evalGradients(const Element& element,
56 const typename Element::Geometry& geometry,
57 const typename FVElementGeometry::GridGeometry& gridGeometry,
59 const typename Element::Geometry::GlobalCoordinate& globalPos,
60 bool ignoreState = false)
61{
62 // determine if all states are the same at all vertices
63 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
64 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
65
67 {
68 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
69
70 // evaluate gradients using the local finite element basis
71 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
72
73 // evaluate the shape function gradients at the scv center
74 using ShapeJacobian = typename std::decay_t< decltype(localBasis) >::Traits::JacobianType;
75 const auto localPos = geometry.local(globalPos);
76 std::vector< ShapeJacobian > shapeJacobian;
77 localBasis.evaluateJacobian(localPos, shapeJacobian);
78
79 // the inverse transposed of the jacobian matrix
80 const auto jacInvT = geometry.jacobianInverseTransposed(localPos);
81
82 // interpolate the gradients
83 Dune::FieldVector<GlobalPosition, PrimaryVariables::dimension> result( GlobalPosition(0.0) );
84 for (int i = 0; i < element.subEntities(Element::Geometry::mydimension); ++i)
85 {
86 // the global shape function gradient
87 GlobalPosition gradN;
88 jacInvT.mv(shapeJacobian[i][0], gradN);
89
90 // add gradient to global privar gradients
91 for (unsigned int pvIdx = 0; pvIdx < PrimaryVariables::dimension; ++pvIdx)
92 {
93 GlobalPosition tmp(gradN);
94 tmp *= elemSol[i][pvIdx];
95 result[pvIdx] += tmp;
96 }
97 }
98
99 return result;
100 }
101 else
102 {
103 DUNE_THROW(Dune::NotImplemented, "Element vertices have different phase states. Enforce calculation by setting ignoreState to true.");
104 }
105}
106
126template<class Element, class FVElementGeometry, class PrimaryVariables>
127Dune::FieldVector<typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension>
128evalGradients(const Element& element,
129 const typename Element::Geometry& geometry,
130 const typename FVElementGeometry::GridGeometry& gridGeometry,
132 const typename Element::Geometry::GlobalCoordinate& globalPos)
133{ DUNE_THROW(Dune::NotImplemented, "General gradient evaluation for cell-centered methods"); }
134
135} // namespace Dumux
136
137#endif
Type traits to be used with matrix types.
A helper function for class member function introspection.
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:55
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
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
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.