version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_EVAL_GRADIENTS_HH
13#define DUMUX_DISCRETIZATION_EVAL_GRADIENTS_HH
14
15#include <dune/common/exceptions.hh>
16
21
22#include "evalsolution.hh"
23
24namespace Dumux {
25
26// some implementation details
27namespace Detail{
43template<class Element, class GridGeometry, class CVFEElemSol>
44auto evalCVFEGradients(const Element& element,
45 const typename Element::Geometry& geometry,
46 const GridGeometry& gridGeometry,
47 const CVFEElemSol& elemSol,
48 const typename Element::Geometry::GlobalCoordinate& globalPos,
49 bool ignoreState = false)
50{
51 // determine if all states are the same at all vertices
52 using HasState = decltype(isValid(Detail::hasState())(elemSol[0]));
53 bool allStatesEqual = ignoreState || Detail::allStatesEqual(elemSol, HasState{});
54
56 {
57 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
58
59 // evaluate gradients using the local finite element basis
60 const auto& localBasis = gridGeometry.feCache().get(geometry.type()).localBasis();
61
62 // evaluate the shape function gradients at the scv center
63 using ShapeJacobian = typename std::decay_t< decltype(localBasis) >::Traits::JacobianType;
64 const auto localPos = geometry.local(globalPos);
65 std::vector< ShapeJacobian > shapeJacobian;
66 localBasis.evaluateJacobian(localPos, shapeJacobian);
67
68 // the inverse transposed of the jacobian matrix
69 const auto jacInvT = geometry.jacobianInverseTransposed(localPos);
70
71 // interpolate the gradients
72 Dune::FieldVector<GlobalPosition, CVFEElemSol::PrimaryVariables::dimension> result( GlobalPosition(0.0) );
73 for (int i = 0; i < shapeJacobian.size(); ++i)
74 {
75 // the global shape function gradient
76 GlobalPosition gradN;
77 jacInvT.mv(shapeJacobian[i][0], gradN);
78
79 // add gradient to global privar gradients
80 for (unsigned int pvIdx = 0; pvIdx < CVFEElemSol::PrimaryVariables::dimension; ++pvIdx)
81 {
82 GlobalPosition tmp(gradN);
83 tmp *= elemSol[i][pvIdx];
84 result[pvIdx] += tmp;
85 }
86 }
87
88 return result;
89 }
90 else
91 {
92 DUNE_THROW(Dune::NotImplemented, "Element dofs have different phase states. Enforce calculation by setting ignoreState to true.");
93 }
94}
95
96} // end namespace Detail
97
109template<class Element, class FVElementGeometry, class PrimaryVariables>
110auto evalGradients(const Element& element,
111 const typename Element::Geometry& geometry,
112 const typename FVElementGeometry::GridGeometry& gridGeometry,
114 const typename Element::Geometry::GlobalCoordinate& globalPos,
115 bool ignoreState = false)
116{
117 return Detail::evalCVFEGradients(element, geometry, gridGeometry, elemSol, globalPos, ignoreState);
118}
119
139template<class Element, class FVElementGeometry, class PrimaryVariables>
140Dune::FieldVector<typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension>
141evalGradients(const Element& element,
142 const typename Element::Geometry& geometry,
143 const typename FVElementGeometry::GridGeometry& gridGeometry,
145 const typename Element::Geometry::GlobalCoordinate& globalPos)
146{ DUNE_THROW(Dune::NotImplemented, "General gradient evaluation for cell-centered methods"); }
147
148} // namespace Dumux
149
150#endif
The local element solution class for cell-centered methods.
The element solution vector.
Definition: cellcentered/elementsolution.hh:28
The element solution vector.
Definition: cvfe/elementsolution.hh:28
The local element solution class for control-volume finite element methods.
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 CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Evaluates the gradient of a given CVFE element solution to a given global position.
Definition: evalgradients.hh:110
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:44
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:81
A helper function for class member function introspection.
bool allStatesEqual(const ElementSolution &elemSol, std::true_type hasState)
returns true if all states in an element solution are the same
Definition: evalsolution.hh:35
Definition: adapt.hh:17
Type traits to be used with matrix types.
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:21