version 3.11-dev
Loading...
Searching...
No Matches
discretization/cvfe/hybrid/fluxvariablescache.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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_HYBRID_CVFE_FLUXVARIABLES_CACHE_HH
13#define DUMUX_DISCRETIZATION_HYBRID_CVFE_FLUXVARIABLES_CACHE_HH
14
15#include <vector>
16
17#include <dune/common/fvector.hh>
18#include <dumux/common/typetraits/localdofs_.hh>
19#include <dumux/common/concepts/ipdata_.hh>
21
22namespace Dumux {
23
30template< class Scalar, class GridGeometry >
32{
33 using GridView = typename GridGeometry::GridView;
34 using Element = typename GridView::template Codim<0>::Entity;
35 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
36 using LocalPosition = typename Element::Geometry::LocalCoordinate;
37
38 using FVElementGeometry = typename GridGeometry::LocalView;
39
40 static const int dim = GridView::dimension;
41 static const int dimWorld = GridView::dimensionworld;
42
43 using CoordScalar = typename GridView::ctype;
44 using FeLocalBasis = typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
45 using ShapeJacobian = typename FeLocalBasis::Traits::JacobianType;
46 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
47 using JacobianInverseTransposed = typename Element::Geometry::JacobianInverseTransposed;
48
49public:
51 static bool constexpr isSolDependent = false;
52
54 template< class Problem, class ElementVolumeVariables >
55 void update(const Problem& problem,
56 const Element& element,
57 const FVElementGeometry& fvGeometry,
58 const ElementVolumeVariables& elemVolVars,
59 const GlobalPosition& globalPos)
60 {
61 update_(fvGeometry, fvGeometry.elementGeometry().local(globalPos), globalPos);
62 }
63
65 template< class Problem, class ElementVolumeVariables, Concept::IpData IpData >
66 void update(const Problem& problem,
67 const Element& element,
68 const FVElementGeometry& fvGeometry,
69 const ElementVolumeVariables& elemVolVars,
70 const IpData& ipData)
71 {
72 update_(fvGeometry, ipData.local(), ipData.global());
73 }
74
76 const GlobalPosition& ipGlobal() const { return ipGlobal_; }
78 const LocalPosition& ipLocal() const { return ipLocal_; }
80 const std::vector<ShapeJacobian>& shapeJacobian() const { return shapeJacobian_; }
82 const std::vector<ShapeValue>& shapeValues() const { return shapeValues_; }
84 const JacobianInverseTransposed& jacInvT() const { return jacInvT_; }
86 const GlobalPosition& gradN(unsigned int scvIdxInElement) const { return gradN_[scvIdxInElement]; }
87
88private:
90 void update_(const FVElementGeometry& fvGeometry,
91 const LocalPosition& localPos,
92 const GlobalPosition& globalPos)
93 {
94 const auto& geometry = fvGeometry.elementGeometry();
95 const auto& localBasis = fvGeometry.feLocalBasis();
96
97 // evaluate shape functions and gradients at the interpolation point
98 ipGlobal_ = globalPos;
99 ipLocal_ = localPos;
100 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal_);
101 localBasis.evaluateJacobian(ipLocal_, shapeJacobian_);
102 localBasis.evaluateFunction(ipLocal_, shapeValues_); // shape values for rho
103
104 // compute the gradN at for every scv/dof
105 gradN_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
106 for (const auto& localDof: localDofs(fvGeometry))
107 jacInvT_.mv(shapeJacobian_[localDof.index()][0], gradN_[localDof.index()]);
108 }
109
110 std::vector<GlobalPosition> gradN_;
111 std::vector<ShapeJacobian> shapeJacobian_;
112 std::vector<ShapeValue> shapeValues_;
113 JacobianInverseTransposed jacInvT_;
114 GlobalPosition ipGlobal_;
115 LocalPosition ipLocal_;
116};
117
118} // end namespace Dumux
119
120#endif
Flux variables cache class for control-volume finite element schemes. For control-volume finite eleme...
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:32
static bool constexpr isSolDependent
whether the cache needs an update when the solution changes
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:51
const LocalPosition & ipLocal() const
returns the local position for which this cache has been updated
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:78
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
returns the shape function gradients in global coordinates at the interpolation point
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:86
const std::vector< ShapeJacobian > & shapeJacobian() const
returns the shape function gradients in local coordinates at the interpolation point
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:80
const JacobianInverseTransposed & jacInvT() const
returns inverse transposed jacobian at the interpolation point
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:84
const std::vector< ShapeValue > & shapeValues() const
returns the shape function values at the interpolation point
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:82
const GlobalPosition & ipGlobal() const
returns the global position for which this cache has been updated
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:76
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const GlobalPosition &globalPos)
update the cache for a given global position
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:55
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const IpData &ipData)
update the cache for interpolation point data
Definition discretization/cvfe/hybrid/fluxvariablescache.hh:66
Class representing dofs on elements for control-volume finite element schemes.
Definition adapt.hh:17
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition localdof.hh:50