ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
SatAdj Class Reference

#include <ERF_SatAdj.H>

Inheritance diagram for SatAdj:
Collaboration diagram for SatAdj:

Public Member Functions

 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) 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
 
void Qmoist_Restart_Vars (const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
 
- Public Member Functions inherited from NullMoist
 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
 

Static Public Member Functions

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonIterSat (int &i, int &j, int &k, const amrex::Real &fac_cond, const amrex::Array4< amrex::Real > &tabs_array, const amrex::Array4< amrex::Real > &pres_array, const amrex::Array4< amrex::Real > &qv_array, const amrex::Array4< amrex::Real > &qc_array)
 

Private Types

using FabPtr = std::shared_ptr< amrex::MultiFab >
 

Private Attributes

int m_qmoist_size = 0
 
int n_qstate_moist_size = 2
 
amrex::Geometry m_geom
 
amrex::Real dt
 
amrex::Real m_fac_cond
 
amrex::Real m_rdOcp
 
bool m_do_cond
 
amrex::Array< FabPtr, MicVar_SatAdj::NumVarsmic_fab_vars
 

Member Typedef Documentation

◆ FabPtr

using SatAdj::FabPtr = std::shared_ptr<amrex::MultiFab>
private

Constructor & Destructor Documentation

◆ SatAdj()

SatAdj::SatAdj ( )
inline
47 {}

◆ ~SatAdj()

virtual SatAdj::~SatAdj ( )
virtualdefault

Member Function Documentation

◆ Advance()

void SatAdj::Advance ( const amrex::Real dt_advance,
const SolverChoice solverChoice 
)
inlineoverridevirtual

Reimplemented from NullMoist.

99  {
100  dt = dt_advance;
101 
102  this->AdvanceSatAdj(solverChoice);
103  }
void AdvanceSatAdj(const SolverChoice &)
Definition: ERF_SatAdj.cpp:8
amrex::Real dt
Definition: ERF_SatAdj.H:203
Here is the call graph for this function:

◆ AdvanceSatAdj()

void SatAdj::AdvanceSatAdj ( const SolverChoice )

Compute Precipitation-related Microphysics quantities.

9 {
10  if (!m_do_cond) { return; }
11 
13 
14  // Expose for GPU
15  Real d_fac_cond = m_fac_cond;
16  Real rdOcp = m_rdOcp;
17 
18  // get the temperature, density, theta, qt and qc from input
19  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
20 
21  const auto& tbx = mfi.tilebox();
22 
23  auto qv_array = mic_fab_vars[MicVar_SatAdj::qv]->array(mfi);
24  auto qc_array = mic_fab_vars[MicVar_SatAdj::qc]->array(mfi);
25  auto tabs_array = mic_fab_vars[MicVar_SatAdj::tabs]->array(mfi);
26  auto theta_array = mic_fab_vars[MicVar_SatAdj::theta]->array(mfi);
27  auto pres_array = mic_fab_vars[MicVar_SatAdj::pres]->array(mfi);
28 
29  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
30  {
31  //qc_array(i,j,k) = std::max(0.0, qc_array(i,j,k));
32 
33  //------- Evaporation/condensation
34  Real qsat;
35  erf_qsatw(tabs_array(i,j,k), pres_array(i,j,k), qsat);
36 
37  // There is enough moisture to drive to equilibrium
38  if ((qv_array(i,j,k)+qc_array(i,j,k)) > qsat) {
39  Real qvprev = qv_array(i,j,k);
40  Real qcprev = qc_array(i,j,k);
41 
42  // clip qc but maintain total water
43  if (qc_array(i,j,k) < 0) {
44  qv_array(i,j,k) += qc_array(i,j,k);
45  qc_array(i,j,k) = 0.0;
46  }
47 
48  // Update temperature
49  tabs_array(i,j,k) = NewtonIterSat(i, j, k ,
50  d_fac_cond, tabs_array, pres_array,
51  qv_array , qc_array );
52 
53  Real qsatnew;
54  erf_qsatw(tabs_array(i,j,k), pres_array(i,j,k), qsatnew);
55 
56  AMREX_ASSERT(std::abs(qv_array(i,j,k)-qsatnew) < 1e-12);
57 
58  amrex::ignore_unused(qvprev);
59  amrex::ignore_unused(qcprev);
60  AMREX_ASSERT(std::abs(qv_array(i,j,k)+qc_array(i,j,k)-qvprev-qcprev) < 1e-14);
61 
62  // Update theta (constant pressure)
63  theta_array(i,j,k) = getThgivenTandP(tabs_array(i,j,k), 100.0*pres_array(i,j,k), rdOcp);
64 
65  //
66  // We cannot blindly relax to qsat, but we can convert qc/qi -> qv.
67  // The concept here is that if we put all the moisture into qv and modify
68  // the temperature, we can then check if qv > qsat occurs (for final T/P/qv).
69  // If the reduction in T/qsat and increase in qv does trigger the
70  // aforementioned condition, we can do Newton iteration to drive qv = qsat.
71  //
72  } else {
73  // Changes in each component
74  Real delta_qc = qc_array(i,j,k);
75 
76  // Partition the change in non-precipitating q
77  qv_array(i,j,k) += qc_array(i,j,k);
78  qc_array(i,j,k) = 0.0;
79 
80  // Update temperature (endothermic since we evap/sublime)
81  tabs_array(i,j,k) -= d_fac_cond * delta_qc;
82 
83  // Update theta
84  theta_array(i,j,k) = getThgivenTandP(tabs_array(i,j,k), 100.0*pres_array(i,j,k), rdOcp);
85 
86  // Verify assumption that qv > qsat does not occur
87  erf_qsatw(tabs_array(i,j,k), pres_array(i,j,k), qsat);
88  if (qv_array(i,j,k) > qsat) {
89  Real qvprev = qv_array(i,j,k);
90  Real qcprev = qc_array(i,j,k);
91  Real Tprev = tabs_array(i,j,k);
92 
93  // Update temperature
94  tabs_array(i,j,k) = NewtonIterSat(i, j, k ,
95  d_fac_cond , tabs_array, pres_array,
96  qv_array , qc_array );
97 
98  Real qsatnew;
99  erf_qsatw(tabs_array(i,j,k), pres_array(i,j,k), qsatnew);
100  amrex::ignore_unused(qvprev);
101  amrex::ignore_unused(qcprev);
102  amrex::ignore_unused(Tprev);
103  AMREX_ASSERT(qv_array(i,j,k) < qvprev);
104  AMREX_ASSERT(qc_array(i,j,k) > qcprev);
105  AMREX_ASSERT(tabs_array(i,j,k) > Tprev);
106  AMREX_ASSERT(std::abs(qv_array(i,j,k)-qsatnew) < 1e-14);
107  AMREX_ASSERT(std::abs(qv_array(i,j,k)+qc_array(i,j,k)-qvprev-qcprev) < 1e-14);
108 
109  // Update theta
110  theta_array(i,j,k) = getThgivenTandP(tabs_array(i,j,k), 100.0*pres_array(i,j,k), rdOcp);
111 
112  }
113  }
114  });
115  }
116 }
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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:163
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Real m_fac_cond
Definition: ERF_SatAdj.H:206
amrex::Array< FabPtr, MicVar_SatAdj::NumVars > mic_fab_vars
Definition: ERF_SatAdj.H:211
amrex::Real m_rdOcp
Definition: ERF_SatAdj.H:207
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonIterSat(int &i, int &j, int &k, const amrex::Real &fac_cond, const amrex::Array4< amrex::Real > &tabs_array, const amrex::Array4< amrex::Real > &pres_array, const amrex::Array4< amrex::Real > &qv_array, const amrex::Array4< amrex::Real > &qc_array)
Definition: ERF_SatAdj.H:130
bool m_do_cond
Definition: ERF_SatAdj.H:208
@ tabs
Definition: ERF_Kessler.H:24
@ theta
Definition: ERF_SatAdj.H:31
@ pres
Definition: ERF_SatAdj.H:33
@ qv
Definition: ERF_SatAdj.H:35
@ tabs
Definition: ERF_SatAdj.H:32
@ qc
Definition: ERF_SatAdj.H:36

