version 3.10-dev
freeflow/navierstokes/iofields.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_NAVIER_STOKES_IO_FIELDS_HH
13#define DUMUX_NAVIER_STOKES_IO_FIELDS_HH
14
16#include <dumux/io/name.hh>
17
18namespace Dumux {
19
25template<class IOFields, class PrimaryVariables, class ModelTraits, class FluidSystem>
26std::function<std::string(int,int)> createCellCenterPVNameFunction(const std::string& paramGroup = "")
27{
28 static constexpr auto offset = ModelTraits::numEq() - PrimaryVariables::dimension;
29
30 if (hasParamInGroup(paramGroup, "LoadSolution.CellCenterPriVarNames"))
31 {
32 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.CellCenterPriVarNames");
33 return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
34 }
35 else
36 // add an offset to the pvIdx so that the velocities are skipped
37 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx + offset, state); };
38}
39
45template<class IOFields, class PrimaryVariables, class ModelTraits, class FluidSystem>
46std::function<std::string(int,int)> createFacePVNameFunction(const std::string& paramGroup = "")
47{
48 if (hasParamInGroup(paramGroup, "LoadSolution.FacePriVarNames"))
49 {
50 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.FacePriVarNames");
51 return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
52 }
53 else
54 // there is only one scalar-type primary variable called "v" living on the faces
55 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx , state); };
56}
57
58// forward declare
59template<class T, class U>
60class StaggeredVtkOutputModule;
61
67{
69 template<class T>
70 struct isStaggered : public std::false_type {};
71
72 template<class... Args>
73 struct isStaggered<StaggeredVtkOutputModule<Args...>>
74 : public std::true_type {};
75
76public:
78 template <class OutputModule>
79 static void initOutputModule(OutputModule& out)
80 {
81 out.addVolumeVariable([](const auto& v){ return v.pressure(); }, IOName::pressure());
82 out.addVolumeVariable([](const auto& v){ return v.density(); }, IOName::density());
83
84 // add discretization-specific fields
85 additionalOutput_(out, isStaggered<OutputModule>());
86 }
87
89 template <class ModelTraits, class FluidSystem = void>
90 static std::string primaryVariableName(int pvIdx = 0, int state = 0)
91 {
92 if (pvIdx < ModelTraits::dim())
93 return "v";
94 else
95 return IOName::pressure();
96 }
97
98private:
99
101 template <class OutputModule>
102 static void additionalOutput_(OutputModule& out, std::false_type)
103 { }
104
106 template <class OutputModule>
107 static void additionalOutput_(OutputModule& out, std::true_type)
108 {
109 const bool writeFaceVars = getParamFromGroup<bool>(out.paramGroup(), "Vtk.WriteFaceData", false);
110 if(writeFaceVars)
111 {
112 auto faceVelocityVector = [](const auto& scvf, const auto& faceVars)
113 {
114 using VelocityVector = std::decay_t<decltype(scvf.unitOuterNormal())>;
115
116 VelocityVector velocity(0.0);
117 velocity[scvf.directionIndex()] = faceVars.velocitySelf();
118 return velocity;
119 };
120
121 out.addFaceVariable(faceVelocityVector, "faceVelocity");
122
123 auto faceNormalVelocity = [](const auto& faceVars)
124 {
125 return faceVars.velocitySelf();
126 };
127
128 out.addFaceVariable(faceNormalVelocity, "v");
129 }
130 }
131};
132
133} // end namespace Dumux
134
135#endif
Adds I/O fields for the Navier-Stokes model.
Definition: freeflow/navierstokes/iofields.hh:67
static void initOutputModule(OutputModule &out)
Initialize the Navier-Stokes specific output fields.
Definition: freeflow/navierstokes/iofields.hh:79
static std::string primaryVariableName(int pvIdx=0, int state=0)
return the names of the primary variables
Definition: freeflow/navierstokes/iofields.hh:90
A VTK output module to simplify writing dumux simulation data to VTK format Specialization for stagge...
Definition: staggeredvtkoutputmodule.hh:38
std::function< std::string(int, int)> createFacePVNameFunction(const std::string &paramGroup="")
helper function to determine the names of primary variables on the cell faces of a model with stagger...
Definition: freeflow/navierstokes/iofields.hh:46
std::function< std::string(int, int)> createCellCenterPVNameFunction(const std::string &paramGroup="")
helper function to determine the names of cell-centered primary variables of a model with staggered g...
Definition: freeflow/navierstokes/iofields.hh:26
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:165
A collection of input/output field names for common physical quantities.
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:22
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.