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,
91  const amrex::MultiFab& /*z_phys_nd*/) override
92  {
93  this->Copy_Micro_to_State(cons_in);
94  }
95 
96  // wrapper to advance micro vars
97  void
98  Advance (const amrex::Real& dt_advance,
99  const SolverChoice& solverChoice) override
100  {
101  dt = dt_advance;
102 
103  this->AdvanceSatAdj(solverChoice);
104  }
105 
106  amrex::MultiFab*
107  Qmoist_Ptr (const int& varIdx) override
108  {
110  return nullptr;
111  }
112 
113  int
114  Qmoist_Size () override { return SatAdj::m_qmoist_size; }
115 
116  int
118 
119  void
121  std::vector<int>& a_idx,
122  std::vector<std::string>& a_names) const override
123  {
124  a_idx.clear();
125  a_names.clear();
126  }
127 
128  void
129  Set_RealWidth (const int real_width) override { m_real_width = real_width; }
130 
131  AMREX_GPU_HOST_DEVICE
132  AMREX_FORCE_INLINE
133  static amrex::Real
134  NewtonIterSat (int& i, int& j, int& k,
135  const amrex::Real& fac_cond,
136  const amrex::Array4<amrex::Real>& tabs_array,
137  const amrex::Array4<amrex::Real>& pres_array,
138  const amrex::Array4<amrex::Real>& qv_array,
139  const amrex::Array4<amrex::Real>& qc_array)
140  {
141  // Solution tolerance
142  amrex::Real tol = amrex::Real(1.0e-8);
143 
144  // Saturation moisture fractions
145  amrex::Real qsat, dqsat;
146 
147  // Newton iteration vars
148  int niter;
149  amrex::Real fff, dfff;
150  amrex::Real dtabs, delta_qv;
151 
152  // Initial guess for temperature & pressure
153  amrex::Real tabs = tabs_array(i,j,k);
154  amrex::Real pres = pres_array(i,j,k);
155 
156  niter = 0;
157  dtabs = 1;
158 
159  //==================================================
160  // Newton iteration to qv=qsat (cloud phase only)
161  //==================================================
162  do {
163  // Saturation moisture fractions
164  erf_qsatw(tabs, pres, qsat);
165  erf_dtqsatw(tabs, pres, dqsat);
166 
167  // Function for root finding:
168  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
169  fff = -tabs + tabs_array(i,j,k) + fac_cond*(qv_array(i,j,k) - qsat);
170 
171  // Derivative of function (T_new iterated on)
172  dfff = -one - fac_cond*dqsat;
173 
174  // Update the temperature
175  dtabs = -fff/dfff;
176  tabs += dtabs;
177 
178  // Update iteration
179  niter = niter+1;
180  } while(std::abs(dtabs) > tol && niter < 20);
181 
182  // Update qsat from last iteration (dq = dq/dt * dt)
183  qsat += dqsat*dtabs;
184 
185  // Changes in each component
186  delta_qv = qv_array(i,j,k) - qsat;
187 
188  // Partition the change in non-precipitating q
189  qv_array(i,j,k) = qsat;
190  qc_array(i,j,k) += delta_qv;
191 
192  // Return to temperature
193  return tabs;
194  }
195 
196 private:
197  // Number of qmoist variables (no precipitating comps to accumulate)
198  int m_qmoist_size = 0;
199 
200  // Number of qstate variables
202 
203  // geometry
204  amrex::Geometry m_geom;
205 
206  // Nudging + set width
207  int m_real_width{0};
208 
209  // timestep
211 
212  // constants
215  bool m_do_cond;
216 
217  // independent variables
218  amrex::Array<FabPtr, MicVar_SatAdj::NumVars> mic_fab_vars;
219 };
220 #endif
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real lcond
Definition: ERF_Constants.H:77
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:186
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:171
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:201
amrex::Real m_fac_cond
Definition: ERF_SatAdj.H:213
int Qmoist_Size() override
Definition: ERF_SatAdj.H:114
void Update_Micro_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_SatAdj.H:83
amrex::Geometry m_geom
Definition: ERF_SatAdj.H:204
void Qmoist_Restart_Vars(const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
Definition: ERF_SatAdj.H:120
amrex::Array< FabPtr, MicVar_SatAdj::NumVars > mic_fab_vars
Definition: ERF_SatAdj.H:218
amrex::Real m_rdOcp
Definition: ERF_SatAdj.H:214
SatAdj()
Definition: ERF_SatAdj.H:47
virtual ~SatAdj()=default
void Update_State_Vars(amrex::MultiFab &cons_in, const amrex::MultiFab &) override
Definition: ERF_SatAdj.H:90
amrex::MultiFab * Qmoist_Ptr(const int &varIdx) override
Definition: ERF_SatAdj.H:107
void Set_RealWidth(const int real_width) override
Definition: ERF_SatAdj.H:129
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:134
int Qstate_Moist_Size() override
Definition: ERF_SatAdj.H:117
std::shared_ptr< amrex::MultiFab > FabPtr
Definition: ERF_SatAdj.H:43
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSatAdj.cpp:12
bool m_do_cond
Definition: ERF_SatAdj.H:215
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:98
int m_qmoist_size
Definition: ERF_SatAdj.H:198
void AdvanceSatAdj(const SolverChoice &)
Definition: ERF_SatAdj.cpp:8
amrex::Real dt
Definition: ERF_SatAdj.H:210
void Define(SolverChoice &sc) override
Definition: ERF_SatAdj.H:57
int m_real_width
Definition: ERF_SatAdj.H:207
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:130
amrex::Real rdOcp
Definition: ERF_DataStruct.H:1169
amrex::Real c_p
Definition: ERF_DataStruct.H:1168
bool use_shoc
Definition: ERF_DataStruct.H:1200