/home/runner/work/amr-wind/amr-wind/amr-wind/transport_models/ConstTransport.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/transport_models/ConstTransport.H Source File
AMR-Wind 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 amr_wind::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 << std::endl;
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 << std::endl;
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 << std::endl;
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 << std::endl;
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 << std::endl;
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 << std::endl;
88 pp_abl.get("reference_temperature", m_reference_temperature);
89 }
90 }
91
92 ~ConstTransport() override = default;
93
94 inline amrex::Real viscosity() const { return m_mu; }
95
96 inline amrex::Real thermal_diffusivity() const { return m_mu / m_Pr; }
97
98 inline amrex::Real laminar_prandtl() const { return m_Pr; }
99
100 inline amrex::Real turbulent_prandtl() const { return m_Prt; }
101
102 static inline amrex::Real laminar_schmidt(const std::string& scalar_name)
103 {
104 amrex::ParmParse pp("transport");
105 const std::string key = scalar_name + "_laminar_schmidt";
106 amrex::Real lam_schmidt = 1.0_rt;
107 pp.query(key, lam_schmidt);
108 return lam_schmidt;
109 }
110
111 static inline amrex::Real turbulent_schmidt(const std::string& scalar_name)
112 {
113 amrex::ParmParse pp("transport");
114 const std::string key = scalar_name + "_turbulent_schmidt";
115 amrex::Real turb_schmidt = 1.0_rt;
116 pp.query(key, turb_schmidt);
117 return turb_schmidt;
118 }
119
121 inline std::unique_ptr<ScratchField> mu() override
122 {
123 auto mu = m_repo.create_scratch_field(1, m_ngrow);
124 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
125 (*mu)(lev).setVal(m_mu);
126 }
127 return mu;
128 }
129
131 inline std::unique_ptr<ScratchField> alpha() override
132 {
133 auto alpha = mu();
134 amrex::Real inv_Pr = 1.0_rt / m_Pr;
135 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
136 (*alpha)(lev).mult(inv_Pr);
137 }
138 return alpha;
139 }
140
141 inline std::unique_ptr<ScratchField>
142 scalar_diffusivity(const std::string& scalar_name) override
143 {
144 amrex::Real lam_schmidt = laminar_schmidt(scalar_name);
145
146 amrex::Real inv_schmidt = 1.0_rt / lam_schmidt;
147 auto diff = mu();
148 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
149 (*diff)(lev).mult(inv_schmidt);
150 }
151
152 return diff;
153 }
154
156 inline std::unique_ptr<ScratchField> beta() const override
157 {
158 auto beta = m_repo.create_scratch_field(1, m_ngrow);
159 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
160#ifdef AMREX_USE_OMP
161#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
162#endif
163 for (amrex::MFIter mfi((*beta)(lev)); mfi.isValid(); ++mfi) {
164 const auto& bx = mfi.tilebox();
165 const auto& beta_arr = (*beta)(lev).array(mfi);
166 beta_impl(lev, mfi, bx, beta_arr);
167 }
168 }
169 return beta;
170 }
171
173 inline void beta_impl(
174 const int lev,
175 const amrex::MFIter& mfi,
176 const amrex::Box& bx,
177 const amrex::Array4<amrex::Real>& beta) const override
178 {
179
180 if (m_constant_beta > 0.0_rt) {
181 const amrex::Real beta_val = m_constant_beta;
182 amrex::ParallelFor(
183 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
184 beta(i, j, k) = beta_val;
185 });
186 } else if (m_repo.field_exists("reference_temperature")) {
187 const auto& temp0 = m_repo.get_field("reference_temperature");
188 const auto& temp0_arr = temp0(lev).const_array(mfi);
189 amrex::ParallelFor(
190 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
191 beta(i, j, k) = 1.0_rt / temp0_arr(i, j, k);
192 });
193 } else {
194 const amrex::Real beta_val = 1.0_rt / m_reference_temperature;
195 amrex::ParallelFor(
196 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
197 beta(i, j, k) = beta_val;
198 });
199 }
200
201 if (m_repo.field_exists("vof")) {
202 const auto& vof = m_repo.get_field("vof");
203 const auto& vof_arr = vof(lev).const_array(mfi);
204 amrex::ParallelFor(
205 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
206 if (vof_arr(i, j, k) > constants::TIGHT_TOL) {
207 beta(i, j, k) = 0.0_rt;
208 }
209 });
210 }
211 }
212
213 inline amrex::Real reference_temperature() const override
214 {
216 }
217
219 inline std::unique_ptr<ScratchField> ref_theta() const override
220 {
221 if (m_reference_temperature < 0.0_rt) {
222 amrex::Abort("Reference temperature was not set");
223 }
224
225 auto ref_theta = m_repo.create_scratch_field(1, m_ngrow);
226 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
227#ifdef AMREX_USE_OMP
228#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
229#endif
230 for (amrex::MFIter mfi((*ref_theta)(lev)); mfi.isValid(); ++mfi) {
231 const auto& bx = mfi.tilebox();
232 const auto& ref_theta_arr = (*ref_theta)(lev).array(mfi);
233 ref_theta_impl(lev, mfi, bx, ref_theta_arr);
234 }
235 }
236 return ref_theta;
237 }
238
240 inline void ref_theta_impl(
241 const int lev,
242 const amrex::MFIter& mfi,
243 const amrex::Box& bx,
244 const amrex::Array4<amrex::Real>& ref_theta) const override
245 {
246 if (m_reference_temperature < 0.0_rt) {
247 amrex::Abort("Reference temperature was not set");
248 }
249
250 if (m_repo.field_exists("reference_temperature")) {
251 auto& temp0 = m_repo.get_field("reference_temperature");
252 const auto& temp0_arr = temp0(lev).const_array(mfi);
253 amrex::ParallelFor(
254 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
255 ref_theta(i, j, k) = temp0_arr(i, j, k);
256 });
257 } else {
258 const amrex::Real ref_theta_val = m_reference_temperature;
259 amrex::ParallelFor(
260 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
261 ref_theta(i, j, k) = ref_theta_val;
262 });
263 }
264 }
265
266private:
269
271 amrex::Real m_mu{1.0e-5_rt};
272
274 amrex::Real m_Pr{1.0_rt};
275
277 amrex::Real m_Prt{1.0_rt};
278
280 amrex::Real m_constant_beta{0.0_rt};
281
283 amrex::Real m_reference_temperature{-1.0_rt};
284};
285
286} // namespace amr_wind::transport
287
288#endif /* CONSTTRANSPORT_H */
Definition CFDSim.H:54
Definition FieldRepo.H:86
void ref_theta_impl(const int lev, const amrex::MFIter &mfi, const amrex::Box &bx, const amrex::Array4< amrex::Real > &ref_theta) const override
Compute the reference temperature.
Definition ConstTransport.H:240
amrex::Real m_mu
(Laminar) dynamic viscosity
Definition ConstTransport.H:271
amrex::Real thermal_diffusivity() const
Definition ConstTransport.H:96
FieldRepo & m_repo
Reference to the field repository (for creating scratch fields)
Definition ConstTransport.H:268
amrex::Real reference_temperature() const override
Definition ConstTransport.H:213
amrex::Real turbulent_prandtl() const
Definition ConstTransport.H:100
void beta_impl(const int lev, const amrex::MFIter &mfi, const amrex::Box &bx, const amrex::Array4< amrex::Real > &beta) const override
Compute the thermal expansion coefficient.
Definition ConstTransport.H:173
amrex::Real m_Prt
Turbulent Prandtl number.
Definition ConstTransport.H:277
static constexpr bool constant_properties
Definition ConstTransport.H:18
amrex::Real viscosity() const
Definition ConstTransport.H:94
static amrex::Real turbulent_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:111
amrex::Real m_Pr
Prandtl number.
Definition ConstTransport.H:274
amrex::Real laminar_prandtl() const
Definition ConstTransport.H:98
static std::string identifier()
Definition ConstTransport.H:20
ConstTransport(const CFDSim &sim)
Definition ConstTransport.H:22
std::unique_ptr< ScratchField > beta() const override
Return the thermal expansion coefficient.
Definition ConstTransport.H:156
amrex::Real m_constant_beta
Constant thermal expansion coefficient.
Definition ConstTransport.H:280
static amrex::Real laminar_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:102
amrex::Real m_reference_temperature
Reference temperature.
Definition ConstTransport.H:283
std::unique_ptr< ScratchField > ref_theta() const override
Return the reference temperature.
Definition ConstTransport.H:219
std::unique_ptr< ScratchField > scalar_diffusivity(const std::string &scalar_name) override
Definition ConstTransport.H:142
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field.
Definition ConstTransport.H:131
std::unique_ptr< ScratchField > mu() override
Return the dynamic visocity field.
Definition ConstTransport.H:121
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:17
Definition CFDSim.H:26