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 TwoPTwoCTests | ||
10 | * \brief Problem where air is injected under a low permeable layer in a depth of 2700m. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_INJECTION_PROBLEM_HH | ||
14 | #define DUMUX_INJECTION_PROBLEM_HH | ||
15 | |||
16 | #include <dumux/common/boundarytypes.hh> | ||
17 | #include <dumux/common/properties.hh> | ||
18 | #include <dumux/common/parameters.hh> | ||
19 | #include <dumux/common/numeqvector.hh> | ||
20 | |||
21 | #include <dumux/porousmediumflow/problem.hh> | ||
22 | #include <dumux/material/fluidsystems/h2on2.hh> | ||
23 | |||
24 | namespace Dumux { | ||
25 | |||
26 | /*! | ||
27 | * \ingroup TwoPTwoCTests | ||
28 | * \brief Problem where air is injected under a low permeable layer in a depth of 2700m. | ||
29 | * | ||
30 | */ | ||
31 | template <class TypeTag> | ||
32 | class InjectionProblem : public PorousMediumFlowProblem<TypeTag> | ||
33 | { | ||
34 | using ParentType = PorousMediumFlowProblem<TypeTag>; | ||
35 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
36 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
37 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
38 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
39 | |||
40 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
41 | using Indices = typename ModelTraits::Indices; | ||
42 | |||
43 | // primary variable indices | ||
44 | enum | ||
45 | { | ||
46 | pressureIdx = Indices::pressureIdx, | ||
47 | switchIdx = Indices::switchIdx | ||
48 | }; | ||
49 | |||
50 | // phase presence | ||
51 | enum { wPhaseOnly = Indices::firstPhaseOnly }; | ||
52 | |||
53 | // equation indices | ||
54 | enum | ||
55 | { | ||
56 | contiH2OEqIdx = Indices::conti0EqIdx + FluidSystem::H2OIdx, | ||
57 | contiN2EqIdx = Indices::conti0EqIdx + FluidSystem::N2Idx, | ||
58 | }; | ||
59 | |||
60 | // phase indices | ||
61 | enum | ||
62 | { | ||
63 | gasPhaseIdx = FluidSystem::N2Idx, | ||
64 | H2OIdx = FluidSystem::H2OIdx, | ||
65 | N2Idx = FluidSystem::N2Idx | ||
66 | }; | ||
67 | |||
68 | using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>; | ||
69 | using NumEqVector = Dumux::NumEqVector<PrimaryVariables>; | ||
70 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
71 | using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView; | ||
72 | using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView; | ||
73 | using BoundaryTypes = Dumux::BoundaryTypes<GetPropType<TypeTag, Properties::ModelTraits>::numEq()>; | ||
74 | using Element = typename GridView::template Codim<0>::Entity; | ||
75 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
76 | using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView; | ||
77 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
78 | |||
79 | //! Property that defines whether mole or mass fractions are used | ||
80 | static constexpr bool useMoles = ModelTraits::useMoles(); | ||
81 | |||
82 | public: | ||
83 | 6 | InjectionProblem(std::shared_ptr<const GridGeometry> gridGeometry) | |
84 |
7/20✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 6 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 6 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 6 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 6 times.
✓ Branch 15 taken 6 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 6 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
18 | : ParentType(gridGeometry) |
85 | { | ||
86 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | nTemperature_ = getParam<int>("Problem.NTemperature"); |
87 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | nPressure_ = getParam<int>("Problem.NPressure"); |
88 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | pressureLow_ = getParam<Scalar>("Problem.PressureLow"); |
89 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | pressureHigh_ = getParam<Scalar>("Problem.PressureHigh"); |
90 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | temperatureLow_ = getParam<Scalar>("Problem.TemperatureLow"); |
91 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | temperatureHigh_ = getParam<Scalar>("Problem.TemperatureHigh"); |
92 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | depthBOR_ = getParam<Scalar>("Problem.DepthBOR"); |
93 |
2/4✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
6 | name_ = getParam<std::string>("Problem.Name"); |
94 | |||
95 | // initialize the tables of the fluid system | ||
96 | 12 | FluidSystem::init(/*Tmin=*/temperatureLow_, | |
97 | /*Tmax=*/temperatureHigh_, | ||
98 | 6 | /*nT=*/nTemperature_, | |
99 | /*pmin=*/pressureLow_, | ||
100 | /*pmax=*/pressureHigh_, | ||
101 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | /*np=*/nPressure_); |
102 | |||
103 | // stating in the console whether mole or mass fractions are used | ||
104 | if(useMoles) | ||
105 |
1/2✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
|
12 | std::cout<<"problem uses mole-fractions"<<std::endl; |
106 | else | ||
107 | std::cout<<"problem uses mass-fractions"<<std::endl; | ||
108 | 6 | } | |
109 | |||
110 | /*! | ||
111 | * \brief Returns the problem name | ||
112 | * | ||
113 | * This is used as a prefix for files generated by the simulation. | ||
114 | */ | ||
115 | const std::string& name() const | ||
116 |
1/2✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
|
6 | { return name_; } |
117 | |||
118 | |||
119 | |||
120 | |||
121 | /*! | ||
122 | * \brief Specifies which kind of boundary condition should be | ||
123 | * used for which equation on a given boundary segment | ||
124 | * | ||
125 | * \param globalPos The global position | ||
126 | */ | ||
127 | ✗ | BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const | |
128 | { | ||
129 | ✗ | BoundaryTypes bcTypes; | |
130 | ✗ | if (globalPos[0] < eps_) | |
131 | bcTypes.setAllDirichlet(); | ||
132 | else | ||
133 | bcTypes.setAllNeumann(); | ||
134 | ✗ | return bcTypes; | |
135 | } | ||
136 | |||
137 | /*! | ||
138 | * \brief Evaluates the boundary conditions for a Dirichlet boundary segment | ||
139 | * | ||
140 | * \param globalPos The global position | ||
141 | */ | ||
142 | PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const | ||
143 | { | ||
144 |
1/2✓ Branch 1 taken 5700 times.
✗ Branch 2 not taken.
|
12500 | return initial_(globalPos); |
145 | } | ||
146 | |||
147 | /*! | ||
148 | * \brief Evaluates the boundary conditions for a Neumann | ||
149 | * boundary segment in dependency on the current solution. | ||
150 | * | ||
151 | * \param element The finite element | ||
152 | * \param fvGeometry The finite volume geometry of the element | ||
153 | * \param elemVolVars All volume variables for the element | ||
154 | * \param elemFluxVarsCache Flux variables caches for all faces in stencil | ||
155 | * \param scvf The sub-control volume face | ||
156 | * | ||
157 | * This method is used for cases, when the Neumann condition depends on the | ||
158 | * solution and requires some quantities that are specific to the fully-implicit method. | ||
159 | * The \a values store the mass flux of each phase normal to the boundary. | ||
160 | * Negative values indicate an inflow. | ||
161 | */ | ||
162 | 115308 | NumEqVector neumann(const Element& element, | |
163 | const FVElementGeometry& fvGeometry, | ||
164 | const ElementVolumeVariables& elemVolVars, | ||
165 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
166 | const SubControlVolumeFace& scvf) const | ||
167 | { | ||
168 | 115308 | NumEqVector values(0.0); | |
169 | |||
170 |
2/2✓ Branch 0 taken 29116 times.
✓ Branch 1 taken 34136 times.
|
115308 | const auto& globalPos = scvf.ipGlobal(); |
171 | |||
172 | 115308 | Scalar injectedPhaseMass = 1e-3; | |
173 |
6/6✓ Branch 0 taken 29116 times.
✓ Branch 1 taken 53348 times.
✓ Branch 2 taken 56112 times.
✓ Branch 3 taken 39984 times.
✓ Branch 4 taken 29116 times.
✓ Branch 5 taken 34136 times.
|
252820 | Scalar moleFracW = elemVolVars[scvf.insideScvIdx()].moleFraction(gasPhaseIdx, H2OIdx); |
174 |
8/8✓ Branch 0 taken 53488 times.
✓ Branch 1 taken 61820 times.
✓ Branch 2 taken 53488 times.
✓ Branch 3 taken 61820 times.
✓ Branch 4 taken 5496 times.
✓ Branch 5 taken 47992 times.
✓ Branch 6 taken 5496 times.
✓ Branch 7 taken 47992 times.
|
230616 | if (globalPos[1] < 14 - eps_ && globalPos[1] > 6.5 - eps_) |
175 | { | ||
176 | 5496 | values[contiN2EqIdx] = -(1-moleFracW)*injectedPhaseMass/FluidSystem::molarMass(N2Idx); //mole/(m^2*s) -> kg/(s*m^2) | |
177 | 10992 | values[contiH2OEqIdx] = -moleFracW*injectedPhaseMass/FluidSystem::molarMass(H2OIdx); //mole/(m^2*s) -> kg/(s*m^2) | |
178 | } | ||
179 | 115308 | return values; | |
180 | } | ||
181 | |||
182 | |||
183 | /*! | ||
184 | * \brief Evaluates the initial values for a control volume. | ||
185 | * | ||
186 | * \param globalPos The global position | ||
187 | */ | ||
188 | PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const | ||
189 | 2386 | { return initial_(globalPos); } | |
190 | |||
191 | // \} | ||
192 | |||
193 | private: | ||
194 | /*! | ||
195 | * \brief Evaluates the initial values for a control volume. | ||
196 | * | ||
197 | * The internal method for the initial condition | ||
198 | * | ||
199 | * \param globalPos The global position | ||
200 | */ | ||
201 | 14886 | PrimaryVariables initial_(const GlobalPosition &globalPos) const | |
202 | { | ||
203 | 14886 | PrimaryVariables priVars(0.0); | |
204 | 14886 | priVars.setState(wPhaseOnly); | |
205 | |||
206 | Scalar densityW = | ||
207 | 29772 | FluidSystem::H2O::liquidDensity(this->spatialParams().temperatureAtPos(globalPos), 1e5); | |
208 | |||
209 | 59544 | Scalar pl = 1e5 - densityW*this->spatialParams().gravity(globalPos)[1]*(depthBOR_ - globalPos[1]); | |
210 | 29772 | Scalar moleFracLiquidN2 = pl*0.95/BinaryCoeff::H2O_N2::henry(this->spatialParams().temperatureAtPos(globalPos)); | |
211 | 14886 | Scalar moleFracLiquidH2O = 1.0 - moleFracLiquidN2; | |
212 | |||
213 | 14886 | Scalar meanM = | |
214 | 14886 | FluidSystem::molarMass(H2OIdx)*moleFracLiquidH2O + | |
215 | 14886 | FluidSystem::molarMass(N2Idx)*moleFracLiquidN2; | |
216 | if(useMoles) | ||
217 | { | ||
218 | //mole-fraction formulation | ||
219 | 29772 | priVars[switchIdx] = moleFracLiquidN2; | |
220 | } | ||
221 | else | ||
222 | { | ||
223 | //mass fraction formulation | ||
224 | Scalar massFracLiquidN2 = moleFracLiquidN2*FluidSystem::molarMass(N2Idx)/meanM; | ||
225 | priVars[switchIdx] = massFracLiquidN2; | ||
226 | } | ||
227 | 29772 | priVars[pressureIdx] = pl; | |
228 | 14886 | return priVars; | |
229 | } | ||
230 | |||
231 | Scalar depthBOR_; | ||
232 | static constexpr Scalar eps_ = 1e-6; | ||
233 | |||
234 | int nTemperature_; | ||
235 | int nPressure_; | ||
236 | Scalar pressureLow_, pressureHigh_; | ||
237 | Scalar temperatureLow_, temperatureHigh_; | ||
238 | |||
239 | std::string name_; | ||
240 | }; | ||
241 | |||
242 | } // end namespace Dumux | ||
243 | |||
244 | #endif | ||
245 |