/home/runner/work/amr-wind/amr-wind/amr-wind/transport_models/ConstTransport.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/transport_models/ConstTransport.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
ConstTransport.H
Go to the documentation of this file.
1#ifndef CONSTTRANSPORT_H
2#define CONSTTRANSPORT_H
3
5#include "AMReX_ParmParse.H"
6
8
13{
14public:
15 static constexpr bool constant_properties = true;
16
17 static std::string identifier() { return "ConstTransport"; }
18
19 explicit ConstTransport(CFDSim& sim) : m_repo(sim.repo())
20 {
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);
25 }
26
27 ~ConstTransport() override = default;
28
29 inline amrex::Real viscosity() const { return m_mu; }
30
31 inline amrex::Real thermal_diffusivity() const { return m_mu / m_Pr; }
32
33 inline amrex::Real laminar_prandtl() const { return m_Pr; }
34
35 inline amrex::Real turbulent_prandtl() const { return m_Prt; }
36
37 static inline amrex::Real laminar_schmidt(const std::string& scalar_name)
38 {
39 amrex::ParmParse pp("transport");
40 const std::string key = scalar_name + "_laminar_schmidt";
41 amrex::Real lam_schmidt = 1.0;
42 pp.query(key.c_str(), lam_schmidt);
43 return lam_schmidt;
44 }
45
46 static inline amrex::Real turbulent_schmidt(const std::string& scalar_name)
47 {
48 amrex::ParmParse pp("transport");
49 const std::string key = scalar_name + "_turbulent_schmidt";
50 amrex::Real turb_schmidt = 1.0;
51 pp.query(key.c_str(), turb_schmidt);
52 return turb_schmidt;
53 }
54
56 inline std::unique_ptr<ScratchField> mu() override
57 {
58 auto mu = m_repo.create_scratch_field(1, 1);
59 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
60 (*mu)(lev).setVal(m_mu);
61 }
62 return mu;
63 }
64
66 inline std::unique_ptr<ScratchField> alpha() override
67 {
68 auto alpha = mu();
69 amrex::Real inv_Pr = 1.0 / m_Pr;
70 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
71 (*alpha)(lev).mult(inv_Pr);
72 }
73 return alpha;
74 }
75
76 inline std::unique_ptr<ScratchField>
77 scalar_diffusivity(const std::string& scalar_name) override
78 {
79 amrex::Real lam_schmidt = laminar_schmidt(scalar_name);
80
81 amrex::Real inv_schmidt = 1.0 / lam_schmidt;
82 auto diff = mu();
83 for (int lev = 0; lev < m_repo.num_active_levels(); ++lev) {
84 (*diff)(lev).mult(inv_schmidt);
85 }
86
87 return diff;
88 }
89
90private:
93
95 amrex::Real m_mu{1.0e-5};
96
98 amrex::Real m_Pr{1.0};
99
101 amrex::Real m_Prt{1.0};
102};
103
104} // namespace amr_wind::transport
105
106#endif /* CONSTTRANSPORT_H */
Definition CFDSim.H:47
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
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
Definition ConstTransport.H:13
amrex::Real m_mu
(Laminar) dynamic viscosity
Definition ConstTransport.H:95
amrex::Real thermal_diffusivity() const
Definition ConstTransport.H:31
FieldRepo & m_repo
Reference to the field repository (for creating scratch fields)
Definition ConstTransport.H:92
amrex::Real turbulent_prandtl() const
Definition ConstTransport.H:35
amrex::Real m_Prt
Turbulent Prandtl number.
Definition ConstTransport.H:101
static constexpr bool constant_properties
Definition ConstTransport.H:15
amrex::Real viscosity() const
Definition ConstTransport.H:29
ConstTransport(CFDSim &sim)
Definition ConstTransport.H:19
static amrex::Real turbulent_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:46
amrex::Real m_Pr
Prandtl number.
Definition ConstTransport.H:98
amrex::Real laminar_prandtl() const
Definition ConstTransport.H:33
static std::string identifier()
Definition ConstTransport.H:17
static amrex::Real laminar_schmidt(const std::string &scalar_name)
Definition ConstTransport.H:37
std::unique_ptr< ScratchField > scalar_diffusivity(const std::string &scalar_name) override
Scalar diffusivity based on Schmidt number.
Definition ConstTransport.H:77
std::unique_ptr< ScratchField > alpha() override
Return the thermal diffusivity field.
Definition ConstTransport.H:66
std::unique_ptr< ScratchField > mu() override
Return the dynamic visocity field.
Definition ConstTransport.H:56
Definition TransportModel.H:29
Definition ConstTransport.H:7