1#ifndef TWOPHASETRANSPORT_H
2#define TWOPHASETRANSPORT_H
6#include "AMReX_ParmParse.H"
9using namespace amrex::literals;
21 static std::string
identifier() {
return "TwoPhaseTransport"; }
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);
34 if (pp.contains(
"viscosity")) {
36 <<
"WARNING: single-phase viscosity has been specified but "
37 "will not be used! (TwoPhaseTransport)\n";
39 if (pp.contains(
"laminar_prandtl")) {
41 <<
"WARNING: single-phase laminar_prandtl has been specified "
42 "but will not be used! (TwoPhaseTransport)\n";
46 amrex::ParmParse pp_boussinesq_buoyancy(
"BoussinesqBuoyancy");
47 amrex::ParmParse pp_abl(
"ABL");
48 if (pp.contains(
"thermal_expansion_coefficient")) {
50 if (pp_boussinesq_buoyancy.contains(
"thermal_expansion_coeff")) {
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 "
59 }
else if (pp_boussinesq_buoyancy.contains(
"thermal_expansion_coeff")) {
61 <<
"WARNING: BoussinesqBuoyancy.thermal_expansion_coeff option "
62 "has been deprecated in favor of "
63 "transport.thermal_expansion_coefficient. Please replace "
66 pp_boussinesq_buoyancy.get(
70 if (pp.contains(
"reference_temperature")) {
72 if (pp_boussinesq_buoyancy.contains(
"reference_temperature")) {
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 "
80 }
else if (pp_abl.contains(
"reference_temperature")) {
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 "
89 }
else if (pp_boussinesq_buoyancy.contains(
"reference_temperature")) {
91 <<
"WARNING: BoussinesqBuoyancy.reference_temperature option "
92 "has been deprecated in favor of "
93 "transport.reference_temperature. Please replace "
96 pp_boussinesq_buoyancy.get(
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 "
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);
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);
134 inline std::unique_ptr<ScratchField>
mu()
override
137 auto mu =
m_repo.create_scratch_field(1, m_ngrow);
142 const auto& vof =
m_repo.get_field(
"vof");
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;
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));
158 amrex::Gpu::streamSynchronize();
162 const auto& levelset =
m_repo.get_field(
"levelset");
163 const auto& geom =
m_repo.mesh().Geom();
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;
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;
186 (1.0_rt + phi_arrs[nbx](i, j, k) / eps +
187 1.0_rt /
static_cast<amrex::Real
>(M_PI) *
189 phi_arrs[nbx](i, j, k) *
190 static_cast<amrex::Real
>(M_PI) / eps));
192 visc_arrs[nbx](i, j, k) =
193 mu1 * smooth_heaviside +
194 mu2 * (1.0_rt - smooth_heaviside);
197 amrex::Gpu::streamSynchronize();
203 inline std::unique_ptr<ScratchField>
alpha()
override
206 auto alpha =
m_repo.create_scratch_field(1, m_ngrow);
211 const auto& vof =
m_repo.get_field(
"vof");
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;
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));
231 const auto& levelset =
m_repo.get_field(
"levelset");
232 const auto& geom =
m_repo.mesh().Geom();
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;
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;
256 (1.0_rt + phi_arrs[nbx](i, j, k) / eps +
257 1.0_rt /
static_cast<amrex::Real
>(M_PI) *
259 phi_arrs[nbx](i, j, k) *
260 static_cast<amrex::Real
>(M_PI) / eps));
262 visc_arrs[nbx](i, j, k) =
263 mu1 / Pr1 * smooth_heaviside +
264 mu2 / Pr2 * (1.0_rt - smooth_heaviside);
268 amrex::Gpu::streamSynchronize();
272 inline std::unique_ptr<ScratchField>
277 amrex::Real inv_schmidt = 1.0_rt / lam_schmidt;
279 for (
int lev = 0; lev <
m_repo.num_active_levels(); ++lev) {
280 (*diff)(lev).mult(inv_schmidt);
287 inline std::unique_ptr<ScratchField>
beta()
const override
289 auto beta =
m_repo.create_scratch_field(1, m_ngrow);
290 for (
int lev = 0; lev <
m_repo.num_active_levels(); ++lev) {
292#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
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);
306 const amrex::MFIter& mfi,
307 const amrex::Box& bx,
308 const amrex::Array4<amrex::Real>&
beta)
const override
314 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
315 beta(i, j, k) = beta_val;
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);
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);
327 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
328 beta(i, j, k) = beta_val;
334 const auto& vof =
m_repo.get_field(
"vof");
335 const auto& vof_arr = vof(lev).const_array(mfi);
337 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
339 beta(i, j, k) = 0.0_rt;
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]);
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;
357 (1.0_rt + phi_arr(i, j, k) / eps +
358 1.0_rt /
static_cast<amrex::Real
>(M_PI) *
361 static_cast<amrex::Real
>(M_PI) / eps));
363 beta(i, j, k) *= (1.0_rt - smooth_heaviside);
374 inline std::unique_ptr<ScratchField>
ref_theta()
const override
377 amrex::Abort(
"Reference temperature was not set");
381 for (
int lev = 0; lev <
m_repo.num_active_levels(); ++lev) {
383#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
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);
397 const amrex::MFIter& mfi,
398 const amrex::Box& bx,
399 const amrex::Array4<amrex::Real>&
ref_theta)
const override
402 amrex::Abort(
"Reference temperature was not set");
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);
409 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
415 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
452 amrex::Abort(
"TwoPhaseTransport requires MultiPhase physics");
455 return multiphase.interface_capturing_method();
Definition FieldRepo.H:86
Definition MultiPhase.H:30
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
~TwoPhaseTransport() override=default
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
InterfaceCapturingMethod
Definition MultiPhase.H:24
@ VOF
Volume of fluid.
Definition MultiPhase.H:25
@ LS
Levelset.
Definition MultiPhase.H:26