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_CAHN_HILLIARD_MODEL_HH | ||
9 | #define DUMUX_EXAMPLES_CAHN_HILLIARD_MODEL_HH | ||
10 | |||
11 | // # Cahn-Hilliard equation model definition | ||
12 | // | ||
13 | // In the file `model.hh`, we define the model equations and | ||
14 | // set all default model properties. The setup consist of four steps: | ||
15 | // 1. Create a model type tag (used to specialize properties) | ||
16 | // 2. Define the volume variables class computing and storing variables for a control volume | ||
17 | // 3. Define the local residual class implementing the discrete equation | ||
18 | // 4. Specialize important properties of the model such that Dumux knows how to assemble the system matrix | ||
19 | // | ||
20 | // We implement the classes | ||
21 | // | ||
22 | // * `CahnHilliardModel` (step 1), | ||
23 | // * `CahnHilliardModelVolumeVariables` (step 2) and | ||
24 | // * `CahnHilliardModelLocalResidual` (step 3). | ||
25 | // | ||
26 | // __Table of contents__ | ||
27 | // | ||
28 | // [TOC] | ||
29 | // | ||
30 | // We start in `model.hh` with the necessary header includes: | ||
31 | // [[details]] includes | ||
32 | #include <dune/common/fvector.hh> | ||
33 | #include <dumux/common/math.hh> | ||
34 | #include <dumux/common/properties.hh> | ||
35 | #include <dumux/common/numeqvector.hh> | ||
36 | #include <dumux/discretization/method.hh> | ||
37 | // [[/details]] | ||
38 | // | ||
39 | // ## 1. Property Tag | ||
40 | // | ||
41 | // The property tag is simply an empty struct with the name `CahnHilliardModel`. | ||
42 | // | ||
43 | // [[content]] | ||
44 | // [[codeblock]] | ||
45 | namespace Dumux::Properties::TTag { | ||
46 | struct CahnHilliardModel {}; | ||
47 | } // end namespace Dumux::Properties::TTag | ||
48 | // [[/codeblock]] | ||
49 | // [[/content]] | ||
50 | // | ||
51 | // ## 2. The volume variables | ||
52 | // | ||
53 | // The volume variables store the control volume variables, both primary and secondary. | ||
54 | // Let's have a look at the class implementation. | ||
55 | // | ||
56 | // [[content]] | ||
57 | // | ||
58 | // We write out a full volume variable class compatible with the Dumux assembly. | ||
59 | // We have to fulfill the same interface as the [`BasicVolumeVariables`](https://git.iws.uni-stuttgart.de/dumux-repositories/dumux/-/blob/master/dumux/common/volumevariables.hh). | ||
60 | // To shorten the code, we could also inherit from `BasicVolumeVariables`. However, | ||
61 | // we want to show the entire interface here as a basis for the implementation of more complex variables. | ||
62 | // Note that in `update` we could, for example, compute some secondary variable that is a | ||
63 | // possibly nonlinear function of the primary variables. That secondary variable could be exposed | ||
64 | // via an interface and then used below in the `CahnHilliardModelLocalResidual` class implementation. | ||
65 | // | ||
66 | // [[codeblock]] | ||
67 | namespace Dumux { | ||
68 | template <class Traits> | ||
69 | 1920 | class CahnHilliardModelVolumeVariables | |
70 | { | ||
71 | using Scalar = typename Traits::PrimaryVariables::value_type; | ||
72 | static_assert(Traits::PrimaryVariables::dimension == Traits::ModelTraits::numEq()); | ||
73 | public: | ||
74 | //! export the type used for the primary variables | ||
75 | using PrimaryVariables = typename Traits::PrimaryVariables; | ||
76 | //! export the indices type | ||
77 | using Indices = typename Traits::ModelTraits::Indices; | ||
78 | // [[/codeblock]] | ||
79 | // | ||
80 | // **Update variables:** The `update` function stores the local primary variables of the current solution and | ||
81 | // potentially computes secondary variables. Secondary variables can be nonlinear functions | ||
82 | // of the primary variables. | ||
83 | // | ||
84 | // [[codeblock]] | ||
85 | //! Update all quantities for a given control volume | ||
86 | template<class ElementSolution, class Problem, class Element, class SubControlVolume> | ||
87 | ✗ | void update(const ElementSolution& elemSol, | |
88 | const Problem& problem, | ||
89 | const Element& element, | ||
90 | const SubControlVolume& scv) | ||
91 | { | ||
92 | 2460160 | priVars_ = elemSol[scv.indexInElement()]; | |
93 | ✗ | } | |
94 | // [[/codeblock]] | ||
95 | // | ||
96 | // **Access functions:** Named and generic functions to access different primary variables. | ||
97 | // The volumevariables are also expected to return an extrusion factor. | ||
98 | // Here we don't extrude our domain and therefore return $1.0$. | ||
99 | // [[codeblock]] | ||
100 | Scalar concentration() const | ||
101 | 51584000 | { return priVars_[Indices::concentrationIdx]; } | |
102 | |||
103 | Scalar chemicalPotential() const | ||
104 | 36892800 | { return priVars_[Indices::chemicalPotentialIdx]; } | |
105 | |||
106 | Scalar priVar(const int pvIdx) const | ||
107 | { return priVars_[pvIdx]; } | ||
108 | |||
109 | const PrimaryVariables& priVars() const | ||
110 | ✗ | { return priVars_; } | |
111 | |||
112 | ✗ | Scalar extrusionFactor() const | |
113 | ✗ | { return 1.0; } | |
114 | |||
115 | private: | ||
116 | PrimaryVariables priVars_; | ||
117 | }; | ||
118 | } // end namespace Dumux | ||
119 | // [[/codeblock]] | ||
120 | // [[/content]] | ||
121 | // | ||
122 | // ## 3. The local residual | ||
123 | // | ||
124 | // The local residual assembles the contribution to the residual for | ||
125 | // all degrees of freedom associated with an element. | ||
126 | // See [examples/diffusion/doc/model.hh](https://git.iws.uni-stuttgart.de/dumux-repositories/dumux/-/blob/master/examples/diffusion/doc/main.md) | ||
127 | // for a more detailed explanation for control-volume finite element method local residuals. | ||
128 | // Let's have a look at the class implementation. | ||
129 | // | ||
130 | // [[content]] | ||
131 | // | ||
132 | // The class `CahnHilliardModelLocalResidual` inherits from a base class set in | ||
133 | // the model properties, depending on the discretization scheme. | ||
134 | // See [examples/diffusion/doc/model.hh](https://git.iws.uni-stuttgart.de/dumux-repositories/dumux/-/blob/master/examples/diffusion/doc/main.md) | ||
135 | // for details on the `BaseLocalResidual`. | ||
136 | namespace Dumux { | ||
137 | template<class TypeTag> | ||
138 | class CahnHilliardModelLocalResidual | ||
139 | : public GetPropType<TypeTag, Properties::BaseLocalResidual> | ||
140 | { | ||
141 | // [[exclude]] | ||
142 | // the base local residual is selected depending on the chosen discretization scheme | ||
143 | using ParentType = GetPropType<TypeTag, Properties::BaseLocalResidual>; | ||
144 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
145 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
146 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
147 | |||
148 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
149 | using VolumeVariables = typename GridVariables::GridVolumeVariables::VolumeVariables; | ||
150 | using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView; | ||
151 | using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView; | ||
152 | |||
153 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
154 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
155 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
156 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
157 | using GridView = typename GridGeometry::GridView; | ||
158 | using Element = typename GridView::template Codim<0>::Entity; | ||
159 | |||
160 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
161 | using Indices = typename ModelTraits::Indices; | ||
162 | static constexpr int dimWorld = GridView::dimensionworld; | ||
163 | public: | ||
164 | 207280 | using ParentType::ParentType; | |
165 | // [[/exclude]] | ||
166 | // | ||
167 | // **Storage term:** The function `computeStorage` receives the volume variables | ||
168 | // at the previous or current time step and computes the value of the storage terms. | ||
169 | // In this case the mass balance equation is a conservation equation of the concentration and | ||
170 | // the equation for the chemical potential does not have a storage term. | ||
171 | // | ||
172 | // [[codeblock]] | ||
173 | ✗ | NumEqVector computeStorage(const Problem& problem, | |
174 | const SubControlVolume& scv, | ||
175 | const VolumeVariables& volVars) const | ||
176 | { | ||
177 | 7345600 | NumEqVector storage; | |
178 | 14691200 | storage[Indices::massBalanceEqIdx] = volVars.concentration(); | |
179 | 14691200 | storage[Indices::chemicalPotentialEqIdx] = 0.0; | |
180 | ✗ | return storage; | |
181 | } | ||
182 | // [[/codeblock]] | ||
183 | |||
184 | // **Flux term:** The function `computeFlux` computes the integrated | ||
185 | // fluxes over a sub control volume face. | ||
186 | // | ||
187 | // ```math | ||
188 | // \begin{aligned} | ||
189 | // F_{K,\sigma,0} &= -M \vert \sigma \vert \sum_{B \in \mathcal{B}_K} \mu_{h,B} \nabla N_B \cdot\boldsymbol{n} \cr | ||
190 | // F_{K,\sigma,1} &= -\gamma \vert \sigma \vert \sum_{B \in \mathcal{B}_K} c_{h,B} \nabla N_B \cdot\boldsymbol{n} | ||
191 | // \end{aligned} | ||
192 | // ```` | ||
193 | // | ||
194 | // [[codeblock]] | ||
195 | 3672800 | NumEqVector computeFlux(const Problem& problem, | |
196 | const Element& element, | ||
197 | const FVElementGeometry& fvGeometry, | ||
198 | const ElementVolumeVariables& elemVolVars, | ||
199 | const SubControlVolumeFace& scvf, | ||
200 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
201 | { | ||
202 | static_assert(DiscretizationMethods::isCVFE<typename GridGeometry::DiscretizationMethod>, | ||
203 | "This local residual is hard-coded to control-volume finite element schemes"); | ||
204 | |||
205 | 3672800 | const auto& fluxVarCache = elemFluxVarsCache[scvf]; | |
206 | 3672800 | Dune::FieldVector<Scalar, dimWorld> gradConcentration(0.0); | |
207 | 3672800 | Dune::FieldVector<Scalar, dimWorld> gradChemicalPotential(0.0); | |
208 |
4/4✓ Branch 0 taken 14691200 times.
✓ Branch 1 taken 3672800 times.
✓ Branch 2 taken 14691200 times.
✓ Branch 3 taken 3672800 times.
|
36728000 | for (const auto& scv : scvs(fvGeometry)) |
209 | { | ||
210 | 14691200 | const auto& volVars = elemVolVars[scv]; | |
211 | // v.axpy(a, w) means v += a*w | ||
212 | 44073600 | gradConcentration.axpy( | |
213 | volVars.concentration(), | ||
214 | fluxVarCache.gradN(scv.indexInElement()) | ||
215 | ); | ||
216 | 58764800 | gradChemicalPotential.axpy( | |
217 | volVars.chemicalPotential(), | ||
218 | fluxVarCache.gradN(scv.indexInElement()) | ||
219 | ); | ||
220 | } | ||
221 | |||
222 | 3672800 | NumEqVector flux; | |
223 | // Compute the flux with `vtmv` (vector transposed times matrix times vector). | ||
224 | // The mobility and surface tension coefficients comes from the `problem` (see Part 2 of the example). | ||
225 | 7345600 | flux[Indices::massBalanceEqIdx] = -1.0*vtmv( | |
226 | scvf.unitOuterNormal(), problem.mobility(), gradChemicalPotential | ||
227 | 3672800 | )*scvf.area(); | |
228 | 7345600 | flux[Indices::chemicalPotentialEqIdx] = -1.0*vtmv( | |
229 | scvf.unitOuterNormal(), problem.surfaceTension(), gradConcentration | ||
230 | 3672800 | )*scvf.area(); | |
231 | 3672800 | return flux; | |
232 | } | ||
233 | // [[/codeblock]] | ||
234 | |||
235 | // **Source term:** The function `computeSource` computes the source terms for a sub control volume. | ||
236 | // We implement a model-specific source term for the chemical potential equation before | ||
237 | // deferring further implementation to the problem where we add the derivative of the free | ||
238 | // energy. | ||
239 | // | ||
240 | // [[codeblock]] | ||
241 | 3672800 | NumEqVector computeSource(const Problem& problem, | |
242 | const Element& element, | ||
243 | const FVElementGeometry& fvGeometry, | ||
244 | const ElementVolumeVariables& elemVolVars, | ||
245 | const SubControlVolume &scv) const | ||
246 | { | ||
247 | 3672800 | NumEqVector source(0.0); | |
248 | 3672800 | source[Indices::massBalanceEqIdx] = 0.0; | |
249 | 11018400 | source[Indices::chemicalPotentialEqIdx] = elemVolVars[scv].chemicalPotential(); | |
250 | // add contributions from problem (e.g. double well potential) | ||
251 | 3672800 | source += problem.source(element, fvGeometry, elemVolVars, scv); | |
252 | 3672800 | return source; | |
253 | } | ||
254 | }; | ||
255 | } // end namespace Dumux | ||
256 | // [[/codeblock]] | ||
257 | // [[/content]] | ||
258 | // | ||
259 | // ## 4. The model properties/traits | ||
260 | // | ||
261 | // In the `Dumux::Properties` namespace, we specialize properties for | ||
262 | // the created type tag `CahnHilliardModel`. | ||
263 | // | ||
264 | // [[content]] | ||
265 | namespace Dumux::Properties { | ||
266 | |||
267 | // The type of the local residual is the class defined above. | ||
268 | template<class TypeTag> | ||
269 | struct LocalResidual<TypeTag, TTag::CahnHilliardModel> | ||
270 | { using type = CahnHilliardModelLocalResidual<TypeTag>; }; | ||
271 | |||
272 | // We compute with double precision floating point numbers. | ||
273 | template<class TypeTag> | ||
274 | struct Scalar<TypeTag, TTag::CahnHilliardModel> | ||
275 | { using type = double; }; | ||
276 | |||
277 | // The model traits specify some information about our equation system. | ||
278 | // Here we have two equations. The indices allow to access primary variables | ||
279 | // and equations with named indices. | ||
280 | template<class TypeTag> | ||
281 | struct ModelTraits<TypeTag, TTag::CahnHilliardModel> | ||
282 | { | ||
283 | struct type | ||
284 | { | ||
285 | struct Indices | ||
286 | { | ||
287 | static constexpr int concentrationIdx = 0; | ||
288 | static constexpr int chemicalPotentialIdx = 1; | ||
289 | |||
290 | static constexpr int massBalanceEqIdx = 0; | ||
291 | static constexpr int chemicalPotentialEqIdx = 1; | ||
292 | }; | ||
293 | |||
294 | static constexpr int numEq() { return 2; } | ||
295 | }; | ||
296 | }; | ||
297 | |||
298 | // The primary variable vector has entries of type `Scalar` and is | ||
299 | // as large as the number of equations (here 2) but we keep it general | ||
300 | // here by obtaining the number of equations from the `ModelTraits`. | ||
301 | template<class TypeTag> | ||
302 | struct PrimaryVariables<TypeTag, TTag::CahnHilliardModel> | ||
303 | { | ||
304 | using type = Dune::FieldVector< | ||
305 | GetPropType<TypeTag, Properties::Scalar>, | ||
306 | GetPropType<TypeTag, Properties::ModelTraits>::numEq() | ||
307 | >; | ||
308 | }; | ||
309 | |||
310 | // Finally, the type of the volume variables is the class defined above. | ||
311 | template<class TypeTag> | ||
312 | struct VolumeVariables<TypeTag, TTag::CahnHilliardModel> | ||
313 | { | ||
314 | struct Traits | ||
315 | { | ||
316 | using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>; | ||
317 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
318 | }; | ||
319 | |||
320 | using type = CahnHilliardModelVolumeVariables<Traits>; | ||
321 | }; | ||
322 | |||
323 | } // end namespace Dumux::Properties | ||
324 | // [[/content]] | ||
325 | // [[exclude]] | ||
326 | #endif | ||
327 | // [[/exclude]] | ||
328 |