1#ifndef TWOPHASETRANSPORT_H
2#define TWOPHASETRANSPORT_H
6#include "AMReX_ParmParse.H"
18 static std::string
identifier() {
return "TwoPhaseTransport"; }
23 amrex::ParmParse pp(
"transport");
24 pp.query(
"viscosity_fluid1",
m_mu1);
25 pp.query(
"viscosity_fluid2",
m_mu2);
26 pp.query(
"laminar_prandtl_fluid1",
m_Pr1);
27 pp.query(
"laminar_prandtl_fluid2",
m_Pr2);
28 pp.query(
"turbulent_prandtl",
m_Prt);
31 if (pp.contains(
"viscosity")) {
33 <<
"WARNING: single-phase viscosity has been specified but "
34 "will not be used! (TwoPhaseTransport)\n";
36 if (pp.contains(
"laminar_prandtl")) {
38 <<
"WARNING: single-phase laminar_prandtl has been specified "
39 "but will not be used! (TwoPhaseTransport)\n";
43 amrex::ParmParse pp_boussinesq_buoyancy(
"BoussinesqBuoyancy");
44 amrex::ParmParse pp_abl(
"ABL");
45 if (pp.contains(
"thermal_expansion_coefficient")) {
47 if (pp_boussinesq_buoyancy.contains(
"thermal_expansion_coeff")) {
49 <<
"WARNING: BoussinesqBuoyancy.thermal_expansion_coeff "
50 "option has been deprecated in favor of "
51 "transport.thermal_expansion_coefficient. Ignoring the "
52 "BoussinesqBuoyancy option in favor of the transport "
56 }
else if (pp_boussinesq_buoyancy.contains(
"thermal_expansion_coeff")) {
58 <<
"WARNING: BoussinesqBuoyancy.thermal_expansion_coeff option "
59 "has been deprecated in favor of "
60 "transport.thermal_expansion_coefficient. Please replace "
63 pp_boussinesq_buoyancy.get(
67 if (pp.contains(
"reference_temperature")) {
69 if (pp_boussinesq_buoyancy.contains(
"reference_temperature")) {
71 <<
"WARNING: BoussinesqBuoyancy.reference_temperature "
72 "option has been deprecated in favor of "
73 "transport.reference_temperature. Ignoring the "
74 "BoussinesqBuoyancy option in favor of the transport "
77 }
else if (pp_abl.contains(
"reference_temperature")) {
79 <<
"WARNING: ABL.reference_temperature "
80 "option has been deprecated in favor of "
81 "transport.reference_temperature. Ignoring the "
82 "ABL option in favor of the transport "
86 }
else if (pp_boussinesq_buoyancy.contains(
"reference_temperature")) {
88 <<
"WARNING: BoussinesqBuoyancy.reference_temperature option "
89 "has been deprecated in favor of "
90 "transport.reference_temperature. Please replace "
93 pp_boussinesq_buoyancy.get(
95 }
else if (pp_abl.contains(
"reference_temperature")) {
96 amrex::Print() <<
"WARNING: ABL.reference_temperature option "
97 "has been deprecated in favor of "
98 "transport.reference_temperature. Please replace "
114 amrex::ParmParse pp(
"transport");
115 const std::string key = scalar_name +
"_laminar_schmidt";
116 amrex::Real lam_schmidt = 1.0;
117 pp.query(key.c_str(), lam_schmidt);
123 amrex::ParmParse pp(
"transport");
124 const std::string key = scalar_name +
"_turbulent_schmidt";
125 amrex::Real turb_schmidt = 1.0;
126 pp.query(key.c_str(), turb_schmidt);
131 inline std::unique_ptr<ScratchField>
mu()
override
142 for (amrex::MFIter mfi((*
mu)(lev)); mfi.isValid(); ++mfi) {
143 const auto& vbx = mfi.growntilebox();
144 const auto& volfrac = vof(lev).const_array(mfi);
145 const auto& visc = (*mu)(lev).array(mfi);
146 const amrex::Real mu1 =
m_mu1;
147 const amrex::Real mu2 =
m_mu2;
150 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
151 visc(i, j, k) = mu1 * volfrac(i, j, k) +
152 mu2 * (1.0 - volfrac(i, j, k));
163 for (amrex::MFIter mfi((*
mu)(lev)); mfi.isValid(); ++mfi) {
164 const auto& vbx = mfi.growntilebox();
165 const auto& dx = geom[lev].CellSizeArray();
166 const auto& visc = (*mu)(lev).array(mfi);
167 const auto& phi = levelset(lev).const_array(mfi);
168 const amrex::Real eps =
169 std::cbrt(2. * dx[0] * dx[1] * dx[2]);
170 const amrex::Real mu1 =
m_mu1;
171 const amrex::Real mu2 =
m_mu2;
174 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
175 amrex::Real smooth_heaviside;
176 if (phi(i, j, k) > eps) {
177 smooth_heaviside = 1.0;
178 }
else if (phi(i, j, k) < -eps) {
179 smooth_heaviside = 0.;
183 (1.0 + phi(i, j, k) / eps +
185 std::sin(phi(i, j, k) * M_PI / eps));
187 visc(i, j, k) = mu1 * smooth_heaviside +
188 mu2 * (1.0 - smooth_heaviside);
197 inline std::unique_ptr<ScratchField>
alpha()
override
208 for (amrex::MFIter mfi((*
alpha)(lev)); mfi.isValid(); ++mfi) {
209 const auto& vbx = mfi.growntilebox();
210 const auto& volfrac = vof(lev).const_array(mfi);
211 const auto& thdiff = (*alpha)(lev).array(mfi);
212 const amrex::Real mu1 =
m_mu1;
213 const amrex::Real mu2 =
m_mu2;
214 const amrex::Real Pr1 =
m_Pr1;
215 const amrex::Real Pr2 =
m_Pr2;
218 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
220 mu1 / Pr1 * volfrac(i, j, k) +
221 mu2 / Pr2 * (1.0 - volfrac(i, j, k));
232 for (amrex::MFIter mfi((*
alpha)(lev)); mfi.isValid(); ++mfi) {
233 const auto& vbx = mfi.growntilebox();
234 const auto& dx = geom[lev].CellSizeArray();
235 const auto& visc = (*alpha)(lev).array(mfi);
236 const auto& phi = levelset(lev).const_array(mfi);
237 const amrex::Real eps =
238 std::cbrt(2. * 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;
245 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
246 amrex::Real smooth_heaviside;
247 if (phi(i, j, k) > eps) {
248 smooth_heaviside = 1.0;
249 }
else if (phi(i, j, k) < -eps) {
250 smooth_heaviside = 0.;
254 (1.0 + phi(i, j, k) / eps +
256 std::sin(phi(i, j, k) * M_PI / eps));
259 mu1 / Pr1 * smooth_heaviside +
260 mu2 / Pr2 * (1.0 - smooth_heaviside);
268 inline std::unique_ptr<ScratchField>
273 amrex::Real inv_schmidt = 1.0 / lam_schmidt;
276 (*diff)(lev).mult(inv_schmidt);
283 inline std::unique_ptr<ScratchField>
beta()
const override
288#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
290 for (amrex::MFIter mfi((*
beta)(lev)); mfi.isValid(); ++mfi) {
291 const auto& bx = mfi.tilebox();
292 const auto& beta_arr = (*beta)(lev).array(mfi);
302 const amrex::MFIter& mfi,
303 const amrex::Box& bx,
304 const amrex::Array4<amrex::Real>&
beta)
const override
310 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
311 beta(i, j, k) = beta_val;
315 const auto& temp0_arr = temp0(lev).const_array(mfi);
317 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
318 beta(i, j, k) = 1.0 / temp0_arr(i, j, k);
323 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
324 beta(i, j, k) = beta_val;
331 const auto& vof_arr = vof(lev).const_array(mfi);
333 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
340 const auto& dx =
m_repo.
mesh().Geom(lev).CellSizeArray();
341 const auto& phi_arr = levelset(lev).const_array(mfi);
342 const amrex::Real eps = std::cbrt(2. * dx[0] * dx[1] * dx[2]);
344 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
345 amrex::Real smooth_heaviside;
346 if (phi_arr(i, j, k) > eps) {
347 smooth_heaviside = 1.0;
348 }
else if (phi_arr(i, j, k) < -eps) {
349 smooth_heaviside = 0.;
352 0.5 * (1.0 + phi_arr(i, j, k) / eps +
354 std::sin(phi_arr(i, j, k) * M_PI / eps));
356 beta(i, j, k) *= (1 - smooth_heaviside);
367 inline std::unique_ptr<ScratchField>
ref_theta()
const override
370 amrex::Abort(
"Reference temperature was not set");
376#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
378 for (amrex::MFIter mfi((*
ref_theta)(lev)); mfi.isValid(); ++mfi) {
379 const auto& bx = mfi.tilebox();
380 const auto& ref_theta_arr = (*ref_theta)(lev).array(mfi);
390 const amrex::MFIter& mfi,
391 const amrex::Box& bx,
392 const amrex::Array4<amrex::Real>&
ref_theta)
const override
395 amrex::Abort(
"Reference temperature was not set");
400 const auto& temp0_arr = temp0(lev).const_array(mfi);
402 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
408 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
445 amrex::Abort(
"TwoPhaseTransport requires MultiPhase physics");
bool contains(const std::string &key) const
Query if an object exists using the lookup key.
Definition CollMgr.H:55
Definition FieldRepo.H:86
std::unique_ptr< ScratchField > create_scratch_field(const std::string &name, const int ncomp=1, const int nghost=0, const FieldLoc floc=FieldLoc::CELL) const
Definition FieldRepo.cpp:302
Field & get_field(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:149
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
bool field_exists(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:215
const amrex::AmrCore & mesh() const
Return a reference to the underlying AMR mesh instance.
Definition FieldRepo.H:358
Definition MultiPhase.H:27
InterfaceCapturingMethod interface_capturing_method() const
Definition MultiPhase.cpp:80
T & get()
Return a concrete physics instance.
Definition Physics.H:104
Definition TwoPhaseTransport.H:14
amrex::Real laminar_prandtl1() const
Definition TwoPhaseTransport.H:107
std::unique_ptr< ScratchField > ref_theta() const override
Return the reference temperature.
Definition TwoPhaseTransport.H:367
amrex::Real reference_temperature() const override
Definition TwoPhaseTransport.H:361
amrex::Real m_Prt
Turbulent Prandtl number.
Definition TwoPhaseTransport.H:434
std::unique_ptr< ScratchField > beta() const override
Return the thermal expansion coefficient.
Definition TwoPhaseTransport.H:283
std::unique_ptr< ScratchField > mu() override
Return the dynamic visocity field.
Definition TwoPhaseTransport.H:131
std::unique_ptr< ScratchField > scalar_diffusivity(const std::string &scalar_name) override
Definition TwoPhaseTransport.H:269
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:388
amrex::Real turbulent_prandtl() const
Definition TwoPhaseTransport.H:110
InterfaceCapturingMethod get_iface_method() const
Definition TwoPhaseTransport.H:442
amrex::Real m_reference_temperature
Reference temperature.
Definition TwoPhaseTransport.H:440
amrex::Real m_Pr1
Phase 1 (liquid) Prandtl number.
Definition TwoPhaseTransport.H:428
~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:300
static constexpr bool constant_properties
Definition TwoPhaseTransport.H:16
TwoPhaseTransport(const CFDSim &sim)
Definition TwoPhaseTransport.H:20
amrex::Real m_constant_beta
Constant thermal expansion coefficient.
Definition TwoPhaseTransport.H:437
amrex::Real laminar_prandtl2() const
Definition TwoPhaseTransport.H:108
amrex::Real m_mu1
Phase 1 (liquid) dynamic molecular viscosity.
Definition TwoPhaseTransport.H:422
FieldRepo & m_repo
Reference to the field repository (for creating scratch fields)
Definition TwoPhaseTransport.H:416
static std::string identifier()
Definition TwoPhaseTransport.H:18
amrex::Real m_mu2
Phase 2 (gas) dynamic molecular viscosity.
Definition TwoPhaseTransport.H:425
static amrex::Real turbulent_schmidt(const std::string &scalar_name)
Definition TwoPhaseTransport.H:121
amrex::Real m_Pr2
Phase 2 (gas) Prandtl number.
Definition TwoPhaseTransport.H:431
static amrex::Real laminar_schmidt(const std::string &scalar_name)
Definition TwoPhaseTransport.H:112
const PhysicsMgr & m_physics_mgr
Reference to the physics manager.
Definition TwoPhaseTransport.H:419
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field (later divided by density, though)
Definition TwoPhaseTransport.H:197
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:19
InterfaceCapturingMethod
Definition MultiPhase.H:21