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_SOLIDSYSTEMS_MICP_SOLID_PHASE_HH | ||
9 | #define DUMUX_SOLIDSYSTEMS_MICP_SOLID_PHASE_HH | ||
10 | |||
11 | // ## The solid system (`biominsolids.hh`) | ||
12 | // | ||
13 | // This file contains the __solidystem class__ which defines all functions needed to describe the solid and their properties, | ||
14 | // which are needed to model biomineralization. | ||
15 | // | ||
16 | // [[content]] | ||
17 | // | ||
18 | // ### Include files | ||
19 | // [[codeblock]] | ||
20 | #include <string> | ||
21 | #include <dune/common/exceptions.hh> | ||
22 | |||
23 | // we include all necessary solid components | ||
24 | #include <examples/biomineralization/material/components/biofilm.hh> | ||
25 | #include <dumux/material/components/calcite.hh> | ||
26 | #include <dumux/material/components/granite.hh> | ||
27 | // [[/codeblock]] | ||
28 | |||
29 | // ### The solidsystem class | ||
30 | // In the BioMinSolidPhase solid system, we define all functions needed to describe the solids and their properties accounted for in our simulation. | ||
31 | // We enter the namespace Dumux. All Dumux functions and classes are in a namespace Dumux, to make sure they don`t clash with symbols from other libraries you may want to use in conjunction with Dumux. | ||
32 | // [[codeblock]] | ||
33 | namespace Dumux::SolidSystems { | ||
34 | |||
35 | template <class Scalar> | ||
36 | class BioMinSolidPhase | ||
37 | { | ||
38 | public: | ||
39 | // [[/codeblock]] | ||
40 | |||
41 | // #### Component definitions | ||
42 | // With the following function we define what solid components will be used by the solid system and define the indices used to distinguish solid components in the course of the simulation. | ||
43 | // [[codeblock]] | ||
44 | // We use convenient declarations that we derive from the property system. | ||
45 | using Biofilm = Components::Biofilm<Scalar>; | ||
46 | using Calcite = Components::Calcite<Scalar>; | ||
47 | using Granite = Components::Granite<Scalar>; | ||
48 | |||
49 | /**************************************** | ||
50 | * Solid phase related static parameters | ||
51 | ****************************************/ | ||
52 | static constexpr int numComponents = 3; | ||
53 | static constexpr int numInertComponents = 1; | ||
54 | static constexpr int BiofilmIdx = 0; | ||
55 | static constexpr int CalciteIdx = 1; | ||
56 | static constexpr int GraniteIdx = 2; | ||
57 | |||
58 | // The component names | ||
59 | 2 | const static std::string componentName(int compIdx) | |
60 | { | ||
61 |
2/4✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2 | switch (compIdx) |
62 | { | ||
63 | 1 | case BiofilmIdx: return Biofilm::name(); | |
64 | 1 | case CalciteIdx: return Calcite::name(); | |
65 | ✗ | case GraniteIdx: return Granite::name(); | |
66 | ✗ | default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx); | |
67 | } | ||
68 | } | ||
69 | |||
70 | // The solid system's name | ||
71 | static std::string name() | ||
72 | { return "BiomineralizationMinSolidPhase"; } | ||
73 | |||
74 | // We assume incompressible solids | ||
75 | static constexpr bool isCompressible(int compIdx) | ||
76 | { return false; } | ||
77 | |||
78 | // we have one inert component, the others may change | ||
79 | static constexpr bool isInert() | ||
80 | { return (numComponents == numInertComponents); } | ||
81 | // [[/codeblock]] | ||
82 | |||
83 | // #### Component properties | ||
84 | // This simply forwards the component properties (name, molar mass, density). | ||
85 | // [[codeblock]] | ||
86 | // The component molar masses | ||
87 | 8119552 | static Scalar molarMass(int compIdx) | |
88 | { | ||
89 |
2/4✓ Branch 0 taken 4059776 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 4059776 times.
|
8119552 | switch (compIdx) |
90 | { | ||
91 | 4059776 | case BiofilmIdx: return Biofilm::molarMass(); | |
92 | case CalciteIdx: return Calcite::molarMass(); | ||
93 | ✗ | case GraniteIdx: return Granite::molarMass(); | |
94 | ✗ | default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx); | |
95 | } | ||
96 | } | ||
97 | |||
98 | // The component densities | ||
99 | template <class SolidState> | ||
100 | 8119552 | static Scalar density(const SolidState& solidState, const int compIdx) | |
101 | { | ||
102 |
2/4✓ Branch 0 taken 4059776 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 4059776 times.
|
8119552 | switch (compIdx) |
103 | { | ||
104 | 4059776 | case BiofilmIdx: return Biofilm::solidDensity(solidState.temperature()); | |
105 | case CalciteIdx: return Calcite::solidDensity(solidState.temperature()); | ||
106 | ✗ | case GraniteIdx: return Granite::solidDensity(solidState.temperature()); | |
107 | ✗ | default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx); | |
108 | } | ||
109 | } | ||
110 | |||
111 | // The component molar densities | ||
112 | template <class SolidState> | ||
113 | 16239104 | static Scalar molarDensity(const SolidState& solidState, const int compIdx) | |
114 | { | ||
115 |
2/4✓ Branch 0 taken 8119552 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 8119552 times.
|
16239104 | switch (compIdx) |
116 | { | ||
117 | 8119552 | case BiofilmIdx: return Biofilm::solidDensity(solidState.temperature())/Biofilm::molarMass(); | |
118 | case CalciteIdx: return Calcite::solidDensity(solidState.temperature())/Calcite::molarMass(); | ||
119 | ✗ | case GraniteIdx: return Granite::solidDensity(solidState.temperature())/Calcite::molarMass(); | |
120 | ✗ | default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx); | |
121 | } | ||
122 | } | ||
123 | // [[/codeblock]] | ||
124 | |||
125 | // #### Averaged solid properties | ||
126 | // The overall solid properties are calculated by a volume-weigthed average of the pure component's properties. | ||
127 | // [[codeblock]] | ||
128 | // The average density | ||
129 | template <class SolidState> | ||
130 | static Scalar density(const SolidState& solidState) | ||
131 | { | ||
132 | const Scalar rho1 = Biofilm::solidDensity(solidState.temperature()); | ||
133 | const Scalar rho2 = Calcite::solidDensity(solidState.temperature()); | ||
134 | const Scalar rho3 = Granite::solidDensity(solidState.temperature()); | ||
135 | const Scalar volFrac1 = solidState.volumeFraction(BiofilmIdx); | ||
136 | const Scalar volFrac2 = solidState.volumeFraction(CalciteIdx); | ||
137 | const Scalar volFrac3 = solidState.volumeFraction(GraniteIdx); | ||
138 | |||
139 | return (rho1*volFrac1+ | ||
140 | rho2*volFrac2+ | ||
141 | rho3*volFrac3) | ||
142 | /(volFrac1+volFrac2+volFrac3); | ||
143 | } | ||
144 | // [[/codeblock]] | ||
145 | // [[exclude]] | ||
146 | // The average thermal conductivity | ||
147 | template <class SolidState> | ||
148 | static Scalar thermalConductivity(const SolidState &solidState) | ||
149 | { | ||
150 | const Scalar lambda1 = Biofilm::solidThermalConductivity(solidState.temperature()); | ||
151 | const Scalar lambda2 = Calcite::solidThermalConductivity(solidState.temperature()); | ||
152 | const Scalar lambda3 = Granite::solidThermalConductivity(solidState.temperature()); | ||
153 | const Scalar volFrac1 = solidState.volumeFraction(BiofilmIdx); | ||
154 | const Scalar volFrac2 = solidState.volumeFraction(CalciteIdx); | ||
155 | const Scalar volFrac3 = solidState.volumeFraction(GraniteIdx); | ||
156 | |||
157 | return (lambda1*volFrac1+ | ||
158 | lambda2*volFrac2+ | ||
159 | lambda3*volFrac3) | ||
160 | /(volFrac1+volFrac2+volFrac3); | ||
161 | } | ||
162 | |||
163 | // The average heat capacity | ||
164 | template <class SolidState> | ||
165 | static Scalar heatCapacity(const SolidState &solidState) | ||
166 | { | ||
167 | const Scalar c1 = Biofilm::solidHeatCapacity(solidState.temperature()); | ||
168 | const Scalar c2 = Calcite::solidHeatCapacity(solidState.temperature()); | ||
169 | const Scalar c3 = Granite::solidHeatCapacity(solidState.temperature()); | ||
170 | const Scalar volFrac1 = solidState.volumeFraction(BiofilmIdx); | ||
171 | const Scalar volFrac2 = solidState.volumeFraction(CalciteIdx); | ||
172 | const Scalar volFrac3 = solidState.volumeFraction(GraniteIdx); | ||
173 | |||
174 | return (c1*volFrac1+ | ||
175 | c2*volFrac2+ | ||
176 | c3*volFrac3) | ||
177 | /(volFrac1+volFrac2+volFrac3); | ||
178 | } | ||
179 | |||
180 | }; | ||
181 | |||
182 | } // end namespace Dumux::SolidSystems | ||
183 | // [[/exclude]] | ||
184 | // [[/content]] | ||
185 | #endif | ||
186 |