GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/examples/diffusion/model.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 13 15 86.7%
Functions: 1 2 50.0%
Branches: 4 4 100.0%

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