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  * SatAdj is a cell-local, fixed-pressure saturation-adjustment source step.
6  * It is not a flux-form finite-volume update.
7  *
8  * The scheme reads rho, rhoTheta, rhoQv, and rhoQc from the conserved state,
9  * adjusts theta/qv/qc locally, and writes rhoTheta/rhoQv/rhoQc back.
10  * Density is copied for diagnostics and conservative reconstruction, but rho is
11  * not a SatAdj prognostic variable and is not modified by this module.
12  */
13 
14 #include <string>
15 #include <vector>
16 #include <memory>
17 
18 #include <AMReX_FArrayBox.H>
19 #include <AMReX_Geometry.H>
20 #include <AMReX_MultiFabUtil.H>
21 #include <AMReX_GpuContainers.H>
22 
23 #include "ERF_EOS.H"
24 #include "ERF_Constants.H"
25 #include "ERF_MicrophysicsUtils.H"
26 #include "ERF_IndexDefines.H"
27 #include "ERF_DataStruct.H"
28 #include "ERF_NullMoist.H"
29 #include "ERF_TileNoZ.H"
30 
31 namespace MicVar_SatAdj {
32  enum {
33  // independent variables
34  rho=0, // density copied from conserved state; not modified by SatAdj
35  theta, // dry potential temperature
36  tabs, // absolute temperature [K]
37  pres, // pressure [mbar/hPa] for saturation helpers
38  // non-precipitating vars
39  qv, // water-vapor mixing ratio
40  qc, // cloud-water mixing ratio
41  NumVars
42  };
43 }
44 
45 class SatAdj : public NullMoist {
46 
47  using FabPtr = std::shared_ptr<amrex::MultiFab>;
48 
49 public:
50  // constructor
51  SatAdj () {}
52 
53  // destructor
54  virtual ~SatAdj () = default;
55 
56  // cloud physics
57  void AdvanceSatAdj (const SolverChoice& /*solverChoice*/);
58 
59  // Set up for first time
60  void
61  Define (SolverChoice& sc) override
62  {
63  m_fac_cond = lcond / sc.c_p;
64  m_rdOcp = sc.rdOcp;
65  m_do_cond = (!sc.use_shoc);
66  }
67 
68  // init
69  void
70  Init (const amrex::MultiFab& cons_in,
71  const amrex::BoxArray& /*grids*/,
72  const amrex::Geometry& geom,
73  const amrex::Real& dt_advance,
74  std::unique_ptr<amrex::MultiFab>& /*z_phys_nd*/,
75  std::unique_ptr<amrex::MultiFab>& /*detJ_cc*/) override;
76 
77  // Copy state into micro vars
78  void
79  Copy_State_to_Micro (const amrex::MultiFab& cons_in) override;
80 
81  // Copy state into micro vars
82  void
83  Copy_Micro_to_State (amrex::MultiFab& cons_in) override;
84 
85  // update micro vars
86  void
87  Update_Micro_Vars (amrex::MultiFab& cons_in) override
88  {
89  this->Copy_State_to_Micro(cons_in);
90  }
91 
92  // update state vars
93  void
94  Update_State_Vars (amrex::MultiFab& cons_in,
95  const amrex::MultiFab& /*z_phys_nd*/) override
96  {
97  this->Copy_Micro_to_State(cons_in);
98  }
99 
100  // wrapper to advance micro vars
101  void
102  Advance (const amrex::Real& dt_advance,
103  const SolverChoice& solverChoice) override
104  {
105  dt = dt_advance;
106 
107  this->AdvanceSatAdj(solverChoice);
108  }
109 
110  amrex::MultiFab*
111  Qmoist_Ptr (const int& varIdx) override
112  {
114  return nullptr;
115  }
116 
117  int
118  Qmoist_Size () override { return SatAdj::m_qmoist_size; }
119 
120  int
122 
123  int
125 
126  void
128  std::vector<int>& a_idx,
129  std::vector<std::string>& a_names) const override
130  {
131  a_idx.clear();
132  a_names.clear();
133  }
134 
135  AMREX_GPU_HOST_DEVICE
136  AMREX_FORCE_INLINE
137  static amrex::Real
138  // Newton solve for the saturated final temperature at fixed pressure.
139  // The invariant is T + (L/cp) * qv for the local cell. The solve finds
140  // T_new such that qv_new = qsat(T_new, p). The pressure argument is held
141  // fixed in mbar/hPa throughout the solve.
142  //
143  // Solves the fixed-pressure saturation equation
144  // 0 = -T_new + T_old + (L/cp) * (qv_old - qsat(T_new, p))
145  // for T_new. On return, qsat is evaluated exactly at the returned
146  // temperature. The caller uses that qsat to set qv and qc while
147  // conserving qv + qc. Density does not enter this solve.
148  // If the Newton iteration reaches max_niter before the tolerance, the
149  // routine returns the last iterate and recomputes qsat at that temperature.
151  const amrex::Real tabs_old,
152  const amrex::Real pres_mbar,
153  const amrex::Real qv_old,
154  amrex::Real& qsat)
155  {
156  constexpr amrex::Real tol = amrex::Real(1.0e-8);
157  constexpr int max_niter = 20;
158 
159  amrex::Real dqsat;
160 
161  int niter = 0;
162  amrex::Real fff, dfff;
163  amrex::Real dtabs = one;
164 
165  amrex::Real tabs = tabs_old;
166 
167  //==================================================
168  // Newton iteration to qv=qsat (cloud phase only)
169  //==================================================
170  do {
171  // Saturation moisture fractions
172  erf_qsatw(tabs, pres_mbar, qsat);
173  erf_dtqsatw(tabs, pres_mbar, dqsat);
174 
175  // Function for root finding:
176  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
177  fff = -tabs + tabs_old + fac_cond*(qv_old - qsat);
178 
179  // Derivative of function (T_new iterated on)
180  dfff = -one - fac_cond*dqsat;
181 
182  // Update the temperature
183  dtabs = -fff/dfff;
184  tabs += dtabs;
185 
186  // Update iteration
187  niter = niter+1;
188  } while(std::abs(dtabs) > tol && niter < max_niter);
189 
190  // Return qsat evaluated at the final temperature, not a tangent extrapolation.
191  erf_qsatw(tabs, pres_mbar, qsat);
192 
193  return tabs;
194  }
195 
196  AMREX_GPU_HOST_DEVICE
197  AMREX_FORCE_INLINE
198  static void
199  // Cell-local fixed-pressure saturation adjustment.
200  // Inputs/outputs are specific quantities for one cell. The density is not
201  // used or modified. The update conserves qv + qc and, at the pressure
202  // passed in, conserves T + (L/cp) * qv. The final state should satisfy
203  // qc >= 0 and either qc == 0 or qv == qsat(T,p).
204  AdjustSatAdjCell (const amrex::Real fac_cond,
205  const amrex::Real rdOcp,
206  amrex::Real& tabs,
207  const amrex::Real pres_mbar,
209  amrex::Real& qv,
210  amrex::Real& qc)
211  {
212  amrex::Real qsat;
213  erf_qsatw(tabs, pres_mbar, qsat);
214 
215  if ((qv + qc) > qsat) {
216  // Total water exceeds saturation at the initial temperature, so the final state
217  // is cloudy and saturated. Solve directly for the final saturated temperature.
218  // This is algebraically equivalent to evaporating existing cloud water first
219  // and then recondensing to the fixed-pressure equilibrium state.
220 #if defined(AMREX_DEBUG)
221  const amrex::Real qvprev = qv;
222  const amrex::Real qcprev = qc;
223 #endif
224 
225  if (qc < 0) {
226  // Repair negative cloud water by transferring the deficit back to vapor
227  // before solving the saturated equilibrium and clamping qc to zero.
228  qv += qc;
229  qc = zero;
230  }
231 
232  tabs = NewtonSolveSatTemperature(fac_cond, tabs, pres_mbar, qv, qsat);
233 
234  const amrex::Real delta_qv = qv - qsat;
235  qv = qsat;
236  qc += delta_qv;
237 
238 #if defined(AMREX_DEBUG)
239  amrex::Real qsatnew;
240  erf_qsatw(tabs, pres_mbar, qsatnew);
241  AMREX_ASSERT(std::abs(qv-qsatnew) < 1e-12);
242  AMREX_ASSERT(std::abs(qv+qc-qvprev-qcprev) < 1e-14);
243 #endif
244 
245  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
246  } else {
247  // Total water does not exceed initial saturation. Evaporate all cloud water,
248  // cooling the cell by latent heat absorption. The cooled cell can become
249  // supersaturated, so recondense if needed.
250  const amrex::Real delta_qc = qc;
251 
252  qv += qc;
253  qc = zero;
254 
255  tabs -= fac_cond * delta_qc;
256  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
257 
258  erf_qsatw(tabs, pres_mbar, qsat);
259  if (qv > qsat) {
260 #if defined(AMREX_DEBUG)
261  const amrex::Real qvprev = qv;
262  const amrex::Real qcprev = qc;
263  const amrex::Real Tprev = tabs;
264 #endif
265 
266  tabs = NewtonSolveSatTemperature(fac_cond, tabs, pres_mbar, qv, qsat);
267 
268  const amrex::Real delta_qv = qv - qsat;
269  qv = qsat;
270  qc += delta_qv;
271 
272 #if defined(AMREX_DEBUG)
273  amrex::Real qsatnew;
274  erf_qsatw(tabs, pres_mbar, qsatnew);
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);
280 #endif
281 
282  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
283  }
284  }
285  }
286 
287 private:
288  // Number of qmoist variables (no precipitating comps to accumulate)
289  int m_qmoist_size = 0;
290 
291  // Number of qstate variables
293 
294  // Number of qstate variables that are number concentrations
296 
297  // geometry
298  amrex::Geometry m_geom;
299 
300  // timestep
302 
303  // constants
306  bool m_do_cond;
307 
308  // independent variables
309  amrex::Array<FabPtr, MicVar_SatAdj::NumVars> mic_fab_vars;
310 };
311 #endif
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real lcond
Definition: ERF_Constants.H:85
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_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:244
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
Definition: ERF_NullMoist.H:8
Definition: ERF_SatAdj.H:45
int n_qstate_moist_size
Definition: ERF_SatAdj.H:292
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
int Qmoist_Size() override
Definition: ERF_SatAdj.H:118
void Update_Micro_Vars(amrex::MultiFab &cons_in) override
Definition: ERF_SatAdj.H:87
amrex::Geometry m_geom
Definition: ERF_SatAdj.H:298
void Qmoist_Restart_Vars(const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
Definition: ERF_SatAdj.H:127
amrex::Array< FabPtr, MicVar_SatAdj::NumVars > mic_fab_vars
Definition: ERF_SatAdj.H:309
amrex::Real m_rdOcp
Definition: ERF_SatAdj.H:305
SatAdj()
Definition: ERF_SatAdj.H:51
virtual ~SatAdj()=default
void Update_State_Vars(amrex::MultiFab &cons_in, const amrex::MultiFab &) override
Definition: ERF_SatAdj.H:94
amrex::MultiFab * Qmoist_Ptr(const int &varIdx) override
Definition: ERF_SatAdj.H:111
int n_qstate_moist_numconc_size
Definition: ERF_SatAdj.H:295
int Qstate_Moist_NumConc_Size() override
Definition: ERF_SatAdj.H:124
int Qstate_Moist_Size() override
Definition: ERF_SatAdj.H:121
std::shared_ptr< amrex::MultiFab > FabPtr
Definition: ERF_SatAdj.H:47
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSatAdj.cpp:13
bool m_do_cond
Definition: ERF_SatAdj.H:306
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:102
int m_qmoist_size
Definition: ERF_SatAdj.H:289
void AdvanceSatAdj(const SolverChoice &)
Definition: ERF_SatAdj.cpp:11
amrex::Real dt
Definition: ERF_SatAdj.H:301
void Define(SolverChoice &sc) override
Definition: ERF_SatAdj.H:61
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
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSatAdj.cpp:44
Definition: ERF_SatAdj.H:31
@ theta
Definition: ERF_SatAdj.H:35
@ rho
Definition: ERF_SatAdj.H:34
@ pres
Definition: ERF_SatAdj.H:37
@ NumVars
Definition: ERF_SatAdj.H:41
@ qv
Definition: ERF_SatAdj.H:39
@ tabs
Definition: ERF_SatAdj.H:36
@ qc
Definition: ERF_SatAdj.H:40
Definition: ERF_DataStruct.H:141
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