Line | Branch | Exec | Source |
---|---|---|---|
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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup NavierStokesModel | ||
10 | * \copydoc Dumux::NavierStokesIOFields | ||
11 | */ | ||
12 | #ifndef DUMUX_NAVIER_STOKES_IO_FIELDS_HH | ||
13 | #define DUMUX_NAVIER_STOKES_IO_FIELDS_HH | ||
14 | |||
15 | #include <dumux/common/parameters.hh> | ||
16 | #include <dumux/io/name.hh> | ||
17 | |||
18 | namespace Dumux { | ||
19 | |||
20 | /*! | ||
21 | * \ingroup NavierStokesModel | ||
22 | * \brief helper function to determine the names of cell-centered primary variables of a model with staggered grid discretization | ||
23 | * \note use this as input for the load solution function | ||
24 | */ | ||
25 | template<class IOFields, class PrimaryVariables, class ModelTraits, class FluidSystem> | ||
26 | 49 | std::function<std::string(int,int)> createCellCenterPVNameFunction(const std::string& paramGroup = "") | |
27 | { | ||
28 | static constexpr auto offset = ModelTraits::numEq() - PrimaryVariables::dimension; | ||
29 | |||
30 |
1/2✗ Branch 2 not taken.
✓ Branch 3 taken 25 times.
|
98 | 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 | 211 | return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx + offset, state); }; | |
38 | } | ||
39 | |||
40 | /*! | ||
41 | * \ingroup NavierStokesModel | ||
42 | * \brief helper function to determine the names of primary variables on the cell faces of a model with staggered grid discretization | ||
43 | * \note use this as input for the load solution function | ||
44 | */ | ||
45 | template<class IOFields, class PrimaryVariables, class ModelTraits, class FluidSystem> | ||
46 | 49 | std::function<std::string(int,int)> createFacePVNameFunction(const std::string& paramGroup = "") | |
47 | { | ||
48 |
1/2✗ Branch 2 not taken.
✓ Branch 3 taken 25 times.
|
98 | 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 | 99 | return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx , state); }; | |
56 | } | ||
57 | |||
58 | // forward declare | ||
59 | template<class T, class U> | ||
60 | class StaggeredVtkOutputModule; | ||
61 | |||
62 | /*! | ||
63 | * \ingroup NavierStokesModel | ||
64 | * \brief Adds I/O fields for the Navier-Stokes model | ||
65 | */ | ||
66 | class NavierStokesIOFields | ||
67 | { | ||
68 | //! Helper strcuts to determine whether a staggered grid discretization is used | ||
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 | |||
76 | public: | ||
77 | //! Initialize the Navier-Stokes specific output fields. | ||
78 | template <class OutputModule> | ||
79 | 177 | static void initOutputModule(OutputModule& out) | |
80 | { | ||
81 |
2/3✓ Branch 2 taken 75 times.
✗ Branch 3 not taken.
✓ Branch 1 taken 78 times.
|
954 | out.addVolumeVariable([](const auto& v){ return v.pressure(); }, IOName::pressure()); |
82 |
1/2✓ Branch 1 taken 153 times.
✗ Branch 2 not taken.
|
954 | out.addVolumeVariable([](const auto& v){ return v.density(); }, IOName::density()); |
83 | |||
84 | // add discretization-specific fields | ||
85 | 99 | additionalOutput_(out, isStaggered<OutputModule>()); | |
86 | 177 | } | |
87 | |||
88 | //! return the names of the primary variables | ||
89 | template <class ModelTraits, class FluidSystem = void> | ||
90 | 196 | static std::string primaryVariableName(int pvIdx = 0, int state = 0) | |
91 | { | ||
92 |
2/2✓ Branch 0 taken 50 times.
✓ Branch 1 taken 50 times.
|
196 | if (pvIdx < ModelTraits::dim()) |
93 | 98 | return "v"; | |
94 | else | ||
95 | return IOName::pressure(); | ||
96 | } | ||
97 | |||
98 | private: | ||
99 | |||
100 | //! Adds discretization-specific fields (nothing by default). | ||
101 | template <class OutputModule> | ||
102 | static void additionalOutput_(OutputModule& out, std::false_type) | ||
103 | { } | ||
104 | |||
105 | //! Adds discretization-specific fields (velocity vectors on the faces for the staggered discretization). | ||
106 | template <class OutputModule> | ||
107 | 99 | static void additionalOutput_(OutputModule& out, std::true_type) | |
108 | { | ||
109 | 99 | const bool writeFaceVars = getParamFromGroup<bool>(out.paramGroup(), "Vtk.WriteFaceData", false); | |
110 |
2/2✓ Branch 0 taken 26 times.
✓ Branch 1 taken 49 times.
|
99 | if(writeFaceVars) |
111 | { | ||
112 | 50055 | auto faceVelocityVector = [](const auto& scvf, const auto& faceVars) | |
113 | { | ||
114 | using VelocityVector = std::decay_t<decltype(scvf.unitOuterNormal())>; | ||
115 | |||
116 | 50055 | VelocityVector velocity(0.0); | |
117 | 50055 | velocity[scvf.directionIndex()] = faceVars.velocitySelf(); | |
118 | 50055 | return velocity; | |
119 | }; | ||
120 | |||
121 |
2/4✓ Branch 2 taken 26 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 26 times.
✗ Branch 7 not taken.
|
100 | out.addFaceVariable(faceVelocityVector, "faceVelocity"); |
122 | |||
123 | 50055 | auto faceNormalVelocity = [](const auto& faceVars) | |
124 | { | ||
125 | 50055 | return faceVars.velocitySelf(); | |
126 | }; | ||
127 | |||
128 |
2/4✓ Branch 1 taken 26 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 26 times.
✗ Branch 5 not taken.
|
100 | out.addFaceVariable(faceNormalVelocity, "v"); |
129 | } | ||
130 | 99 | } | |
131 | }; | ||
132 | |||
133 | } // end namespace Dumux | ||
134 | |||
135 | #endif | ||
136 |