ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_SatAdj.H
Go to the documentation of this file.
1 #ifndef ERF_SATADJ_H
2 #define ERF_SATADJ_H
3 
4 /*
5  * Implementation saturation adjustment microphysics model
6  * The model transports qv and qc and does Newton iterations
7  * to complete condensation/evaporation.
8  */
9 
10 #include <string>
11 #include <vector>
12 #include <memory>
13 
14 #include <AMReX_FArrayBox.H>
15 #include <AMReX_Geometry.H>
16 #include <AMReX_MultiFabUtil.H>
17 #include <AMReX_GpuContainers.H>
18 
19 #include "ERF_EOS.H"
20 #include "ERF_Constants.H"
21 #include "ERF_MicrophysicsUtils.H"
22 #include "ERF_IndexDefines.H"
23 #include "ERF_DataStruct.H"
24 #include "ERF_NullMoist.H"
25 #include "ERF_TileNoZ.H"
26 
27 namespace MicVar_SatAdj {
28  enum {
29  // independent variables
30  rho=0, // density
31  theta, // liquid/ice water potential temperature
32  tabs, // temperature
33  pres, // pressure
34  // non-precipitating vars
35  qv, // cloud vapor
36  qc, // cloud water
37  NumVars
38  };
39 }
40 
41 class SatAdj : public NullMoist {
42 
43  using FabPtr = std::shared_ptr<amrex::MultiFab>;
44 
45 public:
46  // constructor
47  SatAdj () {}
48 
49  // destructor
50  virtual ~SatAdj () = default;
51 
52  // cloud physics
53  void AdvanceSatAdj (const SolverChoice& /*solverChoice*/);
54 
55  // Set up for first time
56  void
57  Define (SolverChoice& sc) override
58  {
59  m_fac_cond = lcond / sc.c_p;
60  m_rdOcp = sc.rdOcp;
61  m_do_cond = (!sc.use_shoc);
62  }
63 
64  // init
65  void
66  Init (const amrex::MultiFab& cons_in,
67  const amrex::BoxArray& /*grids*/,
68  const amrex::Geometry& geom,
69  const amrex::Real& dt_advance,
70  std::unique_ptr<amrex::MultiFab>& /*z_phys_nd*/,
71  std::unique_ptr<amrex::MultiFab>& /*detJ_cc*/) override;
72 
73  // Copy state into micro vars
74  void
75  Copy_State_to_Micro (const amrex::MultiFab& cons_in) override;
76 
77  // Copy state into micro vars
78  void
79  Copy_Micro_to_State (amrex::MultiFab& cons_in) override;
80 
81  // update micro vars
82  void
83  Update_Micro_Vars (amrex::MultiFab& cons_in) override
84  {
85  this->Copy_State_to_Micro(cons_in);
86  }
87 
88  // update state vars
89  void
90  Update_State_Vars (amrex::MultiFab& cons_in) override
91  {
92  this->Copy_Micro_to_State(cons_in);
93  }
94 
95  // wrapper to advance micro vars
96  void
97  Advance (const amrex::Real& dt_advance,
98  const SolverChoice& solverChoice) override
99  {
100  dt = dt_advance;
101 
102  this->AdvanceSatAdj(solverChoice);
103  }
104 
105  amrex::MultiFab*
106  Qmoist_Ptr (const int& varIdx) override
107  {
108  AMREX_ALWAYS_ASSERT(varIdx < m_qmoist_size);
109  return nullptr;
110  }
111 
112  int
113  Qmoist_Size () override { return SatAdj::m_qmoist_size; }
114 
115  int
117 
118  void
120  std::vector<int>& a_idx,
121  std::vector<std::string>& a_names) const override
122  {
123  a_idx.clear();
124  a_names.clear();
125  }
126 
127  AMREX_GPU_HOST_DEVICE
128  AMREX_FORCE_INLINE
129  static amrex::Real
130  NewtonIterSat (int& i, int& j, int& k,
131  const amrex::Real& fac_cond,
132  const amrex::Array4<amrex::Real>& tabs_array,
133  const amrex::Array4<amrex::Real>& pres_array,
134  const amrex::Array4<amrex::Real>& qv_array,
135  const amrex::Array4<amrex::Real>& qc_array)
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  }
191 
192 private:
193  // Number of qmoist variables (no precipitating comps to accumulate)
194  int m_qmoist_size = 0;
195 
196  // Number of qstate variables
198 
199  // geometry
200  amrex::Geometry m_geom;
201 
202  // timestep
204 
205  // constants
208  bool m_do_cond;
209 
210  // independent variables
211  amrex::Array<FabPtr, MicVar_SatAdj::NumVars> mic_fab_vars;
212 };
213 #endif
constexpr amrex::Real lcond
Definition: ERF_Constants.H:66
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:181
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:166
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_NullMoist.H:8
Definition: ERF_SatAdj.H:41
int n_qstate_moist_size
Definition: ERF_SatAdj.H:197
amrex::Real m_fac_cond
Definition: ERF_SatAdj.H:206
int Qmoist_Size() override
Definition: ERF_SatAdj.H:113
void Update_Micro_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_SatAdj.H:83
amrex::Geometry m_geom
Definition: ERF_SatAdj.H:200
void Qmoist_Restart_Vars(const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
Definition: ERF_SatAdj.H:119
amrex::Array< FabPtr, MicVar_SatAdj::NumVars > mic_fab_vars
Definition: ERF_SatAdj.H:211
amrex::Real m_rdOcp
Definition: ERF_SatAdj.H:207
SatAdj()
Definition: ERF_SatAdj.H:47
virtual ~SatAdj()=default
amrex::MultiFab * Qmoist_Ptr(const int &varIdx) override
Definition: ERF_SatAdj.H:106
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
int Qstate_Moist_Size() override
Definition: ERF_SatAdj.H:116
std::shared_ptr< amrex::MultiFab > FabPtr
Definition: ERF_SatAdj.H:43
void Update_State_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_SatAdj.H:90
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSatAdj.cpp:12
bool m_do_cond
Definition: ERF_SatAdj.H:208
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
Definition: ERF_InitSatAdj.cpp:17
void Advance(const amrex::Real &dt_advance, const SolverChoice &solverChoice) override
Definition: ERF_SatAdj.H:97
int m_qmoist_size
Definition: ERF_SatAdj.H:194
void AdvanceSatAdj(const SolverChoice &)
Definition: ERF_SatAdj.cpp:8
amrex::Real dt
Definition: ERF_SatAdj.H:203
void Define(SolverChoice &sc) override
Definition: ERF_SatAdj.H:57
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSatAdj.cpp:40
Definition: ERF_SatAdj.H:27
@ theta
Definition: ERF_SatAdj.H:31
@ rho
Definition: ERF_SatAdj.H:30
@ pres
Definition: ERF_SatAdj.H:33
@ NumVars
Definition: ERF_SatAdj.H:37
@ qv
Definition: ERF_SatAdj.H:35
@ tabs
Definition: ERF_SatAdj.H:32
@ qc
Definition: ERF_SatAdj.H:36
Definition: ERF_DataStruct.H:123
amrex::Real rdOcp
Definition: ERF_DataStruct.H:913
amrex::Real c_p
Definition: ERF_DataStruct.H:912
bool use_shoc
Definition: ERF_DataStruct.H:942