free functions for the evaluation of primary variable gradients inside elements. More...
#include <dune/localfunctions/lagrange/pqkfactory.hh>
#include <dune/common/exceptions.hh>
#include <dumux/common/typetraits/state.hh>
#include <dumux/common/typetraits/isvalid.hh>
#include <dumux/discretization/box/elementsolution.hh>
#include <dumux/discretization/cellcentered/elementsolution.hh>
#include <dumux/discretization/facecentered/diamond/elementsolution.hh>
#include <dumux/discretization/pq1bubble/elementsolution.hh>
#include "evalsolution.hh"
Go to the source code of this file.
free functions for the evaluation of primary variable gradients inside elements.
Namespaces | |
namespace | Dumux |
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2. | |
namespace | Dumux::Detail |
Distance implementation details. | |
Functions | |
template<class Element , class GridGeometry , class CVFEElemSol > | |
auto | Dumux::Detail::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. More... | |
template<class Element , class FVElementGeometry , class PrimaryVariables > | |
auto | Dumux::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. More... | |
template<class Element , class FVElementGeometry , class PrimaryVariables > | |
auto | Dumux::evalGradients (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const FaceCenteredDiamondElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false) |
Evaluates the gradient of a given diamond scheme element solution to a given global position. More... | |
template<class Element , class FVElementGeometry , class PrimaryVariables > | |
auto | Dumux::evalGradients (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const PQ1BubbleElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false) |
Evaluates the gradient of a given pq1bubble scheme element solution to a given global position. More... | |
template<class Element , class FVElementGeometry , class PrimaryVariables > | |
Dune::FieldVector< typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension > | Dumux::evalGradients (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CCElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos) |
Evaluates the gradient of a given CCElementSolution to a given global position. This function is only here for (compilation) compatibility reasons with the box scheme. The solution within the control volumes is constant and thus gradients are zero. One can compute gradients towards the sub-control volume faces after reconstructing the solution on the faces. More... | |