30#ifndef DUMUX_STOKES_SUBPROBLEM_DIFFUSION_COMPARISON_HH
31#define DUMUX_STOKES_SUBPROBLEM_DIFFUSION_COMPARISON_HH
33#include <dune/grid/yaspgrid.hh>
47#define DIFFUSIONTYPE FicksLaw<TypeTag>
51template <
class TypeTag>
52class StokesSubProblem;
61template<
class TypeTag>
65 static constexpr auto phaseIdx = H2OAir::liquidPhaseIdx;
70template<
class TypeTag>
71struct Grid<TypeTag, TTag::StokesOnePTwoC> {
using type = Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
74template<
class TypeTag>
77template<
class TypeTag>
79template<
class TypeTag>
81template<
class TypeTag>
85template<
class TypeTag>
86struct UseMoles<TypeTag, TTag::StokesOnePTwoC> {
static constexpr bool value =
true; };
89template<
class TypeTag>
93template<
class TypeTag>
94struct ReplaceCompEqIdx<TypeTag, TTag::StokesOnePTwoC> {
static constexpr int value = 3; };
103template <
class TypeTag>
113 using FVElementGeometry =
typename GridGeometry::LocalView;
114 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
118 using Element =
typename GridView::template Codim<0>::Entity;
119 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
126 : ParentType(gridGeometry,
"Stokes"), eps_(1e-6), injectionState_(false), couplingManager_(
couplingManager)
128 inletVelocity_ = getParamFromGroup<Scalar>(this->paramGroup(),
"Problem.Velocity");
129 pressure_ = getParamFromGroup<Scalar>(this->paramGroup(),
"Problem.Pressure");
130 inletMoleFraction_ = getParamFromGroup<Scalar>(this->paramGroup(),
"Problem.InletMoleFraction");
147 {
return 273.15 + 10; }
155 {
return NumEqVector(0.0); }
171 const SubControlVolumeFace& scvf)
const
173 BoundaryTypes values;
175 const auto& globalPos = scvf.dofPosition();
177 if(onLeftBoundary_(globalPos))
179 values.setDirichlet(Indices::conti0EqIdx + 1);
180 values.setDirichlet(Indices::velocityXIdx);
181 values.setDirichlet(Indices::velocityYIdx);
183 else if(onRightBoundary_(globalPos))
185 values.setDirichlet(Indices::pressureIdx);
186 values.setOutflow(Indices::conti0EqIdx + 1);
190 values.setDirichlet(Indices::velocityXIdx);
191 values.setDirichlet(Indices::velocityYIdx);
192 values.setNeumann(Indices::conti0EqIdx);
193 values.setNeumann(Indices::conti0EqIdx + 1);
196 if(
couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
198 values.setCouplingNeumann(Indices::conti0EqIdx);
199 values.setCouplingNeumann(Indices::conti0EqIdx + 1);
200 values.setCouplingNeumann(Indices::momentumYBalanceIdx);
201 values.setBJS(Indices::momentumXBalanceIdx);
212 PrimaryVariables values(0.0);
216 if(globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_ &&
isInjectionPeriod())
217 values[Indices::conti0EqIdx + 1] = inletMoleFraction_;
231 template<
class ElementVolumeVariables,
class ElementFaceVariables>
233 const FVElementGeometry& fvGeometry,
234 const ElementVolumeVariables& elemVolVars,
235 const ElementFaceVariables& elemFaceVars,
236 const SubControlVolumeFace& scvf)
const
238 NumEqVector values(0.0);
240 if(
couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
242 values[Indices::momentumYBalanceIdx] =
couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
244 const auto tmp =
couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf, DiffusionCoefficientAveragingType::harmonic);
245 values[Indices::conti0EqIdx] = tmp[0];
246 values[Indices::conti0EqIdx + 1] = tmp[1];
255 {
return *couplingManager_; }
269 PrimaryVariables values(0.0);
270 values[Indices::pressureIdx] = pressure_;
271 values[Indices::velocityXIdx] = inletVelocity_ * (globalPos[1] - this->gridGeometry().bBoxMin()[1])
272 * (this->gridGeometry().bBoxMax()[1] - globalPos[1])
273 / (0.25 * (this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1])
274 * (this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1]));
283 Scalar
permeability(
const Element& element,
const SubControlVolumeFace& scvf)
const
285 return couplingManager().couplingData().darcyPermeability(element, scvf);
292 Scalar
alphaBJ(
const SubControlVolumeFace& scvf)
const
299 injectionState_ = yesNo;
304 return injectionState_;
310 bool onLeftBoundary_(
const GlobalPosition &globalPos)
const
311 {
return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
313 bool onRightBoundary_(
const GlobalPosition &globalPos)
const
314 {
return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
316 bool onLowerBoundary_(
const GlobalPosition &globalPos)
const
317 {
return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_; }
319 bool onUpperBoundary_(
const GlobalPosition &globalPos)
const
320 {
return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
323 Scalar inletVelocity_;
325 Scalar inletMoleFraction_;
326 bool injectionState_;
328 std::shared_ptr<CouplingManager> couplingManager_;
This file contains the data which is required to calculate diffusive mass fluxes due to molecular dif...
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
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
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:38
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
Navier-Stokes problem base class.
Definition: dumux/freeflow/navierstokes/problem.hh:63
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
Definition: 1padapter.hh:46
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
Definition: h2oair.hh:75
DiffusionCoefficientAveragingType
Defines which kind of averanging of diffusion coefficiencients (moleculat diffusion or thermal conduc...
Definition: couplingdata.hh:52
Definition: multidomain/couplingmanager.hh:46
Test problem for the 1pnc (Navier-) Stokes problem.
Definition: 1p_2p/problem_stokes.hh:68
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:146
NumEqVector neumann(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Evaluates the boundary conditions for a Neumann control volume.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:232
bool shouldWriteRestartFile() const
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:138
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:292
BoundaryTypes boundaryTypes(const Element &element, const SubControlVolumeFace &scvf) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:170
const CouplingManager & couplingManager() const
Get the coupling manager.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:254
StokesSubProblem(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< CouplingManager > couplingManager)
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:125
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:154
void setInjectionState(const bool yesNo)
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:297
bool isInjectionPeriod() const
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:302
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:210
Scalar permeability(const Element &element, const SubControlVolumeFace &scvf) const
Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boun...
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:283
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:267
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:57
std::tuple< NavierStokesNC, StaggeredFreeFlowModel > InheritsFrom
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:57
Dune::YaspGrid< 2, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:71
DIFFUSIONTYPE type
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:90
Defines a type tag and some properties for ree-flow models using the staggered scheme....
#define DIFFUSIONTYPE
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:47