3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/porousmediumflow/1pnc/implicit/1p2c/isothermal/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 *****************************************************************************/
26#ifndef DUMUX_1P2C_TEST_PROBLEM_HH
27#define DUMUX_1P2C_TEST_PROBLEM_HH
28
29#if HAVE_UG
30#include <dune/grid/uggrid.hh>
31#endif
32#include <dune/grid/yaspgrid.hh>
33
34#include <dumux/common/math.hh>
42
45
46#include "../spatialparams.hh"
47
48namespace Dumux {
49
50template <class TypeTag>
51class OnePTwoCTestProblem;
52
53namespace Properties {
54
55// Create new type tags
56namespace TTag {
57struct OnePTwoCTest { using InheritsFrom = std::tuple<OnePNC>; };
58struct OnePTwoCTestBox { using InheritsFrom = std::tuple<OnePTwoCTest, BoxModel>; };
59struct OnePTwoCTestCCTpfa { using InheritsFrom = std::tuple<OnePTwoCTest, CCTpfaModel>; };
60struct OnePTwoCTestCCMpfa { using InheritsFrom = std::tuple<OnePTwoCTest, CCMpfaModel>; };
61} // end namespace TTag
62
63// Set the grid type
64#if HAVE_UG
65template<class TypeTag>
66struct Grid<TypeTag, TTag::OnePTwoCTest> { using type = Dune::UGGrid<2>; };
67#else
68template<class TypeTag>
69struct Grid<TypeTag, TTag::OnePTwoCTest> { using type = Dune::YaspGrid<2>; };
70#endif
71
72// Set the problem property
73template<class TypeTag>
74struct Problem<TypeTag, TTag::OnePTwoCTest> { using type = OnePTwoCTestProblem<TypeTag>; };
75
76// Set fluid configuration
77template<class TypeTag>
78struct FluidSystem<TypeTag, TTag::OnePTwoCTest>
79{
83};
84
85// Set the spatial parameters
86template<class TypeTag>
87struct SpatialParams<TypeTag, TTag::OnePTwoCTest>
88{
92};
93
94// Define whether mole(true) or mass (false) fractions are used
95template<class TypeTag>
96struct UseMoles<TypeTag, TTag::OnePTwoCTest> { static constexpr bool value = true; };
97}
98
124template <class TypeTag>
126{
128
135 using Element = typename GridView::template Codim<0>::Entity;
136
138 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
139 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
140
142 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
143 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
144
145 // copy some indices for convenience
147
148 enum
149 {
150 // indices of the primary variables
151 pressureIdx = Indices::pressureIdx,
152
153 // component indices
154 H2OIdx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::H2OIdx),
155 N2Idx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::N2Idx),
156
157 // indices of the equations
158 contiH2OEqIdx = Indices::conti0EqIdx + H2OIdx,
159 contiN2EqIdx = Indices::conti0EqIdx + N2Idx
160 };
161
163 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
165
166 static const int dimWorld = GridView::dimensionworld;
167 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
168
169public:
170 OnePTwoCTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
171 : ParentType(gridGeometry), useNitscheTypeBc_(getParam<bool>("Problem.UseNitscheTypeBc", false))
172 {
173 //initialize fluid system
174 FluidSystem::init();
175
176 // stating in the console whether mole or mass fractions are used
177 if(useMoles)
178 std::cout<<"problem uses mole fractions"<<std::endl;
179 else
180 std::cout<<"problem uses mass fractions"<<std::endl;
181 }
182
186 // \{
187
193 Scalar temperature() const
194 { return 273.15 + 20; } // in [K]
195
196 // \}
197
201 // \{
202
209 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
210 {
211 BoundaryTypes values;
212
213 if(globalPos[0] < eps_)
214 values.setAllDirichlet();
215 else
216 values.setAllNeumann();
217
218 return values;
219 }
220
226 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
227 {
228 PrimaryVariables values = initial_(globalPos);
229
230 // condition for the N2 molefraction at left boundary
231 if (globalPos[0] < eps_ )
232 {
233 values[N2Idx] = 2.0e-5;
234 }
235
236 return values;
237 }
238
257 template<bool useBox = isBox, std::enable_if_t<useBox, int> = 0>
258 NumEqVector neumann(const Element& element,
259 const FVElementGeometry& fvGeometry,
260 const ElementVolumeVariables& elemVolVars,
261 const ElementFluxVariablesCache& elemFluxVarsCache,
262 const SubControlVolumeFace& scvf) const
263 {
264 // no-flow everywhere except at the right boundary
265 NumEqVector values(0.0);
266 const auto xMax = this->gridGeometry().bBoxMax()[0];
267 const auto& ipGlobal = scvf.ipGlobal();
268 if (ipGlobal[0] < xMax - eps_)
269 return values;
270
271 // set a fixed pressure on the right side of the domain
272 static const Scalar dirichletPressure = 1e5;
273 const auto& volVars = elemVolVars[scvf.insideScvIdx()];
274 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
275
276 // if specified in the input file, use a Nitsche type boundary condition for the box model.
277 if (useNitscheTypeBc_)
278 {
279 values[contiH2OEqIdx] = (volVars.pressure() - dirichletPressure)*1e7;
280 values[contiN2EqIdx] = values[contiH2OEqIdx] * (useMoles ? volVars.moleFraction(0, N2Idx)
281 : volVars.massFraction(0, N2Idx));
282 return values;
283 }
284
285 // evaluate the pressure gradient
286 GlobalPosition gradP(0.0);
287 for (const auto& scv : scvs(fvGeometry))
288 {
289 const auto xIp = scv.dofPosition()[0];
290 auto tmp = fluxVarsCache.gradN(scv.localDofIndex());
291 tmp *= xIp > xMax - eps_ ? dirichletPressure
292 : elemVolVars[scv].pressure();
293 gradP += tmp;
294 }
295
296 // compute flux
297 auto phaseFlux = vtmv(scvf.unitOuterNormal(), volVars.permeability(), gradP);
298 phaseFlux *= useMoles ? volVars.molarDensity() : volVars.density();
299 phaseFlux *= volVars.mobility();
300
301 // set Neumann bc values
302 values[contiH2OEqIdx] = phaseFlux;
303 // emulate an outflow condition for the component transport on the right side
304 values[contiN2EqIdx] = phaseFlux * ( useMoles ? volVars.moleFraction(0, N2Idx) : volVars.massFraction(0, N2Idx) );
305
306 return values;
307 }
308
313 template<bool useBox = isBox, std::enable_if_t<!useBox, int> = 0>
314 NumEqVector neumann(const Element& element,
315 const FVElementGeometry& fvGeometry,
316 const ElementVolumeVariables& elemVolVars,
317 const ElementFluxVariablesCache& elemFluxVarsCache,
318 const SubControlVolumeFace& scvf) const
319 {
320 // no-flow everywhere except at the right boundary
321 NumEqVector values(0.0);
322 const auto xMax = this->gridGeometry().bBoxMax()[0];
323 const auto& ipGlobal = scvf.ipGlobal();
324 if (ipGlobal[0] < xMax - eps_)
325 return values;
326
327 // set a fixed pressure on the right side of the domain
328 static const Scalar dirichletPressure = 1e5;
329 const auto& volVars = elemVolVars[scvf.insideScvIdx()];
330
331 auto d = ipGlobal - element.geometry().center();
332 d /= d.two_norm2();
333
334 auto upwindTerm = useMoles ? volVars.molarDensity() : volVars.density();
335 upwindTerm *= volVars.mobility();
336
337 const auto tij = vtmv(scvf.unitOuterNormal(), volVars.permeability(), d);
338 const auto phaseFlux = -1.0*upwindTerm*tij*(dirichletPressure - volVars.pressure());
339
340 // set Neumann bc values
341 values[contiH2OEqIdx] = phaseFlux;
342 // emulate an outflow condition for the component transport on the right side
343 values[contiN2EqIdx] = phaseFlux * (useMoles ? volVars.moleFraction(0, N2Idx) : volVars.massFraction(0, N2Idx));
344
345 return values;
346 }
347
348 // \}
349
353 // \{
354
366 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
367 { return NumEqVector(0.0); }
368
377 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
378 { return initial_(globalPos); }
379
380 // \}
381
382private:
383 // the internal method for the initial condition
384 PrimaryVariables initial_(const GlobalPosition &globalPos) const
385 {
386 PrimaryVariables priVars;
387 priVars[pressureIdx] = 2e5 - 1e5*globalPos[0]; // initial condition for the pressure
388 priVars[N2Idx] = 0.0; // initial condition for the N2 molefraction
389 return priVars;
390 }
391 static constexpr Scalar eps_ = 1e-6;
392 bool useNitscheTypeBc_;
393 };
394
395} // end namespace Dumux
396
397#endif
Defines a type tag and some properties for models using the box scheme.
free functions for the evaluation of primary variables inside elements.
Properties for all models using cell-centered finite volume scheme with TPFA.
Properties for all models using cell-centered finite volume scheme with mpfa.
free functions for the evaluation of primary variable gradients inside elements.
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
Define some often used mathematical functions.
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:840
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition: parameters.hh:428
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
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
The DUNE grid type.
Definition: common/properties.hh:57
UndefinedProperty 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 type of the spatial parameters object.
Definition: common/properties.hh:221
The type of the fluid system to use.
Definition: common/properties.hh:223
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
Definition: 1padapter.hh:46
Policy for the H2O-N2 fluid system.
Definition: h2on2.hh:52
A two-phase fluid system with two components water Nitrogen for non-equilibrium models.
Definition: h2on2.hh:69
Base class for all fully implicit porous media problems.
Definition: dumux/porousmediumflow/problem.hh:39
Definition of a problem for the 1pnc problem: Component transport of nitrogen dissolved in the water ...
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:126
NumEqVector 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 (implementation for the box scheme).
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:258
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:377
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet boundary segment.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:226
Scalar temperature() const
Returns the temperature within the domain [K].
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:193
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:209
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Evaluates the source term for all phases within a given sub-control volume.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:366
OnePTwoCTestProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:170
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:57
std::tuple< OnePNC > InheritsFrom
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:57
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:58
std::tuple< OnePTwoCTest, BoxModel > InheritsFrom
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:58
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:59
std::tuple< OnePTwoCTest, CCTpfaModel > InheritsFrom
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:59
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:60
std::tuple< OnePTwoCTest, CCMpfaModel > InheritsFrom
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:60
Dune::YaspGrid< 2 > type
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:69
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:80
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:90
GetPropType< TypeTag, Properties::GridGeometry > GridGeometry
Definition: test/porousmediumflow/1pnc/implicit/1p2c/isothermal/problem.hh:89
Definition of the spatial parameters for the 1pnc test problems.
Definition: porousmediumflow/1pnc/implicit/1p2c/spatialparams.hh:41
Base class for all porous media problems.
Adaption of the fully implicit model to the one-phase n-component flow model.