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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/transport_models/TwoPhaseTransport.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
TwoPhaseTransport.H
Go to the documentation of this file.
1#ifndef TWOPHASETRANSPORT_H
2#define TWOPHASETRANSPORT_H
3
6#include "AMReX_ParmParse.H"
7#include "AMReX_REAL.H"
8
9using namespace amrex::literals;
10
11namespace amr_wind::transport {
12
16class TwoPhaseTransport : public TransportModel::Register<TwoPhaseTransport>
17{
18public:
19 static constexpr bool constant_properties = false;
20
21 static std::string identifier() { return "TwoPhaseTransport"; }
22
23 explicit TwoPhaseTransport(const CFDSim& sim)
24 : m_repo(sim.repo()), m_physics_mgr(sim.physics_manager())
25 {
26 amrex::ParmParse pp("transport");
27 pp.query("viscosity_fluid1", m_mu1);
28 pp.query("viscosity_fluid2", m_mu2);
29 pp.query("laminar_prandtl_fluid1", m_Pr1);
30 pp.query("laminar_prandtl_fluid2", m_Pr2);
31 pp.query("turbulent_prandtl", m_Prt);
32
33 // Check for single-phase quantities and warn
34 if (pp.contains("viscosity")) {
35 amrex::Print()
36 << "WARNING: single-phase viscosity has been specified but "
37 "will not be used! (TwoPhaseTransport)\n";
38 }
39 if (pp.contains("laminar_prandtl")) {
40 amrex::Print()
41 << "WARNING: single-phase laminar_prandtl has been specified "
42 "but will not be used! (TwoPhaseTransport)\n";
43 }
44
45 // Backwards compatibility
46 amrex::ParmParse pp_boussinesq_buoyancy("BoussinesqBuoyancy");
47 amrex::ParmParse pp_abl("ABL");
48 if (pp.contains("thermal_expansion_coefficient")) {
49 pp.get("thermal_expansion_coefficient", m_constant_beta);
50 if (pp_boussinesq_buoyancy.contains("thermal_expansion_coeff")) {
51 amrex::Print()
52 << "WARNING: BoussinesqBuoyancy.thermal_expansion_coeff "
53 "option has been deprecated in favor of "
54 "transport.thermal_expansion_coefficient. Ignoring the "
55 "BoussinesqBuoyancy option in favor of the transport "
56 "option."
57 << std::endl;
58 }
59 } else if (pp_boussinesq_buoyancy.contains("thermal_expansion_coeff")) {
60 amrex::Print()
61 << "WARNING: BoussinesqBuoyancy.thermal_expansion_coeff option "
62 "has been deprecated in favor of "
63 "transport.thermal_expansion_coefficient. Please replace "
64 "this option."
65 << std::endl;
66 pp_boussinesq_buoyancy.get(
67 "thermal_expansion_coeff", m_constant_beta);
68 }
69
70 if (pp.contains("reference_temperature")) {
71 pp.get("reference_temperature", m_reference_temperature);
72 if (pp_boussinesq_buoyancy.contains("reference_temperature")) {
73 amrex::Print()
74 << "WARNING: BoussinesqBuoyancy.reference_temperature "
75 "option has been deprecated in favor of "
76 "transport.reference_temperature. Ignoring the "
77 "BoussinesqBuoyancy option in favor of the transport "
78 "option."
79 << std::endl;
80 } else if (pp_abl.contains("reference_temperature")) {
81 amrex::Print()
82 << "WARNING: ABL.reference_temperature "
83 "option has been deprecated in favor of "
84 "transport.reference_temperature. Ignoring the "
85 "ABL option in favor of the transport "
86 "option."
87 << std::endl;
88 }
89 } else if (pp_boussinesq_buoyancy.contains("reference_temperature")) {
90 amrex::Print()
91 << "WARNING: BoussinesqBuoyancy.reference_temperature option "
92 "has been deprecated in favor of "
93 "transport.reference_temperature. Please replace "
94 "this option."
95 << std::endl;
96 pp_boussinesq_buoyancy.get(
97 "reference_temperature", m_reference_temperature);
98 } else if (pp_abl.contains("reference_temperature")) {
99 amrex::Print() << "WARNING: ABL.reference_temperature option "
100 "has been deprecated in favor of "
101 "transport.reference_temperature. Please replace "
102 "this option."
103 << std::endl;
104 pp_abl.get("reference_temperature", m_reference_temperature);
105 }
106 }
107
108 ~TwoPhaseTransport() override = default;
109
110 inline amrex::Real laminar_prandtl1() const { return m_Pr1; }
111 inline amrex::Real laminar_prandtl2() const { return m_Pr2; }
112
113 inline amrex::Real turbulent_prandtl() const { return m_Prt; }
114
115 static inline amrex::Real laminar_schmidt(const std::string& scalar_name)
116 {
117 amrex::ParmParse pp("transport");
118 const std::string key = scalar_name + "_laminar_schmidt";
119 amrex::Real lam_schmidt = 1.0_rt;
120 pp.query(key, lam_schmidt);
121 return lam_schmidt;
122 }
123
124 static inline amrex::Real turbulent_schmidt(const std::string& scalar_name)
125 {
126 amrex::ParmParse pp("transport");
127 const std::string key = scalar_name + "_turbulent_schmidt";
128 amrex::Real turb_schmidt = 1.0_rt;
129 pp.query(key, turb_schmidt);
130 return turb_schmidt;
131 }
132
134 inline std::unique_ptr<ScratchField> mu() override
135 {
136 // Select the interface capturing method
137 auto mu = m_repo.create_scratch_field(1, m_ngrow);
138
139 auto ifacetype = get_iface_method();
140 if (ifacetype == InterfaceCapturingMethod::VOF) {
141
142 const auto& vof = m_repo.get_field("vof");
143
144 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
145 const auto& volfrac_arrs = vof(lev).const_arrays();
146 const auto& visc_arrs = (*mu)(lev).arrays();
147 const amrex::Real mu1 = m_mu1;
148 const amrex::Real mu2 = m_mu2;
149 amrex::ParallelFor(
150 (*mu)(lev), mu->num_grow(),
151 [=] AMREX_GPU_DEVICE(
152 int nbx, int i, int j, int k) noexcept {
153 visc_arrs[nbx](i, j, k) =
154 mu1 * volfrac_arrs[nbx](i, j, k) +
155 mu2 * (1.0_rt - volfrac_arrs[nbx](i, j, k));
156 });
157 }
158 amrex::Gpu::streamSynchronize();
159
160 } else if (ifacetype == InterfaceCapturingMethod::LS) {
161
162 const auto& levelset = m_repo.get_field("levelset");
163 const auto& geom = m_repo.mesh().Geom();
164
165 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
166 const auto& dx = geom[lev].CellSizeArray();
167 const auto& visc_arrs = (*mu)(lev).arrays();
168 const auto& phi_arrs = levelset(lev).const_arrays();
169 const amrex::Real eps =
170 std::cbrt(2.0_rt * dx[0] * dx[1] * dx[2]);
171 const amrex::Real mu1 = m_mu1;
172 const amrex::Real mu2 = m_mu2;
173
174 amrex::ParallelFor(
175 (*mu)(lev), mu->num_grow(),
176 [=] AMREX_GPU_DEVICE(
177 int nbx, int i, int j, int k) noexcept {
178 amrex::Real smooth_heaviside;
179 if (phi_arrs[nbx](i, j, k) > eps) {
180 smooth_heaviside = 1.0_rt;
181 } else if (phi_arrs[nbx](i, j, k) < -eps) {
182 smooth_heaviside = 0.0_rt;
183 } else {
184 smooth_heaviside =
185 0.5_rt *
186 (1.0_rt + phi_arrs[nbx](i, j, k) / eps +
187 1.0_rt / static_cast<amrex::Real>(M_PI) *
188 std::sin(
189 phi_arrs[nbx](i, j, k) *
190 static_cast<amrex::Real>(M_PI) / eps));
191 }
192 visc_arrs[nbx](i, j, k) =
193 mu1 * smooth_heaviside +
194 mu2 * (1.0_rt - smooth_heaviside);
195 });
196 }
197 amrex::Gpu::streamSynchronize();
198 }
199 return mu;
200 }
201
203 inline std::unique_ptr<ScratchField> alpha() override
204 {
205 // Select the interface capturing method
206 auto alpha = m_repo.create_scratch_field(1, m_ngrow);
207
208 auto ifacetype = get_iface_method();
209 if (ifacetype == InterfaceCapturingMethod::VOF) {
210
211 const auto& vof = m_repo.get_field("vof");
212
213 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
214 const auto& volfrac_arrs = vof(lev).const_arrays();
215 const auto& thdiff_arrs = (*alpha)(lev).arrays();
216 const amrex::Real mu1 = m_mu1;
217 const amrex::Real mu2 = m_mu2;
218 const amrex::Real Pr1 = m_Pr1;
219 const amrex::Real Pr2 = m_Pr2;
220 amrex::ParallelFor(
221 (*alpha)(lev), alpha->num_grow(),
222 [=] AMREX_GPU_DEVICE(
223 int nbx, int i, int j, int k) noexcept {
224 thdiff_arrs[nbx](i, j, k) =
225 mu1 / Pr1 * volfrac_arrs[nbx](i, j, k) +
226 mu2 / Pr2 * (1.0_rt - volfrac_arrs[nbx](i, j, k));
227 });
228 }
229
230 } else if (ifacetype == InterfaceCapturingMethod::LS) {
231 const auto& levelset = m_repo.get_field("levelset");
232 const auto& geom = m_repo.mesh().Geom();
233
234 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
235 const auto& dx = geom[lev].CellSizeArray();
236 const auto& visc_arrs = (*alpha)(lev).arrays();
237 const auto& phi_arrs = levelset(lev).const_arrays();
238 const amrex::Real eps =
239 std::cbrt(2.0_rt * dx[0] * dx[1] * dx[2]);
240 const amrex::Real mu1 = m_mu1;
241 const amrex::Real mu2 = m_mu2;
242 const amrex::Real Pr1 = m_Pr1;
243 const amrex::Real Pr2 = m_Pr2;
244 amrex::ParallelFor(
245 (*alpha)(lev), alpha->num_grow(),
246 [=] AMREX_GPU_DEVICE(
247 int nbx, int i, int j, int k) noexcept {
248 amrex::Real smooth_heaviside;
249 if (phi_arrs[nbx](i, j, k) > eps) {
250 smooth_heaviside = 1.0_rt;
251 } else if (phi_arrs[nbx](i, j, k) < -eps) {
252 smooth_heaviside = 0.0_rt;
253 } else {
254 smooth_heaviside =
255 0.5_rt *
256 (1.0_rt + phi_arrs[nbx](i, j, k) / eps +
257 1.0_rt / static_cast<amrex::Real>(M_PI) *
258 std::sin(
259 phi_arrs[nbx](i, j, k) *
260 static_cast<amrex::Real>(M_PI) / eps));
261 }
262 visc_arrs[nbx](i, j, k) =
263 mu1 / Pr1 * smooth_heaviside +
264 mu2 / Pr2 * (1.0_rt - smooth_heaviside);
265 });
266 }
267 }
268 amrex::Gpu::streamSynchronize();
269 return alpha;
270 }
271
272 inline std::unique_ptr<ScratchField>
273 scalar_diffusivity(const std::string& scalar_name) override
274 {
275 amrex::Real lam_schmidt = laminar_schmidt(scalar_name);
276
277 amrex::Real inv_schmidt = 1.0_rt / lam_schmidt;
278 auto diff = mu();
279 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
280 (*diff)(lev).mult(inv_schmidt);
281 }
282
283 return diff;
284 }
285
287 inline std::unique_ptr<ScratchField> beta() const override
288 {
289 auto beta = m_repo.create_scratch_field(1, m_ngrow);
290 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
291#ifdef AMREX_USE_OMP
292#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
293#endif
294 for (amrex::MFIter mfi((*beta)(lev)); mfi.isValid(); ++mfi) {
295 const auto& bx = mfi.tilebox();
296 const auto& beta_arr = (*beta)(lev).array(mfi);
297 beta_impl(lev, mfi, bx, beta_arr);
298 }
299 }
300 return beta;
301 }
302
304 inline void beta_impl(
305 const int lev,
306 const amrex::MFIter& mfi,
307 const amrex::Box& bx,
308 const amrex::Array4<amrex::Real>& beta) const override
309 {
310
311 if (m_constant_beta > 0.0_rt) {
312 const amrex::Real beta_val = m_constant_beta;
313 amrex::ParallelFor(
314 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
315 beta(i, j, k) = beta_val;
316 });
317 } else if (m_repo.field_exists("reference_temperature")) {
318 const auto& temp0 = m_repo.get_field("reference_temperature");
319 const auto& temp0_arr = temp0(lev).const_array(mfi);
320 amrex::ParallelFor(
321 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
322 beta(i, j, k) = 1.0_rt / temp0_arr(i, j, k);
323 });
324 } else {
325 const amrex::Real beta_val = 1.0_rt / m_reference_temperature;
326 amrex::ParallelFor(
327 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
328 beta(i, j, k) = beta_val;
329 });
330 }
331
332 auto ifacetype = get_iface_method();
333 if (ifacetype == InterfaceCapturingMethod::VOF) {
334 const auto& vof = m_repo.get_field("vof");
335 const auto& vof_arr = vof(lev).const_array(mfi);
336 amrex::ParallelFor(
337 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
338 if (vof_arr(i, j, k) > constants::TIGHT_TOL) {
339 beta(i, j, k) = 0.0_rt;
340 }
341 });
342 } else if (ifacetype == InterfaceCapturingMethod::LS) {
343 const auto& levelset = m_repo.get_field("levelset");
344 const auto& dx = m_repo.mesh().Geom(lev).CellSizeArray();
345 const auto& phi_arr = levelset(lev).const_array(mfi);
346 const amrex::Real eps = std::cbrt(2.0_rt * dx[0] * dx[1] * dx[2]);
347 amrex::ParallelFor(
348 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
349 amrex::Real smooth_heaviside;
350 if (phi_arr(i, j, k) > eps) {
351 smooth_heaviside = 1.0_rt;
352 } else if (phi_arr(i, j, k) < -eps) {
353 smooth_heaviside = 0.0_rt;
354 } else {
355 smooth_heaviside =
356 0.5_rt *
357 (1.0_rt + phi_arr(i, j, k) / eps +
358 1.0_rt / static_cast<amrex::Real>(M_PI) *
359 std::sin(
360 phi_arr(i, j, k) *
361 static_cast<amrex::Real>(M_PI) / eps));
362 }
363 beta(i, j, k) *= (1.0_rt - smooth_heaviside);
364 });
365 }
366 }
367
368 inline amrex::Real reference_temperature() const override
369 {
371 }
372
374 inline std::unique_ptr<ScratchField> ref_theta() const override
375 {
376 if (m_reference_temperature < 0.0_rt) {
377 amrex::Abort("Reference temperature was not set");
378 }
379
380 auto ref_theta = m_repo.create_scratch_field(1, m_ngrow);
381 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
382#ifdef AMREX_USE_OMP
383#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
384#endif
385 for (amrex::MFIter mfi((*ref_theta)(lev)); mfi.isValid(); ++mfi) {
386 const auto& bx = mfi.tilebox();
387 const auto& ref_theta_arr = (*ref_theta)(lev).array(mfi);
388 ref_theta_impl(lev, mfi, bx, ref_theta_arr);
389 }
390 }
391 return ref_theta;
392 }
393
395 inline void ref_theta_impl(
396 const int lev,
397 const amrex::MFIter& mfi,
398 const amrex::Box& bx,
399 const amrex::Array4<amrex::Real>& ref_theta) const override
400 {
401 if (m_reference_temperature < 0.0_rt) {
402 amrex::Abort("Reference temperature was not set");
403 }
404
405 if (m_repo.field_exists("reference_temperature")) {
406 auto& temp0 = m_repo.get_field("reference_temperature");
407 const auto& temp0_arr = temp0(lev).const_array(mfi);
408 amrex::ParallelFor(
409 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
410 ref_theta(i, j, k) = temp0_arr(i, j, k);
411 });
412 } else {
413 const amrex::Real ref_theta_val = m_reference_temperature;
414 amrex::ParallelFor(
415 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
416 ref_theta(i, j, k) = ref_theta_val;
417 });
418 }
419 }
420
421private:
424
427
429 amrex::Real m_mu1{1.0e-3_rt};
430
432 amrex::Real m_mu2{1.0e-5_rt};
433
435 amrex::Real m_Pr1{7.2_rt};
436
438 amrex::Real m_Pr2{0.7_rt};
439
441 amrex::Real m_Prt{1.0_rt};
442
444 amrex::Real m_constant_beta{0.0_rt};
445
447 amrex::Real m_reference_temperature{-1.0_rt};
448
450 {
451 if (!m_physics_mgr.contains("MultiPhase")) {
452 amrex::Abort("TwoPhaseTransport requires MultiPhase physics");
453 }
454 const auto& multiphase = m_physics_mgr.get<MultiPhase>();
455 return multiphase.interface_capturing_method();
456 }
457};
458
459} // namespace amr_wind::transport
460
461#endif /* TWOPHASETRANSPORT_H */
Definition CFDSim.H:54
Definition FieldRepo.H:86
Definition MultiPhase.H:30
Definition Physics.H:100
amrex::Real laminar_prandtl1() const
Definition TwoPhaseTransport.H:110
std::unique_ptr< ScratchField > ref_theta() const override
Return the reference temperature.
Definition TwoPhaseTransport.H:374
amrex::Real reference_temperature() const override
Definition TwoPhaseTransport.H:368
amrex::Real m_Prt
Turbulent Prandtl number.
Definition TwoPhaseTransport.H:441
std::unique_ptr< ScratchField > beta() const override
Return the thermal expansion coefficient.
Definition TwoPhaseTransport.H:287
std::unique_ptr< ScratchField > mu() override
Return the dynamic visocity field.
Definition TwoPhaseTransport.H:134
std::unique_ptr< ScratchField > scalar_diffusivity(const std::string &scalar_name) override
Definition TwoPhaseTransport.H:273
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 TwoPhaseTransport.H:395
amrex::Real turbulent_prandtl() const
Definition TwoPhaseTransport.H:113
InterfaceCapturingMethod get_iface_method() const
Definition TwoPhaseTransport.H:449
amrex::Real m_reference_temperature
Reference temperature.
Definition TwoPhaseTransport.H:447
amrex::Real m_Pr1
Phase 1 (liquid) Prandtl number.
Definition TwoPhaseTransport.H:435
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 TwoPhaseTransport.H:304
static constexpr bool constant_properties
Definition TwoPhaseTransport.H:19
TwoPhaseTransport(const CFDSim &sim)
Definition TwoPhaseTransport.H:23
amrex::Real m_constant_beta
Constant thermal expansion coefficient.
Definition TwoPhaseTransport.H:444
amrex::Real laminar_prandtl2() const
Definition TwoPhaseTransport.H:111
amrex::Real m_mu1
Phase 1 (liquid) dynamic molecular viscosity.
Definition TwoPhaseTransport.H:429
FieldRepo & m_repo
Reference to the field repository (for creating scratch fields)
Definition TwoPhaseTransport.H:423
static std::string identifier()
Definition TwoPhaseTransport.H:21
amrex::Real m_mu2
Phase 2 (gas) dynamic molecular viscosity.
Definition TwoPhaseTransport.H:432
static amrex::Real turbulent_schmidt(const std::string &scalar_name)
Definition TwoPhaseTransport.H:124
amrex::Real m_Pr2
Phase 2 (gas) Prandtl number.
Definition TwoPhaseTransport.H:438
static amrex::Real laminar_schmidt(const std::string &scalar_name)
Definition TwoPhaseTransport.H:115
const PhysicsMgr & m_physics_mgr
Reference to the physics manager.
Definition TwoPhaseTransport.H:426
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field (later divided by density, though)
Definition TwoPhaseTransport.H:203
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:17
Definition height_functions.H:8
Definition CFDSim.H:26
InterfaceCapturingMethod
Definition MultiPhase.H:24
@ VOF
Volume of fluid.
Definition MultiPhase.H:25
@ LS
Levelset.
Definition MultiPhase.H:26