3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/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 *****************************************************************************/
27#ifndef DUMUX_TWOPTWOC_NONEQUILIBRIUM_PROBLEM_HH
28#define DUMUX_TWOPTWOC_NONEQUILIBRIUM_PROBLEM_HH
29
30#include <dune/grid/yaspgrid.hh>
31
32#include <dune/common/parametertreeparser.hh>
33
36
39
42
43#include "spatialparams.hh"
44
45namespace Dumux {
46
47template <class TypeTag>
48class TwoPTwoCChemicalNonequilibriumProblem;
49
50namespace Properties {
51// Create new type tags
52namespace TTag {
53struct TwoPTwoCChemicalNonequilibrium { using InheritsFrom = std::tuple<TwoPTwoCNINonEquil>; };
54struct TwoPTwoCChemicalNonequilibriumBox { using InheritsFrom = std::tuple<TwoPTwoCChemicalNonequilibrium, BoxModel>; };
55struct TwoPTwoCChemicalNonequilibriumCC { using InheritsFrom = std::tuple<TwoPTwoCChemicalNonequilibrium, CCTpfaModel>; };
56} // end namespace TTag
57
58// Set the grid type
59template<class TypeTag>
60struct Grid<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium> { using type = Dune::YaspGrid<2>; };
61
62// Set the problem property
63template<class TypeTag>
64struct Problem<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium>
66
67// Set fluid configuration
68template<class TypeTag>
69struct FluidSystem<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium>
70{
71private:
73public:
74
75 using type = FluidSystems::H2OAir<Scalar,
77 FluidSystems::H2OAirDefaultPolicy</*fastButSimplifiedRelations=*/true>,
78 true /*useKelvinEquation*/>;
79};
80
81// Set the spatial parameters
82template<class TypeTag>
83struct SpatialParams<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium>
84{
88};
89
90// decide which type to use for floating values (double / quad)
91template<class TypeTag>
92struct Scalar<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium> { using type = double; };
93template<class TypeTag>
94struct Formulation<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium>
95{
96public:
98};
99
100template<class TypeTag>
101struct UseMoles<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium>
102{ static constexpr bool value = true; };
103
104template<class TypeTag>
105struct EnableThermalNonEquilibrium<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium>
106{ static constexpr bool value = false; };
107
108template<class TypeTag>
109struct HeatConductionType<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium>
111
112template<class TypeTag>
113struct EnergyLocalResidual<TypeTag, TTag::TwoPTwoCChemicalNonequilibrium>
115
116} // end namespace Properties
117
124template <class TypeTag>
126{
134 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
135 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
136 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
137 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
139 using Element = typename GridView::template Codim<0>::Entity;
140 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
144
146 using Indices = typename ModelTraits::Indices;
147
148 static constexpr int dim = GridView::dimension;
149 static constexpr int dimWorld = GridView::dimensionworld;
150 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethod::box;
151 enum { dofCodim = isBox ? dim : 0 };
152public:
153 TwoPTwoCChemicalNonequilibriumProblem(std::shared_ptr<const GridGeometry> gridGeometry)
155 {
156 temperature_ = 273.15 + 25; // -> 25°C
157
158 // initialize the tables of the fluid system
159 Scalar Tmin = temperature_ - 1.0;
160 Scalar Tmax = temperature_ + 1.0;
161 int nT = 3;
162
163 Scalar pmin = 1.0e5 * 0.75;
164 Scalar pmax = 2.0e5 * 1.25;
165 int np = 1000;
166
167 FluidSystem::init(Tmin, Tmax, nT, pmin, pmax, np);
168 name_ = getParam<std::string>("Problem.Name");
169 }
170
171 void setGridVariables(std::shared_ptr<GridVariables> gridVariables)
172 { gridVariables_ = gridVariables; }
173
174 const GridVariables& gridVariables() const
175 { return *gridVariables_; }
176
180 // \{
181
187 const std::string name() const
188 { return name_; }
189
194 Scalar temperature() const
195 { return temperature_; }
196
204 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
205 {
206 BoundaryTypes bcTypes;
207 if (globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_)
208 bcTypes.setAllDirichlet();
209 else
210 bcTypes.setAllNeumann();
211 return bcTypes;
212 }
213
219 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
220 {
221 PrimaryVariables values;
222 values = initial_(globalPos);
223 return values;
224 }
225
237 NeumannFluxes neumann(const Element& element,
238 const FVElementGeometry& fvGeometry,
239 const ElementVolumeVariables& elemVolVars,
240 const ElementFluxVariablesCache& elemFluxVarsCache,
241 const SubControlVolumeFace& scvf) const
242 {
243 NeumannFluxes values(0.0);
244 const auto& globalPos = scvf.ipGlobal();
245 const auto& volVars = elemVolVars[scvf.insideScvIdx()];
246 Scalar boundaryLayerThickness = 0.0016;
247 //right side
248 if (globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_)
249 {
250 Scalar moleFracH2OInside = volVars.moleFraction(FluidSystem::gasPhaseIdx, FluidSystem::H2OIdx);
251 Scalar moleFracRefH2O = 0.0;
252 Scalar evaporationRate = volVars.diffusionCoefficient(FluidSystem::gasPhaseIdx, FluidSystem::H2OIdx)
253 * (moleFracH2OInside - moleFracRefH2O)
254 / boundaryLayerThickness
255 * volVars.molarDensity(FluidSystem::gasPhaseIdx);
256 values[Indices::conti0EqIdx+2] = evaporationRate;
257
258 values[Indices::energyEqIdx] = FluidSystem::enthalpy(volVars.fluidState(), FluidSystem::gasPhaseIdx) * evaporationRate;
259 values[Indices::energyEqIdx] += FluidSystem::thermalConductivity(volVars.fluidState(), FluidSystem::gasPhaseIdx)
260 * (volVars.temperature() - temperature_)/boundaryLayerThickness;
261 }
262 return values;
263 }
264
265 // \}
266
275 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
276 {
277 return initial_(globalPos);
278 }
279
285 template<class VTKWriter>
286 void addVtkFields(VTKWriter& vtk)
287 {
288 vtk.addField(xEquilxwn_, "xEquil^Air_liq");
289 vtk.addField(xEquilxnw_, "xEquil^H2O_gas");
290 }
291
292 void updateVtkFields(const SolutionVector& curSol)
293 {
294 const auto& gridView = this->gridGeometry().gridView();
295 xEquilxwn_.resize(gridView.size(dofCodim));
296 xEquilxnw_.resize(gridView.size(dofCodim));
297 for (const auto& element : elements(this->gridGeometry().gridView()))
298 {
299 auto elemSol = elementSolution(element, curSol, this->gridGeometry());
300
301 auto fvGeometry = localView(this->gridGeometry());
302 fvGeometry.bindElement(element);
303
304 for (auto&& scv : scvs(fvGeometry))
305 {
306 VolumeVariables volVars;
307 volVars.update(elemSol, *this, element, scv);
308 const auto dofIdxGlobal = scv.dofIndex();
309 xEquilxwn_[dofIdxGlobal] = volVars.xEquil(0,1);
310 xEquilxnw_[dofIdxGlobal] = volVars.xEquil(1,0);
311 }
312 }
313 }
314
315 // \}
316
317private:
318 // the internal method for the initial condition
319 PrimaryVariables initial_(const GlobalPosition &globalPos) const
320 {
321 PrimaryVariables values(0.0);
322 values[Indices::pressureIdx] = 1e5; // water pressure
323 values[Indices::switchIdx] = 0.8; // gas saturation
324 values[2] = 5e-4; // xwn higher than equil, equil is 3.4e-5
325 values[3] = 1e-2; // xnw lower than 1.3e-2
326 values[Indices::temperatureIdx] = temperature_;
327 values.setState(Indices::bothPhases);
328
329 return values;
330 }
331
332 Scalar temperature_;
333 static constexpr Scalar eps_ = 1e-6;
334 std::string name_;
335
336 std::shared_ptr<GridVariables> gridVariables_;
337 std::vector<double> xEquilxnw_;
338 std::vector<double> xEquilxwn_;
339};
340} // end namespace Dumux
341
342#endif
Defines a type tag and some properties for models using the box scheme.
Properties for all models using cell-centered finite volume scheme with TPFA.
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
TwoPFormulation
Enumerates the formulations which the two-phase model accepts.
Definition: formulation.hh:35
@ p0s1
first phase pressure and second phase saturation as primary variables
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
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
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:50
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:588
Property to specify the type of scalar values.
Definition: common/properties.hh:53
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 local residual of the energy equation.
Definition: common/properties.hh:206
The type for the calculation of the heat conduction fluxes.
Definition: common/properties.hh:216
The type of the spatial parameters object.
Definition: common/properties.hh:221
The type of the fluid system to use.
Definition: common/properties.hh:223
The formulation of the model.
Definition: common/properties.hh:237
Definition: common/properties.hh:328
forward declaration of the method-specific implementation
Definition: fourierslaw.hh:37
Tabulates all thermodynamic properties of a given untabulated chemical species.
Definition: tabulatedcomponent.hh:82
Policy for the H2O-air fluid system.
Definition: h2oair.hh:52
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
Definition: h2oair.hh:75
Definition: porousmediumflow/nonisothermal/localresidual.hh:35
Base class for all fully implicit porous media problems.
Definition: dumux/porousmediumflow/problem.hh:39
Problem where air is injected in a unsaturated porous medium.
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:126
void setGridVariables(std::shared_ptr< GridVariables > gridVariables)
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:171
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:275
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:204
void updateVtkFields(const SolutionVector &curSol)
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:292
const GridVariables & gridVariables() const
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:174
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet oundary segment.
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:219
TwoPTwoCChemicalNonequilibriumProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:153
void addVtkFields(VTKWriter &vtk)
Adds additional VTK output data to the VTKWriter.
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:286
Scalar temperature() const
Returns the temperature .
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:194
const std::string name() const
Returns the problem name.
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:187
NeumannFluxes neumann(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluates the boundary conditions for a Neumann boundary segment.
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:237
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:53
std::tuple< TwoPTwoCNINonEquil > InheritsFrom
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:53
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:54
std::tuple< TwoPTwoCChemicalNonequilibrium, BoxModel > InheritsFrom
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:54
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:55
std::tuple< TwoPTwoCChemicalNonequilibrium, CCTpfaModel > InheritsFrom
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:55
Dune::YaspGrid< 2 > type
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:60
GetPropType< TypeTag, Properties::GridGeometry > GridGeometry
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:85
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:86
double type
Definition: test/porousmediumflow/2p2c/implicit/chemicalnonequilibrium/problem.hh:92
The spatial parameters for the 2p2c chemical nonequilibrium problem.
Definition: porousmediumflow/2p2c/implicit/chemicalnonequilibrium/spatialparams.hh:51
Adaption of the fully implicit scheme to the two-phase two-component fully implicit model.
Base class for all porous media problems.
Definition of the spatial parameters for the MaxwellStefan problem.