/home/runner/work/kynema-sgf/kynema-sgf/src/transport_models/ConstTransport.H Source File

Kynema-SGF API: /home/runner/work/kynema-sgf/kynema-sgf/src/transport_models/ConstTransport.H Source File
Kynema-SGF API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
ConstTransport.H
Go to the documentation of this file.
1#ifndef CONSTTRANSPORT_H
2#define CONSTTRANSPORT_H
3
5#include "AMReX_ParmParse.H"
6#include "AMReX_REAL.H"
7
8using namespace amrex::literals;
9
10namespace kynema_sgf::transport {
11
15class ConstTransport : public TransportModel::Register<ConstTransport>
16{
17public:
18 static constexpr bool constant_properties = true;
19
20 static std::string identifier() { return "ConstTransport"; }
21
22 explicit ConstTransport(const CFDSim& sim) : m_repo(sim.repo())
23 {
24 amrex::ParmParse pp("transport");
25 pp.query("viscosity", m_mu);
26 pp.query("laminar_prandtl", m_Pr);
27 pp.query("turbulent_prandtl", m_Prt);
28
29 // Backwards compatibility
30 amrex::ParmParse pp_boussinesq_buoyancy("BoussinesqBuoyancy");
31 amrex::ParmParse pp_abl("ABL");
32 if (pp.contains("thermal_expansion_coefficient")) {
33 pp.get("thermal_expansion_coefficient", m_constant_beta);
34 if (pp_boussinesq_buoyancy.contains("thermal_expansion_coeff")) {
35 amrex::Print()
36 << "WARNING: BoussinesqBuoyancy.thermal_expansion_coeff "
37 "option has been deprecated in favor of "
38 "transport.thermal_expansion_coefficient. Ignoring the "
39 "BoussinesqBuoyancy option in favor of the transport "
40 "option."
41 << '\n';
42 }
43 } else if (pp_boussinesq_buoyancy.contains("thermal_expansion_coeff")) {
44 amrex::Print()
45 << "WARNING: BoussinesqBuoyancy.thermal_expansion_coeff option "
46 "has been deprecated in favor of "
47 "transport.thermal_expansion_coefficient. Please replace "
48 "this option."
49 << '\n';
50 pp_boussinesq_buoyancy.get(
51 "thermal_expansion_coeff", m_constant_beta);
52 }
53
54 if (pp.contains("reference_temperature")) {
55 pp.get("reference_temperature", m_reference_temperature);
56 if (pp_boussinesq_buoyancy.contains("reference_temperature")) {
57 amrex::Print()
58 << "WARNING: BoussinesqBuoyancy.reference_temperature "
59 "option has been deprecated in favor of "
60 "transport.reference_temperature. Ignoring the "
61 "BoussinesqBuoyancy option in favor of the transport "
62 "option."
63 << '\n';
64 } else if (pp_abl.contains("reference_temperature")) {
65 amrex::Print()
66 << "WARNING: ABL.reference_temperature "
67 "option has been deprecated in favor of "
68 "transport.reference_temperature. Ignoring the "
69 "ABL option in favor of the transport "
70 "option."
71 << '\n';
72 }
73 } else if (pp_boussinesq_buoyancy.contains("reference_temperature")) {
74 amrex::Print()
75 << "WARNING: BoussinesqBuoyancy.reference_temperature option "
76 "has been deprecated in favor of "
77 "transport.reference_temperature. Please replace "
78 "this option."
79 << '\n';
80 pp_boussinesq_buoyancy.get(
81 "reference_temperature", m_reference_temperature);
82 } else if (pp_abl.contains("reference_temperature")) {
83 amrex::Print() << "WARNING: ABL.reference_temperature option "
84 "has been deprecated in favor of "
85 "transport.reference_temperature. Please replace "
86 "this option."
87 << '\n';
88 pp_abl.get("reference_temperature", m_reference_temperature);
89 }
90 }
91
92 ~ConstTransport() override = default;
93
94 [[nodiscard]] amrex::Real viscosity() const { return m_mu; }
95
96 [[nodiscard]] amrex::Real thermal_diffusivity() const
97 {
98 return m_mu / m_Pr;
99 }
100
101 [[nodiscard]] amrex::Real laminar_prandtl() const { return m_Pr; }
102
103 [[nodiscard]] amrex::Real turbulent_prandtl() const { return m_Prt; }
104
105 static amrex::Real laminar_schmidt(const std::string& scalar_name)
106 {
107 amrex::ParmParse pp("transport");
108 const std::string key = scalar_name + "_laminar_schmidt";
109 amrex::Real lam_schmidt = 1.0_rt;
110 pp.query(key, lam_schmidt);
111 return lam_schmidt;
112 }
113
114 static amrex::Real turbulent_schmidt(const std::string& scalar_name)
115 {
116 amrex::ParmParse pp("transport");
117 const std::string key = scalar_name + "_turbulent_schmidt";
118 amrex::Real turb_schmidt = 1.0_rt;
119 pp.query(key, turb_schmidt);
120 return turb_schmidt;
121 }
122
124 std::unique_ptr<ScratchField> mu() override
125 {
126 auto mu = m_repo.create_scratch_field(1, m_ngrow);
127 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
128 (*mu)(lev).setVal(m_mu);
129 }
130 return mu;
131 }
132
134 std::unique_ptr<ScratchField> alpha() override
135 {
136 auto alpha = mu();
137 amrex::Real inv_Pr = 1.0_rt / m_Pr;
138 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
139 (*alpha)(lev).mult(inv_Pr);
140 }
141 return alpha;
142 }
143
144 std::unique_ptr<ScratchField>
145 scalar_diffusivity(const std::string& scalar_name) override
146 {
147 amrex::Real lam_schmidt = laminar_schmidt(scalar_name);
148
149 amrex::Real inv_schmidt = 1.0_rt / lam_schmidt;
150 auto diff = mu();
151 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
152 (*diff)(lev).mult(inv_schmidt);
153 }
154
155 return diff;
156 }
157
159 [[nodiscard]] std::unique_ptr<ScratchField> beta() const override
160 {
161 auto beta = m_repo.create_scratch_field(1, m_ngrow);
162 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
163 beta_fill(lev, (*beta)(lev));
164 }
165 return beta;
166 }
167
168 [[nodiscard]] amrex::Real reference_temperature() const override
169 {
171 }
172
173 void beta_fill(const int lev, amrex::MultiFab& beta_mf) const override
174 {
175 if (m_constant_beta > 0.0_rt) {
176 beta_mf.setVal(m_constant_beta);
177 } else if (m_repo.field_exists("reference_temperature")) {
178 const auto& temp0 = m_repo.get_field("reference_temperature");
179 auto const& beta_arrs = beta_mf.arrays();
180 auto const& temp0_arrs = temp0(lev).const_arrays();
181 amrex::ParallelFor(
182 beta_mf, amrex::IntVect(0), 1,
183 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int) {
184 beta_arrs[nbx](i, j, k) = 1.0_rt / temp0_arrs[nbx](i, j, k);
185 });
186 amrex::Gpu::streamSynchronize();
187 } else {
188 beta_mf.setVal(1.0_rt / m_reference_temperature);
189 }
190
191 if (m_repo.field_exists("vof")) {
192 const auto& vof = m_repo.get_field("vof");
193 auto const& beta_arrs = beta_mf.arrays();
194 auto const& vof_arrs = vof(lev).const_arrays();
195 amrex::ParallelFor(
196 beta_mf, amrex::IntVect(0), 1,
197 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int) {
198 if (vof_arrs[nbx](i, j, k) > constants::TIGHT_TOL) {
199 beta_arrs[nbx](i, j, k) = 0.0_rt;
200 }
201 });
202 amrex::Gpu::streamSynchronize();
203 }
204 }
205
206 void
207 ref_theta_fill(const int lev, amrex::MultiFab& ref_theta_mf) const override
208 {
209 if (m_reference_temperature < 0.0_rt) {
210 amrex::Abort("Reference temperature was not set");
211 }
212
213 if (m_repo.field_exists("reference_temperature")) {
214 const auto& temp0 = m_repo.get_field("reference_temperature");
215 amrex::MultiFab::Copy(ref_theta_mf, temp0(lev), 0, 0, 1, 0);
216 } else {
217 ref_theta_mf.setVal(m_reference_temperature);
218 }
219 }
220
222 [[nodiscard]] std::unique_ptr<ScratchField> ref_theta() const override
223 {
224 if (m_reference_temperature < 0.0_rt) {
225 amrex::Abort("Reference temperature was not set");
226 }
227
228 auto ref_theta = m_repo.create_scratch_field(1, m_ngrow);
229 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
230 ref_theta_fill(lev, (*ref_theta)(lev));
231 }
232 return ref_theta;
233 }
234
235private:
238
240 amrex::Real m_mu{1.0e-5_rt};
241
243 amrex::Real m_Pr{1.0_rt};
244
246 amrex::Real m_Prt{1.0_rt};
247
249 amrex::Real m_constant_beta{0.0_rt};
250
252 amrex::Real m_reference_temperature{-1.0_rt};
253};
254
255} // namespace kynema_sgf::transport
256
257#endif /* CONSTTRANSPORT_H */
Definition CFDSim.H:55
Definition FieldRepo.H:86
void ref_theta_fill(const int lev, amrex::MultiFab &ref_theta_mf) const override
Definition ConstTransport.H:207
amrex::Real laminar_prandtl() const
Definition ConstTransport.H:101
static constexpr bool constant_properties
Definition ConstTransport.H:18
static amrex::Real turbulent_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:114
std::unique_ptr< ScratchField > ref_theta() const override
Return the reference temperature.
Definition ConstTransport.H:222
std::unique_ptr< ScratchField > mu() override
Return the dynamic visocity field.
Definition ConstTransport.H:124
static std::string identifier()
Definition ConstTransport.H:20
static amrex::Real laminar_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:105
amrex::Real m_mu
(Laminar) dynamic viscosity
Definition ConstTransport.H:240
amrex::Real m_constant_beta
Constant thermal expansion coefficient.
Definition ConstTransport.H:249
void beta_fill(const int lev, amrex::MultiFab &beta_mf) const override
Definition ConstTransport.H:173
ConstTransport(const CFDSim &sim)
Definition ConstTransport.H:22
std::unique_ptr< ScratchField > scalar_diffusivity(const std::string &scalar_name) override
Definition ConstTransport.H:145
FieldRepo & m_repo
Reference to the field repository (for creating scratch fields)
Definition ConstTransport.H:237
amrex::Real m_Prt
Turbulent Prandtl number.
Definition ConstTransport.H:246
amrex::Real m_reference_temperature
Reference temperature.
Definition ConstTransport.H:252
amrex::Real thermal_diffusivity() const
Definition ConstTransport.H:96
std::unique_ptr< ScratchField > beta() const override
Return the thermal expansion coefficient.
Definition ConstTransport.H:159
amrex::Real m_Pr
Prandtl number.
Definition ConstTransport.H:243
amrex::Real reference_temperature() const override
Definition ConstTransport.H:168
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field.
Definition ConstTransport.H:134
amrex::Real turbulent_prandtl() const
Definition ConstTransport.H:103
amrex::Real viscosity() const
Definition ConstTransport.H:94
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:17
Definition CFDSim.H:27