25#ifndef DUMUX_STOKES_SUBPROBLEM_HH
26#define DUMUX_STOKES_SUBPROBLEM_HH
28#include <dune/grid/yaspgrid.hh>
38template <
class TypeTag>
48template<
class TypeTag>
57template<
class TypeTag>
58struct Grid<TypeTag,
TTag::StokesOnePTwoC> {
using type = Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
61template<
class TypeTag>
64template<
class TypeTag>
66template<
class TypeTag>
67struct EnableGridFluxVariablesCache<TypeTag, TTag::StokesOnePTwoC> {
static constexpr bool value =
true; };
68template<
class TypeTag>
69struct EnableGridVolumeVariablesCache<TypeTag, TTag::StokesOnePTwoC> {
static constexpr bool value =
true; };
72template<
class TypeTag>
73struct UseMoles<TypeTag,
TTag::StokesOnePTwoC> {
static constexpr bool value =
true; };
76template<
class TypeTag>
86template <
class TypeTag>
92 using Indices =
typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
95 using FVElementGeometry =
typename GridGeometry::LocalView;
96 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
100 using Element =
typename GridView::template Codim<0>::Entity;
101 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
104 using TimeLoopPtr = std::shared_ptr<TimeLoop<Scalar>>;
108 : ParentType(gridGeometry,
"Stokes"), eps_(1e-6), couplingManager_(
couplingManager)
113 verticalFlow_ = problemName_.find(
"vertical") != std::string::npos;
119 const std::string&
name()
const
135 {
return 273.15 + 10; }
143 {
return NumEqVector(0.0); }
159 const SubControlVolumeFace& scvf)
const
161 BoundaryTypes values;
163 const auto& globalPos = scvf.dofPosition();
168 if(onUpperBoundary_(globalPos))
170 values.setDirichlet(Indices::velocityXIdx);
171 values.setDirichlet(Indices::velocityYIdx);
172 values.setDirichlet(Indices::conti0EqIdx + 1);
176 if (onRightBoundary_(globalPos) || (onLeftBoundary_(globalPos)))
178 values.setDirichlet(Indices::velocityXIdx);
179 values.setDirichlet(Indices::velocityYIdx);
180 values.setNeumann(Indices::conti0EqIdx);
181 values.setNeumann(Indices::conti0EqIdx + 1);
187 if (onLeftBoundary_(globalPos))
189 values.setDirichlet(Indices::conti0EqIdx + 1);
190 values.setDirichlet(Indices::velocityXIdx);
191 values.setDirichlet(Indices::velocityYIdx);
193 else if (onRightBoundary_(globalPos))
195 values.setDirichlet(Indices::pressureIdx);
196 values.setOutflow(Indices::conti0EqIdx + 1);
200 values.setDirichlet(Indices::velocityXIdx);
201 values.setDirichlet(Indices::velocityYIdx);
202 values.setNeumann(Indices::conti0EqIdx);
203 values.setNeumann(Indices::conti0EqIdx + 1);
207 if(
couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
209 values.setCouplingNeumann(Indices::conti0EqIdx);
210 values.setCouplingNeumann(Indices::conti0EqIdx + 1);
211 values.setCouplingNeumann(Indices::momentumYBalanceIdx);
212 values.setBJS(Indices::momentumXBalanceIdx);
223 PrimaryVariables values(0.0);
229 static const bool isDiffusionProblem = problemName_.find(
"diffusion") != std::string::npos;
231 Scalar topMoleFraction = 1e-3;
233 if (isDiffusionProblem)
238 topMoleFraction = 0.0;
244 values[Indices::velocityYIdx] *= -1.0;
247 if(globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_)
248 values[Indices::conti0EqIdx + 1] = topMoleFraction;
253 if(globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_)
254 values[Indices::conti0EqIdx + 1] = inletMoleFraction;
270 template<
class ElementVolumeVariables,
class ElementFaceVariables>
272 const FVElementGeometry& fvGeometry,
273 const ElementVolumeVariables& elemVolVars,
274 const ElementFaceVariables& elemFaceVars,
275 const SubControlVolumeFace& scvf)
const
277 NumEqVector values(0.0);
279 if(
couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
281 values[Indices::momentumYBalanceIdx] =
couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
283 const auto tmp =
couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
284 values[Indices::conti0EqIdx] = tmp[0];
285 values[Indices::conti0EqIdx + 1] = tmp[1];
294 {
return *couplingManager_; }
308 PrimaryVariables values(0.0);
309 values[Indices::pressureIdx] = 1e5;
313 auto parabolicProfile = [&](
const GlobalPosition& globalPos,
int coord)
315 return vMax * (globalPos[coord] - this->gridGeometry().bBoxMin()[coord])
316 * (this->gridGeometry().bBoxMax()[coord] - globalPos[coord])
317 / (0.25 * (this->gridGeometry().bBoxMax()[coord] - this->gridGeometry().bBoxMin()[coord])
318 * (this->gridGeometry().bBoxMax()[coord] - this->gridGeometry().bBoxMin()[coord]));
322 values[Indices::velocityYIdx] = parabolicProfile(globalPos, 0);
324 values[Indices::velocityXIdx] = parabolicProfile(globalPos, 1);
333 Scalar
permeability(
const Element& element,
const SubControlVolumeFace& scvf)
const
335 return couplingManager().couplingData().darcyPermeability(element, scvf);
342 Scalar
alphaBJ(
const SubControlVolumeFace& scvf)
const
351 { timeLoop_ = timeLoop; }
357 {
return timeLoop_->time(); }
362 bool onLeftBoundary_(
const GlobalPosition &globalPos)
const
363 {
return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
365 bool onRightBoundary_(
const GlobalPosition &globalPos)
const
366 {
return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
368 bool onLowerBoundary_(
const GlobalPosition &globalPos)
const
369 {
return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_; }
371 bool onUpperBoundary_(
const GlobalPosition &globalPos)
const
372 {
return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
376 std::string problemName_;
377 std::shared_ptr<CouplingManager> couplingManager_;
378 TimeLoopPtr timeLoop_;
A single-phase, isothermal Navier-Stokes model.
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...
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:438
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
Definition common/properties.hh:47
Type tag for numeric models.
Definition grid.hh:35
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
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
NavierStokesProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The constructor.
Definition dumux/freeflow/navierstokes/problem.hh:99
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
static constexpr int liquidPhaseIdx
Definition h2oair.hh:87
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/problem_stokes.hh:134
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/problem_stokes.hh:271
void setTimeLoop(TimeLoopPtr timeLoop)
Sets the time loop pointer.
Definition 1p2c_1p2c/problem_stokes.hh:350
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition 1p2c_1p2c/problem_stokes.hh:342
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/problem_stokes.hh:158
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/problem_stokes.hh:107
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition 1p2c_1p2c/problem_stokes.hh:142
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition 1p2c_1p2c/problem_stokes.hh:221
Scalar time() const
Returns the time.
Definition 1p2c_1p2c/problem_stokes.hh:356
const std::string & name() const
The problem name.
Definition 1p2c_1p2c/problem_stokes.hh:119
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/problem_stokes.hh:333
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
static constexpr auto phaseIdx
Definition 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:65
FluidSystems::H2OAir< GetPropType< TypeTag, Properties::Scalar > > H2OAir
Definition 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:64
FluidSystems::OnePAdapter< H2OAir, phaseIdx > type
Definition 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:66
Dune::YaspGrid< 2, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:71
Dumux::StokesSubProblem< TypeTag > type
Definition 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:75
Defines a type tag and some properties for ree-flow models using the staggered scheme....
Navier-Stokes problem base class.