3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/freeflow/navierstokesnc/maxwellstefan/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_MAXWELL_STEFAN_TEST_PROBLEM_HH
26#define DUMUX_CHANNEL_MAXWELL_STEFAN_TEST_PROBLEM_HH
27
28#include <dune/grid/yaspgrid.hh>
29
32
35
38
40
41namespace Dumux {
42template <class TypeTag>
43class MaxwellStefanNCTestProblem;
44
45namespace Properties {
46
47// Create new type tags
48namespace TTag {
49struct MaxwellStefanNCTest { using InheritsFrom = std::tuple<NavierStokesNC, StaggeredFreeFlowModel>; };
50} // end namespace TTag
51
52template<class TypeTag>
53struct ReplaceCompEqIdx<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr int value = 0; };
54
55// Set the grid type
56template<class TypeTag>
57struct Grid<TypeTag, TTag::MaxwellStefanNCTest> { using type = Dune::YaspGrid<2>; };
58
59// Set the problem property
60template<class TypeTag>
61struct Problem<TypeTag, TTag::MaxwellStefanNCTest> { using type = Dumux::MaxwellStefanNCTestProblem<TypeTag> ; };
62
63template<class TypeTag>
64struct EnableGridGeometryCache<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr bool value = true; };
65
66template<class TypeTag>
67struct EnableGridFluxVariablesCache<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr bool value = true; };
68template<class TypeTag>
69struct EnableGridVolumeVariablesCache<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr bool value = true; };
70
71template<class TypeTag>
72struct UseMoles<TypeTag, TTag::MaxwellStefanNCTest> { static constexpr bool value = true; };
73
74
76template<class TypeTag>
77struct MolecularDiffusionType<TypeTag, TTag::MaxwellStefanNCTest> { using type = MaxwellStefansLaw<TypeTag>; };
78
79
84template<class TypeTag>
86: public FluidSystems::Base<GetPropType<TypeTag, Properties::Scalar>, MaxwellStefanFluidSystem<TypeTag>>
87
88{
92
93public:
95 static constexpr int numPhases = 1;
96 static constexpr int numComponents = 3;
97
98 static constexpr int H2Idx = 0;//first major component
99 static constexpr int N2Idx = 1;//second major component
100 static constexpr int CO2Idx = 2;//secondary component
101
103 static std::string componentName(int compIdx)
104 { return "MaxwellStefan_" + std::to_string(compIdx); }
105
107 static std::string phaseName(int phaseIdx = 0)
108 { return "Gas"; }
109
111 static Scalar molarMass(unsigned int compIdx)
112 { return 0.02896; }
113
114
126 template <class FluidState>
127 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
128 int phaseIdx,
129 int compIIdx,
130 int compJIdx)
131 {
132 if (compIIdx > compJIdx)
133 {
134 using std::swap;
135 swap(compIIdx, compJIdx);
136 }
137
138 if (compIIdx == H2Idx && compJIdx == N2Idx)
139 return 83.3e-6;
140 if (compIIdx == H2Idx && compJIdx == CO2Idx)
141 return 68.0e-6;
142 if (compIIdx == N2Idx && compJIdx == CO2Idx)
143 return 16.8e-6;
144 DUNE_THROW(Dune::InvalidStateException,
145 "Binary diffusion coefficient of components "
146 << compIIdx << " and " << compJIdx << " is undefined!\n");
147 }
148 using Base::density;
158 template <class FluidState>
159 static Scalar density(const FluidState &fluidState,
160 const int phaseIdx)
161 {
162 return 1;
163 }
164
165 using Base::viscosity;
172 template <class FluidState>
173 static Scalar viscosity(const FluidState &fluidState,
174 int phaseIdx)
175 {
176 return 1e-6;
177 }
178
179 using Base::molarDensity;
189 template <class FluidState>
190 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
191 {
192 return density(fluidState, phaseIdx)/molarMass(0);
193 }
194};
195
196template<class TypeTag>
197struct FluidSystem<TypeTag, TTag::MaxwellStefanNCTest> { using type = MaxwellStefanFluidSystem<TypeTag>; };
198
199} // end namespace Properties
200
205template <class TypeTag>
207{
208 using ParentType = NavierStokesProblem<TypeTag>;
209
217
218 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
219 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
220
221 using TimeLoopPtr = std::shared_ptr<CheckPointTimeLoop<Scalar>>;
222
223 enum {
224 compTwoIdx = Indices::conti0EqIdx + FluidSystem::N2Idx,
225 compThreeIdx = Indices::conti0EqIdx + FluidSystem::CO2Idx,
226 };
227
228public:
229 MaxwellStefanNCTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
230 : ParentType(gridGeometry)
231 {
232 plotOutput_ = getParam<bool>("Problem.PlotOutput", false);
233 }
234
238 // \{
239
249 template<class SolutionVector, class GridVariables>
250 void plotComponentsOverTime(const SolutionVector& curSol,
251 const GridVariables& gridVariables,
252 const Scalar time)
253 {
254
255 if (plotOutput_)
256 {
257 Scalar x_co2_left = 0.0;
258 Scalar x_n2_left = 0.0;
259 Scalar x_co2_right = 0.0;
260 Scalar x_n2_right = 0.0;
261 Scalar x_h2_left = 0.0;
262 Scalar x_h2_right = 0.0;
263 Scalar i = 0.0;
264 Scalar j = 0.0;
265 for (const auto& element : elements(this->gridGeometry().gridView()))
266 {
267 auto fvGeometry = localView(this->gridGeometry());
268 fvGeometry.bindElement(element);
269
270 auto elemVolVars = localView(gridVariables.curGridVolVars());
271 elemVolVars.bind(element, fvGeometry, curSol);
272 for (auto&& scv : scvs(fvGeometry))
273 {
274 const auto globalPos = scv.dofPosition();
275
276 if (globalPos[0] < 0.5)
277 {
278 x_co2_left += elemVolVars[scv].moleFraction(FluidSystem::CO2Idx);
279 x_n2_left += elemVolVars[scv].moleFraction(FluidSystem::N2Idx);
280 x_h2_left += elemVolVars[scv].moleFraction(FluidSystem::H2Idx);
281 i +=1;
282 }
283 else
284 {
285 x_co2_right += elemVolVars[scv].moleFraction(FluidSystem::CO2Idx);
286 x_n2_right += elemVolVars[scv].moleFraction(FluidSystem::N2Idx);
287 x_h2_right += elemVolVars[scv].moleFraction(FluidSystem::H2Idx);
288 j +=1;
289 }
290
291 }
292 }
293 x_co2_left /= i;
294 x_n2_left /= i;
295 x_h2_left /= i;
296 x_co2_right /= j;
297 x_n2_right /= j;
298 x_h2_right /= j;
299
300 //do a gnuplot
301 x_.push_back(time); // in seconds
302 y1_.push_back(x_n2_left);
303 y2_.push_back(x_n2_right);
304 y3_.push_back(x_co2_left);
305 y4_.push_back(x_co2_right);
306 y5_.push_back(x_h2_left);
307 y6_.push_back(x_h2_right);
308
309 gnuplot_.resetPlot();
310 gnuplot_.setXRange(0, std::max(time, 72000.0));
311 gnuplot_.setYRange(0.4, 0.6);
312 gnuplot_.setXlabel("time [s]");
313 gnuplot_.setYlabel("mole fraction mol/mol");
314 gnuplot_.addDataSetToPlot(x_, y1_, "N2_left.dat", "w l t 'N_2 left'");
315 gnuplot_.addDataSetToPlot(x_, y2_, "N2_right.dat", "w l t 'N_2 right'");
316 gnuplot_.plot("mole_fraction_N2");
317
318 gnuplot2_.resetPlot();
319 gnuplot2_.setXRange(0, std::max(time, 72000.0));
320 gnuplot2_.setYRange(0.0, 0.6);
321 gnuplot2_.setXlabel("time [s]");
322 gnuplot2_.setYlabel("mole fraction mol/mol");
323 gnuplot2_.addDataSetToPlot(x_, y3_, "CO2_left.dat", "w l t 'CO_2 left'");
324 gnuplot2_.addDataSetToPlot(x_, y4_, "CO2_right.dat", "w l t 'CO_2 right'");
325 gnuplot2_.plot("mole_fraction_C02");
326
327 gnuplot3_.resetPlot();
328 gnuplot3_.setXRange(0, std::max(time, 72000.0));
329 gnuplot3_.setYRange(0.0, 0.6);
330 gnuplot3_.setXlabel("time [s]");
331 gnuplot3_.setYlabel("mole fraction mol/mol");
332 gnuplot3_.addDataSetToPlot(x_, y5_, "H2_left.dat", "w l t 'H_2 left'");
333 gnuplot3_.addDataSetToPlot(x_, y6_, "H2_right.dat", "w l t 'H_2 right'");
334 gnuplot3_.plot("mole_fraction_H2");
335 }
336
337 }
338
344 Scalar temperature() const
345 { return 273.15 + 10; } // 10C
346
352 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
353 {
354 return NumEqVector(0.0);
355 }
356
357 // \}
361 // \{
362
369 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
370 {
371 BoundaryTypes values;
372 // set Dirichlet values for the velocity everywhere
373 values.setDirichlet(Indices::velocityXIdx);
374 values.setDirichlet(Indices::velocityYIdx);
375 values.setOutflow(compTwoIdx);
376 values.setOutflow(compThreeIdx);
377 return values;
378 }
379
385 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
386 {
387 return initialAtPos(globalPos);
388 }
389
393 // \{
394
400 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
401 {
402 PrimaryVariables initialValues(0.0);
403 if (globalPos[0] < 0.5)
404 {
405 initialValues[compTwoIdx] = 0.50086;
406 initialValues[compThreeIdx] = 0.49914;
407 }
408 else
409 {
410 initialValues[compTwoIdx] = 0.49879;
411 initialValues[compThreeIdx] = 0.0;
412 }
413
414 initialValues[Indices::pressureIdx] = 1.1e+5;
415 initialValues[Indices::velocityXIdx] = 0.0;
416 initialValues[Indices::velocityYIdx] = 0.0;
417
418 return initialValues;
419 }
420
421private:
422
423 bool plotOutput_;
424
425 const Scalar eps_{1e-6};
426
430
431 std::vector<Scalar> x_;
432 std::vector<Scalar> y1_;
433 std::vector<Scalar> y2_;
434 std::vector<Scalar> y3_;
435 std::vector<Scalar> y4_;
436 std::vector<Scalar> y5_;
437 std::vector<Scalar> y6_;
438};
439} // end namespace Dumux
440
441#endif
This file contains the data which is required to calculate diffusive mass fluxes due to molecular dif...
Interface for passing data sets to a file and plotting them, if gnuplot is installed.
A much simpler (and thus potentially less buggy) version of pure water.
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 for the calculation of the molecular diffusion fluxes.
Definition: common/properties.hh:212
The type of the fluid system to use.
Definition: common/properties.hh:223
The type of the fluid state to use.
Definition: common/properties.hh:225
Definition: maxwellstefanslaw.hh:37
Navier-Stokes problem base class.
Definition: dumux/freeflow/navierstokes/problem.hh:63
Interface for passing data sets to a file and plotting them, if gnuplot is installed.
Definition: gnuplotinterface.hh:57
Fluid system base class.
Definition: fluidsystems/base.hh:45
static Scalar density(const FluidState &fluidState, int phaseIdx)
Calculate the density of a fluid phase.
Definition: fluidsystems/base.hh:134
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: fluidsystems/base.hh:326
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
Calculate the molar density of a fluid phase.
Definition: fluidsystems/base.hh:160
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: fluidsystems/base.hh:236
Test problem for the Maxwell-Stefan model.
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:207
MaxwellStefanNCTestProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:229
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:400
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:352
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:344
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/maxwellstefan/problem.hh:369
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:385
void plotComponentsOverTime(const SolutionVector &curSol, const GridVariables &gridVariables, const Scalar time)
Writes out the diffusion rates from left to right.
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:250
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:49
std::tuple< NavierStokesNC, StaggeredFreeFlowModel > InheritsFrom
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:49
Dune::YaspGrid< 2 > type
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:57
A simple fluid system with three components for testing the multi-component diffusion with the Maxwel...
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:88
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, returns the binary diffusion coefficient for ...
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:127
static constexpr int numComponents
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:96
static constexpr int numPhases
The number of phases.
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:95
static constexpr int CO2Idx
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:100
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculates the dynamic viscosity of a fluid phase .
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:173
static constexpr int H2Idx
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:98
static Scalar molarMass(unsigned int compIdx)
Molar mass in kg/mol of the component with index compIdx.
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:111
static std::string componentName(int compIdx)
Human readable component name (index compIdx) (for vtk output)
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:103
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:190
static std::string phaseName(int phaseIdx=0)
Human readable phase name (index phaseIdx) (for velocity vtk output)
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:107
static constexpr int N2Idx
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:99
static Scalar density(const FluidState &fluidState, const int phaseIdx)
Given a phase's composition, temperature, pressure, and the partial pressures of all components,...
Definition: test/freeflow/navierstokesnc/maxwellstefan/problem.hh:159
Defines a type tag and some properties for ree-flow models using the staggered scheme....