24#ifndef DUMUX_NAVIER_STOKES_IO_FIELDS_HH
25#define DUMUX_NAVIER_STOKES_IO_FIELDS_HH
37template<
class IOFields,
class PrimaryVariables,
class ModelTraits,
class Flu
idSystem>
40 static constexpr auto offset = ModelTraits::numEq() - PrimaryVariables::dimension;
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]; };
49 return [](
int pvIdx,
int state = 0){
return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx + offset, state); };
57template<
class IOFields,
class PrimaryVariables,
class ModelTraits,
class Flu
idSystem>
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]; };
67 return [](
int pvIdx,
int state = 0){
return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx , state); };
71template<
class T,
class U>
72class StaggeredVtkOutputModule;
82 struct isStaggered :
public std::false_type {};
84 template<
class... Args>
86 :
public std::true_type {};
90 template <
class OutputModule>
93 out.addVolumeVariable([](
const auto& v){
return v.pressure(); },
IOName::pressure());
94 out.addVolumeVariable([](
const auto& v){
return v.density(); },
IOName::density());
97 additionalOutput_(out, isStaggered<OutputModule>());
101 template <
class ModelTraits,
class Flu
idSystem =
void>
104 if (pvIdx < ModelTraits::dim())
113 template <
class OutputModule>
114 static void additionalOutput_(OutputModule& out, std::false_type)
118 template <
class OutputModule>
119 static void additionalOutput_(OutputModule& out, std::true_type)
121 const bool writeFaceVars = getParamFromGroup<bool>(out.paramGroup(),
"Vtk.WriteFaceData",
false);
124 auto faceVelocityVector = [](
const auto& scvf,
const auto& faceVars)
126 using VelocityVector = std::decay_t<
decltype(scvf.unitOuterNormal())>;
128 VelocityVector velocity(0.0);
129 velocity[scvf.directionIndex()] = faceVars.velocitySelf();
133 out.addFaceVariable(faceVelocityVector,
"faceVelocity");
135 auto faceNormalVelocity = [](
const auto& faceVars)
137 return faceVars.velocitySelf();
140 out.addFaceVariable(faceNormalVelocity,
"v");
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 ¶mGroup="")
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 ¶mGroup="")
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 ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:177
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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