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