3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/cvfe/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 * 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_CVFE_FLUXVARIABLES_CACHE_HH
25#define DUMUX_DISCRETIZATION_CVFE_FLUXVARIABLES_CACHE_HH
26
27#include <dune/common/fvector.hh>
28
29namespace Dumux {
30
37template< class Scalar, class GridGeometry >
39{
40 using GridView = typename GridGeometry::GridView;
41 using Element = typename GridView::template Codim<0>::Entity;
42 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
43
44 using FVElementGeometry = typename GridGeometry::LocalView;
45 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
46
47 static const int dim = GridView::dimension;
48 static const int dimWorld = GridView::dimensionworld;
49
50 using CoordScalar = typename GridView::ctype;
51 using FeLocalBasis = typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
52 using ShapeJacobian = typename FeLocalBasis::Traits::JacobianType;
53 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
54 using JacobianInverseTransposed = typename Element::Geometry::JacobianInverseTransposed;
55
56public:
57
59 template< class Problem, class ElementVolumeVariables >
60 void update(const Problem& problem,
61 const Element& element,
62 const FVElementGeometry& fvGeometry,
63 const ElementVolumeVariables& elemVolVars,
64 const SubControlVolumeFace& scvf)
65 {
66 update(problem, element, fvGeometry, elemVolVars, scvf.ipGlobal());
67 }
68
70 template< class Problem, class ElementVolumeVariables >
71 void update(const Problem& problem,
72 const Element& element,
73 const FVElementGeometry& fvGeometry,
74 const ElementVolumeVariables& elemVolVars,
75 const GlobalPosition& globalPos)
76 {
77 const auto geometry = element.geometry();
78 const auto& localBasis = fvGeometry.feLocalBasis();
79
80 // evaluate shape functions and gradients at the integration point
81 ipGlobal_ = globalPos;
82 const auto ipLocal = geometry.local(ipGlobal_);
83 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal);
84 localBasis.evaluateJacobian(ipLocal, shapeJacobian_);
85 localBasis.evaluateFunction(ipLocal, shapeValues_); // shape values for rho
86
87 // compute the gradN at for every scv/dof
88 gradN_.resize(fvGeometry.numScv());
89 for (const auto& scv: scvs(fvGeometry))
90 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
91 }
92
94 const GlobalPosition& ipGlobal() const { return ipGlobal_; }
96 const std::vector<ShapeJacobian>& shapeJacobian() const { return shapeJacobian_; }
98 const std::vector<ShapeValue>& shapeValues() const { return shapeValues_; }
100 const JacobianInverseTransposed& jacInvT() const { return jacInvT_; }
102 const GlobalPosition& gradN(unsigned int scvIdxInElement) const { return gradN_[scvIdxInElement]; }
103
104private:
105 GlobalPosition ipGlobal_;
106 std::vector<GlobalPosition> gradN_;
107 std::vector<ShapeJacobian> shapeJacobian_;
108 std::vector<ShapeValue> shapeValues_;
109 JacobianInverseTransposed jacInvT_;
110};
111
112} // end namespace Dumux
113
114#endif
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Flux variables cache class for control-volume finite element schemes. For control-volume finite eleme...
Definition: discretization/cvfe/fluxvariablescache.hh:39
const JacobianInverseTransposed & jacInvT() const
returns inverse transposed jacobian at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:100
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
update the cache for an scvf
Definition: discretization/cvfe/fluxvariablescache.hh:60
const std::vector< ShapeValue > & shapeValues() const
returns the shape function values at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:98
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
returns the shape function gradients in global coordinates at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:102
const GlobalPosition & ipGlobal() const
returns the global position for which this cache has been updated
Definition: discretization/cvfe/fluxvariablescache.hh:94
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/fluxvariablescache.hh:71
const std::vector< ShapeJacobian > & shapeJacobian() const
returns the shape function gradients in local coordinates at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:96