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

Kynema-SGF API: /home/runner/work/kynema-sgf/kynema-sgf/src/transport_models/TwoPhaseTransport.H Source File
Kynema-SGF 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
4#include <numbers>
7#include "AMReX_ParmParse.H"
8#include "AMReX_REAL.H"
9
10using namespace amrex::literals;
11
12namespace kynema_sgf::transport {
13
17class TwoPhaseTransport : public TransportModel::Register<TwoPhaseTransport>
18{
19public:
20 static constexpr bool constant_properties = false;
21
22 static std::string identifier() { return "TwoPhaseTransport"; }
23
24 explicit TwoPhaseTransport(const CFDSim& sim)
25 : m_repo(sim.repo()), m_physics_mgr(sim.physics_manager())
26 {
27 amrex::ParmParse pp("transport");
28 pp.query("viscosity_fluid1", m_mu1);
29 pp.query("viscosity_fluid2", m_mu2);
30 pp.query("laminar_prandtl_fluid1", m_Pr1);
31 pp.query("laminar_prandtl_fluid2", m_Pr2);
32 pp.query("turbulent_prandtl", m_Prt);
33
34 // Check for single-phase quantities and warn
35 if (pp.contains("viscosity")) {
36 amrex::Print()
37 << "WARNING: single-phase viscosity has been specified but "
38 "will not be used! (TwoPhaseTransport)\n";
39 }
40 if (pp.contains("laminar_prandtl")) {
41 amrex::Print()
42 << "WARNING: single-phase laminar_prandtl has been specified "
43 "but will not be used! (TwoPhaseTransport)\n";
44 }
45
46 // Backwards compatibility
47 amrex::ParmParse pp_boussinesq_buoyancy("BoussinesqBuoyancy");
48 amrex::ParmParse pp_abl("ABL");
49 if (pp.contains("thermal_expansion_coefficient")) {
50 pp.get("thermal_expansion_coefficient", m_constant_beta);
51 if (pp_boussinesq_buoyancy.contains("thermal_expansion_coeff")) {
52 amrex::Print()
53 << "WARNING: BoussinesqBuoyancy.thermal_expansion_coeff "
54 "option has been deprecated in favor of "
55 "transport.thermal_expansion_coefficient. Ignoring the "
56 "BoussinesqBuoyancy option in favor of the transport "
57 "option."
58 << '\n';
59 }
60 } else if (pp_boussinesq_buoyancy.contains("thermal_expansion_coeff")) {
61 amrex::Print()
62 << "WARNING: BoussinesqBuoyancy.thermal_expansion_coeff option "
63 "has been deprecated in favor of "
64 "transport.thermal_expansion_coefficient. Please replace "
65 "this option."
66 << '\n';
67 pp_boussinesq_buoyancy.get(
68 "thermal_expansion_coeff", m_constant_beta);
69 }
70
71 if (pp.contains("reference_temperature")) {
72 pp.get("reference_temperature", m_reference_temperature);
73 if (pp_boussinesq_buoyancy.contains("reference_temperature")) {
74 amrex::Print()
75 << "WARNING: BoussinesqBuoyancy.reference_temperature "
76 "option has been deprecated in favor of "
77 "transport.reference_temperature. Ignoring the "
78 "BoussinesqBuoyancy option in favor of the transport "
79 "option."
80 << '\n';
81 } else if (pp_abl.contains("reference_temperature")) {
82 amrex::Print()
83 << "WARNING: ABL.reference_temperature "
84 "option has been deprecated in favor of "
85 "transport.reference_temperature. Ignoring the "
86 "ABL option in favor of the transport "
87 "option."
88 << '\n';
89 }
90 } else if (pp_boussinesq_buoyancy.contains("reference_temperature")) {
91 amrex::Print()
92 << "WARNING: BoussinesqBuoyancy.reference_temperature option "
93 "has been deprecated in favor of "
94 "transport.reference_temperature. Please replace "
95 "this option."
96 << '\n';
97 pp_boussinesq_buoyancy.get(
98 "reference_temperature", m_reference_temperature);
99 } else if (pp_abl.contains("reference_temperature")) {
100 amrex::Print() << "WARNING: ABL.reference_temperature option "
101 "has been deprecated in favor of "
102 "transport.reference_temperature. Please replace "
103 "this option."
104 << '\n';
105 pp_abl.get("reference_temperature", m_reference_temperature);
106 }
107 }
108
109 ~TwoPhaseTransport() override = default;
110
111 [[nodiscard]] amrex::Real laminar_prandtl1() const { return m_Pr1; }
112 [[nodiscard]] amrex::Real laminar_prandtl2() const { return m_Pr2; }
113
114 [[nodiscard]] amrex::Real turbulent_prandtl() const { return m_Prt; }
115
116 static amrex::Real laminar_schmidt(const std::string& scalar_name)
117 {
118 amrex::ParmParse pp("transport");
119 const std::string key = scalar_name + "_laminar_schmidt";
120 amrex::Real lam_schmidt = 1.0_rt;
121 pp.query(key, lam_schmidt);
122 return lam_schmidt;
123 }
124
125 static amrex::Real turbulent_schmidt(const std::string& scalar_name)
126 {
127 amrex::ParmParse pp("transport");
128 const std::string key = scalar_name + "_turbulent_schmidt";
129 amrex::Real turb_schmidt = 1.0_rt;
130 pp.query(key, turb_schmidt);
131 return turb_schmidt;
132 }
133
135 std::unique_ptr<ScratchField> mu() override
136 {
137 // Select the interface capturing method
138 auto mu = m_repo.create_scratch_field(1, m_ngrow);
139
140 auto ifacetype = get_iface_method();
141 if (ifacetype == InterfaceCapturingMethod::VOF) {
142
143 const auto& vof = m_repo.get_field("vof");
144
145 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
146 const auto& volfrac_arrs = vof(lev).const_arrays();
147 const auto& visc_arrs = (*mu)(lev).arrays();
148 const amrex::Real mu1 = m_mu1;
149 const amrex::Real mu2 = m_mu2;
150 amrex::ParallelFor(
151 (*mu)(lev), mu->num_grow(),
152 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) {
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(int nbx, int i, int j, int k) {
177 amrex::Real smooth_heaviside;
178 if (phi_arrs[nbx](i, j, k) > eps) {
179 smooth_heaviside = 1.0_rt;
180 } else if (phi_arrs[nbx](i, j, k) < -eps) {
181 smooth_heaviside = 0.0_rt;
182 } else {
183 smooth_heaviside =
184 0.5_rt *
185 (1.0_rt + phi_arrs[nbx](i, j, k) / eps +
186 1.0_rt / std::numbers::pi_v<amrex::Real> *
187 std::sin(
188 phi_arrs[nbx](i, j, k) *
189 std::numbers::pi_v<amrex::Real> /
190 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 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(int nbx, int i, int j, int k) {
223 thdiff_arrs[nbx](i, j, k) =
224 (mu1 / Pr1 * volfrac_arrs[nbx](i, j, k)) +
225 (mu2 / Pr2 * (1.0_rt - volfrac_arrs[nbx](i, j, k)));
226 });
227 }
228
229 } else if (ifacetype == InterfaceCapturingMethod::LS) {
230 const auto& levelset = m_repo.get_field("levelset");
231 const auto& geom = m_repo.mesh().Geom();
232
233 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
234 const auto& dx = geom[lev].CellSizeArray();
235 const auto& visc_arrs = (*alpha)(lev).arrays();
236 const auto& phi_arrs = levelset(lev).const_arrays();
237 const amrex::Real eps =
238 std::cbrt(2.0_rt * dx[0] * dx[1] * dx[2]);
239 const amrex::Real mu1 = m_mu1;
240 const amrex::Real mu2 = m_mu2;
241 const amrex::Real Pr1 = m_Pr1;
242 const amrex::Real Pr2 = m_Pr2;
243 amrex::ParallelFor(
244 (*alpha)(lev), alpha->num_grow(),
245 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) {
246 amrex::Real smooth_heaviside;
247 if (phi_arrs[nbx](i, j, k) > eps) {
248 smooth_heaviside = 1.0_rt;
249 } else if (phi_arrs[nbx](i, j, k) < -eps) {
250 smooth_heaviside = 0.0_rt;
251 } else {
252 smooth_heaviside =
253 0.5_rt *
254 (1.0_rt + phi_arrs[nbx](i, j, k) / eps +
255 1.0_rt / std::numbers::pi_v<amrex::Real> *
256 std::sin(
257 phi_arrs[nbx](i, j, k) *
258 std::numbers::pi_v<amrex::Real> /
259 eps));
260 }
261 visc_arrs[nbx](i, j, k) =
262 (mu1 / Pr1 * smooth_heaviside) +
263 (mu2 / Pr2 * (1.0_rt - smooth_heaviside));
264 });
265 }
266 }
267 amrex::Gpu::streamSynchronize();
268 return alpha;
269 }
270
271 std::unique_ptr<ScratchField>
272 scalar_diffusivity(const std::string& scalar_name) override
273 {
274 amrex::Real lam_schmidt = laminar_schmidt(scalar_name);
275
276 amrex::Real inv_schmidt = 1.0_rt / lam_schmidt;
277 auto diff = mu();
278 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
279 (*diff)(lev).mult(inv_schmidt);
280 }
281
282 return diff;
283 }
284
286 [[nodiscard]] std::unique_ptr<ScratchField> beta() const override
287 {
288 auto beta = m_repo.create_scratch_field(1, m_ngrow);
289 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
290 beta_fill(lev, (*beta)(lev));
291 }
292 return beta;
293 }
294
295 [[nodiscard]] amrex::Real reference_temperature() const override
296 {
298 }
299
300 void beta_fill(const int lev, amrex::MultiFab& beta_mf) const override
301 {
302 if (m_constant_beta > 0.0_rt) {
303 beta_mf.setVal(m_constant_beta);
304 } else if (m_repo.field_exists("reference_temperature")) {
305 const auto& temp0 = m_repo.get_field("reference_temperature");
306 auto const& beta_arrs = beta_mf.arrays();
307 auto const& temp0_arrs = temp0(lev).const_arrays();
308 amrex::ParallelFor(
309 beta_mf, amrex::IntVect(0), 1,
310 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int) {
311 beta_arrs[nbx](i, j, k) = 1.0_rt / temp0_arrs[nbx](i, j, k);
312 });
313 amrex::Gpu::streamSynchronize();
314 } else {
315 beta_mf.setVal(1.0_rt / m_reference_temperature);
316 }
317
318 auto ifacetype = get_iface_method();
319 if (ifacetype == InterfaceCapturingMethod::VOF) {
320 const auto& vof = m_repo.get_field("vof");
321 auto const& beta_arrs = beta_mf.arrays();
322 auto const& vof_arrs = vof(lev).const_arrays();
323 amrex::ParallelFor(
324 beta_mf, amrex::IntVect(0), 1,
325 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int) {
326 if (vof_arrs[nbx](i, j, k) > constants::TIGHT_TOL) {
327 beta_arrs[nbx](i, j, k) = 0.0_rt;
328 }
329 });
330 amrex::Gpu::streamSynchronize();
331 } else if (ifacetype == InterfaceCapturingMethod::LS) {
332 const auto& levelset = m_repo.get_field("levelset");
333 const auto& dx = m_repo.mesh().Geom(lev).CellSizeArray();
334 const amrex::Real eps = std::cbrt(2.0_rt * dx[0] * dx[1] * dx[2]);
335 auto const& beta_arrs = beta_mf.arrays();
336 auto const& phi_arrs = levelset(lev).const_arrays();
337 amrex::ParallelFor(
338 beta_mf, amrex::IntVect(0), 1,
339 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int) {
340 amrex::Real smooth_heaviside;
341 if (phi_arrs[nbx](i, j, k) > eps) {
342 smooth_heaviside = 1.0_rt;
343 } else if (phi_arrs[nbx](i, j, k) < -eps) {
344 smooth_heaviside = 0.0_rt;
345 } else {
346 smooth_heaviside =
347 0.5_rt *
348 (1.0_rt + phi_arrs[nbx](i, j, k) / eps +
349 1.0_rt / std::numbers::pi_v<amrex::Real> *
350 std::sin(
351 phi_arrs[nbx](i, j, k) *
352 std::numbers::pi_v<amrex::Real> / eps));
353 }
354 beta_arrs[nbx](i, j, k) *= (1.0_rt - smooth_heaviside);
355 });
356 amrex::Gpu::streamSynchronize();
357 }
358 }
359
360 void
361 ref_theta_fill(const int lev, amrex::MultiFab& ref_theta_mf) const override
362 {
363 if (m_reference_temperature < 0.0_rt) {
364 amrex::Abort("Reference temperature was not set");
365 }
366
367 if (m_repo.field_exists("reference_temperature")) {
368 const auto& temp0 = m_repo.get_field("reference_temperature");
369 amrex::MultiFab::Copy(ref_theta_mf, temp0(lev), 0, 0, 1, 0);
370 } else {
371 ref_theta_mf.setVal(m_reference_temperature);
372 }
373 }
374
376 [[nodiscard]] std::unique_ptr<ScratchField> ref_theta() const override
377 {
378 if (m_reference_temperature < 0.0_rt) {
379 amrex::Abort("Reference temperature was not set");
380 }
381
382 auto ref_theta = m_repo.create_scratch_field(1, m_ngrow);
383 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
384 ref_theta_fill(lev, (*ref_theta)(lev));
385 }
386 return ref_theta;
387 }
388
389private:
392
395
397 amrex::Real m_mu1{1.0e-3_rt};
398
400 amrex::Real m_mu2{1.0e-5_rt};
401
403 amrex::Real m_Pr1{7.2_rt};
404
406 amrex::Real m_Pr2{0.7_rt};
407
409 amrex::Real m_Prt{1.0_rt};
410
412 amrex::Real m_constant_beta{0.0_rt};
413
415 amrex::Real m_reference_temperature{-1.0_rt};
416
418 {
419 if (!m_physics_mgr.contains("MultiPhase")) {
420 amrex::Abort("TwoPhaseTransport requires MultiPhase physics");
421 }
422 const auto& multiphase = m_physics_mgr.get<MultiPhase>();
423 return multiphase.interface_capturing_method();
424 }
425};
426
427} // namespace kynema_sgf::transport
428
429#endif /* TWOPHASETRANSPORT_H */
Definition CFDSim.H:55
Definition FieldRepo.H:86
Definition MultiPhase.H:31
Definition Physics.H:100
static amrex::Real turbulent_schmidt(const std::string &scalar_name)
Definition TwoPhaseTransport.H:125
InterfaceCapturingMethod get_iface_method() const
Definition TwoPhaseTransport.H:417
amrex::Real reference_temperature() const override
Definition TwoPhaseTransport.H:295
std::unique_ptr< ScratchField > beta() const override
Return the thermal expansion coefficient.
Definition TwoPhaseTransport.H:286
amrex::Real laminar_prandtl1() const
Definition TwoPhaseTransport.H:111
amrex::Real m_Pr1
Phase 1 (liquid) Prandtl number.
Definition TwoPhaseTransport.H:403
void beta_fill(const int lev, amrex::MultiFab &beta_mf) const override
Definition TwoPhaseTransport.H:300
std::unique_ptr< ScratchField > mu() override
Return the dynamic visocity field.
Definition TwoPhaseTransport.H:135
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field (later divided by density, though)
Definition TwoPhaseTransport.H:203
FieldRepo & m_repo
Reference to the field repository (for creating scratch fields)
Definition TwoPhaseTransport.H:391
TwoPhaseTransport(const CFDSim &sim)
Definition TwoPhaseTransport.H:24
const PhysicsMgr & m_physics_mgr
Reference to the physics manager.
Definition TwoPhaseTransport.H:394
amrex::Real m_constant_beta
Constant thermal expansion coefficient.
Definition TwoPhaseTransport.H:412
amrex::Real turbulent_prandtl() const
Definition TwoPhaseTransport.H:114
amrex::Real m_reference_temperature
Reference temperature.
Definition TwoPhaseTransport.H:415
std::unique_ptr< ScratchField > ref_theta() const override
Return the reference temperature.
Definition TwoPhaseTransport.H:376
static std::string identifier()
Definition TwoPhaseTransport.H:22
amrex::Real m_mu2
Phase 2 (gas) dynamic molecular viscosity.
Definition TwoPhaseTransport.H:400
static amrex::Real laminar_schmidt(const std::string &scalar_name)
Definition TwoPhaseTransport.H:116
void ref_theta_fill(const int lev, amrex::MultiFab &ref_theta_mf) const override
Definition TwoPhaseTransport.H:361
amrex::Real m_mu1
Phase 1 (liquid) dynamic molecular viscosity.
Definition TwoPhaseTransport.H:397
std::unique_ptr< ScratchField > scalar_diffusivity(const std::string &scalar_name) override
Definition TwoPhaseTransport.H:272
static constexpr bool constant_properties
Definition TwoPhaseTransport.H:20
amrex::Real m_Pr2
Phase 2 (gas) Prandtl number.
Definition TwoPhaseTransport.H:406
amrex::Real laminar_prandtl2() const
Definition TwoPhaseTransport.H:112
amrex::Real m_Prt
Turbulent Prandtl number.
Definition TwoPhaseTransport.H:409
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:17
Definition height_functions.H:8
Definition CFDSim.H:27
InterfaceCapturingMethod
Definition MultiPhase.H:25
@ VOF
Volume of fluid.
Definition MultiPhase.H:26
@ LS
Levelset.
Definition MultiPhase.H:27