version 3.10-dev
porousmediumflow/boxdfm/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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
15#define DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
16
17#include <dune/common/fvector.hh>
18
22
23namespace Dumux {
24
31template<class TypeTag>
33{
37 using GridView = typename GridGeometry::GridView;
39 using FVElementGeometry = typename GridGeometry::LocalView;
40 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
41 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
42 using Element = typename GridView::template Codim<0>::Entity;
43 using GridIndexType = typename GridView::IndexSet::IndexType;
44 using Stencil = std::vector<GridIndexType>;
45 using TransmissibilityVector = std::vector<GridIndexType>;
46
47 using CoordScalar = typename GridView::ctype;
48 static const int dim = GridView::dimension;
49 static const int dimWorld = GridView::dimensionworld;
50
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 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
56
57public:
58 static constexpr bool isSolDependent = false;
59
60 void update(const Problem& problem,
61 const Element& element,
62 const FVElementGeometry& fvGeometry,
63 const ElementVolumeVariables& elemVolVars,
64 const SubControlVolumeFace& scvf)
65 {
66 const auto geometry = element.geometry();
67 const auto& localBasis = fvGeometry.feLocalBasis();
68
69 // evaluate shape functions and gradients at the integration point
70 std::vector<ShapeValue> shapeVals;
71 const auto ipLocal = geometry.local(scvf.ipGlobal());
72 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal);
73 localBasis.evaluateJacobian(ipLocal, shapeJacobian_);
74 localBasis.evaluateFunction(ipLocal, shapeVals);
75
76 // set the shape values
77 shapeValues_.resize(fvGeometry.numScv(), 0.0);
78 if (!scvf.isOnFracture())
79 std::copy(shapeVals.begin(), shapeVals.end(), shapeValues_.begin());
80 else
81 {
82 const auto thisFacetIdx = scvf.facetIndexInElement();
83 for (const auto& scv: scvs(fvGeometry))
84 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
85 shapeValues_[scv.indexInElement()] = shapeVals[scv.localDofIndex()];
86 }
87
88 // set the shape value gradients
89 gradN_.resize(fvGeometry.numScv(), GlobalPosition(0.0));
90 if (!scvf.isOnFracture())
91 {
92 for (const auto& scv: scvs(fvGeometry))
93 if (!scv.isOnFracture())
94 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
95 }
96 else
97 {
98 const auto thisFacetIdx = scvf.facetIndexInElement();
99
100 // first, find all local dofs on this facet
101 std::vector<unsigned int> facetLocalDofs;
102 for (const auto& scv : scvs(fvGeometry))
103 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
104 facetLocalDofs.push_back(scv.localDofIndex());
105
106 for (const auto& scv: scvs(fvGeometry))
107 {
108 // now, create entries for all fracture scvs on this same facet ...
109 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
110 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
111
112 // ... and those non-fracture scvs that are not on this facet
113 else if (!scv.isOnFracture()
114 && std::find( facetLocalDofs.begin(),
115 facetLocalDofs.end(),
116 scv.localDofIndex() ) == facetLocalDofs.end())
117 {
118 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
119 }
120 }
121 }
122 }
123
125 const std::vector<ShapeJacobian>& shapeJacobian() const { return shapeJacobian_; }
127 const std::vector<ShapeValue>& shapeValues() const { return shapeValues_; }
129 const JacobianInverseTransposed& jacInvT() const { return jacInvT_; }
131 const GlobalPosition& gradN(unsigned int scvIdxInElement) const { return gradN_[scvIdxInElement]; }
132
133private:
134 std::vector<GlobalPosition> gradN_;
135 std::vector<ShapeJacobian> shapeJacobian_;
136 std::vector<ShapeValue> shapeValues_;
137 JacobianInverseTransposed jacInvT_;
138};
139
140} // end namespace Dumux
141
142#endif
We only store discretization-related quantities for the box method. However, we cannot reuse the cach...
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:33
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
Returns the shape value gradients corresponding to an scv.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:131
static constexpr bool isSolDependent
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:58
const std::vector< ShapeValue > & shapeValues() const
Returns the shape values for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:127
const std::vector< ShapeJacobian > & shapeJacobian() const
Returns the Jacobian of the shape functions at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:125
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:60
const JacobianInverseTransposed & jacInvT() const
Returns the shape value gradients for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:129
Defines all properties used in Dumux.
Classes related to flux variables caching.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
Definition: adapt.hh:17