1#ifndef CONSTTRANSPORT_H
2#define CONSTTRANSPORT_H
5#include "AMReX_ParmParse.H"
17 static std::string
identifier() {
return "ConstTransport"; }
21 amrex::ParmParse pp(
"transport");
22 pp.query(
"viscosity",
m_mu);
23 pp.query(
"laminar_prandtl",
m_Pr);
24 pp.query(
"turbulent_prandtl",
m_Prt);
27 amrex::ParmParse pp_boussinesq_buoyancy(
"BoussinesqBuoyancy");
28 amrex::ParmParse pp_abl(
"ABL");
29 if (pp.contains(
"thermal_expansion_coefficient")) {
31 if (pp_boussinesq_buoyancy.contains(
"thermal_expansion_coeff")) {
33 <<
"WARNING: BoussinesqBuoyancy.thermal_expansion_coeff "
34 "option has been deprecated in favor of "
35 "transport.thermal_expansion_coefficient. Ignoring the "
36 "BoussinesqBuoyancy option in favor of the transport "
40 }
else if (pp_boussinesq_buoyancy.contains(
"thermal_expansion_coeff")) {
42 <<
"WARNING: BoussinesqBuoyancy.thermal_expansion_coeff option "
43 "has been deprecated in favor of "
44 "transport.thermal_expansion_coefficient. Please replace "
47 pp_boussinesq_buoyancy.get(
51 if (pp.contains(
"reference_temperature")) {
53 if (pp_boussinesq_buoyancy.contains(
"reference_temperature")) {
55 <<
"WARNING: BoussinesqBuoyancy.reference_temperature "
56 "option has been deprecated in favor of "
57 "transport.reference_temperature. Ignoring the "
58 "BoussinesqBuoyancy option in favor of the transport "
61 }
else if (pp_abl.contains(
"reference_temperature")) {
63 <<
"WARNING: ABL.reference_temperature "
64 "option has been deprecated in favor of "
65 "transport.reference_temperature. Ignoring the "
66 "ABL option in favor of the transport "
70 }
else if (pp_boussinesq_buoyancy.contains(
"reference_temperature")) {
72 <<
"WARNING: BoussinesqBuoyancy.reference_temperature option "
73 "has been deprecated in favor of "
74 "transport.reference_temperature. Please replace "
77 pp_boussinesq_buoyancy.get(
79 }
else if (pp_abl.contains(
"reference_temperature")) {
80 amrex::Print() <<
"WARNING: ABL.reference_temperature option "
81 "has been deprecated in favor of "
82 "transport.reference_temperature. Please replace "
101 amrex::ParmParse pp(
"transport");
102 const std::string key = scalar_name +
"_laminar_schmidt";
103 amrex::Real lam_schmidt = 1.0;
104 pp.query(key.c_str(), lam_schmidt);
110 amrex::ParmParse pp(
"transport");
111 const std::string key = scalar_name +
"_turbulent_schmidt";
112 amrex::Real turb_schmidt = 1.0;
113 pp.query(key.c_str(), turb_schmidt);
118 inline std::unique_ptr<ScratchField>
mu()
override
122 (*mu)(lev).setVal(
m_mu);
118 inline std::unique_ptr<ScratchField>
mu()
override {
…}
128 inline std::unique_ptr<ScratchField>
alpha()
override
131 amrex::Real inv_Pr = 1.0 /
m_Pr;
133 (*alpha)(lev).mult(inv_Pr);
128 inline std::unique_ptr<ScratchField>
alpha()
override {
…}
138 inline std::unique_ptr<ScratchField>
143 amrex::Real inv_schmidt = 1.0 / lam_schmidt;
146 (*diff)(lev).mult(inv_schmidt);
153 inline std::unique_ptr<ScratchField>
beta()
const override
158#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
160 for (amrex::MFIter mfi((*
beta)(lev)); mfi.isValid(); ++mfi) {
161 const auto& bx = mfi.tilebox();
162 const auto& beta_arr = (*beta)(lev).array(mfi);
153 inline std::unique_ptr<ScratchField>
beta()
const override {
…}
172 const amrex::MFIter& mfi,
173 const amrex::Box& bx,
174 const amrex::Array4<amrex::Real>&
beta)
const override
180 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
181 beta(i, j, k) = beta_val;
185 const auto& temp0_arr = temp0(lev).const_array(mfi);
187 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
188 beta(i, j, k) = 1.0 / temp0_arr(i, j, k);
193 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
194 beta(i, j, k) = beta_val;
200 const auto& vof_arr = vof(lev).const_array(mfi);
202 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
216 inline std::unique_ptr<ScratchField>
ref_theta()
const override
219 amrex::Abort(
"Reference temperature was not set");
225#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
227 for (amrex::MFIter mfi((*
ref_theta)(lev)); mfi.isValid(); ++mfi) {
228 const auto& bx = mfi.tilebox();
229 const auto& ref_theta_arr = (*ref_theta)(lev).array(mfi);
216 inline std::unique_ptr<ScratchField>
ref_theta()
const override {
…}
239 const amrex::MFIter& mfi,
240 const amrex::Box& bx,
241 const amrex::Array4<amrex::Real>&
ref_theta)
const override
244 amrex::Abort(
"Reference temperature was not set");
249 const auto& temp0_arr = temp0(lev).const_array(mfi);
251 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
257 bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
noexcept {
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
Definition ConstTransport.H:13
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 ConstTransport.H:237
amrex::Real m_mu
(Laminar) dynamic viscosity
Definition ConstTransport.H:268
amrex::Real thermal_diffusivity() const
Definition ConstTransport.H:93
FieldRepo & m_repo
Reference to the field repository (for creating scratch fields)
Definition ConstTransport.H:265
amrex::Real reference_temperature() const override
Definition ConstTransport.H:210
amrex::Real turbulent_prandtl() const
Definition ConstTransport.H:97
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 ConstTransport.H:170
~ConstTransport() override=default
amrex::Real m_Prt
Turbulent Prandtl number.
Definition ConstTransport.H:274
static constexpr bool constant_properties
Definition ConstTransport.H:15
amrex::Real viscosity() const
Definition ConstTransport.H:91
static amrex::Real turbulent_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:108
amrex::Real m_Pr
Prandtl number.
Definition ConstTransport.H:271
amrex::Real laminar_prandtl() const
Definition ConstTransport.H:95
static std::string identifier()
Definition ConstTransport.H:17
ConstTransport(const CFDSim &sim)
Definition ConstTransport.H:19
std::unique_ptr< ScratchField > beta() const override
Return the thermal expansion coefficient.
Definition ConstTransport.H:153
amrex::Real m_constant_beta
Constant thermal expansion coefficient.
Definition ConstTransport.H:277
static amrex::Real laminar_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:99
amrex::Real m_reference_temperature
Reference temperature.
Definition ConstTransport.H:280
std::unique_ptr< ScratchField > ref_theta() const override
Return the reference temperature.
Definition ConstTransport.H:216
std::unique_ptr< ScratchField > scalar_diffusivity(const std::string &scalar_name) override
Definition ConstTransport.H:139
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field.
Definition ConstTransport.H:128
std::unique_ptr< ScratchField > mu() override
Return the dynamic visocity field.
Definition ConstTransport.H:118
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:19