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 const auto& volfrac_arrs = vof(lev).const_arrays();
143 const auto& visc_arrs = (*mu)(lev).arrays();
144 const amrex::Real mu1 =
m_mu1;
145 const amrex::Real mu2 =
m_mu2;
147 (*
mu)(lev),
mu->num_grow(),
148 [=] AMREX_GPU_DEVICE(
149 int nbx,
int i,
int j,
int k)
noexcept {
150 visc_arrs[nbx](i, j, k) =
151 mu1 * volfrac_arrs[nbx](i, j, k) +
152 mu2 * (1.0 - volfrac_arrs[nbx](i, j, k));
155 amrex::Gpu::streamSynchronize();
163 const auto& dx = geom[lev].CellSizeArray();
164 const auto& visc_arrs = (*mu)(lev).arrays();
165 const auto& phi_arrs = levelset(lev).const_arrays();
166 const amrex::Real eps = std::cbrt(2. * dx[0] * dx[1] * dx[2]);
167 const amrex::Real mu1 =
m_mu1;
168 const amrex::Real mu2 =
m_mu2;
171 (*
mu)(lev),
mu->num_grow(),
172 [=] AMREX_GPU_DEVICE(
173 int nbx,
int i,
int j,
int k)
noexcept {
174 amrex::Real smooth_heaviside;
175 if (phi_arrs[nbx](i, j, k) > eps) {
176 smooth_heaviside = 1.0;
177 }
else if (phi_arrs[nbx](i, j, k) < -eps) {
178 smooth_heaviside = 0.;
182 (1.0 + phi_arrs[nbx](i, j, k) / eps +
185 phi_arrs[nbx](i, j, k) * M_PI / eps));
187 visc_arrs[nbx](i, j, k) =
188 mu1 * smooth_heaviside +
189 mu2 * (1.0 - smooth_heaviside);
192 amrex::Gpu::streamSynchronize();
131 inline std::unique_ptr<ScratchField>
mu()
override {
…}
198 inline std::unique_ptr<ScratchField>
alpha()
override
209 const auto& volfrac_arrs = vof(lev).const_arrays();
210 const auto& thdiff_arrs = (*alpha)(lev).arrays();
211 const amrex::Real mu1 =
m_mu1;
212 const amrex::Real mu2 =
m_mu2;
213 const amrex::Real Pr1 =
m_Pr1;
214 const amrex::Real Pr2 =
m_Pr2;
217 [=] AMREX_GPU_DEVICE(
218 int nbx,
int i,
int j,
int k)
noexcept {
219 thdiff_arrs[nbx](i, j, k) =
220 mu1 / Pr1 * volfrac_arrs[nbx](i, j, k) +
221 mu2 / Pr2 * (1.0 - volfrac_arrs[nbx](i, j, k));
230 const auto& dx = geom[lev].CellSizeArray();
231 const auto& visc_arrs = (*alpha)(lev).arrays();
232 const auto& phi_arrs = levelset(lev).const_arrays();
233 const amrex::Real eps = std::cbrt(2. * dx[0] * dx[1] * dx[2]);
234 const amrex::Real mu1 =
m_mu1;
235 const amrex::Real mu2 =
m_mu2;
236 const amrex::Real Pr1 =
m_Pr1;
237 const amrex::Real Pr2 =
m_Pr2;
240 [=] AMREX_GPU_DEVICE(
241 int nbx,
int i,
int j,
int k)
noexcept {
242 amrex::Real smooth_heaviside;
243 if (phi_arrs[nbx](i, j, k) > eps) {
244 smooth_heaviside = 1.0;
245 }
else if (phi_arrs[nbx](i, j, k) < -eps) {
246 smooth_heaviside = 0.;
250 (1.0 + phi_arrs[nbx](i, j, k) / eps +
253 phi_arrs[nbx](i, j, k) * M_PI / eps));
255 visc_arrs[nbx](i, j, k) =
256 mu1 / Pr1 * smooth_heaviside +
257 mu2 / Pr2 * (1.0 - smooth_heaviside);
261 amrex::Gpu::streamSynchronize();
198 inline std::unique_ptr<ScratchField>
alpha()
override {
…}
265 inline std::unique_ptr<ScratchField>
270 amrex::Real inv_schmidt = 1.0 / lam_schmidt;
273 (*diff)(lev).mult(inv_schmidt);
280 inline std::unique_ptr<ScratchField>
beta()
const override
285#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
287 for (amrex::MFIter mfi((*
beta)(lev)); mfi.isValid(); ++mfi) {
288 const auto& bx = mfi.tilebox();
289 const auto& beta_arr = (*beta)(lev).array(mfi);
280 inline std::unique_ptr<ScratchField>
beta()
const override {
…}
299 const amrex::MFIter& mfi,
300 const amrex::Box& bx,
301 const amrex::Array4<amrex::Real>&
beta)
const override
307 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
308 beta(i, j, k) = beta_val;
312 const auto& temp0_arr = temp0(lev).const_array(mfi);
314 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
315 beta(i, j, k) = 1.0 / temp0_arr(i, j, k);
320 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
321 beta(i, j, k) = beta_val;
328 const auto& vof_arr = vof(lev).const_array(mfi);
330 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
337 const auto& dx =
m_repo.
mesh().Geom(lev).CellSizeArray();
338 const auto& phi_arr = levelset(lev).const_array(mfi);
339 const amrex::Real eps = std::cbrt(2. * dx[0] * dx[1] * dx[2]);
341 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
342 amrex::Real smooth_heaviside;
343 if (phi_arr(i, j, k) > eps) {
344 smooth_heaviside = 1.0;
345 }
else if (phi_arr(i, j, k) < -eps) {
346 smooth_heaviside = 0.;
349 0.5 * (1.0 + phi_arr(i, j, k) / eps +
351 std::sin(phi_arr(i, j, k) * M_PI / eps));
353 beta(i, j, k) *= (1 - smooth_heaviside);
364 inline std::unique_ptr<ScratchField>
ref_theta()
const override
367 amrex::Abort(
"Reference temperature was not set");
373#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
375 for (amrex::MFIter mfi((*
ref_theta)(lev)); mfi.isValid(); ++mfi) {
376 const auto& bx = mfi.tilebox();
377 const auto& ref_theta_arr = (*ref_theta)(lev).array(mfi);
364 inline std::unique_ptr<ScratchField>
ref_theta()
const override {
…}
387 const amrex::MFIter& mfi,
388 const amrex::Box& bx,
389 const amrex::Array4<amrex::Real>&
ref_theta)
const override
392 amrex::Abort(
"Reference temperature was not set");
397 const auto& temp0_arr = temp0(lev).const_array(mfi);
399 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
405 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
442 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:364
amrex::Real reference_temperature() const override
Definition TwoPhaseTransport.H:358
amrex::Real m_Prt
Turbulent Prandtl number.
Definition TwoPhaseTransport.H:431
std::unique_ptr< ScratchField > beta() const override
Return the thermal expansion coefficient.
Definition TwoPhaseTransport.H:280
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:266
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:385
amrex::Real turbulent_prandtl() const
Definition TwoPhaseTransport.H:110
InterfaceCapturingMethod get_iface_method() const
Definition TwoPhaseTransport.H:439
amrex::Real m_reference_temperature
Reference temperature.
Definition TwoPhaseTransport.H:437
amrex::Real m_Pr1
Phase 1 (liquid) Prandtl number.
Definition TwoPhaseTransport.H:425
~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:297
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:434
amrex::Real laminar_prandtl2() const
Definition TwoPhaseTransport.H:108
amrex::Real m_mu1
Phase 1 (liquid) dynamic molecular viscosity.
Definition TwoPhaseTransport.H:419
FieldRepo & m_repo
Reference to the field repository (for creating scratch fields)
Definition TwoPhaseTransport.H:413
static std::string identifier()
Definition TwoPhaseTransport.H:18
amrex::Real m_mu2
Phase 2 (gas) dynamic molecular viscosity.
Definition TwoPhaseTransport.H:422
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:428
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:416
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field (later divided by density, though)
Definition TwoPhaseTransport.H:198
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:19
InterfaceCapturingMethod
Definition MultiPhase.H:21