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 | * \file | ||
9 | * \ingroup NavierStokesNCTests | ||
10 | * \brief The properties of the channel flow test for the multi-component staggered grid (Navier-)Stokes model. | ||
11 | */ | ||
12 | #ifndef DUMUX_CHANNEL_MAXWELL_STEPHAN_TEST_PROPERTIES_HH | ||
13 | #define DUMUX_CHANNEL_MAXWELL_STEPHAN_TEST_PROPERTIES_HH | ||
14 | |||
15 | #include <dune/grid/yaspgrid.hh> | ||
16 | |||
17 | #include <dumux/freeflow/navierstokes/mass/1pnc/model.hh> | ||
18 | #include <dumux/freeflow/navierstokes/mass/problem.hh> | ||
19 | |||
20 | #include <dumux/freeflow/navierstokes/momentum/problem.hh> | ||
21 | #include <dumux/freeflow/navierstokes/momentum/model.hh> | ||
22 | |||
23 | #include <dumux/flux/maxwellstefanslaw.hh> | ||
24 | #include <dumux/material/fluidsystems/base.hh> | ||
25 | |||
26 | #include <dumux/discretization/fcstaggered.hh> | ||
27 | #include <dumux/discretization/cctpfa.hh> | ||
28 | |||
29 | #include <dumux/multidomain/traits.hh> | ||
30 | #include <dumux/multidomain/freeflow/couplingmanager.hh> | ||
31 | |||
32 | #include "problem.hh" | ||
33 | |||
34 | namespace Dumux::Properties { | ||
35 | |||
36 | // Create new type tags | ||
37 | namespace TTag { | ||
38 | struct MaxwellStefanNCTest {}; | ||
39 | struct MaxwellStefanTestMomentum { using InheritsFrom = std::tuple<MaxwellStefanNCTest, NavierStokesMomentum, FaceCenteredStaggeredModel>; }; | ||
40 | struct MaxwellStefanTestMass { using InheritsFrom = std::tuple<MaxwellStefanNCTest, NavierStokesMassOnePNC, CCTpfaModel>; }; | ||
41 | } // end namespace TTag | ||
42 | |||
43 | template<class TypeTag> | ||
44 | struct ReplaceCompEqIdx<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr int value = 0; }; | ||
45 | |||
46 | // Set the grid type | ||
47 | template<class TypeTag> | ||
48 | struct Grid<TypeTag, TTag::MaxwellStefanNCTest> | ||
49 | { using type = Dune::YaspGrid<2>; }; | ||
50 | |||
51 | // Set the problem property | ||
52 | template<class TypeTag> | ||
53 | struct Problem<TypeTag, TTag::MaxwellStefanTestMomentum> | ||
54 | { using type = MaxwellStefanNCTestProblem<TypeTag, Dumux::NavierStokesMomentumProblem<TypeTag>>; }; | ||
55 | |||
56 | template<class TypeTag> | ||
57 | struct Problem<TypeTag, TTag::MaxwellStefanTestMass> | ||
58 | { using type = MaxwellStefanNCTestProblem<TypeTag, Dumux::NavierStokesMassProblem<TypeTag>>; }; | ||
59 | |||
60 | template<class TypeTag> | ||
61 | struct EnableGridGeometryCache<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr bool value = true; }; | ||
62 | template<class TypeTag> | ||
63 | struct EnableGridFluxVariablesCache<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr bool value = true; }; | ||
64 | template<class TypeTag> | ||
65 | struct EnableGridVolumeVariablesCache<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr bool value = true; }; | ||
66 | |||
67 | template<class TypeTag> | ||
68 | struct UseMoles<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr bool value = true; }; | ||
69 | |||
70 | //! Here we set FicksLaw or MaxwellStefansLaw | ||
71 | template<class TypeTag> | ||
72 | struct MolecularDiffusionType<TypeTag, TTag::MaxwellStefanNCTest> { using type = MaxwellStefansLaw<TypeTag>; }; | ||
73 | |||
74 | /*! | ||
75 | * \ingroup NavierStokesNCTests | ||
76 | * \brief A simple fluid system with three components for testing the multi-component diffusion with the Maxwell-Stefan formulation. | ||
77 | */ | ||
78 | template<class TypeTag> | ||
79 | class MaxwellStefanFluidSystem | ||
80 | : public FluidSystems::Base<GetPropType<TypeTag, Properties::Scalar>, MaxwellStefanFluidSystem<TypeTag>> | ||
81 | |||
82 | { | ||
83 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
84 | using ThisType = MaxwellStefanFluidSystem<TypeTag>; | ||
85 | using Base = FluidSystems::Base<Scalar, ThisType>; | ||
86 | |||
87 | public: | ||
88 | //! The number of phases | ||
89 | static constexpr int numPhases = 1; | ||
90 | static constexpr int numComponents = 3; | ||
91 | |||
92 | static constexpr int H2Idx = 0;//first major component | ||
93 | static constexpr int N2Idx = 1;//second major component | ||
94 | static constexpr int CO2Idx = 2;//secondary component | ||
95 | |||
96 | //! Human readable component name (index compIdx) (for vtk output) | ||
97 | 8 | static std::string componentName(int compIdx) | |
98 |
2/6✓ Branch 2 taken 8 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
8 | { return "MaxwellStefan_" + std::to_string(compIdx); } |
99 | |||
100 | //! Human readable phase name (index phaseIdx) (for velocity vtk output) | ||
101 | ✗ | static std::string phaseName(int phaseIdx = 0) | |
102 |
2/8✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 19 taken 2 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
|
16 | { return "Gas"; } |
103 | |||
104 | //! Molar mass in kg/mol of the component with index compIdx | ||
105 | ✗ | static Scalar molarMass(unsigned int compIdx) | |
106 | ✗ | { return 0.02896; } | |
107 | |||
108 | //! Returns whether the fluids are miscible | ||
109 | static constexpr bool isMiscible() | ||
110 | { return false; } | ||
111 | |||
112 | //! Returns whether the fluids are compressible | ||
113 | static constexpr bool isCompressible(int phaseIdx = 0) | ||
114 | { return false; } | ||
115 | |||
116 | using Base::binaryDiffusionCoefficient; | ||
117 | /*! | ||
118 | * \brief Given a phase's composition, temperature and pressure, | ||
119 | * returns the binary diffusion coefficient \f$\mathrm{[m^2/s]}\f$ for components | ||
120 | * \f$i\f$ and \f$j\f$ in this phase. | ||
121 | * | ||
122 | * \param fluidState An arbitrary fluid state | ||
123 | * \param phaseIdx The index of the fluid phase to consider | ||
124 | * \param compIIdx The index of the first component to consider | ||
125 | * \param compJIdx The index of the second component to consider | ||
126 | */ | ||
127 | template <class FluidState> | ||
128 | 1287900 | static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, | |
129 | int phaseIdx, | ||
130 | int compIIdx, | ||
131 | int compJIdx) | ||
132 | { | ||
133 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1287900 times.
|
1287900 | if (compIIdx > compJIdx) |
134 | { | ||
135 | using std::swap; | ||
136 | ✗ | swap(compIIdx, compJIdx); | |
137 | } | ||
138 | |||
139 |
4/4✓ Branch 0 taken 858600 times.
✓ Branch 1 taken 429300 times.
✓ Branch 2 taken 429300 times.
✓ Branch 3 taken 429300 times.
|
1287900 | if (compIIdx == H2Idx && compJIdx == N2Idx) |
140 | return 83.3e-6; | ||
141 |
3/4✓ Branch 0 taken 429300 times.
✓ Branch 1 taken 429300 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 429300 times.
|
858600 | if (compIIdx == H2Idx && compJIdx == CO2Idx) |
142 | return 68.0e-6; | ||
143 |
2/4✓ Branch 0 taken 429300 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 429300 times.
|
429300 | if (compIIdx == N2Idx && compJIdx == CO2Idx) |
144 | return 16.8e-6; | ||
145 | ✗ | DUNE_THROW(Dune::InvalidStateException, | |
146 | "Binary diffusion coefficient of components " | ||
147 | << compIIdx << " and " << compJIdx << " is undefined!\n"); | ||
148 | } | ||
149 | using Base::density; | ||
150 | /*! | ||
151 | * \brief Given a phase's composition, temperature, pressure, and | ||
152 | * the partial pressures of all components, returns its | ||
153 | * density \f$\mathrm{[kg/m^3]}\f$. | ||
154 | * | ||
155 | * \param phaseIdx index of the phase | ||
156 | * \param fluidState the fluid state | ||
157 | * | ||
158 | */ | ||
159 | template <class FluidState> | ||
160 | ✗ | static Scalar density(const FluidState &fluidState, | |
161 | const int phaseIdx) | ||
162 | { | ||
163 | ✗ | return 1; | |
164 | } | ||
165 | |||
166 | using Base::viscosity; | ||
167 | /*! | ||
168 | * \brief Calculates the dynamic viscosity of a fluid phase \f$\mathrm{[Pa*s]}\f$ | ||
169 | * | ||
170 | * \param fluidState An arbitrary fluid state | ||
171 | * \param phaseIdx The index of the fluid phase to consider | ||
172 | */ | ||
173 | template <class FluidState> | ||
174 | ✗ | static Scalar viscosity(const FluidState &fluidState, | |
175 | int phaseIdx) | ||
176 | { | ||
177 | ✗ | return 1e-6; | |
178 | } | ||
179 | |||
180 | using Base::molarDensity; | ||
181 | /*! | ||
182 | * \brief The molar density \f$\rho_{mol,\alpha}\f$ | ||
183 | * of a fluid phase \f$\alpha\f$ in \f$\mathrm{[mol/m^3]}\f$ | ||
184 | * | ||
185 | * The molar density for the simple relation is defined by the | ||
186 | * mass density \f$\rho_\alpha\f$ and the molar mass of the main component \f$M_\kappa\f$: | ||
187 | * | ||
188 | * \f[\rho_{mol,\alpha} = \frac{\rho_\alpha}{M_\kappa} \;.\f] | ||
189 | */ | ||
190 | template <class FluidState> | ||
191 | ✗ | static Scalar molarDensity(const FluidState &fluidState, int phaseIdx) | |
192 | { | ||
193 | 429300 | return density(fluidState, phaseIdx)/molarMass(0); | |
194 | } | ||
195 | }; | ||
196 | |||
197 | template<class TypeTag> | ||
198 | struct FluidSystem<TypeTag, TTag::MaxwellStefanNCTest> { using type = MaxwellStefanFluidSystem<TypeTag>; }; | ||
199 | |||
200 | // Set the coupling manager | ||
201 | template<class TypeTag> | ||
202 | struct CouplingManager<TypeTag, TTag::MaxwellStefanNCTest> | ||
203 | { | ||
204 | using Traits = MultiDomainTraits<TTag::MaxwellStefanTestMomentum, TTag::MaxwellStefanTestMass>; | ||
205 | using type = FreeFlowCouplingManager<Traits>; | ||
206 | }; | ||
207 | |||
208 | |||
209 | } // end namespace Dumux::Properties | ||
210 | |||
211 | #endif | ||
212 |