version 3.8
freeflow/navierstokes/velocityoutput.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_FREEFLOW_NAVIERSTOKES_VELOCITYOUTPUT_HH
13#define DUMUX_FREEFLOW_NAVIERSTOKES_VELOCITYOUTPUT_HH
14
15#include <type_traits>
16#include <dune/common/exceptions.hh>
21
22namespace Dumux {
23
28template<class GridVariables>
29class NavierStokesVelocityOutput : public VelocityOutput<GridVariables>
30{
32 using GridGeometry = typename GridVariables::GridGeometry;
33 using FVElementGeometry = typename GridGeometry::LocalView;
34 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
35 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
36 using ElementFluxVarsCache = typename GridVariables::GridFluxVariablesCache::LocalView;
37 using VolumeVariables = typename GridVariables::VolumeVariables;
38 using FluidSystem = typename VolumeVariables::FluidSystem;
39 using GridView = typename GridGeometry::GridView;
40 using Element = typename GridView::template Codim<0>::Entity;
41 using FieldType = typename ParentType::FieldType;
42
43public:
45
46 NavierStokesVelocityOutput(const std::string& paramGroup = "")
47 {
48 enableOutput_ = getParamFromGroup<bool>(paramGroup, "Vtk.AddVelocity", true);
49 }
50
52 bool enableOutput() const override { return enableOutput_; }
53
55 std::string phaseName(int phaseIdx) const override { return FluidSystem::phaseName(phaseIdx); }
56
58 int numFluidPhases() const override { return VolumeVariables::numFluidPhases(); }
59
61 FieldType fieldType() const override { return FieldType::element; }
62
66 const Element& element,
67 const FVElementGeometry& fvGeometry,
68 const ElementVolumeVariables& elemVolVars,
69 const ElementFluxVarsCache& elemFluxVarsCache,
70 int phaseIdx) const override
71 {
72 using CouplingManager = std::decay_t<decltype(elemVolVars.gridVolVars().problem().couplingManager())>;
73 using MomGG = std::decay_t<decltype(std::declval<CouplingManager>().problem(CouplingManager::freeFlowMomentumIndex).gridGeometry())>;
74 if constexpr (MomGG::discMethod == DiscretizationMethods::fcstaggered)
75 calculateVelocityForStaggeredGrid_(velocity, element, fvGeometry, elemVolVars);
76 else if constexpr (DiscretizationMethods::isCVFE<typename MomGG::DiscretizationMethod>)
77 calculateVelocityForCVFESchemes_(velocity, element, fvGeometry, elemVolVars);
78 else
79 DUNE_THROW(Dune::NotImplemented, "Navier-Stokes velocity output for scheme " << MomGG::discMethod);
80 }
81
82private:
83 void calculateVelocityForStaggeredGrid_(VelocityVector& velocity,
84 const Element& element,
85 const FVElementGeometry& fvGeometry,
86 const ElementVolumeVariables& elemVolVars) const
87 {
88 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
89 const auto getFaceVelocity = [&](const FVElementGeometry& fvG, const auto& scvf)
90 {
91 return elemVolVars.gridVolVars().problem().faceVelocity(element, fvGeometry, scvf);
92 };
93
94 velocity[eIdx] = StaggeredVelocityReconstruction::cellCenterVelocity(getFaceVelocity, fvGeometry);
95 }
96
97 void calculateVelocityForCVFESchemes_(VelocityVector& velocity,
98 const Element& element,
99 const FVElementGeometry& fvGeometry,
100 const ElementVolumeVariables& elemVolVars) const
101 {
102 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
103 velocity[eIdx] = elemVolVars.gridVolVars().problem().elementVelocity(fvGeometry);
104 }
105
106
107 bool enableOutput_;
108};
109
110} // end namespace Dumux
111
112#endif
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:309
Velocity output for staggered free-flow models.
Definition: freeflow/navierstokes/velocityoutput.hh:30
NavierStokesVelocityOutput(const std::string &paramGroup="")
Definition: freeflow/navierstokes/velocityoutput.hh:46
void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const override
Definition: freeflow/navierstokes/velocityoutput.hh:65
std::string phaseName(int phaseIdx) const override
returns the phase name of a given phase index
Definition: freeflow/navierstokes/velocityoutput.hh:55
typename ParentType::VelocityVector VelocityVector
Definition: freeflow/navierstokes/velocityoutput.hh:44
bool enableOutput() const override
Returns whether to enable the velocity output or not.
Definition: freeflow/navierstokes/velocityoutput.hh:52
FieldType fieldType() const override
returns the field type
Definition: freeflow/navierstokes/velocityoutput.hh:61
int numFluidPhases() const override
returns the number of phases
Definition: freeflow/navierstokes/velocityoutput.hh:58
Velocity output for implicit (porous media) models.
Definition: io/velocityoutput.hh:29
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: io/velocityoutput.hh:38
FieldType
A container for possible velocity data types.
Definition: io/velocityoutput.hh:44
Default velocity output policy for porous media models.
The available discretization methods in Dumux.
constexpr FCStaggered fcstaggered
Definition: method.hh:151
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporenetwork/couplingmanager.hh:34
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
static auto cellCenterVelocity(const VelocityHelper &getFaceVelocity, const FVElementGeometry &fvGeometry)
Return the velocity vector at the center of the primal grid.
Definition: velocityreconstruction.hh:27