Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup BoxDFMModel | ||
10 | * \brief Cache class for the flux variables to be used | ||
11 | * in conjunction with the box discrete fracture scheme. | ||
12 | */ | ||
13 | |||
14 | #ifndef DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH | ||
15 | #define DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH | ||
16 | |||
17 | #include <dune/common/fvector.hh> | ||
18 | |||
19 | #include <dumux/common/properties.hh> | ||
20 | #include <dumux/discretization/method.hh> | ||
21 | #include <dumux/flux/fluxvariablescaching.hh> | ||
22 | |||
23 | namespace Dumux { | ||
24 | |||
25 | /*! | ||
26 | * \ingroup BoxDFMModel | ||
27 | * \brief We only store discretization-related quantities for the box method. | ||
28 | * However, we cannot reuse the cache of the standard box method as we have | ||
29 | * to take into account the scvs that lie on fracture facets. | ||
30 | */ | ||
31 | template<class TypeTag> | ||
32 | 43761720 | class BoxDfmFluxVariablesCache | |
33 | { | ||
34 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
35 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
36 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
37 | using GridView = typename GridGeometry::GridView; | ||
38 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
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 | |||
57 | public: | ||
58 | static constexpr bool isSolDependent = false; | ||
59 | |||
60 | 8752344 | void update(const Problem& problem, | |
61 | const Element& element, | ||
62 | const FVElementGeometry& fvGeometry, | ||
63 | const ElementVolumeVariables& elemVolVars, | ||
64 | const SubControlVolumeFace& scvf) | ||
65 | { | ||
66 | 8752344 | const auto geometry = element.geometry(); | |
67 |
1/2✓ Branch 1 taken 4376172 times.
✗ Branch 2 not taken.
|
8752344 | const auto& localBasis = fvGeometry.feLocalBasis(); |
68 | |||
69 | // evaluate shape functions and gradients at the integration point | ||
70 |
2/4✓ Branch 1 taken 8752344 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 4376172 times.
✗ Branch 4 not taken.
|
21880860 | std::vector<ShapeValue> shapeVals; |
71 |
2/4✓ Branch 1 taken 4376172 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 4376172 times.
✗ Branch 5 not taken.
|
17504688 | const auto ipLocal = geometry.local(scvf.ipGlobal()); |
72 |
1/2✓ Branch 1 taken 4376172 times.
✗ Branch 2 not taken.
|
8752344 | jacInvT_ = geometry.jacobianInverseTransposed(ipLocal); |
73 |
1/2✓ Branch 1 taken 8752344 times.
✗ Branch 2 not taken.
|
8752344 | localBasis.evaluateJacobian(ipLocal, shapeJacobian_); |
74 |
1/2✓ Branch 1 taken 8752344 times.
✗ Branch 2 not taken.
|
8752344 | localBasis.evaluateFunction(ipLocal, shapeVals); |
75 | |||
76 | // set the shape values | ||
77 |
3/6✓ Branch 1 taken 8752344 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 8752344 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 8752344 times.
✗ Branch 8 not taken.
|
26257032 | shapeValues_.resize(fvGeometry.numScv(), 0.0); |
78 |
2/2✓ Branch 0 taken 8459356 times.
✓ Branch 1 taken 292988 times.
|
8752344 | if (!scvf.isOnFracture()) |
79 | 33837424 | std::copy(shapeVals.begin(), shapeVals.end(), shapeValues_.begin()); | |
80 | else | ||
81 | { | ||
82 | 292988 | const auto thisFacetIdx = scvf.facetIndexInElement(); | |
83 |
4/4✓ Branch 0 taken 1606024 times.
✓ Branch 1 taken 292988 times.
✓ Branch 2 taken 1606024 times.
✓ Branch 3 taken 292988 times.
|
2192000 | for (const auto& scv: scvs(fvGeometry)) |
84 |
6/6✓ Branch 0 taken 661672 times.
✓ Branch 1 taken 944352 times.
✓ Branch 2 taken 648016 times.
✓ Branch 3 taken 13656 times.
✓ Branch 4 taken 648016 times.
✓ Branch 5 taken 13656 times.
|
1606024 | if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx) |
85 | 1944048 | shapeValues_[scv.indexInElement()] = shapeVals[scv.localDofIndex()]; | |
86 | } | ||
87 | |||
88 | // set the shape value gradients | ||
89 |
2/6✓ Branch 1 taken 8752344 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 8752344 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
26257032 | gradN_.resize(fvGeometry.numScv(), GlobalPosition(0.0)); |
90 |
2/2✓ Branch 0 taken 8459356 times.
✓ Branch 1 taken 292988 times.
|
8752344 | if (!scvf.isOnFracture()) |
91 | { | ||
92 |
4/4✓ Branch 0 taken 28912380 times.
✓ Branch 1 taken 8459356 times.
✓ Branch 2 taken 28912380 times.
✓ Branch 3 taken 8459356 times.
|
45831092 | for (const auto& scv: scvs(fvGeometry)) |
93 |
2/2✓ Branch 0 taken 27100464 times.
✓ Branch 1 taken 1811916 times.
|
28912380 | if (!scv.isOnFracture()) |
94 | 108401856 | jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]); | |
95 | } | ||
96 | else | ||
97 | { | ||
98 | 292988 | const auto thisFacetIdx = scvf.facetIndexInElement(); | |
99 | |||
100 | // first, find all local dofs on this facet | ||
101 |
1/2✓ Branch 0 taken 292988 times.
✗ Branch 1 not taken.
|
878964 | std::vector<unsigned int> facetLocalDofs; |
102 |
4/4✓ Branch 0 taken 1606024 times.
✓ Branch 1 taken 292988 times.
✓ Branch 2 taken 1606024 times.
✓ Branch 3 taken 292988 times.
|
2192000 | for (const auto& scv : scvs(fvGeometry)) |
103 |
6/6✓ Branch 0 taken 661672 times.
✓ Branch 1 taken 944352 times.
✓ Branch 2 taken 648016 times.
✓ Branch 3 taken 13656 times.
✓ Branch 4 taken 648016 times.
✓ Branch 5 taken 13656 times.
|
1606024 | if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx) |
104 |
1/4✓ Branch 1 taken 648016 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
648016 | facetLocalDofs.push_back(scv.localDofIndex()); |
105 | |||
106 |
4/4✓ Branch 0 taken 1606024 times.
✓ Branch 1 taken 292988 times.
✓ Branch 2 taken 1606024 times.
✓ Branch 3 taken 292988 times.
|
2192000 | for (const auto& scv: scvs(fvGeometry)) |
107 | { | ||
108 | // now, create entries for all fracture scvs on this same facet ... | ||
109 |
6/6✓ Branch 0 taken 661672 times.
✓ Branch 1 taken 944352 times.
✓ Branch 2 taken 13656 times.
✓ Branch 3 taken 648016 times.
✓ Branch 4 taken 13656 times.
✓ Branch 5 taken 648016 times.
|
1606024 | if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx) |
110 | 2592064 | jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]); | |
111 | |||
112 | // ... and those non-fracture scvs that are not on this facet | ||
113 | 958008 | else if (!scv.isOnFracture() | |
114 |
2/2✓ Branch 0 taken 944352 times.
✓ Branch 1 taken 13656 times.
|
1606024 | && std::find( facetLocalDofs.begin(), |
115 | facetLocalDofs.end(), | ||
116 |
4/4✓ Branch 2 taken 648016 times.
✓ Branch 3 taken 296336 times.
✓ Branch 4 taken 648016 times.
✓ Branch 5 taken 296336 times.
|
1888704 | scv.localDofIndex() ) == facetLocalDofs.end()) |
117 | { | ||
118 | 1185344 | jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]); | |
119 | } | ||
120 | } | ||
121 | } | ||
122 | 8752344 | } | |
123 | |||
124 | //! Returns the Jacobian of the shape functions at the integration point. | ||
125 | const std::vector<ShapeJacobian>& shapeJacobian() const { return shapeJacobian_; } | ||
126 | //! Returns the shape values for all scvs at the integration point. | ||
127 | 132628712 | const std::vector<ShapeValue>& shapeValues() const { return shapeValues_; } | |
128 | //! Returns the shape value gradients for all scvs at the integration point. | ||
129 | const JacobianInverseTransposed& jacInvT() const { return jacInvT_; } | ||
130 | //! Returns the shape value gradients corresponding to an scv. | ||
131 | 972730976 | const GlobalPosition& gradN(unsigned int scvIdxInElement) const { return gradN_[scvIdxInElement]; } | |
132 | |||
133 | private: | ||
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 | ||
143 |