Referenced by Advance().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Copy_Micro_to_State()

void SatAdj::Copy_Micro_to_State ( amrex::MultiFab &  cons_in)
overridevirtual

Updates conserved and microphysics variables in the provided MultiFabs from the internal MultiFabs that store Microphysics module data.

Parameters
[out]consConserved variables
[out]qmoistqv, qc

Reimplemented from NullMoist.

13 {
14  // Get the temperature, density, theta, qt and qp from input
15  for (MFIter mfi(cons,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
16  const auto& tbx = mfi.tilebox();
17 
18  auto states_arr = cons.array(mfi);
19 
20  auto rho_arr = mic_fab_vars[MicVar_SatAdj::rho]->array(mfi);
21  auto theta_arr = mic_fab_vars[MicVar_SatAdj::theta]->array(mfi);
22  auto qv_arr = mic_fab_vars[MicVar_SatAdj::qv]->array(mfi);
23  auto qc_arr = mic_fab_vars[MicVar_SatAdj::qc]->array(mfi);
24 
25  // get potential total density, temperature, qt, qp
26  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
27  {
28  states_arr(i,j,k,RhoTheta_comp) = rho_arr(i,j,k)*theta_arr(i,j,k);
29  states_arr(i,j,k,RhoQ1_comp) = rho_arr(i,j,k)*qv_arr(i,j,k);
30  states_arr(i,j,k,RhoQ2_comp) = rho_arr(i,j,k)*qc_arr(i,j,k);
31  });
32  }
33 
34  // Fill interior ghost cells and periodic boundaries
35  cons.FillBoundary(m_geom.periodicity());
36 }
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
amrex::Geometry m_geom
Definition: ERF_SatAdj.H:200
@ rho
Definition: ERF_SatAdj.H:30
@ cons
Definition: ERF_IndexDefines.H:140

