version 3.8
freeflow/navierstokes/momentum/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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_PROBLEM_HH
13#define DUMUX_NAVIERSTOKES_MOMENTUM_PROBLEM_HH
14
15#include <dune/common/exceptions.hh>
16#include <dune/common/typetraits.hh>
22
23namespace Dumux {
24
25template<class TypeTag, class DiscretizationMethod>
27
28template<class TypeTag>
29class NavierStokesMomentumProblemImpl<TypeTag, DiscretizationMethods::FCStaggered>
30: public FVProblemWithSpatialParams<TypeTag>
31{
33 using Implementation = GetPropType<TypeTag, Properties::Problem>;
34
36 using GridView = typename GridGeometry::GridView;
37 using Element = typename GridView::template Codim<0>::Entity;
38
40 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
41 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
43
44 using FVElementGeometry = typename GridGeometry::LocalView;
45 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
46 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
48
49 enum {
50 dim = GridView::dimension,
51 dimWorld = GridView::dimensionworld
52 };
53
54 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
55 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
56 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
59
60 static constexpr bool isCoupled_ = !std::is_empty_v<CouplingManager>;
61
62
63public:
67 using InitialValues = Dune::FieldVector<Scalar, dimWorld>;
68 using Sources = Dune::FieldVector<Scalar, dimWorld>;
69 using DirichletValues = Dune::FieldVector<Scalar, dimWorld>;
70 using BoundaryFluxes = Dune::FieldVector<Scalar, dimWorld>;
71
72 using MomentumFluxType = Dune::FieldVector<Scalar, dimWorld>;
73
76
78 static constexpr bool isMomentumProblem() { return true; }
79
86 NavierStokesMomentumProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
87 std::shared_ptr<CouplingManager> couplingManager,
88 const std::string& paramGroup = "")
89 : ParentType(gridGeometry, paramGroup)
90 , gravity_(0.0)
91 , couplingManager_(couplingManager)
92 {
93 if (getParamFromGroup<bool>(paramGroup, "Problem.EnableGravity"))
94 gravity_[dim-1] = -9.81;
95
96 enableInertiaTerms_ = getParamFromGroup<bool>(paramGroup, "Problem.EnableInertiaTerms");
97 }
98
104 NavierStokesMomentumProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
105 const std::string& paramGroup = "")
106 : NavierStokesMomentumProblemImpl(gridGeometry, {}, paramGroup)
107 {}
108
127 Sources source(const Element& element,
128 const FVElementGeometry& fvGeometry,
129 const ElementVolumeVariables& elemVolVars,
130 const SubControlVolume& scv) const
131 {
132 // forward to solution independent, fully-implicit specific interface
133 return asImp_().sourceAtPos(scv.dofPosition());
134 }
135
149 Sources sourceAtPos(const GlobalPosition& globalPos) const
150 {
153 return Sources(0.0);
154 }
155
163 auto boundaryTypes(const Element& element,
164 const SubControlVolumeFace& scvf) const
165 {
166 // Forward it to the method which only takes the global coordinate.
167 // We evaluate the boundary type at the center of the sub control volume face
168 // in order to avoid ambiguities at domain corners.
169 return asImp_().boundaryTypesAtPos(scvf.center());
170 }
171
180 DirichletValues dirichlet(const Element& element, const SubControlVolumeFace& scvf) const
181 {
182 return asImp_().dirichletAtPos(scvf.ipGlobal());
183 }
184
194 template<class ElementFluxVariablesCache>
195 BoundaryFluxes neumann(const Element& element,
196 const FVElementGeometry& fvGeometry,
197 const ElementVolumeVariables& elemVolVars,
198 const ElementFluxVariablesCache& elemFluxVarsCache,
199 const SubControlVolumeFace& scvf) const
200 { return asImp_().neumannAtPos(scvf.ipGlobal()); }
201
205 BoundaryFluxes neumannAtPos(const GlobalPosition& globalPos) const
206 { return BoundaryFluxes(0.0); }
207
214 const GravityVector& gravity() const
215 { return gravity_; }
216
221 { return enableInertiaTerms_; }
222
227 Scalar pressure(const Element& element,
228 const FVElementGeometry& fvGeometry,
229 const SubControlVolumeFace& scvf) const
230 {
231 if constexpr (isCoupled_)
232 return couplingManager_->pressure(element, fvGeometry, scvf);
233 else
234 return asImp_().pressureAtPos(scvf.ipGlobal());
235 }
236
240 Scalar pressureAtPos(const GlobalPosition&) const
241 {
242 DUNE_THROW(Dune::NotImplemented, "pressureAtPos not implemented");
243 }
244
251 Scalar referencePressure(const Element& element,
252 const FVElementGeometry& fvGeometry,
253 const SubControlVolumeFace& scvf) const
254 { return 0.0; }
255
260 Scalar density(const Element& element,
261 const FVElementGeometry& fvGeometry,
262 const SubControlVolumeFace& scvf) const
263 {
264 if constexpr (isCoupled_)
265 return couplingManager_->density(element, fvGeometry, scvf);
266 else
267 return asImp_().densityAtPos(scvf.ipGlobal());
268 }
269
274 Scalar density(const Element& element,
275 const SubControlVolume& scv,
276 const bool isPreviousTimeStep = false) const
277 {
278 if constexpr (isCoupled_)
279 return couplingManager_->density(element, scv, isPreviousTimeStep);
280 else
281 return asImp_().densityAtPos(scv.dofPosition());
282 }
283
284 auto insideAndOutsideDensity(const Element& element,
285 const FVElementGeometry& fvGeometry,
286 const SubControlVolumeFace& scvf,
287 const bool isPreviousTimeStep = false) const
288 {
289 if constexpr (isCoupled_)
290 return couplingManager_->insideAndOutsideDensity(element, fvGeometry, scvf, isPreviousTimeStep);
291 else
292 {
293 const auto rho = asImp_().densityAtPos(scvf.ipGlobal());
294 return std::make_pair(rho, rho);
295 }
296 }
297
301 Scalar densityAtPos(const GlobalPosition&) const
302 {
303 DUNE_THROW(Dune::NotImplemented, "densityAtPos not implemented");
304 }
305
310 Scalar effectiveViscosity(const Element& element,
311 const FVElementGeometry& fvGeometry,
312 const SubControlVolumeFace& scvf) const
313 {
314 if constexpr (isCoupled_)
315 return couplingManager_->effectiveViscosity(element, fvGeometry, scvf);
316 else
317 return asImp_().effectiveViscosityAtPos(scvf.ipGlobal());
318 }
319
323 Scalar effectiveViscosityAtPos(const GlobalPosition&) const
324 {
325 DUNE_THROW(Dune::NotImplemented, "effectiveViscosityAtPos not implemented");
326 }
327
332 template<class SolutionVector>
333 void applyInitialSolution(SolutionVector& sol) const
334 {
335 sol.resize(this->gridGeometry().numDofs());
336 std::vector<bool> dofHandled(this->gridGeometry().numDofs(), false);
337 auto fvGeometry = localView(this->gridGeometry());
338 for (const auto& element : elements(this->gridGeometry().gridView()))
339 {
340 fvGeometry.bindElement(element);
341 for (const auto& scv : scvs(fvGeometry))
342 {
343 const auto dofIdx = scv.dofIndex();
344 if (!dofHandled[dofIdx])
345 {
346 dofHandled[dofIdx] = true;
347 sol[dofIdx] = asImp_().initial(scv)[scv.dofAxis()];
348 }
349 }
350 }
351 }
352
356 InitialValues initial(const SubControlVolume& scv) const
357 {
358 return asImp_().initialAtPos(scv.dofPosition());
359 }
360
362 Scalar pseudo3DWallFriction(const Element& element,
363 const FVElementGeometry& fvGeometry,
364 const ElementVolumeVariables& elemVolVars,
365 const SubControlVolume& scv,
366 const Scalar height,
367 const Scalar factor = 8.0) const
368 {
369 const Scalar velocity = elemVolVars[scv].velocity();
370 const auto scvf = scvfs(fvGeometry, scv).begin();
371 const Scalar viscosity = effectiveViscosity(element, fvGeometry, *scvf);
372 return pseudo3DWallFriction(velocity, viscosity, height, factor);
373 }
374
388 Scalar pseudo3DWallFriction(const Scalar velocity,
389 const Scalar viscosity,
390 const Scalar height,
391 const Scalar factor = 8.0) const
392 {
393 static_assert(dim == 2, "Pseudo 3D wall friction may only be used in 2D");
394 return -factor * velocity * viscosity / (height*height);
395 }
396
402 bool onSlipBoundary(const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf) const
403 { return asImp_().onSlipBoundaryAtPos(scvf.center()); }
404
409 bool onSlipBoundaryAtPos(const GlobalPosition& pos) const
410 { return false; }
411
416 Scalar permeability(const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf) const
417 {
418 DUNE_THROW(Dune::NotImplemented,
419 "When using the Beavers-Joseph-Saffman boundary condition, "
420 "the permeability must be returned in the actual problem"
421 );
422 }
423
428 Scalar alphaBJ(const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf) const
429 {
430 DUNE_THROW(Dune::NotImplemented,
431 "When using the Beavers-Joseph-Saffman boundary condition, "
432 "the alpha value must be returned in the actual problem"
433 );
434 }
435
439 Scalar betaBJ(const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf, const GlobalPosition& tangentialVector) const
440 {
441 const Scalar interfacePermeability = interfacePermeability_(fvGeometry, scvf, tangentialVector);
442 using std::sqrt;
443 return asImp_().alphaBJ(fvGeometry, scvf) / sqrt(interfacePermeability);
444 }
445
449 VelocityVector porousMediumVelocity(const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf) const
450 {
451 return VelocityVector(0.0);
452 }
453
458 const VelocityVector beaversJosephVelocity(const FVElementGeometry& fvGeometry,
459 const SubControlVolumeFace& scvf,
460 const ElementVolumeVariables& elemVolVars,
461 Scalar tangentialVelocityGradient) const
462 {
463 assert(scvf.isLateral());
464 assert(scvf.boundary());
465
466 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
467
468 // create a unit normal vector oriented in positive coordinate direction
469 GlobalPosition orientation(0.0);
470 orientation[scv.dofAxis()] = 1.0;
471
472 // du/dy + dv/dx = beta * (u_boundary-uPM)
473 // beta = alpha/sqrt(K)
474 const Scalar betaBJ = asImp_().betaBJ(fvGeometry, scvf, orientation);
475 const Scalar distanceNormalToBoundary = (scvf.ipGlobal() - scv.dofPosition()).two_norm();
476
477 static const bool onlyNormalGradient = getParamFromGroup<bool>(this->paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph", false);
478 if (onlyNormalGradient)
479 tangentialVelocityGradient = 0.0;
480
481 const Scalar scalarSlipVelocity = (tangentialVelocityGradient*distanceNormalToBoundary
482 + asImp_().porousMediumVelocity(fvGeometry, scvf) * orientation * betaBJ * distanceNormalToBoundary
483 + elemVolVars[scv].velocity()) / (betaBJ*distanceNormalToBoundary + 1.0);
484
485 orientation *= scalarSlipVelocity;
486 return orientation;
487 }
488
489 const CouplingManager& couplingManager() const
490 {
491 if constexpr (isCoupled_)
492 return *couplingManager_;
493 else
494 DUNE_THROW(Dune::InvalidStateException,
495 "Accessing coupling manager of an uncoupled problem is not possible."
496 );
497 }
498
499private:
501 template<class Scvf>
502 Scalar interfacePermeability_(const FVElementGeometry& fvGeometry, const Scvf& scvf, const GlobalPosition& tangentialVector) const
503 {
504 const auto& K = asImp_().permeability(fvGeometry, scvf);
505
506 // use t*K*t for permeability tensors
507 if constexpr (Dune::IsNumber<std::decay_t<decltype(K)>>::value)
508 return K;
509 else
510 return vtmv(tangentialVector, K, tangentialVector);
511 }
512
514 Implementation& asImp_()
515 { return *static_cast<Implementation *>(this); }
516
518 const Implementation& asImp_() const
519 { return *static_cast<const Implementation *>(this); }
520
521 GravityVector gravity_;
522 bool enableInertiaTerms_;
523 std::shared_ptr<CouplingManager> couplingManager_;
524};
525
529template<class TypeTag, class DM>
530class NavierStokesMomentumProblemImpl<TypeTag, DiscretizationMethods::CVFE<DM>>
531: public FVProblemWithSpatialParams<TypeTag>
532{
534 using Implementation = GetPropType<TypeTag, Properties::Problem>;
535
537 using GridView = typename GridGeometry::GridView;
538 using Element = typename GridView::template Codim<0>::Entity;
539
541 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
542 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
544
545 using FVElementGeometry = typename GridGeometry::LocalView;
546 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
547 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
549
550 static constexpr int dim = GridView::dimension;
551 static constexpr int dimWorld = GridView::dimensionworld;
552
553 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
554 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
555 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
558
559 static constexpr bool isCVFE = DiscretizationMethods::isCVFE<
560 typename GridGeometry::DiscretizationMethod
561 >;
562
563public:
567 using InitialValues = Dune::FieldVector<Scalar, dimWorld>;
568 using Sources = Dune::FieldVector<Scalar, dimWorld>;
569 using DirichletValues = Dune::FieldVector<Scalar, dimWorld>;
570 using BoundaryFluxes = Dune::FieldVector<Scalar, dimWorld>;
571
574
576 static constexpr bool isMomentumProblem() { return true; }
577
584 NavierStokesMomentumProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
585 std::shared_ptr<CouplingManager> couplingManager,
586 const std::string& paramGroup = "")
587 : ParentType(gridGeometry, paramGroup)
588 , gravity_(0.0)
589 , couplingManager_(couplingManager)
590 {
591 if (getParamFromGroup<bool>(paramGroup, "Problem.EnableGravity"))
592 gravity_[dim-1] = -9.81;
593
594 enableInertiaTerms_ = getParamFromGroup<bool>(paramGroup, "Problem.EnableInertiaTerms");
595 }
596
602 NavierStokesMomentumProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
603 const std::string& paramGroup = "")
604 : NavierStokesMomentumProblemImpl(gridGeometry, {}, paramGroup)
605 {}
606
625 template<class ElementVolumeVariables>
626 Sources source(const Element &element,
627 const FVElementGeometry& fvGeometry,
628 const ElementVolumeVariables& elemVolVars,
629 const SubControlVolume &scv) const
630 {
631 // forward to solution independent, fully-implicit specific interface
632 return asImp_().sourceAtPos(scv.center());
633 }
634
648 Sources sourceAtPos(const GlobalPosition &globalPos) const
649 {
652 return Sources(0.0);
653 }
654
662 BoundaryTypes boundaryTypes(const Element& element,
663 const SubControlVolume& scv) const
664 {
665 if (!isCVFE)
666 DUNE_THROW(Dune::InvalidStateException,
667 "boundaryTypes(..., scv) called for for a non CVFE method.");
668
669 // forward it to the method which only takes the global coordinate
670 return asImp_().boundaryTypesAtPos(scv.dofPosition());
671 }
672
680 BoundaryTypes boundaryTypes(const Element& element,
681 const SubControlVolumeFace& scvf) const
682 {
683 if (isCVFE)
684 DUNE_THROW(Dune::InvalidStateException,
685 "boundaryTypes(..., scvf) called for a CVFE method.");
686
687 // forward it to the method which only takes the global coordinate
688 return asImp_().boundaryTypesAtPos(scvf.ipGlobal());
689 }
690
698 DirichletValues dirichlet(const Element& element, const SubControlVolume& scv) const
699 {
700 // forward it to the method which only takes the global coordinate
701 if (!isCVFE)
702 DUNE_THROW(Dune::InvalidStateException, "dirichlet(scv) called for a non CVFE method.");
703 else
704 return asImp_().dirichletAtPos(scv.dofPosition());
705 }
706
714 DirichletValues dirichlet(const Element& element, const SubControlVolumeFace& scvf) const
715 {
716 // forward it to the method which only takes the global coordinate
717 if (isCVFE)
718 DUNE_THROW(Dune::InvalidStateException, "dirichlet(scvf) called for CVFE method.");
719 else
720 return asImp_().dirichletAtPos(scvf.ipGlobal());
721 }
722
729 const GravityVector& gravity() const
730 { return gravity_; }
731
736 { return enableInertiaTerms_; }
737
744 Scalar referencePressure() const
745 { return 0.0; }
746
751 Scalar pressure(const Element& element,
752 const FVElementGeometry& fvGeometry,
753 const SubControlVolumeFace& scvf) const
754 {
755 if constexpr (std::is_empty_v<CouplingManager>)
756 return asImp_().pressureAtPos(scvf.ipGlobal());
757 else
758 return couplingManager_->pressure(element, fvGeometry, scvf);
759 }
760
765 Scalar pressure(const Element& element,
766 const FVElementGeometry& fvGeometry,
767 const SubControlVolume& scv,
768 const bool isPreviousTimeStep = false) const
769 {
770 if constexpr (std::is_empty_v<CouplingManager>)
771 return asImp_().pressureAtPos(scv.dofPosition());
772 else
773 return couplingManager_->pressure(element, fvGeometry, scv, isPreviousTimeStep);
774 }
775
779 Scalar pressureAtPos(const GlobalPosition&) const
780 {
781 DUNE_THROW(Dune::NotImplemented, "pressureAtPos not implemented");
782 }
783
788 Scalar density(const Element& element,
789 const FVElementGeometry& fvGeometry,
790 const SubControlVolumeFace& scvf) const
791 {
792 if constexpr (std::is_empty_v<CouplingManager>)
793 return asImp_().densityAtPos(scvf.ipGlobal());
794 else
795 return couplingManager_->density(element, fvGeometry, scvf);
796 }
797
802 Scalar density(const Element& element,
803 const FVElementGeometry& fvGeometry,
804 const SubControlVolume& scv,
805 const bool isPreviousTimeStep = false) const
806 {
807 if constexpr (std::is_empty_v<CouplingManager>)
808 return asImp_().densityAtPos(scv.dofPosition());
809 else
810 return couplingManager_->density(element, fvGeometry, scv, isPreviousTimeStep);
811 }
812
813
817 Scalar densityAtPos(const GlobalPosition&) const
818 {
819 DUNE_THROW(Dune::NotImplemented, "densityAtPos not implemented");
820 }
821
826 Scalar effectiveViscosity(const Element& element,
827 const FVElementGeometry& fvGeometry,
828 const SubControlVolumeFace& scvf) const
829 {
830 if constexpr (std::is_empty_v<CouplingManager>)
831 return asImp_().effectiveViscosityAtPos(scvf.ipGlobal());
832 else
833 return couplingManager_->effectiveViscosity(element, fvGeometry, scvf);
834 }
835
840 Scalar effectiveViscosity(const Element& element,
841 const FVElementGeometry& fvGeometry,
842 const SubControlVolume& scv,
843 const bool isPreviousTimeStep = false) const
844 {
845 if constexpr (std::is_empty_v<CouplingManager>)
846 return asImp_().effectiveViscosityAtPos(scv.dofPosition());
847 else
848 return couplingManager_->effectiveViscosity(element, fvGeometry, scv, isPreviousTimeStep);
849 }
850
854 Scalar effectiveViscosityAtPos(const GlobalPosition&) const
855 {
856 DUNE_THROW(Dune::NotImplemented, "effectiveViscosityAtPos not implemented");
857 }
858
863 template<class SolutionVector>
864 void applyInitialSolution(SolutionVector& sol) const
865 {
866 static_assert(GridGeometry::discMethod == DiscretizationMethods::CVFE<DM>{});
867 sol.resize(this->gridGeometry().numDofs());
868 std::vector<bool> dofHandled(this->gridGeometry().numDofs(), false);
869 auto fvGeometry = localView(this->gridGeometry());
870 for (const auto& element : elements(this->gridGeometry().gridView()))
871 {
872 fvGeometry.bindElement(element);
873 for (const auto& scv : scvs(fvGeometry))
874 {
875 const auto dofIdx = scv.dofIndex();
876 if (!dofHandled[dofIdx])
877 {
878 dofHandled[dofIdx] = true;
879 sol[dofIdx] = asImp_().initial(scv);
880 }
881 }
882 }
883 }
884
888 InitialValues initial(const SubControlVolume& scv) const
889 {
890 static_assert(GridGeometry::discMethod == DiscretizationMethods::CVFE<DM>{});
891 return asImp_().initialAtPos(scv.dofPosition());
892 }
893
894private:
896 Implementation &asImp_()
897 { return *static_cast<Implementation *>(this); }
898
900 const Implementation &asImp_() const
901 { return *static_cast<const Implementation *>(this); }
902
903 GravityVector gravity_;
904 bool enableInertiaTerms_;
905 std::shared_ptr<CouplingManager> couplingManager_ = {};
906};
907
914template<class TypeTag>
917>;
918
919} // end namespace Dumux
920
921#endif
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:26
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:43
Base class for all finite-volume problems using spatial parameters.
Definition: fvproblemwithspatialparams.hh:29
Class to specify the type of a boundary condition for the Navier-Stokes model.
Definition: freeflow/navierstokes/momentum/boundarytypes.hh:25
Scalar effectiveViscosity(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the effective dynamic viscosity at a given sub control volume face.
Definition: freeflow/navierstokes/momentum/problem.hh:826
Dune::FieldVector< Scalar, dimWorld > Sources
Definition: freeflow/navierstokes/momentum/problem.hh:568
NavierStokesMomentumProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< CouplingManager > couplingManager, const std::string &paramGroup="")
The constructor.
Definition: freeflow/navierstokes/momentum/problem.hh:584
void applyInitialSolution(SolutionVector &sol) const
Applies the initial solution for all degrees of freedom of the grid.
Definition: freeflow/navierstokes/momentum/problem.hh:864
Scalar pressureAtPos(const GlobalPosition &) const
Returns the pressure at a given position.
Definition: freeflow/navierstokes/momentum/problem.hh:779
Scalar pressure(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, const bool isPreviousTimeStep=false) const
Returns the pressure at a given sub control volume.
Definition: freeflow/navierstokes/momentum/problem.hh:765
static constexpr bool isMomentumProblem()
This problem is used for the momentum balance model.
Definition: freeflow/navierstokes/momentum/problem.hh:576
Scalar densityAtPos(const GlobalPosition &) const
Returns the density at a given position.
Definition: freeflow/navierstokes/momentum/problem.hh:817
InitialValues initial(const SubControlVolume &scv) const
Evaluate the initial value at an sub control volume.
Definition: freeflow/navierstokes/momentum/problem.hh:888
Sources sourceAtPos(const GlobalPosition &globalPos) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: freeflow/navierstokes/momentum/problem.hh:648
Sources source(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: freeflow/navierstokes/momentum/problem.hh:626
Scalar effectiveViscosityAtPos(const GlobalPosition &) const
Returns the effective dynamic viscosity at a given position.
Definition: freeflow/navierstokes/momentum/problem.hh:854
Dune::FieldVector< Scalar, dimWorld > InitialValues
Definition: freeflow/navierstokes/momentum/problem.hh:567
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: freeflow/navierstokes/momentum/problem.hh:680
const GravityVector & gravity() const
Returns the acceleration due to gravity.
Definition: freeflow/navierstokes/momentum/problem.hh:729
Scalar density(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the density at a given sub control volume face.
Definition: freeflow/navierstokes/momentum/problem.hh:788
Dune::FieldVector< Scalar, dimWorld > DirichletValues
Definition: freeflow/navierstokes/momentum/problem.hh:569
DirichletValues dirichlet(const Element &element, const SubControlVolume &scv) const
Evaluate the boundary conditions for a Dirichlet control volume.
Definition: freeflow/navierstokes/momentum/problem.hh:698
bool enableInertiaTerms() const
Returns whether inertia terms should be considered.
Definition: freeflow/navierstokes/momentum/problem.hh:735
NavierStokesMomentumProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor for usage without a coupling manager.
Definition: freeflow/navierstokes/momentum/problem.hh:602
BoundaryTypes boundaryTypes(const Element &element, const SubControlVolume &scv) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: freeflow/navierstokes/momentum/problem.hh:662
Scalar density(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, const bool isPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: freeflow/navierstokes/momentum/problem.hh:802
Dune::FieldVector< Scalar, dimWorld > BoundaryFluxes
Definition: freeflow/navierstokes/momentum/problem.hh:570
Scalar referencePressure() const
Returns a reference pressure This pressure is subtracted from the actual pressure for the momentum ba...
Definition: freeflow/navierstokes/momentum/problem.hh:744
Scalar effectiveViscosity(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, const bool isPreviousTimeStep=false) const
Returns the effective dynamic viscosity at a given sub control volume.
Definition: freeflow/navierstokes/momentum/problem.hh:840
Scalar pressure(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the pressure at a given sub control volume face.
Definition: freeflow/navierstokes/momentum/problem.hh:751
DirichletValues dirichlet(const Element &element, const SubControlVolumeFace &scvf) const
Evaluate the boundary conditions for a Dirichlet control volume face.
Definition: freeflow/navierstokes/momentum/problem.hh:714
Sources source(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: freeflow/navierstokes/momentum/problem.hh:127
NavierStokesMomentumProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor for usage without a coupling manager.
Definition: freeflow/navierstokes/momentum/problem.hh:104
Scalar pressure(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the pressure at a given sub control volume face.
Definition: freeflow/navierstokes/momentum/problem.hh:227
auto 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: freeflow/navierstokes/momentum/problem.hh:163
Scalar pseudo3DWallFriction(const Scalar velocity, const Scalar viscosity, const Scalar height, const Scalar factor=8.0) const
An additional drag term can be included as source term for the momentum balance to mimic 3D flow beha...
Definition: freeflow/navierstokes/momentum/problem.hh:388
const VelocityVector beaversJosephVelocity(const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, Scalar tangentialVelocityGradient) const
Returns the slip velocity at a porous boundary based on the Beavers-Joseph(-Saffman) condition.
Definition: freeflow/navierstokes/momentum/problem.hh:458
Scalar effectiveViscosity(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the effective dynamic viscosity at a given sub control volume face.
Definition: freeflow/navierstokes/momentum/problem.hh:310
bool enableInertiaTerms() const
Returns whether inertia terms should be considered.
Definition: freeflow/navierstokes/momentum/problem.hh:220
Scalar permeability(const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boun...
Definition: freeflow/navierstokes/momentum/problem.hh:416
NavierStokesMomentumProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< CouplingManager > couplingManager, const std::string &paramGroup="")
The constructor.
Definition: freeflow/navierstokes/momentum/problem.hh:86
Scalar betaBJ(const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const GlobalPosition &tangentialVector) const
Returns the beta value which is the alpha value divided by the square root of the (scalar-valued) int...
Definition: freeflow/navierstokes/momentum/problem.hh:439
VelocityVector porousMediumVelocity(const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the velocity in the porous medium (which is 0 by default according to Saffmann).
Definition: freeflow/navierstokes/momentum/problem.hh:449
Scalar pseudo3DWallFriction(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv, const Scalar height, const Scalar factor=8.0) const
Convenience function for staggered grid implementation.
Definition: freeflow/navierstokes/momentum/problem.hh:362
auto insideAndOutsideDensity(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const bool isPreviousTimeStep=false) const
Definition: freeflow/navierstokes/momentum/problem.hh:284
Scalar effectiveViscosityAtPos(const GlobalPosition &) const
Returns the effective dynamic viscosity at a given position.
Definition: freeflow/navierstokes/momentum/problem.hh:323
const CouplingManager & couplingManager() const
Definition: freeflow/navierstokes/momentum/problem.hh:489
Dune::FieldVector< Scalar, dimWorld > MomentumFluxType
Definition: freeflow/navierstokes/momentum/problem.hh:72
static constexpr bool isMomentumProblem()
This problem is used for the momentum balance model.
Definition: freeflow/navierstokes/momentum/problem.hh:78
Scalar referencePressure(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns a reference pressure at a given sub control volume face. This pressure is subtracted from the...
Definition: freeflow/navierstokes/momentum/problem.hh:251
void applyInitialSolution(SolutionVector &sol) const
Applies the initial solution for all degrees of freedom of the grid.
Definition: freeflow/navierstokes/momentum/problem.hh:333
Scalar density(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the density at a given sub control volume face.
Definition: freeflow/navierstokes/momentum/problem.hh:260
const GravityVector & gravity() const
A default, i.e. if the user's does not overload any neumann method.
Definition: freeflow/navierstokes/momentum/problem.hh:214
bool onSlipBoundary(const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns true if the scvf is located on a boundary with a slip condition.
Definition: freeflow/navierstokes/momentum/problem.hh:402
Dune::FieldVector< Scalar, dimWorld > InitialValues
Definition: freeflow/navierstokes/momentum/problem.hh:67
Scalar density(const Element &element, const SubControlVolume &scv, const bool isPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: freeflow/navierstokes/momentum/problem.hh:274
BoundaryFluxes neumannAtPos(const GlobalPosition &globalPos) const
Returns the neumann flux at a given position.
Definition: freeflow/navierstokes/momentum/problem.hh:205
InitialValues initial(const SubControlVolume &scv) const
Evaluate the initial value at an sub control volume.
Definition: freeflow/navierstokes/momentum/problem.hh:356
Scalar alphaBJ(const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: freeflow/navierstokes/momentum/problem.hh:428
Scalar densityAtPos(const GlobalPosition &) const
Returns the density at a given position.
Definition: freeflow/navierstokes/momentum/problem.hh:301
Dune::FieldVector< Scalar, dimWorld > BoundaryFluxes
Definition: freeflow/navierstokes/momentum/problem.hh:70
DirichletValues dirichlet(const Element &element, const SubControlVolumeFace &scvf) const
Evaluate the boundary conditions for a dirichlet control volume face.
Definition: freeflow/navierstokes/momentum/problem.hh:180
Dune::FieldVector< Scalar, dimWorld > DirichletValues
Definition: freeflow/navierstokes/momentum/problem.hh:69
Dune::FieldVector< Scalar, dimWorld > Sources
Definition: freeflow/navierstokes/momentum/problem.hh:68
BoundaryFluxes neumann(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluates the boundary conditions for a Neumann control volume.
Definition: freeflow/navierstokes/momentum/problem.hh:195
bool onSlipBoundaryAtPos(const GlobalPosition &pos) const
Returns true if the scvf is located on a boundary with a slip condition.
Definition: freeflow/navierstokes/momentum/problem.hh:409
Scalar pressureAtPos(const GlobalPosition &) const
Returns the pressure at a given position.
Definition: freeflow/navierstokes/momentum/problem.hh:240
Sources sourceAtPos(const GlobalPosition &globalPos) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: freeflow/navierstokes/momentum/problem.hh:149
Definition: freeflow/navierstokes/momentum/problem.hh:26
Defines all properties used in Dumux.
Base class for all finite volume problems that are parameterized.
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:851
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
constexpr bool isCVFE
Definition: method.hh:67
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:62
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
Definition: method.hh:46