3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_NAVIER_STOKES_IO_FIELDS_HH
25#define DUMUX_NAVIER_STOKES_IO_FIELDS_HH
26
28#include <dumux/io/name.hh>
29
30namespace Dumux {
31
37template<class IOFields, class PrimaryVariables, class ModelTraits, class FluidSystem>
38std::function<std::string(int,int)> createCellCenterPVNameFunction(const std::string& paramGroup = "")
39{
40 static constexpr auto offset = ModelTraits::numEq() - PrimaryVariables::dimension;
41
42 if (hasParamInGroup(paramGroup, "LoadSolution.CellCenterPriVarNames"))
43 {
44 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.CellCenterPriVarNames");
45 return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
46 }
47 else
48 // add an offset to the pvIdx so that the velocities are skipped
49 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx + offset, state); };
50}
51
57template<class IOFields, class PrimaryVariables, class ModelTraits, class FluidSystem>
58std::function<std::string(int,int)> createFacePVNameFunction(const std::string& paramGroup = "")
59{
60 if (hasParamInGroup(paramGroup, "LoadSolution.FacePriVarNames"))
61 {
62 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.FacePriVarNames");
63 return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
64 }
65 else
66 // there is only one scalar-type primary variable called "v" living on the faces
67 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx , state); };
68}
69
70// forward declare
71template<class T, class U>
72class StaggeredVtkOutputModule;
73
79{
81 template<class T>
82 struct isStaggered : public std::false_type {};
83
84 template<class... Args>
85 struct isStaggered<StaggeredVtkOutputModule<Args...>>
86 : public std::true_type {};
87
88public:
90 template <class OutputModule>
91 static void initOutputModule(OutputModule& out)
92 {
93 out.addVolumeVariable([](const auto& v){ return v.pressure(); }, IOName::pressure());
94 out.addVolumeVariable([](const auto& v){ return v.density(); }, IOName::density());
95
96 // add discretization-specific fields
97 additionalOutput_(out, isStaggered<OutputModule>());
98 }
99
101 template <class ModelTraits, class FluidSystem = void>
102 static std::string primaryVariableName(int pvIdx = 0, int state = 0)
103 {
104 if (pvIdx < ModelTraits::dim())
105 return "v";
106 else
107 return IOName::pressure();
108 }
109
110private:
111
113 template <class OutputModule>
114 static void additionalOutput_(OutputModule& out)
115 { }
116
118 template <class OutputModule>
119 static void additionalOutput_(OutputModule& out, std::true_type)
120 {
121 const bool writeFaceVars = getParamFromGroup<bool>(out.paramGroup(), "Vtk.WriteFaceData", false);
122 if(writeFaceVars)
123 {
124 auto faceVelocityVector = [](const auto& scvf, const auto& faceVars)
125 {
126 using VelocityVector = std::decay_t<decltype(scvf.unitOuterNormal())>;
127
128 VelocityVector velocity(0.0);
129 velocity[scvf.directionIndex()] = faceVars.velocitySelf();
130 return velocity;
131 };
132
133 out.addFaceVariable(faceVelocityVector, "faceVelocity");
134
135 auto faceNormalVelocity = [](const auto& faceVars)
136 {
137 return faceVars.velocitySelf();
138 };
139
140 out.addFaceVariable(faceNormalVelocity, "v");
141 }
142 }
143};
144
145} // end namespace Dumux
146
147#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A collection of input/output field names for common physical quantities.
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:58
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:38
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:391
Definition: adapt.hh:29
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
A VTK output module to simplify writing dumux simulation data to VTK format Specialization for stagge...
Definition: staggeredvtkoutputmodule.hh:50
Adds I/O fields for the Navier-Stokes model.
Definition: freeflow/navierstokes/iofields.hh:79
static void initOutputModule(OutputModule &out)
Initialize the Navier-Stokes specific output fields.
Definition: freeflow/navierstokes/iofields.hh:91
static std::string primaryVariableName(int pvIdx=0, int state=0)
return the names of the primary variables
Definition: freeflow/navierstokes/iofields.hh:102