Referenced by Update_State_Vars().

Here is the caller graph for this function:

◆ Copy_State_to_Micro()

void SatAdj::Copy_State_to_Micro ( const amrex::MultiFab &  cons_in)
overridevirtual

Initializes the Microphysics module.

Parameters
[in]cons_inConserved variables input

Reimplemented from NullMoist.

41 {
42  // Get the temperature, density, theta, qt and qp from input
43  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
44  const auto& tbx = mfi.tilebox();
45 
46  auto states_array = cons_in.array(mfi);
47 
48  auto qv_array = mic_fab_vars[MicVar_SatAdj::qv]->array(mfi);
49  auto qc_array = mic_fab_vars[MicVar_SatAdj::qc]->array(mfi);
50 
51  auto rho_array = mic_fab_vars[MicVar_SatAdj::rho]->array(mfi);
52  auto theta_array = mic_fab_vars[MicVar_SatAdj::theta]->array(mfi);
53  auto tabs_array = mic_fab_vars[MicVar_SatAdj::tabs]->array(mfi);
54  auto pres_array = mic_fab_vars[MicVar_SatAdj::pres]->array(mfi);
55 
56  // Get pressure, theta, temperature, density, and qt, qp
57  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
58  {
59  rho_array(i,j,k) = states_array(i,j,k,Rho_comp);
60  theta_array(i,j,k) = states_array(i,j,k,RhoTheta_comp)/states_array(i,j,k,Rho_comp);
61  qv_array(i,j,k) = states_array(i,j,k,RhoQ1_comp)/states_array(i,j,k,Rho_comp);
62  qc_array(i,j,k) = states_array(i,j,k,RhoQ2_comp)/states_array(i,j,k,Rho_comp);
63 
64  tabs_array(i,j,k) = getTgivenRandRTh(states_array(i,j,k,Rho_comp),
65  states_array(i,j,k,RhoTheta_comp),
66  qv_array(i,j,k));
67 
68  // Pressure in [mbar] for qsat evaluation
69  pres_array(i,j,k) = getPgivenRTh(states_array(i,j,k,RhoTheta_comp), qv_array(i,j,k)) * 0.01;
70  });
71  }
72 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:46
#define Rho_comp
Definition: ERF_IndexDefines.H:36

Referenced by Update_Micro_Vars().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Define()

void SatAdj::Define ( SolverChoice sc)
inlineoverridevirtual

Reimplemented from NullMoist.

58  {
59  m_fac_cond = lcond / sc.c_p;
60  m_rdOcp = sc.rdOcp;
61  m_do_cond = (!sc.use_shoc);
62  }
constexpr amrex::Real lcond
Definition: ERF_Constants.H:66
amrex::Real rdOcp
Definition: ERF_DataStruct.H:866
amrex::Real c_p
Definition: ERF_DataStruct.H:865
bool use_shoc
Definition: ERF_DataStruct.H:895

◆ 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_inConserved variables input
[in]qc_inCloud variables input
[in,out]qv_inVapor variables input
[in]qi_inIce variables input
[in]gridsThe boxes on which we will evolve the solution
[in]geomGeometry associated with these MultiFabs and grids
[in]dt_advanceTimestep for the advance

Reimplemented from NullMoist.

23 {
24  dt = dt_advance;
25  m_geom = geom;
26 
27  // initialize microphysics variables
28  for (auto ivar = 0; ivar < MicVar_SatAdj::NumVars; ++ivar) {
29  mic_fab_vars[ivar] = std::make_shared<MultiFab>(cons_in.boxArray(), cons_in.DistributionMap(),
30  1, cons_in.nGrowVect());
31  mic_fab_vars[ivar]->setVal(0.);
32  }
33 }
@ NumVars
Definition: ERF_SatAdj.H:37

