3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/freeflow/navierstokesnc/channel/problem.hh
Go to the documentation of this file.
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_CHANNEL_NC_TEST_PROBLEM_HH
26#define DUMUX_CHANNEL_NC_TEST_PROBLEM_HH
27
28#ifndef ENABLECACHING
29#define ENABLECACHING 1
30#endif
31
32#include <dune/grid/yaspgrid.hh>
33
37
40
42
43namespace Dumux {
44
45template <class TypeTag>
46class ChannelNCTestProblem;
47
48namespace Properties {
49
50// Create new type tags
51namespace TTag {
52#if !NONISOTHERMAL
53struct ChannelNCTest { using InheritsFrom = std::tuple<NavierStokesNC, StaggeredFreeFlowModel>; };
54#else
55struct ChannelNCTest { using InheritsFrom = std::tuple<NavierStokesNCNI, StaggeredFreeFlowModel>; };
56#endif
57} // end namespace TTag
58
59// Select the fluid system
60template<class TypeTag>
61struct FluidSystem<TypeTag, TTag::ChannelNCTest>
62{
64 static constexpr int phaseIdx = H2OAir::liquidPhaseIdx;
66};
67
68template<class TypeTag>
69struct ReplaceCompEqIdx<TypeTag, TTag::ChannelNCTest> { static constexpr int value = 0; };
70
71// Set the grid type
72template<class TypeTag>
73struct Grid<TypeTag, TTag::ChannelNCTest> { using type = Dune::YaspGrid<2>; };
74
75// Set the problem property
76template<class TypeTag>
77struct Problem<TypeTag, TTag::ChannelNCTest> { using type = Dumux::ChannelNCTestProblem<TypeTag> ; };
78
79template<class TypeTag>
80struct EnableGridGeometryCache<TypeTag, TTag::ChannelNCTest> { static constexpr bool value = ENABLECACHING; };
81template<class TypeTag>
82struct EnableGridFluxVariablesCache<TypeTag, TTag::ChannelNCTest> { static constexpr bool value = ENABLECACHING; };
83template<class TypeTag>
84struct EnableGridVolumeVariablesCache<TypeTag, TTag::ChannelNCTest> { static constexpr bool value = ENABLECACHING; };
85template<class TypeTag>
86struct EnableGridFaceVariablesCache<TypeTag, TTag::ChannelNCTest> { static constexpr bool value = ENABLECACHING; };
87
88// Use mole fraction formulation
89#if USE_MASS
90template<class TypeTag>
91struct UseMoles<TypeTag, TTag::ChannelNCTest> { static constexpr bool value = false; };
92#else
93template<class TypeTag>
94struct UseMoles<TypeTag, TTag::ChannelNCTest> { static constexpr bool value = true; };
95#endif
96} // end namespace Properties
97
107template <class TypeTag>
109{
110 using ParentType = NavierStokesProblem<TypeTag>;
111
119
120 static constexpr auto dimWorld = GetPropType<TypeTag, Properties::GridView>::dimensionworld;
121 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
122
123 using TimeLoopPtr = std::shared_ptr<CheckPointTimeLoop<Scalar>>;
124
125 static constexpr auto compIdx = 1;
126 static constexpr auto transportCompIdx = Indices::conti0EqIdx + compIdx;
127 static constexpr auto transportEqIdx = Indices::conti0EqIdx + compIdx;
128
129public:
130 ChannelNCTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
131 : ParentType(gridGeometry), eps_(1e-6)
132 {
133 inletVelocity_ = getParam<Scalar>("Problem.InletVelocity");
134 FluidSystem::init();
135 deltaP_.resize(this->gridGeometry().numCellCenterDofs());
136 }
137
141 // \{
142
144 {
145 return false;
146 }
147
153 Scalar temperature() const
154 { return 273.15 + 10; } // 10C
155
161 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
162 {
163 return NumEqVector(0.0);
164 }
165
166 // \}
170 // \{
171
178 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
179 {
180 BoundaryTypes values;
181
182 if(isInlet_(globalPos))
183 {
184 values.setDirichlet(Indices::velocityXIdx);
185 values.setDirichlet(Indices::velocityYIdx);
186 values.setDirichlet(transportCompIdx);
187#if NONISOTHERMAL
188 values.setDirichlet(Indices::temperatureIdx);
189#endif
190 }
191 else if(isOutlet_(globalPos))
192 {
193 values.setDirichlet(Indices::pressureIdx);
194 values.setOutflow(transportEqIdx);
195#if NONISOTHERMAL
196 values.setOutflow(Indices::energyEqIdx);
197#endif
198 }
199 else
200 {
201 // set Dirichlet values for the velocity everywhere
202 values.setDirichlet(Indices::velocityXIdx);
203 values.setDirichlet(Indices::velocityYIdx);
204 values.setNeumann(Indices::conti0EqIdx);
205 values.setNeumann(transportEqIdx);
206#if NONISOTHERMAL
207 values.setNeumann(Indices::energyEqIdx);
208#endif
209 }
210
211 return values;
212 }
213
219 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
220 {
221 PrimaryVariables values = initialAtPos(globalPos);
222
223 // give the system some time so that the pressure can equilibrate, then start the injection of the tracer
224 if(isInlet_(globalPos))
225 {
226 if(time() >= 10.0 || inletVelocity_ < eps_)
227 {
228 Scalar moleFracTransportedComp = 1e-3;
229#if USE_MASS
230 Scalar averageMolarMassPhase = moleFracTransportedComp * FluidSystem::molarMass(compIdx)
231 + (1. - moleFracTransportedComp) * FluidSystem::molarMass(1-compIdx);
232 values[transportCompIdx] = moleFracTransportedComp * FluidSystem::molarMass(compIdx)
233 / averageMolarMassPhase;
234#else
235 values[transportCompIdx] = moleFracTransportedComp;
236#endif
237#if NONISOTHERMAL
238 values[Indices::temperatureIdx] = 293.15;
239#endif
240 }
241 }
242
243 return values;
244 }
245
246 // \}
247
251 // \{
252
258 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
259 {
260 PrimaryVariables values;
261 values[Indices::pressureIdx] = 1.1e+5;
262 values[transportCompIdx] = 0.0;
263#if NONISOTHERMAL
264 values[Indices::temperatureIdx] = 283.15;
265#endif
266
267 // parabolic velocity profile
268 values[Indices::velocityXIdx] = inletVelocity_*(globalPos[1] - this->gridGeometry().bBoxMin()[1])*(this->gridGeometry().bBoxMax()[1] - globalPos[1])
269 / (0.25*(this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1])*(this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1]));
270
271 values[Indices::velocityYIdx] = 0.0;
272
273 return values;
274 }
275
284 template<class GridVariables, class SolutionVector>
285 void calculateDeltaP(const GridVariables& gridVariables, const SolutionVector& sol)
286 {
287 for (const auto& element : elements(this->gridGeometry().gridView()))
288 {
289 auto fvGeometry = localView(this->gridGeometry());
290 fvGeometry.bindElement(element);
291 for (auto&& scv : scvs(fvGeometry))
292 {
293 auto ccDofIdx = scv.dofIndex();
294
295 auto elemVolVars = localView(gridVariables.curGridVolVars());
296 elemVolVars.bindElement(element, fvGeometry, sol);
297
298 deltaP_[ccDofIdx] = elemVolVars[scv].pressure() - 1.1e5;
299 }
300 }
301 }
302
303 auto& getDeltaP() const
304 { return deltaP_; }
305
306 // \}
307
308 void setTimeLoop(TimeLoopPtr timeLoop)
309 {
310 timeLoop_ = timeLoop;
311 }
312
313 Scalar time() const
314 {
315 return timeLoop_->time();
316 }
317
318private:
319
320 bool isInlet_(const GlobalPosition& globalPos) const
321 {
322 return globalPos[0] < eps_;
323 }
324
325 bool isOutlet_(const GlobalPosition& globalPos) const
326 {
327 return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_;
328 }
329
330 const Scalar eps_;
331 Scalar inletVelocity_;
332 TimeLoopPtr timeLoop_;
333 std::vector<Scalar> deltaP_;
334};
335} // end namespace Dumux
336
337#endif
A much simpler (and thus potentially less buggy) version of pure water.
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
The DUNE grid type.
Definition: common/properties.hh:57
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Property whether to use moles or kg as amount unit for balance equations.
Definition: common/properties.hh:102
The component balance index that should be replaced by the total mass/mole balance.
Definition: common/properties.hh:104
Definition: common/properties.hh:169
If disabled, the volume variables are not stored (reduces memory, but is slower)
Definition: common/properties.hh:178
specifies if data on flux vars should be saved (faster, but more memory consuming)
Definition: common/properties.hh:188
The type of the fluid system to use.
Definition: common/properties.hh:223
Switch on/off caching of face variables.
Definition: common/properties.hh:303
Navier-Stokes problem base class.
Definition: dumux/freeflow/navierstokes/problem.hh:63
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
Definition: 1padapter.hh:46
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
Definition: h2oair.hh:75
Test problem for the one-phase (Navier-)Stokes model.
Definition: test/freeflow/navierstokesnc/channel/problem.hh:109
Scalar time() const
Definition: test/freeflow/navierstokesnc/channel/problem.hh:313
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Specifies which kind of boundary condition should be used for which equation on a given boundary cont...
Definition: test/freeflow/navierstokesnc/channel/problem.hh:178
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: test/freeflow/navierstokesnc/channel/problem.hh:161
auto & getDeltaP() const
Definition: test/freeflow/navierstokesnc/channel/problem.hh:303
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/freeflow/navierstokesnc/channel/problem.hh:258
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: test/freeflow/navierstokesnc/channel/problem.hh:153
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: test/freeflow/navierstokesnc/channel/problem.hh:219
void calculateDeltaP(const GridVariables &gridVariables, const SolutionVector &sol)
Adds additional VTK output data to the VTKWriter.
Definition: test/freeflow/navierstokesnc/channel/problem.hh:285
void setTimeLoop(TimeLoopPtr timeLoop)
Definition: test/freeflow/navierstokesnc/channel/problem.hh:308
ChannelNCTestProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/navierstokesnc/channel/problem.hh:130
bool shouldWriteRestartFile() const
Definition: test/freeflow/navierstokesnc/channel/problem.hh:143
Definition: test/freeflow/navierstokesnc/channel/problem.hh:53
std::tuple< NavierStokesNC, StaggeredFreeFlowModel > InheritsFrom
Definition: test/freeflow/navierstokesnc/channel/problem.hh:53
Dune::YaspGrid< 2 > type
Definition: test/freeflow/navierstokesnc/channel/problem.hh:73
Defines a type tag and some properties for ree-flow models using the staggered scheme....
#define ENABLECACHING
Definition: test/freeflow/navierstokesnc/channel/problem.hh:29