GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/test/porousmediumflow/1pnc/1p2c/isothermal/problem.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 42 65 64.6%
Functions: 6 21 28.6%
Branches: 52 110 47.3%

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 OnePNCTests
10 * \brief Definition of a problem for the 1pnc problem:
11 * Component transport of nitrogen dissolved in the water phase.
12 */
13
14 #ifndef DUMUX_1P2C_TEST_PROBLEM_HH
15 #define DUMUX_1P2C_TEST_PROBLEM_HH
16
17 #include <dumux/common/properties.hh>
18 #include <dumux/common/parameters.hh>
19
20 #include <dumux/common/boundarytypes.hh>
21 #include <dumux/common/numeqvector.hh>
22 #include <dumux/porousmediumflow/problem.hh>
23
24
25 namespace Dumux {
26
27 /*!
28 * \ingroup OnePNCTests
29 * \brief Definition of a problem for the 1pnc problem:
30 * Component transport of nitrogen dissolved in the water phase.
31 *
32 * Nitrogen is dissolved in the water phase and is transported with the
33 * water flow from the left side to the right.
34 *
35 * The model domain is specified in the input file and
36 * we use homogeneous soil properties.
37 * Initially, the domain is filled with pure water.
38 *
39 * At the left side, a Dirichlet condition defines the nitrogen mole fraction.
40 * The water phase flows from the left side to the right if the applied pressure
41 * gradient is >0. The nitrogen is transported with the water flow
42 * and leaves the domain at the boundary, where again Dirichlet boundary
43 * conditions are applied.
44 *
45 * This problem uses the \ref OnePNCModel model.
46 */
47 template <class TypeTag>
48 3 class OnePTwoCTestProblem : public PorousMediumFlowProblem<TypeTag>
49 {
50 using ParentType = PorousMediumFlowProblem<TypeTag>;
51
52 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
53 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
54 using BoundaryTypes = Dumux::BoundaryTypes<GetPropType<TypeTag, Properties::ModelTraits>::numEq()>;
55 using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
56 using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
57 using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
58 using Element = typename GridView::template Codim<0>::Entity;
59
60 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
61 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
62 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
63
64 using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
65 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
66 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
67
68 // copy some indices for convenience
69 using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
70
71 enum
72 {
73 // indices of the primary variables
74 pressureIdx = Indices::pressureIdx,
75
76 // component indices
77 H2OIdx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::H2OIdx),
78 N2Idx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::N2Idx),
79
80 // indices of the equations
81 contiH2OEqIdx = Indices::conti0EqIdx + H2OIdx,
82 contiN2EqIdx = Indices::conti0EqIdx + N2Idx
83 };
84
85 //! Property that defines whether mole or mass fractions are used
86 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
87 static const bool isBox = GetPropType<TypeTag, Properties::GridGeometry>::discMethod == DiscretizationMethods::box;
88
89 static const int dimWorld = GridView::dimensionworld;
90 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
91
92 public:
93 3 OnePTwoCTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
94
7/18
✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 3 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 3 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 3 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 3 times.
✓ Branch 15 taken 3 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 3 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
9 : ParentType(gridGeometry), useNitscheTypeBc_(getParam<bool>("Problem.UseNitscheTypeBc", false))
95 {
96 //initialize fluid system
97
1/2
✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
3 FluidSystem::init();
98
99 // stating in the console whether mole or mass fractions are used
100 if(useMoles)
101
1/2
✓ Branch 2 taken 3 times.
✗ Branch 3 not taken.
6 std::cout<<"problem uses mole fractions"<<std::endl;
102 else
103 std::cout<<"problem uses mass fractions"<<std::endl;
104 3 }
105
106 /*!
107 * \name Problem parameters
108 */
109 // \{
110
111 // \}
112
113 /*!
114 * \name Boundary conditions
115 */
116 // \{
117
118 /*!
119 * \brief Specifies which kind of boundary condition should be
120 * used for which equation on a given boundary segment.
121 *
122 * \param globalPos The position for which the bc type should be evaluated
123 */
124 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
125 {
126 BoundaryTypes values;
127
128 if(globalPos[0] < eps_)
129 values.setAllDirichlet();
130 else
131 values.setAllNeumann();
132
133 return values;
134 }
135
136 /*!
137 * \brief Evaluates the boundary conditions for a Dirichlet boundary segment.
138 *
139 * \param globalPos The position for which the bc type should be evaluated
140 */
141 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
142 {
143 1344 PrimaryVariables values = initial_(globalPos);
144
145 // condition for the N2 molefraction at left boundary
146
3/8
✗ Branch 0 not taken.
✗ Branch 1 not taken.
✓ Branch 2 taken 184 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 168 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 320 times.
✗ Branch 7 not taken.
672 if (globalPos[0] < eps_ )
147 {
148 1344 values[N2Idx] = 2.0e-5;
149 }
150
151 return values;
152 }
153
154 /*!
155 * \brief Evaluates the boundary conditions for a neumann
156 * boundary segment (implementation for the box scheme).
157 *
158 * This is the method for the case where the Neumann condition is
159 * potentially solution dependent and requires some quantities that
160 * are specific to the fully-implicit method.
161 *
162 * \param element The finite element
163 * \param fvGeometry The finite-volume geometry
164 * \param elemVolVars All volume variables for the element
165 * \param elemFluxVarsCache Flux variables caches for all faces in stencil
166 * \param scvf The sub-control volume face
167 *
168 * For this method, the \a values parameter stores the flux
169 * in normal direction of each phase. Negative values mean influx.
170 * E.g. for the mass balance that would be the mass flux in \f$ [ kg / (m^2 \cdot s)] \f$.
171 */
172 template<bool useBox = isBox, std::enable_if_t<useBox, int> = 0>
173 16810 NumEqVector neumann(const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const ElementFluxVariablesCache& elemFluxVarsCache,
177 const SubControlVolumeFace& scvf) const
178 {
179 // no-flow everywhere except at the right boundary
180 16810 NumEqVector values(0.0);
181
6/6
✓ Branch 0 taken 15990 times.
✓ Branch 1 taken 820 times.
✓ Branch 2 taken 15990 times.
✓ Branch 3 taken 820 times.
✓ Branch 4 taken 15990 times.
✓ Branch 5 taken 820 times.
50430 const auto xMax = this->gridGeometry().bBoxMax()[0];
182
2/2
✓ Branch 0 taken 15990 times.
✓ Branch 1 taken 820 times.
16810 const auto& ipGlobal = scvf.ipGlobal();
183
4/4
✓ Branch 0 taken 15990 times.
✓ Branch 1 taken 820 times.
✓ Branch 2 taken 15990 times.
✓ Branch 3 taken 820 times.
33620 if (ipGlobal[0] < xMax - eps_)
184 15990 return values;
185
186 // set a fixed pressure on the right side of the domain
187 static const Scalar dirichletPressure = 1e5;
188
2/4
✓ Branch 0 taken 820 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 820 times.
✗ Branch 3 not taken.
1640 const auto& volVars = elemVolVars[scvf.insideScvIdx()];
189
1/2
✓ Branch 0 taken 820 times.
✗ Branch 1 not taken.
820 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
190
191 // if specified in the input file, use a Nitsche type boundary condition for the box model.
192
1/2
✓ Branch 0 taken 820 times.
✗ Branch 1 not taken.
820 if (useNitscheTypeBc_)
193 {
194 1640 values[contiH2OEqIdx] = (volVars.pressure() - dirichletPressure)*1e7;
195 2460 values[contiN2EqIdx] = values[contiH2OEqIdx] * (useMoles ? volVars.moleFraction(0, N2Idx)
196 : volVars.massFraction(0, N2Idx));
197 820 return values;
198 }
199
200 // evaluate the pressure gradient
201 GlobalPosition gradP(0.0);
202 for (const auto& scv : scvs(fvGeometry))
203 {
204 const auto xIp = scv.dofPosition()[0];
205 auto tmp = fluxVarsCache.gradN(scv.localDofIndex());
206 tmp *= xIp > xMax - eps_ ? dirichletPressure
207 : elemVolVars[scv].pressure();
208 gradP += tmp;
209 }
210
211 // compute flux
212 auto phaseFlux = vtmv(scvf.unitOuterNormal(), volVars.permeability(), gradP);
213 phaseFlux *= useMoles ? volVars.molarDensity() : volVars.density();
214 phaseFlux *= volVars.mobility();
215
216 // set Neumann bc values
217 values[contiH2OEqIdx] = phaseFlux;
218 // emulate an outflow condition for the component transport on the right side
219 values[contiN2EqIdx] = phaseFlux * ( useMoles ? volVars.moleFraction(0, N2Idx) : volVars.massFraction(0, N2Idx) );
220
221 return values;
222 }
223
224 /*!
225 * \brief Evaluates the boundary conditions for a neumann
226 * boundary segment (implementation for the tpfa scheme).
227 */
228 template<bool useBox = isBox, std::enable_if_t<!useBox, int> = 0>
229 14544 NumEqVector neumann(const Element& element,
230 const FVElementGeometry& fvGeometry,
231 const ElementVolumeVariables& elemVolVars,
232 const ElementFluxVariablesCache& elemFluxVarsCache,
233 const SubControlVolumeFace& scvf) const
234 {
235 // no-flow everywhere except at the right boundary
236 14544 NumEqVector values(0.0);
237
6/6
✓ Branch 0 taken 13944 times.
✓ Branch 1 taken 600 times.
✓ Branch 2 taken 13944 times.
✓ Branch 3 taken 600 times.
✓ Branch 4 taken 13944 times.
✓ Branch 5 taken 600 times.
43632 const auto xMax = this->gridGeometry().bBoxMax()[0];
238
2/2
✓ Branch 0 taken 13944 times.
✓ Branch 1 taken 600 times.
14544 const auto& ipGlobal = scvf.ipGlobal();
239
4/4
✓ Branch 0 taken 13944 times.
✓ Branch 1 taken 600 times.
✓ Branch 2 taken 13944 times.
✓ Branch 3 taken 600 times.
29088 if (ipGlobal[0] < xMax - eps_)
240 13944 return values;
241
242 // set a fixed pressure on the right side of the domain
243 static const Scalar dirichletPressure = 1e5;
244 764 const auto& volVars = elemVolVars[scvf.insideScvIdx()];
245
246 1200 auto d = ipGlobal - element.geometry().center();
247 600 d /= d.two_norm2();
248
249 600 auto upwindTerm = useMoles ? volVars.molarDensity() : volVars.density();
250 600 upwindTerm *= volVars.mobility();
251
252 1200 const auto tij = vtmv(scvf.unitOuterNormal(), volVars.permeability(), d);
253 600 const auto phaseFlux = -1.0*upwindTerm*tij*(dirichletPressure - volVars.pressure());
254
255 // set Neumann bc values
256 600 values[contiH2OEqIdx] = phaseFlux;
257 // emulate an outflow condition for the component transport on the right side
258 1200 values[contiN2EqIdx] = phaseFlux * (useMoles ? volVars.moleFraction(0, N2Idx) : volVars.massFraction(0, N2Idx));
259
260 600 return values;
261 }
262
263 // \}
264
265 /*!
266 * \name Volume terms
267 */
268 // \{
269
270 /*!
271 * \brief Evaluates the source term for all phases within a given
272 * sub-control volume.
273 *
274 * For this method, the \a priVars parameter stores the rate mass
275 * of a component is generated or annihilated per volume
276 * unit. Positive values mean that mass is created, negative ones
277 * mean that it vanishes.
278 *
279 * The units must be according to either using mole or mass fractions (mole/(m^3*s) or kg/(m^3*s)).
280 */
281 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
282 77120 { return NumEqVector(0.0); }
283
284 /*!
285 * \brief Evaluates the initial value for a control volume.
286 *
287 * \param globalPos The position for which the initial condition should be evaluated
288 *
289 * For this method, the \a values parameter stores primary
290 * variables.
291 */
292 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
293 160 { return initial_(globalPos); }
294
295 // \}
296
297 private:
298 // the internal method for the initial condition
299 PrimaryVariables initial_(const GlobalPosition &globalPos) const
300 {
301 752 PrimaryVariables priVars;
302
6/16
✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 184 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 184 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 168 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 168 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 320 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 320 times.
✗ Branch 15 not taken.
1504 priVars[pressureIdx] = 2e5 - 1e5*globalPos[0]; // initial condition for the pressure
303
6/16
✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 184 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 184 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 168 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 168 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 320 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 320 times.
✗ Branch 15 not taken.
1504 priVars[N2Idx] = 0.0; // initial condition for the N2 molefraction
304 return priVars;
305 }
306 static constexpr Scalar eps_ = 1e-6;
307 bool useNitscheTypeBc_;
308 };
309
310 } // end namespace Dumux
311
312 #endif
313