GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: dumux/test/porousmediumflow/1pnc/dispersion/problem.hh
Date: 2025-04-12 19:19:20
Exec Total Coverage
Lines: 58 58 100.0%
Functions: 12 12 100.0%
Branches: 40 52 76.9%

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 * \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 a 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 theboundary, where again Dirichlet boundary
43 * conditions are applied.
44 *
45 * This problem uses the \ref OnePNCModel model.
46 */
47 template <class TypeTag>
48 7 class OnePNCDispersionProblem : 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 #if NONISOTHERMAL
74 energyEqIdx = Indices::energyEqIdx,
75 #endif
76 // indices of the primary variables
77 pressureIdx = Indices::pressureIdx,
78
79 // component indices
80 H2OIdx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::H2OIdx),
81 N2Idx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::N2Idx),
82 // indices of the equations
83 contiH2OEqIdx = Indices::conti0EqIdx + H2OIdx,
84 contiN2EqIdx = Indices::conti0EqIdx + N2Idx
85 };
86
87 //! Property that defines whether mole or mass fractions are used
88 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
89 static constexpr int numFluidComps = FluidSystem::numComponents;
90 static const int dimWorld = GridView::dimensionworld;
91 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
92
93 public:
94 7 OnePNCDispersionProblem(std::shared_ptr<const GridGeometry> gridGeometry)
95
3/6
✓ Branch 2 taken 7 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 7 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 7 times.
✗ Branch 9 not taken.
21 : ParentType(gridGeometry)
96 {
97 //initialize fluid system
98 FluidSystem::init();
99
100 // stating in the console whether mole or mass fractions are used
101 if constexpr (useMoles)
102
1/2
✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
7 std::cout<<"problem uses mole fractions" << '\n';
103 else
104 std::cout<<"problem uses mass fractions" << '\n';
105
106
1/2
✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
7 boundaryConcentration_ = getParam<Scalar>("Problem.BoundaryConcentration");
107
1/2
✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
7 temperatureDifference_ = getParam<Scalar>("Problem.BoundaryTemperatureDifference");
108
1/2
✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
7 counterFlowRate_ = getParam<Scalar>("Problem.CounterFlowRate");
109
1/2
✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
7 pressure_ = getParam<Scalar>("Problem.Pressure");
110
1/2
✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
7 problemName_ = getParam<std::string>("Problem.Name");
111 7 }
112
113 /*!
114 * \brief The problem name.
115 */
116 7 const std::string& name() const
117
1/2
✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
7 { return problemName_; }
118
119 /*!
120 * \brief Specifies which kind of boundary condition should be
121 * used for which equation on a given boundary segment.
122 */
123 170196 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
124 {
125
2/2
✓ Branch 0 taken 418704 times.
✓ Branch 1 taken 170196 times.
588900 BoundaryTypes values;
126
127 170196 values.setAllDirichlet();
128
129 296752 if (isRightBoundary_(globalPos) ||
130
4/4
✓ Branch 0 taken 126556 times.
✓ Branch 1 taken 43640 times.
✓ Branch 2 taken 42549 times.
✓ Branch 3 taken 84007 times.
170196 isLowerBoundary_(globalPos) ||
131
2/2
✓ Branch 0 taken 42549 times.
✓ Branch 1 taken 41458 times.
84007 isUpperBoundary_(globalPos) )
132 170196 values.setAllNeumann();
133
134 170196 return values;
135 }
136
137 /*!
138 * \brief Evaluate the boundary conditions for a dirichlet
139 * boundary segment.
140 *
141 * For this method, the \a values parameter stores primary variables.
142 */
143 41458 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
144 {
145 41458 PrimaryVariables values = initialAtPos(globalPos);
146
147 if (isLeftBoundary_(globalPos))
148 {
149 41458 values[contiN2EqIdx] = boundaryConcentration_;
150 // assumes simplified density
151 if (!useMoles)
152 values[contiN2EqIdx] *= FluidSystem::molarMass(N2Idx) / FluidSystem::molarMass(H2OIdx);
153 #if NONISOTHERMAL
154 19076 values[energyEqIdx] = this->spatialParams().temperatureAtPos(globalPos) + temperatureDifference_;
155 #endif
156 }
157
158 41458 return values;
159 }
160
161 /*!
162 * \brief Evaluates the boundary conditions for a neumann
163 * boundary segment (implementation for the box scheme).
164 *
165 * This is the method for the case where the Neumann condition is
166 * potentially solution dependent and requires some quantities that
167 * are specific to the fully-implicit method.
168 *
169 * \param element The finite element
170 * \param fvGeometry The finite-volume geometry
171 * \param elemVolVars All volume variables for the element
172 * \param elemFluxVarsCache Flux variables caches for all faces in stencil
173 * \param scvf The sub-control volume face
174 *
175 * For this method, the \a values parameter stores the flux
176 * in normal direction of each phase. Negative values mean influx.
177 * E.g. for the mass balance that would be the mass flux in \f$ [ kg / (m^2 \cdot s)] \f$.
178 */
179 783606 NumEqVector neumann(const Element& element,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars,
182 const ElementFluxVariablesCache& elemFluxVarsCache,
183 const SubControlVolumeFace& scvf) const
184 {
185 783606 const auto globalPos = element.geometry().corner(scvf.insideScvIdx());
186 783606 NumEqVector values(0.0);
187
2/2
✓ Branch 0 taken 269766 times.
✓ Branch 1 taken 513840 times.
783606 if (isRightBoundary_(globalPos))
188 269766 values[contiH2OEqIdx] = -1.0 * counterFlowRate_ * (useMoles ? 1.0 : FluidSystem::molarMass(H2OIdx));
189 783606 return values;
190 }
191
192 /*!
193 * \brief Evaluates the source term for all phases within a given
194 * sub-control volume.
195 *
196 * For this method, the \a priVars parameter stores the rate mass
197 * of a component is generated or annihilated per volume
198 * unit. Positive values mean that mass is created, negative ones
199 * mean that it vanishes.
200 *
201 * The units must be according to either using mole or mass fractions (mole/(m^3*s) or kg/(m^3*s)).
202 */
203 10276800 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
204
2/2
✓ Branch 0 taken 16867200 times.
✓ Branch 1 taken 5622400 times.
37420800 { return NumEqVector(0.0); }
205
206 //#### initial value for a control volume
207 44545 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
208 {
209 44545 PrimaryVariables values;
210
211 44545 GetPropType<TypeTag, Properties::FluidState> fluidState;
212 44545 const Scalar depth = this->gridGeometry().bBoxMax()[1] - globalPos[1];
213 44545 const Scalar gravity = this->spatialParams().gravity(globalPos)[1];
214 44545 const Scalar moleFraction = 0.0; //initial condition for the molefraction
215 44545 fluidState.setTemperature(this->spatialParams().temperatureAtPos(globalPos));
216 44545 fluidState.setPressure(0, pressure_);
217 44545 fluidState.setMoleFraction(0, N2Idx, moleFraction);
218 44545 fluidState.setMoleFraction(0, H2OIdx, 1.0-moleFraction);
219 44545 const Scalar density = FluidSystem::density(fluidState, 0);
220 44545 const Scalar molarMass = density / FluidSystem::molarDensity(fluidState, 0);
221 44545 const Scalar massFraction = moleFraction * FluidSystem::molarMass(N2Idx) / molarMass;
222
223 44545 values[Indices::pressureIdx] = pressure_ - (density * gravity * depth); //initial condition for the pressure
224 44545 values[contiN2EqIdx] = useMoles ? moleFraction : massFraction;
225 #if NONISOTHERMAL
226 20399 values[energyEqIdx] = this->spatialParams().temperatureAtPos(globalPos);
227 #endif
228 44545 return values;
229 }
230
231 private:
232
1/2
✓ Branch 0 taken 41458 times.
✗ Branch 1 not taken.
41458 bool isLeftBoundary_(const GlobalPosition& globalPos) const
233
1/2
✓ Branch 0 taken 41458 times.
✗ Branch 1 not taken.
41458 { return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
234
235
4/4
✓ Branch 0 taken 269766 times.
✓ Branch 1 taken 513840 times.
✓ Branch 2 taken 126556 times.
✓ Branch 3 taken 43640 times.
953802 bool isRightBoundary_(const GlobalPosition& globalPos) const
236
4/4
✓ Branch 0 taken 269766 times.
✓ Branch 1 taken 513840 times.
✓ Branch 2 taken 126556 times.
✓ Branch 3 taken 43640 times.
953802 { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
237
238
2/2
✓ Branch 0 taken 42549 times.
✓ Branch 1 taken 84007 times.
126556 bool isLowerBoundary_(const GlobalPosition& globalPos) const
239
2/2
✓ Branch 0 taken 42549 times.
✓ Branch 1 taken 84007 times.
126556 { return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_; }
240
241
2/2
✓ Branch 0 taken 42549 times.
✓ Branch 1 taken 41458 times.
84007 bool isUpperBoundary_(const GlobalPosition& globalPos) const
242
2/2
✓ Branch 0 taken 42549 times.
✓ Branch 1 taken 41458 times.
84007 { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
243
244 static constexpr Scalar eps_ = 1e-6;
245 Scalar pressure_;
246 Scalar counterFlowRate_;
247 Scalar temperatureDifference_;
248 Scalar boundaryConcentration_;
249 std::string problemName_;
250 };
251
252 } // end namespace Dumux
253
254 #endif
255