1#ifndef CONSTTRANSPORT_H
2#define CONSTTRANSPORT_H
5#include "AMReX_ParmParse.H"
8using namespace amrex::literals;
20 static std::string
identifier() {
return "ConstTransport"; }
24 amrex::ParmParse pp(
"transport");
25 pp.query(
"viscosity",
m_mu);
26 pp.query(
"laminar_prandtl",
m_Pr);
27 pp.query(
"turbulent_prandtl",
m_Prt);
30 amrex::ParmParse pp_boussinesq_buoyancy(
"BoussinesqBuoyancy");
31 amrex::ParmParse pp_abl(
"ABL");
32 if (pp.contains(
"thermal_expansion_coefficient")) {
34 if (pp_boussinesq_buoyancy.contains(
"thermal_expansion_coeff")) {
36 <<
"WARNING: BoussinesqBuoyancy.thermal_expansion_coeff "
37 "option has been deprecated in favor of "
38 "transport.thermal_expansion_coefficient. Ignoring the "
39 "BoussinesqBuoyancy option in favor of the transport "
43 }
else if (pp_boussinesq_buoyancy.contains(
"thermal_expansion_coeff")) {
45 <<
"WARNING: BoussinesqBuoyancy.thermal_expansion_coeff option "
46 "has been deprecated in favor of "
47 "transport.thermal_expansion_coefficient. Please replace "
50 pp_boussinesq_buoyancy.get(
54 if (pp.contains(
"reference_temperature")) {
56 if (pp_boussinesq_buoyancy.contains(
"reference_temperature")) {
58 <<
"WARNING: BoussinesqBuoyancy.reference_temperature "
59 "option has been deprecated in favor of "
60 "transport.reference_temperature. Ignoring the "
61 "BoussinesqBuoyancy option in favor of the transport "
64 }
else if (pp_abl.contains(
"reference_temperature")) {
66 <<
"WARNING: ABL.reference_temperature "
67 "option has been deprecated in favor of "
68 "transport.reference_temperature. Ignoring the "
69 "ABL option in favor of the transport "
73 }
else if (pp_boussinesq_buoyancy.contains(
"reference_temperature")) {
75 <<
"WARNING: BoussinesqBuoyancy.reference_temperature option "
76 "has been deprecated in favor of "
77 "transport.reference_temperature. Please replace "
80 pp_boussinesq_buoyancy.get(
82 }
else if (pp_abl.contains(
"reference_temperature")) {
83 amrex::Print() <<
"WARNING: ABL.reference_temperature option "
84 "has been deprecated in favor of "
85 "transport.reference_temperature. Please replace "
107 amrex::ParmParse pp(
"transport");
108 const std::string key = scalar_name +
"_laminar_schmidt";
109 amrex::Real lam_schmidt = 1.0_rt;
110 pp.query(key, lam_schmidt);
116 amrex::ParmParse pp(
"transport");
117 const std::string key = scalar_name +
"_turbulent_schmidt";
118 amrex::Real turb_schmidt = 1.0_rt;
119 pp.query(key, turb_schmidt);
124 std::unique_ptr<ScratchField>
mu()
override
126 auto mu =
m_repo.create_scratch_field(1, m_ngrow);
127 for (
int lev = 0; lev <
m_repo.num_active_levels(); ++lev) {
128 (*mu)(lev).setVal(
m_mu);
134 std::unique_ptr<ScratchField>
alpha()
override
137 amrex::Real inv_Pr = 1.0_rt /
m_Pr;
138 for (
int lev = 0; lev <
m_repo.num_active_levels(); ++lev) {
139 (*alpha)(lev).mult(inv_Pr);
144 std::unique_ptr<ScratchField>
149 amrex::Real inv_schmidt = 1.0_rt / lam_schmidt;
151 for (
int lev = 0; lev <
m_repo.num_active_levels(); ++lev) {
152 (*diff)(lev).mult(inv_schmidt);
159 [[nodiscard]] std::unique_ptr<ScratchField>
beta()
const override
161 auto beta =
m_repo.create_scratch_field(1, m_ngrow);
162 for (
int lev = 0; lev <
m_repo.num_active_levels(); ++lev) {
164#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
166 for (amrex::MFIter mfi((*
beta)(lev)); mfi.isValid(); ++mfi) {
167 const auto& bx = mfi.tilebox();
168 const auto& beta_arr = (*beta)(lev).array(mfi);
178 const amrex::MFIter& mfi,
179 const amrex::Box& bx,
180 const amrex::Array4<amrex::Real>&
beta)
const override
185 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
186 beta(i, j, k) = beta_val;
188 }
else if (
m_repo.field_exists(
"reference_temperature")) {
189 const auto& temp0 =
m_repo.get_field(
"reference_temperature");
190 const auto& temp0_arr = temp0(lev).const_array(mfi);
191 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
192 beta(i, j, k) = 1.0_rt / temp0_arr(i, j, k);
196 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
197 beta(i, j, k) = beta_val;
201 if (
m_repo.field_exists(
"vof")) {
202 const auto& vof =
m_repo.get_field(
"vof");
203 const auto& vof_arr = vof(lev).const_array(mfi);
204 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
206 beta(i, j, k) = 0.0_rt;
218 [[nodiscard]] std::unique_ptr<ScratchField>
ref_theta()
const override
221 amrex::Abort(
"Reference temperature was not set");
225 for (
int lev = 0; lev <
m_repo.num_active_levels(); ++lev) {
227#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
229 for (amrex::MFIter mfi((*
ref_theta)(lev)); mfi.isValid(); ++mfi) {
230 const auto& bx = mfi.tilebox();
231 const auto& ref_theta_arr = (*ref_theta)(lev).array(mfi);
241 const amrex::MFIter& mfi,
242 const amrex::Box& bx,
243 const amrex::Array4<amrex::Real>&
ref_theta)
const override
246 amrex::Abort(
"Reference temperature was not set");
249 if (
m_repo.field_exists(
"reference_temperature")) {
250 auto& temp0 =
m_repo.get_field(
"reference_temperature");
251 const auto& temp0_arr = temp0(lev).const_array(mfi);
252 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
257 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
Definition FieldRepo.H:86
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:239
amrex::Real m_mu
(Laminar) dynamic viscosity
Definition ConstTransport.H:268
amrex::Real thermal_diffusivity() const
Definition ConstTransport.H:96
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:212
amrex::Real turbulent_prandtl() const
Definition ConstTransport.H:103
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:176
~ConstTransport() override=default
amrex::Real m_Prt
Turbulent Prandtl number.
Definition ConstTransport.H:274
static constexpr bool constant_properties
Definition ConstTransport.H:18
amrex::Real viscosity() const
Definition ConstTransport.H:94
static amrex::Real turbulent_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:114
amrex::Real m_Pr
Prandtl number.
Definition ConstTransport.H:271
amrex::Real laminar_prandtl() const
Definition ConstTransport.H:101
static std::string identifier()
Definition ConstTransport.H:20
ConstTransport(const CFDSim &sim)
Definition ConstTransport.H:22
std::unique_ptr< ScratchField > beta() const override
Return the thermal expansion coefficient.
Definition ConstTransport.H:159
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:105
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:218
std::unique_ptr< ScratchField > scalar_diffusivity(const std::string &scalar_name) override
Definition ConstTransport.H:145
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field.
Definition ConstTransport.H:134
std::unique_ptr< ScratchField > mu() override
Return the dynamic visocity field.
Definition ConstTransport.H:124
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:17