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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | |||
8 | #ifndef DUMUX_EXAMPLES_DIFFUSION_MODEL_HH | ||
9 | #define DUMUX_EXAMPLES_DIFFUSION_MODEL_HH | ||
10 | |||
11 | // In the file `model.hh`, we define the model equations and | ||
12 | // set all default model properties. The setup consist of three steps: | ||
13 | // 1. Create a model type tag (used to specialize properties) | ||
14 | // 2. Define the local residual class implementing the discrete equation | ||
15 | // 3. Specialize important properties of the model such that Dumux knows how to assemble the system matrix | ||
16 | // | ||
17 | // __Table of contents__ | ||
18 | // | ||
19 | // [TOC] | ||
20 | // | ||
21 | // We start in `model.hh` with the necessary header includes: | ||
22 | // [[details]] includes | ||
23 | #include <dune/common/fvector.hh> | ||
24 | #include <dumux/common/math.hh> | ||
25 | #include <dumux/common/properties.hh> | ||
26 | #include <dumux/common/numeqvector.hh> | ||
27 | #include <dumux/common/volumevariables.hh> | ||
28 | #include <dumux/discretization/method.hh> | ||
29 | #include <dumux/discretization/defaultlocaloperator.hh> | ||
30 | // [[/details]] | ||
31 | // | ||
32 | // ## 1. Property Tag | ||
33 | // | ||
34 | // The property tag is simply an empty struct with the name `DiffusionModel` | ||
35 | // | ||
36 | // [[content]] | ||
37 | // [[codeblock]] | ||
38 | namespace Dumux::Properties::TTag { | ||
39 | //! The diffusion model tag that we can specialize properties for | ||
40 | struct DiffusionModel {}; | ||
41 | } // end namespace Dumux::Properties::TTag | ||
42 | // [[/codeblock]] | ||
43 | // [[/content]] | ||
44 | // | ||
45 | // ## 2. The local (element-wise) residual | ||
46 | // | ||
47 | // The local residual assembles the contribution to the residual for | ||
48 | // all degrees of freedom associated with an element. Here, we use the | ||
49 | // Box method which is based on $P_1$ basis functions (piece-wise linears) | ||
50 | // and the degrees of freedom are on the nodes. Each node is associated with | ||
51 | // exactly one sub control volume (`scv`) per element and several ($2$ in $\mathbb{R}^2$) | ||
52 | // sub control volume faces (`scvf`). In the local residual, we can implement the | ||
53 | // contribution for one `scv` (storage and source terms) or one `scvf` (flux terms). | ||
54 | // | ||
55 | // Let's have a look at the class implementation. | ||
56 | // | ||
57 | // [[content]] | ||
58 | // | ||
59 | // The class `DiffusionModelLocalResidual` inherits from something called `DiscretizationDefaultLocalOperator`. | ||
60 | // This base class differs depending on the chosen discretization scheme. For the box method | ||
61 | // (which is a control-volume finite element scheme) used in this example, | ||
62 | // `DiscretizationDefaultLocalOperator` is specialized to `CVFELocalResidual` | ||
63 | // in [dumux/discretization/box.hh](https://git.iws.uni-stuttgart.de/dumux-repositories/dumux/-/blob/master/dumux/discretization/box.hh). | ||
64 | // Since the following implementation only works with control-volume finite element schemes due to | ||
65 | // the chosen flux implementation, we could have also chosen to | ||
66 | // directly inherit from `public CVFELocalResidual<TypeTag>`. | ||
67 | namespace Dumux { | ||
68 | template<class TypeTag> | ||
69 | class DiffusionModelLocalResidual | ||
70 | : public DiscretizationDefaultLocalOperator<TypeTag> | ||
71 | { | ||
72 | // [[exclude]] | ||
73 | // the base local operator is selected depending on the chosen discretization scheme | ||
74 | using ParentType = DiscretizationDefaultLocalOperator<TypeTag>; | ||
75 | |||
76 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
77 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
78 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
79 | |||
80 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
81 | using VolumeVariables = typename GridVariables::GridVolumeVariables::VolumeVariables; | ||
82 | using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView; | ||
83 | using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView; | ||
84 | |||
85 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
86 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
87 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
88 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
89 | using GridView = typename GridGeometry::GridView; | ||
90 | using Element = typename GridView::template Codim<0>::Entity; | ||
91 | |||
92 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
93 | using Indices = typename ModelTraits::Indices; | ||
94 | static constexpr int dimWorld = GridView::dimensionworld; | ||
95 | |||
96 | public: | ||
97 | 3672 | using ParentType::ParentType; | |
98 | // [[/exclude]] | ||
99 | // | ||
100 | // **Storage term:** Evaluate the rate of change of all conserved quantities | ||
101 | // [[codeblock]] | ||
102 | 15840 | NumEqVector computeStorage(const Problem& problem, | |
103 | const SubControlVolume& scv, | ||
104 | const VolumeVariables& volVars) const | ||
105 | { | ||
106 | 15840 | NumEqVector storage; | |
107 | 15840 | storage[Indices::massBalanceEqIdx] = volVars.priVar(Indices::concentrationIdx); | |
108 | return storage; | ||
109 | } | ||
110 | // [[/codeblock]] | ||
111 | |||
112 | // **Flux term:** Evaluate the fluxes over a face of a sub control volume. | ||
113 | // Here we evaluate the (integrated) flux | ||
114 | // | ||
115 | // ```math | ||
116 | // F_{K,\sigma} = -D \sum_{B \in \mathcal{B}_K} c_{h,B} \nabla N_B \cdot\boldsymbol{n} \vert \sigma \vert | ||
117 | // ```` | ||
118 | // | ||
119 | // [[codeblock]] | ||
120 | 15840 | NumEqVector computeFlux(const Problem& problem, | |
121 | const Element& element, | ||
122 | const FVElementGeometry& fvGeometry, | ||
123 | const ElementVolumeVariables& elemVolVars, | ||
124 | const SubControlVolumeFace& scvf, | ||
125 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
126 | { | ||
127 | static_assert(DiscretizationMethods::isCVFE<typename GridGeometry::DiscretizationMethod>, | ||
128 | "This local residual is hard-coded to control-volume finite element schemes"); | ||
129 | |||
130 | // Compute ∇c at the integration point of this sub control volume face. | ||
131 | 15840 | const auto& fluxVarCache = elemFluxVarsCache[scvf]; | |
132 | 15840 | Dune::FieldVector<Scalar, dimWorld> gradConcentration(0.0); | |
133 |
2/2✓ Branch 0 taken 63360 times.
✓ Branch 1 taken 15840 times.
|
79200 | for (const auto& scv : scvs(fvGeometry)) |
134 | { | ||
135 | 63360 | const auto& volVars = elemVolVars[scv]; | |
136 | // v.axpy(a, w) means v += a*w | ||
137 | 63360 | gradConcentration.axpy( | |
138 | 63360 | volVars.priVar(Indices::concentrationIdx), | |
139 | 63360 | fluxVarCache.gradN(scv.indexInElement()) | |
140 | ); | ||
141 | } | ||
142 | |||
143 | 15840 | NumEqVector flux; | |
144 | |||
145 | // Compute the flux with `vtmv` (vector transposed times matrix times vector) or -n^T D ∇c A. | ||
146 | // The diffusion coefficient comes from the `problem` (see Part 2 of the example). | ||
147 | 15840 | flux[Indices::massBalanceEqIdx] = -1.0*vtmv( | |
148 | 15840 | scvf.unitOuterNormal(), problem.diffusionCoefficient(), gradConcentration | |
149 | 15840 | )*scvf.area(); | |
150 | |||
151 | 15840 | return flux; | |
152 | } | ||
153 | }; | ||
154 | } // end namespace Dumux | ||
155 | // [[/codeblock]] | ||
156 | // [[/content]] | ||
157 | // | ||
158 | // ## 3. The model properties | ||
159 | // | ||
160 | // By specializing properties for our type tag `DiffusionModel`, | ||
161 | // every other class that knows about the type tag (this will be | ||
162 | // for example the assembler or the problem), can extract the | ||
163 | // type information that we specify here. | ||
164 | // | ||
165 | // Note that these types can be overwritten for specific problem | ||
166 | // definitions if this is needed (we will show this on the next page). | ||
167 | // | ||
168 | // [[content]] | ||
169 | namespace Dumux::Properties { | ||
170 | |||
171 | // The type of the local residual is the class defined above | ||
172 | template<class TypeTag> | ||
173 | struct LocalResidual<TypeTag, TTag::DiffusionModel> | ||
174 | { using type = DiffusionModelLocalResidual<TypeTag>; }; | ||
175 | |||
176 | // The default scalar type is double | ||
177 | // we compute with double precision floating point numbers | ||
178 | template<class TypeTag> | ||
179 | struct Scalar<TypeTag, TTag::DiffusionModel> | ||
180 | { using type = double; }; | ||
181 | |||
182 | // The model traits specify some information about our equation system. | ||
183 | // Here we have just one equation. We still specify indices so in the | ||
184 | // places where we access primary variables, we can do so with a named variable. | ||
185 | template<class TypeTag> | ||
186 | struct ModelTraits<TypeTag, TTag::DiffusionModel> | ||
187 | { | ||
188 | struct type | ||
189 | { | ||
190 | struct Indices | ||
191 | { | ||
192 | static constexpr int concentrationIdx = 0; | ||
193 | static constexpr int massBalanceEqIdx = 0; | ||
194 | }; | ||
195 | |||
196 | static constexpr int numEq() { return 1; } | ||
197 | }; | ||
198 | }; | ||
199 | |||
200 | // The primary variable vector has entries of type `Scalar` and is | ||
201 | // as large as the number of equations (here 1) but we keep it general. | ||
202 | template<class TypeTag> | ||
203 | struct PrimaryVariables<TypeTag, TTag::DiffusionModel> | ||
204 | { | ||
205 | using type = Dune::FieldVector< | ||
206 | GetPropType<TypeTag, Properties::Scalar>, | ||
207 | GetPropType<TypeTag, Properties::ModelTraits>::numEq() | ||
208 | >; | ||
209 | }; | ||
210 | |||
211 | // The `BasicVolumeVariables` are the simplest class of volume variables. | ||
212 | // They only store one instance of `PrimaryVariables` for the | ||
213 | // degree of freedom (here: vertex dof) that they are attached to. | ||
214 | template<class TypeTag> | ||
215 | struct VolumeVariables<TypeTag, TTag::DiffusionModel> | ||
216 | { | ||
217 | struct Traits | ||
218 | { | ||
219 | using PrimaryVariables | ||
220 | = GetPropType<TypeTag, Properties::PrimaryVariables>; | ||
221 | }; | ||
222 | using type = BasicVolumeVariables<Traits>; | ||
223 | }; | ||
224 | |||
225 | } // end namespace Dumux::Properties | ||
226 | // [[/content]] | ||
227 | // [[exclude]] | ||
228 | #endif | ||
229 | // [[/exclude]] | ||
230 |