3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
26#ifndef DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
27#define DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
28
32
33namespace Dumux {
34
41template<class TypeTag>
43{
47 using GridView = typename GridGeometry::GridView;
49 using FVElementGeometry = typename GridGeometry::LocalView;
50 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
51 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
52 using Element = typename GridView::template Codim<0>::Entity;
53 using GridIndexType = typename GridView::IndexSet::IndexType;
54 using Stencil = std::vector<GridIndexType>;
55 using TransmissibilityVector = std::vector<GridIndexType>;
56
57 using CoordScalar = typename GridView::ctype;
58 static const int dim = GridView::dimension;
59 static const int dimWorld = GridView::dimensionworld;
60
61 using FeLocalBasis = typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
62 using ShapeJacobian = typename FeLocalBasis::Traits::JacobianType;
63 using ShapeValue = typename Dune::FieldVector<Scalar, 1>;
64 using JacobianInverseTransposed = typename Element::Geometry::JacobianInverseTransposed;
65 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
66
67public:
68
69 void update(const Problem& problem,
70 const Element& element,
71 const FVElementGeometry& fvGeometry,
72 const ElementVolumeVariables& elemVolVars,
73 const SubControlVolumeFace& scvf)
74 {
75 const auto geometry = element.geometry();
76 const auto& localBasis = fvGeometry.feLocalBasis();
77
78 // evaluate shape functions and gradients at the integration point
79 std::vector<ShapeValue> shapeVals;
80 const auto ipLocal = geometry.local(scvf.ipGlobal());
81 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal);
82 localBasis.evaluateJacobian(ipLocal, shapeJacobian_);
83 localBasis.evaluateFunction(ipLocal, shapeVals);
84
85 // set the shape values
86 shapeValues_.resize(fvGeometry.numScv(), 0.0);
87 if (!scvf.isOnFracture())
88 std::copy(shapeVals.begin(), shapeVals.end(), shapeValues_.begin());
89 else
90 {
91 const auto thisFacetIdx = scvf.facetIndexInElement();
92 for (const auto& scv: scvs(fvGeometry))
93 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
94 shapeValues_[scv.indexInElement()] = shapeVals[scv.localDofIndex()];
95 }
96
97 // set the shape value gradients
98 gradN_.resize(fvGeometry.numScv(), GlobalPosition(0.0));
99 if (!scvf.isOnFracture())
100 {
101 for (const auto& scv: scvs(fvGeometry))
102 if (!scv.isOnFracture())
103 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
104 }
105 else
106 {
107 const auto thisFacetIdx = scvf.facetIndexInElement();
108
109 // first, find all local dofs on this facet
110 std::vector<unsigned int> facetLocalDofs;
111 for (const auto& scv : scvs(fvGeometry))
112 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
113 facetLocalDofs.push_back(scv.localDofIndex());
114
115 for (const auto& scv: scvs(fvGeometry))
116 {
117 // now, create entries for all fracture scvs on this same facet ...
118 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
119 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
120
121 // ... and those non-fracture scvs that are not on this facet
122 else if (!scv.isOnFracture()
123 && std::find( facetLocalDofs.begin(),
124 facetLocalDofs.end(),
125 scv.localDofIndex() ) == facetLocalDofs.end())
126 {
127 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
128 }
129 }
130 }
131 }
132
134 const std::vector<ShapeJacobian>& shapeJacobian() const { return shapeJacobian_; }
136 const std::vector<ShapeValue>& shapeValues() const { return shapeValues_; }
138 const JacobianInverseTransposed& jacInvT() const { return jacInvT_; }
140 const GlobalPosition& gradN(unsigned int scvIdxInElement) const { return gradN_[scvIdxInElement]; }
141
142private:
143 std::vector<GlobalPosition> gradN_;
144 std::vector<ShapeJacobian> shapeJacobian_;
145 std::vector<ShapeValue> shapeValues_;
146 JacobianInverseTransposed jacInvT_;
147};
148
149} // end namespace Dumux
150
151#endif
Classes related to flux variables caching.
The available discretization methods in Dumux.
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
We only store discretization-related quantities for the box method. However, we cannot reuse the cach...
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:43
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
Returns the shape value gradients corresponding to an scv.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:140
const std::vector< ShapeValue > & shapeValues() const
Returns the shape values for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:136
const std::vector< ShapeJacobian > & shapeJacobian() const
Returns the Jacobian of the shape functions at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:134
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:69
const JacobianInverseTransposed & jacInvT() const
Returns the shape value gradients for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:138
Declares all properties used in Dumux.