27#ifndef DUMUX_1P2CNI_CONVECTION_TEST_PROBLEM_HH
28#define DUMUX_1P2CNI_CONVECTION_TEST_PROBLEM_HH
31#include <dune/grid/uggrid.hh>
33#include <dune/grid/yaspgrid.hh>
45#include "../../spatialparams.hh"
49template <
class TypeTag>
50class OnePTwoCNIConvectionProblem;
63template<
class TypeTag>
64struct Grid<TypeTag, TTag::OnePTwoCNIConvection> {
using type = Dune::UGGrid<2>; };
66template<
class TypeTag>
67struct Grid<TypeTag, TTag::OnePTwoCNIConvection> {
using type = Dune::YaspGrid<2>; };
71template<
class TypeTag>
75template<
class TypeTag>
84template<
class TypeTag>
93template<
class TypeTag>
94struct UseMoles<TypeTag, TTag::OnePTwoCNIConvection> {
static constexpr bool value =
true; };
122template <
class TypeTag>
138 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
139 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
141 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
142 using Element =
typename GridView::template Codim<0>::Entity;
151 pressureIdx = Indices::pressureIdx,
152 temperatureIdx = Indices::temperatureIdx,
155 H2OIdx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::H2OIdx),
156 N2Idx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::N2Idx),
159 contiH2OEqIdx = Indices::conti0EqIdx + H2OIdx,
160 contiN2EqIdx = Indices::conti0EqIdx + N2Idx,
161 energyEqIdx = Indices::energyEqIdx
165 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;
182 temperatureExact_.resize(
gridGeometry->numDofs(), 290.0);
184 darcyVelocity_ = getParam<Scalar>(
"Problem.DarcyVelocity");
186 temperatureHigh_ = 291.;
187 temperatureLow_ = 290.;
195 return temperatureExact_;
201 const auto someElement = *(elements(this->
gridGeometry().gridView()).begin());
204 const auto someInitSol =
initialAtPos(someElement.geometry().center());
207 someFvGeometry.bindElement(someElement);
208 const auto someScv = *(scvs(someFvGeometry).begin());
210 VolumeVariables volVars;
211 volVars.update(someElemSol, *
this, someElement, someScv);
214 const auto densityW = volVars.density();
216 const auto storageW = densityW*heatCapacityW*
porosity;
217 const auto densityS = volVars.solidDensity();
218 const auto heatCapacityS = volVars.solidHeatCapacity();
219 const auto storageTotal = storageW + densityS*heatCapacityS*(1 -
porosity);
220 std::cout <<
"storage: " << storageTotal <<
'\n';
223 time = max(time, 1e-10);
224 const Scalar retardedFrontVelocity = darcyVelocity_*storageW/storageTotal/
porosity;
225 std::cout <<
"retarded velocity: " << retardedFrontVelocity <<
'\n';
227 for (
const auto& element : elements(this->
gridGeometry().gridView()))
230 fvGeometry.bindElement(element);
231 for (
auto&& scv : scvs(fvGeometry))
233 auto dofIdxGlobal = scv.dofIndex();
234 auto dofPosition = scv.dofPosition();
235 temperatureExact_[dofIdxGlobal] = (dofPosition[0] < retardedFrontVelocity*time) ? temperatureHigh_ : temperatureLow_;
251 {
return 273.15 + 20; }
268 BoundaryTypes values;
270 if(globalPos[0] > this->
gridGeometry().bBoxMax()[0] - eps_)
272 values.setAllDirichlet();
276 values.setAllNeumann();
289 PrimaryVariables values = initial_(globalPos);
311 const FVElementGeometry& fvGeometry,
312 const ElementVolumeVariables& elemVolVars,
313 const ElementFluxVariablesCache& elemFluxVarsCache,
314 const SubControlVolumeFace& scvf)
const
316 NumEqVector flux(0.0);
317 const auto& globalPos = scvf.ipGlobal();
318 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
320 if(globalPos[0] < eps_)
322 flux[contiH2OEqIdx] = -darcyVelocity_*elemVolVars[scv].molarDensity();
323 flux[contiN2EqIdx] = -darcyVelocity_*elemVolVars[scv].molarDensity()*elemVolVars[scv].moleFraction(0, N2Idx);
324 flux[energyEqIdx] = -darcyVelocity_
325 *elemVolVars[scv].density()
351 {
return NumEqVector(0.0); }
362 {
return initial_(globalPos); }
368 PrimaryVariables initial_(
const GlobalPosition &globalPos)
const
370 PrimaryVariables priVars;
371 priVars[pressureIdx] = pressureLow_;
372 priVars[N2Idx] = 1e-10;
373 priVars[temperatureIdx] = temperatureLow_;
376 static constexpr Scalar eps_ = 1e-6;
377 Scalar temperatureHigh_;
378 Scalar temperatureLow_;
379 Scalar pressureHigh_;
381 Scalar darcyVelocity_;
383 std::vector<Scalar> temperatureExact_;
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
Material properties of pure water .
Properties for all models using cell-centered finite volume scheme with mpfa.
Element solution classes and factory functions.
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.
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
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string porosity() noexcept
I/O name of porosity.
Definition: name.hh:139
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
Material properties of pure water .
Definition: h2o.hh:61
static const Scalar liquidHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of liquid water .
Definition: h2o.hh:281
static const Scalar liquidEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of liquid water .
Definition: h2o.hh:217
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
SpatialParams & spatialParams()
Returns the spatial parameters object.
Definition: dumux/porousmediumflow/problem.hh:146
Test for the OnePTwoCModel in combination with the NI model for a convection problem.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:124
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:361
OnePTwoCNIConvectionProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:170
void updateExactTemperature(const SolutionVector &curSol, Scalar time)
Udpate the analytical temperature.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:199
const std::vector< Scalar > & getExactTemperature()
Get the analytical temperature.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:193
Scalar temperature() const
Returns the temperature within the domain [K].
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:250
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/nonisothermal/convection/problem.hh:350
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.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:310
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet Sboundary segment.
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:287
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/nonisothermal/convection/problem.hh:266
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:55
std::tuple< OnePNCNI > InheritsFrom
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:55
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:56
std::tuple< OnePTwoCNIConvection, CCTpfaModel > InheritsFrom
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:56
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:57
std::tuple< OnePTwoCNIConvection, CCMpfaModel > InheritsFrom
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:57
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:58
std::tuple< OnePTwoCNIConvection, BoxModel > InheritsFrom
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:58
Dune::YaspGrid< 2 > type
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:67
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:78
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:88
GetPropType< TypeTag, Properties::GridGeometry > GridGeometry
Definition: test/porousmediumflow/1pnc/implicit/1p2c/nonisothermal/convection/problem.hh:87
Definition of the spatial parameters for the 1pnc test problems.
Definition: porousmediumflow/1pnc/implicit/1p2c/spatialparams.hh:41
Adaption of the fully implicit model to the one-phase n-component flow model.
Base class for all porous media problems.