26#ifndef DUMUX_1P2C_TEST_PROBLEM_HH
27#define DUMUX_1P2C_TEST_PROBLEM_HH
30#include <dune/grid/uggrid.hh>
32#include <dune/grid/yaspgrid.hh>
46#include "../spatialparams.hh"
50template <
class TypeTag>
51class OnePTwoCTestProblem;
65template<
class TypeTag>
66struct Grid<TypeTag, TTag::OnePTwoCTest> {
using type = Dune::UGGrid<2>; };
68template<
class TypeTag>
69struct Grid<TypeTag, TTag::OnePTwoCTest> {
using type = Dune::YaspGrid<2>; };
73template<
class TypeTag>
77template<
class TypeTag>
86template<
class TypeTag>
95template<
class TypeTag>
96struct UseMoles<TypeTag, TTag::OnePTwoCTest> {
static constexpr bool value =
true; };
124template <
class TypeTag>
135 using Element =
typename GridView::template Codim<0>::Entity;
139 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
142 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
143 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
151 pressureIdx = Indices::pressureIdx,
154 H2OIdx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::H2OIdx),
155 N2Idx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::N2Idx),
158 contiH2OEqIdx = Indices::conti0EqIdx + H2OIdx,
159 contiN2EqIdx = Indices::conti0EqIdx + N2Idx
163 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
166 static const int dimWorld = GridView::dimensionworld;
167 using GlobalPosition =
typename SubControlVolumeFace::GlobalPosition;
178 std::cout<<
"problem uses mole fractions"<<std::endl;
180 std::cout<<
"problem uses mass fractions"<<std::endl;
194 {
return 273.15 + 20; }
211 BoundaryTypes values;
213 if(globalPos[0] < eps_)
214 values.setAllDirichlet();
216 values.setAllNeumann();
228 PrimaryVariables values = initial_(globalPos);
231 if (globalPos[0] < eps_ )
233 values[N2Idx] = 2.0e-5;
257 template<
bool useBox = isBox, std::enable_if_t<useBox,
int> = 0>
259 const FVElementGeometry& fvGeometry,
260 const ElementVolumeVariables& elemVolVars,
261 const ElementFluxVariablesCache& elemFluxVarsCache,
262 const SubControlVolumeFace& scvf)
const
265 NumEqVector values(0.0);
267 const auto& ipGlobal = scvf.ipGlobal();
268 if (ipGlobal[0] < xMax - eps_)
272 static const Scalar dirichletPressure = 1e5;
273 const auto& volVars = elemVolVars[scvf.insideScvIdx()];
274 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
277 if (useNitscheTypeBc_)
279 values[contiH2OEqIdx] = (volVars.pressure() - dirichletPressure)*1e7;
280 values[contiN2EqIdx] = values[contiH2OEqIdx] * (useMoles ? volVars.moleFraction(0, N2Idx)
281 : volVars.massFraction(0, N2Idx));
286 GlobalPosition gradP(0.0);
287 for (
const auto& scv : scvs(fvGeometry))
289 const auto xIp = scv.dofPosition()[0];
290 auto tmp = fluxVarsCache.gradN(scv.localDofIndex());
291 tmp *= xIp > xMax - eps_ ? dirichletPressure
292 : elemVolVars[scv].pressure();
297 auto phaseFlux =
vtmv(scvf.unitOuterNormal(), volVars.permeability(), gradP);
298 phaseFlux *= useMoles ? volVars.molarDensity() : volVars.density();
299 phaseFlux *= volVars.mobility();
302 values[contiH2OEqIdx] = phaseFlux;
304 values[contiN2EqIdx] = phaseFlux * ( useMoles ? volVars.moleFraction(0, N2Idx) : volVars.massFraction(0, N2Idx) );
313 template<
bool useBox = isBox, std::enable_if_t<!useBox,
int> = 0>
315 const FVElementGeometry& fvGeometry,
316 const ElementVolumeVariables& elemVolVars,
317 const ElementFluxVariablesCache& elemFluxVarsCache,
318 const SubControlVolumeFace& scvf)
const
321 NumEqVector values(0.0);
323 const auto& ipGlobal = scvf.ipGlobal();
324 if (ipGlobal[0] < xMax - eps_)
328 static const Scalar dirichletPressure = 1e5;
329 const auto& volVars = elemVolVars[scvf.insideScvIdx()];
331 auto d = ipGlobal - element.geometry().center();
334 auto upwindTerm = useMoles ? volVars.molarDensity() : volVars.density();
335 upwindTerm *= volVars.mobility();
337 const auto tij =
vtmv(scvf.unitOuterNormal(), volVars.permeability(), d);
338 const auto phaseFlux = -1.0*upwindTerm*tij*(dirichletPressure - volVars.pressure());
341 values[contiH2OEqIdx] = phaseFlux;
343 values[contiN2EqIdx] = phaseFlux * (useMoles ? volVars.moleFraction(0, N2Idx) : volVars.massFraction(0, N2Idx));
367 {
return NumEqVector(0.0); }
378 {
return initial_(globalPos); }
384 PrimaryVariables initial_(
const GlobalPosition &globalPos)
const
386 PrimaryVariables priVars;
387 priVars[pressureIdx] = 2e5 - 1e5*globalPos[0];
388 priVars[N2Idx] = 0.0;
391 static constexpr Scalar eps_ = 1e-6;
392 bool useNitscheTypeBc_;
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.