◆ NewtonIterSat()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real SatAdj::NewtonIterSat ( int &  i,
int &  j,
int &  k,
const amrex::Real fac_cond,
const amrex::Array4< amrex::Real > &  tabs_array,
const amrex::Array4< amrex::Real > &  pres_array,
const amrex::Array4< amrex::Real > &  qv_array,
const amrex::Array4< amrex::Real > &  qc_array 
)
inlinestatic
136  {
137  // Solution tolerance
138  amrex::Real tol = 1.0e-8;
139 
140  // Saturation moisture fractions
141  amrex::Real qsat, dqsat;
142 
143  // Newton iteration vars
144  int niter;
145  amrex::Real fff, dfff;
146  amrex::Real dtabs, delta_qv;
147 
148  // Initial guess for temperature & pressure
149  amrex::Real tabs = tabs_array(i,j,k);
150  amrex::Real pres = pres_array(i,j,k);
151 
152  niter = 0;
153  dtabs = 1;
154 
155  //==================================================
156  // Newton iteration to qv=qsat (cloud phase only)
157  //==================================================
158  do {
159  // Saturation moisture fractions
160  erf_qsatw(tabs, pres, qsat);
161  erf_dtqsatw(tabs, pres, dqsat);
162 
163  // Function for root finding:
164  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
165  fff = -tabs + tabs_array(i,j,k) + fac_cond*(qv_array(i,j,k) - qsat);
166 
167  // Derivative of function (T_new iterated on)
168  dfff = -1.0 - fac_cond*dqsat;
169 
170  // Update the temperature
171  dtabs = -fff/dfff;
172  tabs += dtabs;
173 
174  // Update iteration
175  niter = niter+1;
176  } while(std::abs(dtabs) > tol && niter < 20);
177 
178  // Update qsat from last iteration (dq = dq/dt * dt)
179  qsat += dqsat*dtabs;
180 
181  // Changes in each component
182  delta_qv = qv_array(i,j,k) - qsat;
183 
184  // Partition the change in non-precipitating q
185  qv_array(i,j,k) = qsat;
186  qc_array(i,j,k) += delta_qv;
187 
188  // Return to temperature
189  return tabs;
190  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:178
Here is the call graph for this function:

◆ Qmoist_Ptr()

amrex::MultiFab* SatAdj::Qmoist_Ptr ( const int &  varIdx)
inlineoverridevirtual

Reimplemented from NullMoist.

107  {
108  AMREX_ALWAYS_ASSERT(varIdx < m_qmoist_size);
109  return nullptr;
110  }
int m_qmoist_size
Definition: ERF_SatAdj.H:194

◆ Qmoist_Restart_Vars()

void SatAdj::Qmoist_Restart_Vars ( const SolverChoice ,
std::vector< int > &  a_idx,
std::vector< std::string > &  a_names 
) const
inlineoverridevirtual

Reimplemented from NullMoist.

122  {
123  a_idx.clear();
124  a_names.clear();
125  }

◆ Qmoist_Size()

int SatAdj::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

113 { return SatAdj::m_qmoist_size; }

◆ Qstate_Moist_Size()

int SatAdj::Qstate_Moist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

116 { return SatAdj::n_qstate_moist_size; }
int n_qstate_moist_size
Definition: ERF_SatAdj.H:197

◆ Update_Micro_Vars()

void SatAdj::Update_Micro_Vars ( amrex::MultiFab &  cons_in)
inlineoverridevirtual

Reimplemented from NullMoist.

84  {
85  this->Copy_State_to_Micro(cons_in);
86  }
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSatAdj.cpp:40
Here is the call graph for this function:

◆ Update_State_Vars()

void SatAdj::Update_State_Vars ( amrex::MultiFab &  cons_in)
inlineoverridevirtual

Reimplemented from NullMoist.

91  {
92  this->Copy_Micro_to_State(cons_in);
93  }
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSatAdj.cpp:12
Here is the call graph for this function:

Member Data Documentation

◆ dt

amrex::Real SatAdj::dt
private

Referenced by Advance().

◆ m_do_cond

bool SatAdj::m_do_cond
private

Referenced by Define().

◆ m_fac_cond

amrex::Real SatAdj::m_fac_cond
private

Referenced by Define().

◆ m_geom

amrex::Geometry SatAdj::m_geom
private

◆ m_qmoist_size

int SatAdj::m_qmoist_size = 0
private

Referenced by Qmoist_Ptr(), and Qmoist_Size().

◆ m_rdOcp

amrex::Real SatAdj::m_rdOcp
private

Referenced by Define().

◆ mic_fab_vars

amrex::Array<FabPtr, MicVar_SatAdj::NumVars> SatAdj::mic_fab_vars
private

◆ n_qstate_moist_size

int SatAdj::n_qstate_moist_size = 2
private

Referenced by Qstate_Moist_Size().


The documentation for this class was generated from the following files: