#include <ERF_SatAdj.H>
|
| | SatAdj () |
| |
| virtual | ~SatAdj ()=default |
| |
| void | AdvanceSatAdj (const SolverChoice &) |
| |
| void | Define (SolverChoice &sc) override |
| |
| void | Init (const amrex::MultiFab &cons_in, const amrex::BoxArray &, const amrex::Geometry &geom, const amrex::Real &dt_advance, std::unique_ptr< amrex::MultiFab > &, std::unique_ptr< amrex::MultiFab > &) override |
| |
| void | Copy_State_to_Micro (const amrex::MultiFab &cons_in) override |
| |
| void | Copy_Micro_to_State (amrex::MultiFab &cons_in) override |
| |
| void | Update_Micro_Vars (amrex::MultiFab &cons_in) override |
| |
| void | Update_State_Vars (amrex::MultiFab &cons_in, const amrex::MultiFab &) override |
| |
| void | Advance (const amrex::Real &dt_advance, const SolverChoice &solverChoice) override |
| |
| amrex::MultiFab * | Qmoist_Ptr (const int &varIdx) override |
| |
| int | Qmoist_Size () override |
| |
| int | Qstate_Moist_Size () override |
| |
| int | Qstate_Moist_NumConc_Size () override |
| |
| void | Qmoist_Restart_Vars (const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override |
| |
| | NullMoist () |
| |
| virtual | ~NullMoist ()=default |
| |
| virtual int | Qstate_NonMoist_Size () |
| |
| virtual void | GetPlotVarNames (amrex::Vector< std::string > &a_vec) const |
| |
| virtual void | GetPlotVar (const std::string &, amrex::MultiFab &) const |
| |
| virtual void | GetPlotVar (const std::string &a_name, amrex::MultiFab &a_mf, const int) const |
| |
| virtual void | SetCurrentLevel (const int &) |
| |
| virtual void | InitLevel (const int, const amrex::MultiFab &) |
| |
| virtual int | getDiagnosticsInterval () const |
| |
| virtual void | Set_dzmin (const amrex::Real) |
| |
| virtual void | Set_RealWidth (const int) |
| |
|
| AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real | NewtonSolveSatTemperature (const amrex::Real fac_cond, const amrex::Real tabs_old, const amrex::Real pres_mbar, const amrex::Real qv_old, amrex::Real &qsat) |
| |
| AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void | AdjustSatAdjCell (const amrex::Real fac_cond, const amrex::Real rdOcp, amrex::Real &tabs, const amrex::Real pres_mbar, amrex::Real &theta, amrex::Real &qv, amrex::Real &qc) |
| |
|
| using | FabPtr = std::shared_ptr< amrex::MultiFab > |
| |
◆ FabPtr
◆ SatAdj()
◆ ~SatAdj()
| virtual SatAdj::~SatAdj |
( |
| ) |
|
|
virtualdefault |
◆ AdjustSatAdjCell()
215 if ((
qv +
qc) > qsat) {
220 #if defined(AMREX_DEBUG)
238 #if defined(AMREX_DEBUG)
241 AMREX_ASSERT(std::abs(
qv-qsatnew) < 1e-12);
242 AMREX_ASSERT(std::abs(
qv+
qc-qvprev-qcprev) < 1e-14);
255 tabs -= fac_cond * delta_qc;
260 #if defined(AMREX_DEBUG)
272 #if defined(AMREX_DEBUG)
275 AMREX_ASSERT(
qv < qvprev);
276 AMREX_ASSERT(
qc > qcprev);
277 AMREX_ASSERT(
tabs > Tprev);
278 AMREX_ASSERT(std::abs(
qv-qsatnew) < 1e-14);
279 AMREX_ASSERT(std::abs(
qv+
qc-qvprev-qcprev) < 1e-14);
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getThgivenTandP(const amrex::Real T, const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:18
const Real rdOcp
Definition: ERF_InitCustomPert_Bomex.H:16
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:228
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonSolveSatTemperature(const amrex::Real fac_cond, const amrex::Real tabs_old, const amrex::Real pres_mbar, const amrex::Real qv_old, amrex::Real &qsat)
Definition: ERF_SatAdj.H:150
@ theta
Definition: ERF_SatAdj.H:35
@ qv
Definition: ERF_SatAdj.H:39
@ tabs
Definition: ERF_SatAdj.H:36
@ qc
Definition: ERF_SatAdj.H:40
◆ Advance()
Reimplemented from NullMoist.
void AdvanceSatAdj(const SolverChoice &)
Definition: ERF_SatAdj.cpp:11
amrex::Real dt
Definition: ERF_SatAdj.H:301
◆ AdvanceSatAdj()
AdvanceSatAdj is a local valid-cell update over MFIter tileboxes. It uses pressure diagnosed in Copy_State_to_Micro and holds that pressure fixed inside each cell adjustment. There are no stencils, face fluxes, or ghost-cell reads in this kernel.
23 for ( MFIter mfi(*
tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
25 auto tbx = mfi.tilebox();
33 ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
35 Real T = tabs_array(i,j,k);
36 Real p = pres_array(i,j,k);
37 Real th = theta_array(i,j,k);
43 tabs_array(i,j,k) =
T;
44 theta_array(i,j,k) = th;
Real T
Definition: ERF_InitCustomPert_Bubble.H:105
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void AdjustSatAdjCell(const amrex::Real fac_cond, const amrex::Real rdOcp, amrex::Real &tabs, const amrex::Real pres_mbar, amrex::Real &theta, amrex::Real &qv, amrex::Real &qc)
Definition: ERF_SatAdj.H:204
amrex::Real m_fac_cond
Definition: ERF_SatAdj.H:304
amrex::Array< FabPtr, MicVar_SatAdj::NumVars > mic_fab_vars
Definition: ERF_SatAdj.H:309
amrex::Real m_rdOcp
Definition: ERF_SatAdj.H:305
bool m_do_cond
Definition: ERF_SatAdj.H:306
@ tabs
Definition: ERF_Kessler.H:25
@ qv
Definition: ERF_Kessler.H:29
@ pres
Definition: ERF_SatAdj.H:37
@ p
Definition: ERF_WSM6.H:176
Referenced by Advance().
◆ Copy_Micro_to_State()
| void SatAdj::Copy_Micro_to_State |
( |
amrex::MultiFab & |
cons_in | ) |
|
|
overridevirtual |
Copies SatAdj's adjusted specific variables back into ERF's conserved state. The copied rho is used only to reconstruct rhoTheta, rhoQv, and rhoQc. Density itself is intentionally not written back because SatAdj does not own or modify rho.
- Parameters
-
| [out] | cons | Conserved variables |
Reimplemented from NullMoist.
15 for (MFIter mfi(
cons,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
16 const auto& tbx = mfi.tilebox();
18 auto states_arr =
cons.array(mfi);
28 ParallelFor(tbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
auto qv_arr
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:210
auto rho_arr
Definition: ERF_UpdateWSubsidence_SineMassFlux.H:3
amrex::Geometry m_geom
Definition: ERF_SatAdj.H:298
@ rho
Definition: ERF_SatAdj.H:34
@ cons
Definition: ERF_IndexDefines.H:174
Referenced by Update_State_Vars().
◆ Copy_State_to_Micro()
| void SatAdj::Copy_State_to_Micro |
( |
const amrex::MultiFab & |
cons_in | ) |
|
|
overridevirtual |
Copies conserved ERF state into SatAdj's internal specific variables. rho is copied but not modified by SatAdj. theta, qv, and qc are formed by dividing conserved densities by rho. tabs and pressure are diagnosed from the same conserved state snapshot. pressure is converted from Pa to mbar/hPa for saturation helper calls.
- Parameters
-
| [in] | cons_in | Conserved variables input |
Reimplemented from NullMoist.
46 for (MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
47 const auto& tbx = mfi.tilebox();
49 auto states_array = cons_in.array(mfi);
62 ParallelFor(tbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
70 rho_array(i,j,k) =
rho;
71 theta_array(i,j,k) =
theta;
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:46
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:81
#define Rho_comp
Definition: ERF_IndexDefines.H:36
rho
Definition: ERF_InitCustomPert_Bubble.H:106
@ theta
Definition: ERF_MM5.H:20
Referenced by Update_Micro_Vars().
◆ Define()
Reimplemented from NullMoist.
constexpr amrex::Real lcond
Definition: ERF_Constants.H:85
amrex::Real rdOcp
Definition: ERF_DataStruct.H:1222
amrex::Real c_p
Definition: ERF_DataStruct.H:1221
bool use_shoc
Definition: ERF_DataStruct.H:1253
◆ Init()
| void SatAdj::Init |
( |
const amrex::MultiFab & |
cons_in, |
|
|
const amrex::BoxArray & |
, |
|
|
const amrex::Geometry & |
geom, |
|
|
const amrex::Real & |
dt_advance, |
|
|
std::unique_ptr< amrex::MultiFab > & |
, |
|
|
std::unique_ptr< amrex::MultiFab > & |
|
|
) |
| |
|
overridevirtual |
Initializes the Microphysics module.
- Parameters
-
| [in] | cons_in | Conserved variables input |
| [in] | qc_in | Cloud variables input |
| [in,out] | qv_in | Vapor variables input |
| [in] | qi_in | Ice variables input |
| [in] | grids | The boxes on which we will evolve the solution |
| [in] | geom | Geometry associated with these MultiFabs and grids |
| [in] | dt_advance | Timestep for the advance |
Reimplemented from NullMoist.
29 mic_fab_vars[ivar] = std::make_shared<MultiFab>(cons_in.boxArray(), cons_in.DistributionMap(),
30 1, cons_in.nGrowVect());
@ NumVars
Definition: ERF_SatAdj.H:41
◆ NewtonSolveSatTemperature()
157 constexpr
int max_niter = 20;
177 fff = -
tabs + tabs_old + fac_cond*(qv_old - qsat);
180 dfff = -
one - fac_cond*dqsat;
188 }
while(std::abs(dtabs) > tol && niter < max_niter);
constexpr amrex::Real one
Definition: ERF_Constants.H:7
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:244
Referenced by AdjustSatAdjCell().
◆ Qmoist_Ptr()
| amrex::MultiFab* SatAdj::Qmoist_Ptr |
( |
const int & |
varIdx | ) |
|
|
inlineoverridevirtual |
Reimplemented from NullMoist.
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
int m_qmoist_size
Definition: ERF_SatAdj.H:289
◆ Qmoist_Restart_Vars()
| void SatAdj::Qmoist_Restart_Vars |
( |
const SolverChoice & |
, |
|
|
std::vector< int > & |
a_idx, |
|
|
std::vector< std::string > & |
a_names |
|
) |
| const |
|
inlineoverridevirtual |
◆ Qmoist_Size()
| int SatAdj::Qmoist_Size |
( |
| ) |
|
|
inlineoverridevirtual |
◆ Qstate_Moist_NumConc_Size()
| int SatAdj::Qstate_Moist_NumConc_Size |
( |
| ) |
|
|
inlineoverridevirtual |
Reimplemented from NullMoist.
int n_qstate_moist_numconc_size
Definition: ERF_SatAdj.H:295
◆ Qstate_Moist_Size()
| int SatAdj::Qstate_Moist_Size |
( |
| ) |
|
|
inlineoverridevirtual |
Reimplemented from NullMoist.
int n_qstate_moist_size
Definition: ERF_SatAdj.H:292
◆ Update_Micro_Vars()
| void SatAdj::Update_Micro_Vars |
( |
amrex::MultiFab & |
cons_in | ) |
|
|
inlineoverridevirtual |
Reimplemented from NullMoist.
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSatAdj.cpp:44
◆ Update_State_Vars()
| void SatAdj::Update_State_Vars |
( |
amrex::MultiFab & |
cons_in, |
|
|
const amrex::MultiFab & |
|
|
) |
| |
|
inlineoverridevirtual |
Reimplemented from NullMoist.
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSatAdj.cpp:13
◆ dt
◆ m_do_cond
◆ m_fac_cond
◆ m_geom
| amrex::Geometry SatAdj::m_geom |
|
private |
◆ m_qmoist_size
| int SatAdj::m_qmoist_size = 0 |
|
private |
◆ m_rdOcp
◆ mic_fab_vars
◆ n_qstate_moist_numconc_size
| int SatAdj::n_qstate_moist_numconc_size = 0 |
|
private |
◆ n_qstate_moist_size
| int SatAdj::n_qstate_moist_size = 2 |
|
private |
The documentation for this class was generated from the following files: