ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_ABLMost.H
Go to the documentation of this file.
1 #ifndef ERF_ABLMOST_H
2 #define ERF_ABLMOST_H
3 
4 #include <AMReX_Geometry.H>
5 #include <AMReX_ParmParse.H>
6 #include <AMReX_FArrayBox.H>
7 #include <AMReX_MultiFab.H>
8 #include <AMReX_iMultiFab.H>
9 #include <AMReX_Interpolater.H>
10 
11 #include <ERF_IndexDefines.H>
12 #include <ERF_Constants.H>
13 #include <ERF_MOSTAverage.H>
14 #include <ERF_MOSTStress.H>
15 #include <ERF_TerrainMetrics.H>
16 #include <ERF_PBLHeight.H>
17 
18 /** Monin-Obukhov surface layer profile
19  *
20  * van der Laan, P., Kelly, M. C., & Sørensen, N. N. (2017). A new k-epsilon
21  * model consistent with Monin-Obukhov similarity theory. Wind Energy,
22  * 20(3), 479–489. https://doi.org/10.1002/we.2017
23  *
24  * Consistent with Dyer (1974) formulation from page 57, Chapter 2, Modeling
25  * the vertical ABL structure in Modelling of Atmospheric Flow Fields,
26  * Demetri P Lalas and Corrado F Ratto, January 1996,
27  * https://doi.org/10.1142/2975.
28  */
29 class ABLMost
30 {
31 
32 public:
33 
34  // Constructor
35  explicit ABLMost (const amrex::Vector<amrex::Geometry>& geom,
36  bool& use_exp_most,
37  bool& use_rot_most,
38  amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_old,
39  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Theta_prim,
40  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Qv_prim,
41  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Qr_prim,
42  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& z_phys_nd,
43  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>>& sst_lev,
44  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::iMultiFab>>>& lmask_lev,
45  amrex::Vector<amrex::Vector<amrex::MultiFab*>> lsm_data,
46  amrex::Vector<amrex::Vector<amrex::MultiFab*>> lsm_flux,
47  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Hwave,
48  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Lwave,
49  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& eddyDiffs,
50  amrex::Real start_bdy_time = 0.0,
51  amrex::Real bdy_time_interval = 0.0)
52  : m_exp_most(use_exp_most),
53  m_rotate(use_rot_most),
54  m_start_bdy_time(start_bdy_time),
55  m_bdy_time_interval(bdy_time_interval),
56  m_geom(geom),
57  m_ma(geom,vars_old,Theta_prim,Qv_prim,Qr_prim,z_phys_nd)
58  {
59  // We have a moisture model if Qv_prim is a valid pointer
60  use_moisture = (Qv_prim[0].get());
61 
62  // Get roughness
63  amrex::ParmParse pp("erf");
64  pp.query("most.z0", z0_const);
65 
66  // Specify how to compute the flux
67  if (use_rot_most) {
69  } else {
70  std::string flux_string{"moeng"};
71  pp.query("most.flux_type", flux_string);
72  if (flux_string == "donelan") {
74  } else if (flux_string == "moeng") {
76  } else if (flux_string == "custom") {
78  } else {
79  amrex::Abort("Undefined MOST flux type!");
80  }
81  }
82 
83  // Include w* to handle free convection (Beljaars 1995, QJRMS)
84  pp.query("most.include_wstar", m_include_wstar);
85 
86  std::string pblh_string{"none"};
87  pp.query("most.pblh_calc", pblh_string);
88  if (pblh_string == "none") {
90  } else if (pblh_string == "MYNN2.5") {
92  } else if (pblh_string == "YSU") {
94  } else {
95  amrex::Abort("Undefined PBLH calc type!");
96  }
97 
98  // Get surface temperature
99  auto erf_st = pp.query("most.surf_temp", surf_temp);
100 
101  // Custom type user must specify the fluxes
104  pp.query("most.ustar", custom_ustar);
105  pp.query("most.tstar", custom_tstar);
106  pp.query("most.qstar", custom_qstar);
107  if (custom_qstar != 0) {
108  AMREX_ASSERT_WITH_MESSAGE(use_moisture, "Specified custom MOST qv flux without moisture model!");
109  }
110  amrex::Print() << "Using specified ustar, tstar, qstar for MOST = "
111  << custom_ustar << " "
112  << custom_tstar << " "
113  << custom_qstar << std::endl;
114 
115  // Specify surface temperature or surface flux
116  } else {
117  if (erf_st) {
119  pp.query("most.surf_heating_rate", surf_heating_rate); // [K/h]
120  surf_heating_rate = surf_heating_rate / 3600.0; // [K/s]
121  if (pp.query("most.surf_temp_flux", surf_temp_flux)) {
122  amrex::Abort("Can only specify one of surf_temp_flux or surf_heating_rate");
123  }
124  } else {
125  pp.query("most.surf_temp_flux", surf_temp_flux);
126  if (pp.query("most.surf_heating_rate", surf_heating_rate)) {
127  amrex::Abort("Can only specify one of surf_temp_flux or surf_heating_rate");
128  }
129  if (std::abs(surf_temp_flux) > std::numeric_limits<amrex::Real>::epsilon()) {
131  } else {
133  }
134  }
135  }
136 
137  // Specify how to compute the surface flux over land (if there is any)
138  std::string rough_land_string{"constant"};
139  pp.query("most.roughness_type_land", rough_land_string);
140  if (rough_land_string == "constant") {
142  } else {
143  amrex::Abort("Undefined MOST roughness type for land!");
144  }
145 
146  // Specify how to compute the surface flux over sea (if there is any)
147  std::string rough_sea_string{"charnock"};
148  pp.query("most.roughness_type_sea", rough_sea_string);
149  if (rough_sea_string == "charnock") {
151  pp.query("most.charnock_constant",cnk_a);
152  pp.query("most.charnock_viscosity",cnk_visc);
153  if (cnk_a > 0) {
154  amrex::Print() << "If there is water, Charnock relation with C_a=" << cnk_a
155  << (cnk_visc? " and viscosity" : "")
156  << " will be used"
157  << std::endl;
158  } else {
159  amrex::Print() << "If there is water, Charnock relation with variable Charnock parameter (COARE3.0)"
160  << (cnk_visc? " and viscosity" : "")
161  << " will be used"
162  << std::endl;
163  }
164  } else if (rough_sea_string == "coare3.0") {
166  amrex::Print() << "If there is water, Charnock relation with variable Charnock parameter (COARE3.0)"
167  << (cnk_visc? " and viscosity" : "")
168  << " will be used"
169  << std::endl;
170  cnk_a = -1;
171  } else if (rough_sea_string == "donelan") {
173  } else if (rough_sea_string == "modified_charnock") {
175  pp.query("most.modified_charnock_depth",depth);
176  } else if (rough_sea_string == "wave_coupled") {
178  } else if (rough_sea_string == "constant") {
180  } else {
181  amrex::Abort("Undefined MOST roughness type for sea!");
182  }
183 
184  // Check if there is a user-specified roughness file to be read
185  std::string fname;
186  bool read_z0 = false;
187  if ( (flux_type == FluxCalcType::MOENG) ||
189  read_z0 = pp.query("most.roughness_file_name", fname);
190  }
191 
192  // Do we have a time-varying surface roughness that needs to be saved?
193  const int nghost = 0; // ghost cells not included
194  const int level = 0;
195  int lmask_min = lmask_min_reduce(*lmask_lev[level][0].get(), nghost);
196  amrex::ParallelDescriptor::ReduceIntMin(lmask_min);
197 
198  m_var_z0 = (lmask_min < 1) & (rough_type_sea != RoughCalcType::CONSTANT);
199  if (m_var_z0) {
200  amrex::Print() << "Variable sea roughness (type " << rough_sea_string << ")" << std::endl;
201  }
202 
203  // Size the MOST params for all levels
204  int nlevs = m_geom.size();
205  z_0.resize(nlevs);
206  u_star.resize(nlevs);
207  w_star.resize(nlevs);
208  t_star.resize(nlevs);
209  q_star.resize(nlevs);
210  t_surf.resize(nlevs);
211  olen.resize(nlevs);
212  pblh.resize(nlevs);
213 
214  // Get pointers to SST and LANDMASK data
215  m_sst_lev.resize(nlevs);
216  m_lmask_lev.resize(nlevs);
217 
218  for (int lev(0); lev<nlevs; ++lev) {
219  int nt_tot_sst = sst_lev[lev].size();
220  m_sst_lev[lev].resize(nt_tot_sst);
221  for (int nt(0); nt<nt_tot_sst; ++nt) {
222  m_sst_lev[lev][nt] = sst_lev[lev][nt].get();
223  }
224  int nt_tot_lmask = lmask_lev[lev].size();
225  m_lmask_lev[lev].resize(nt_tot_lmask);
226  for (int nt(0); nt<nt_tot_lmask; ++nt) {
227  m_lmask_lev[lev][nt] = lmask_lev[lev][nt].get();
228  }
229  } // lev
230 
231  // Get pointers to LSM data and Fluxes
232  m_lsm_data_lev.resize(nlevs);
233  m_lsm_flux_lev.resize(nlevs);
234  for (int lev(0); lev<nlevs; ++lev) {
235  int nvar = lsm_data[lev].size();
236  m_lsm_data_lev[lev].resize(nvar);
237  m_lsm_flux_lev[lev].resize(nvar);
238  for (int n(0); n<nvar; ++n) {
239  m_lsm_data_lev[lev][n] = lsm_data[lev][n];
240  m_lsm_flux_lev[lev][n] = lsm_flux[lev][n];
241  }
242  } // lev
243 
244  // Get pointers to wave data
245  m_Hwave_lev.resize(nlevs);
246  m_Lwave_lev.resize(nlevs);
247  m_eddyDiffs_lev.resize(nlevs);
248  for (int lev(0); lev<nlevs; ++lev) {
249  m_Hwave_lev[lev] = Hwave[lev].get();
250  m_Lwave_lev[lev] = Lwave[lev].get();
251  m_eddyDiffs_lev[lev] = eddyDiffs[lev].get();
252  }
253 
254  for (int lev = 0; lev < nlevs; lev++) {
255  // Attributes for MFs and FABs
256  //--------------------------------------------------------
257  auto& mf = vars_old[lev][Vars::cons];
258  // Create a 2D ba, dm, & ghost cells
259  amrex::BoxArray ba = mf.boxArray();
260  amrex::BoxList bl2d = ba.boxList();
261  for (auto& b : bl2d) {
262  b.setRange(2,0);
263  }
264  amrex::BoxArray ba2d(std::move(bl2d));
265  const amrex::DistributionMapping& dm = mf.DistributionMap();
266  const int ncomp = 1;
267  amrex::IntVect ng = mf.nGrowVect(); ng[2]=0;
268 
269  // Z0 heights FAB
270  //--------------------------------------------------------
271  amrex::Arena* Arena_Used = amrex::The_Arena();
272 #ifdef AMREX_USE_GPU
273  Arena_Used = amrex::The_Pinned_Arena();
274 #endif
275  amrex::Box bx = amrex::grow(m_geom[lev].Domain(),ng);
276  bx.setSmall(2,0);
277  bx.setBig(2,0);
278  z_0[lev].resize(bx,1,Arena_Used);
279  z_0[lev].setVal<amrex::RunOn::Device>(z0_const);
280  if (read_z0) read_custom_roughness(lev,fname);
281 
282  // 2D MFs for U*, T*, T_surf
283  //--------------------------------------------------------
284  u_star[lev] = std::make_unique<amrex::MultiFab>(ba2d,dm,ncomp,ng);
285  u_star[lev]->setVal(1.E34);
286 
287  w_star[lev] = std::make_unique<amrex::MultiFab>(ba2d,dm,ncomp,ng);
288  w_star[lev]->setVal(1.E34);
289 
290  t_star[lev] = std::make_unique<amrex::MultiFab>(ba2d,dm,ncomp,ng);
291  t_star[lev]->setVal(0.0); // default to neutral
292 
293  q_star[lev] = std::make_unique<amrex::MultiFab>(ba2d,dm,ncomp,ng);
294  q_star[lev]->setVal(0.0); // default to dry
295 
296  olen[lev] = std::make_unique<amrex::MultiFab>(ba2d,dm,ncomp,ng);
297  olen[lev]->setVal(1.E34);
298 
299  pblh[lev] = std::make_unique<amrex::MultiFab>(ba2d,dm,ncomp,ng);
300  pblh[lev]->setVal(1.E34);
301 
302  t_surf[lev] = std::make_unique<amrex::MultiFab>(ba2d,dm,ncomp,ng);
303 
304  // TODO: Do we want lsm_data to have theta at 0 index always?
305  // Do we want an external enum struct for indexing?
306  if (m_sst_lev[lev][0] || m_lsm_data_lev[lev][0]) {
307  // Valid SST or LSM data; t_surf set before computing fluxes (avoids extended lambda capture)
309  } else if (erf_st) {
310  t_surf[lev]->setVal(surf_temp);
311  } else {
312  t_surf[lev]->setVal(0.0);
313  }
314  }// lev
315  }
316 
317  void
318  update_fluxes (const int& lev,
319  const amrex::Real& time,
320  int max_iters = 25);
321 
322  template <typename FluxIter>
323  void
324  compute_fluxes (const int& lev,
325  const int& max_iters,
326  const FluxIter& most_flux,
327  bool is_land);
328 
329  void
330  impose_most_bcs (const int& lev,
331  const amrex::Vector<amrex::MultiFab*>& mfs,
332  amrex::MultiFab* xxmom_flux,
333  amrex::MultiFab* yymom_flux,
334  amrex::MultiFab* zzmom_flux,
335  amrex::MultiFab* xymom_flux, amrex::MultiFab* yxmom_flux,
336  amrex::MultiFab* xzmom_flux, amrex::MultiFab* zxmom_flux,
337  amrex::MultiFab* yzmom_flux, amrex::MultiFab* zymom_flux,
338  amrex::MultiFab* xheat_flux,
339  amrex::MultiFab* yheat_flux,
340  amrex::MultiFab* zheat_flux,
341  amrex::MultiFab* xqv_flux,
342  amrex::MultiFab* yqv_flux,
343  amrex::MultiFab* zqv_flux,
344  amrex::MultiFab* z_phys);
345 
346  template <typename FluxCalc>
347  void
348  compute_most_bcs (const int& lev,
349  const amrex::Vector<amrex::MultiFab*>& mfs,
350  amrex::MultiFab* xxmom_flux,
351  amrex::MultiFab* yymom_flux,
352  amrex::MultiFab* zzmom_flux,
353  amrex::MultiFab* xymom_flux, amrex::MultiFab* yxmom_flux,
354  amrex::MultiFab* xzmom_flux, amrex::MultiFab* zxmom_flux,
355  amrex::MultiFab* yzmom_flux, amrex::MultiFab* zymom_flux,
356  amrex::MultiFab* xheat_flux,
357  amrex::MultiFab* yheat_flux,
358  amrex::MultiFab* zheat_flux,
359  amrex::MultiFab* xqv_flux,
360  amrex::MultiFab* yqv_flux,
361  amrex::MultiFab* zqv_flux,
362  amrex::MultiFab* z_phys,
363  const FluxCalc& flux_comp);
364 
365  void
366  time_interp_sst (const int& lev,
367  const amrex::Real& time);
368 
369  void
370  get_lsm_tsurf (const int& lev);
371 
372  /* wrapper around compute_pblh */
373  void
374  update_pblh (const int& lev,
375  amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars,
376  amrex::MultiFab* z_phys_cc,
377  const int RhoQv_comp,
378  const int RhoQc_comp,
379  const int RhoQr_comp);
380 
381  template <typename PBLHeightEstimator>
382  void
383  compute_pblh (const int& lev,
384  amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars,
385  amrex::MultiFab* z_phys_cc,
386  const PBLHeightEstimator& est,
387  const int RhoQv_comp,
388  const int RhoQc_comp,
389  const int RhoQr_comp);
390 
391  void
392  read_custom_roughness (const int& lev,
393  const std::string& fname);
394 
395  void
396  update_surf_temp (const amrex::Real& time)
397  {
398  if (surf_heating_rate != 0) {
399  int nlevs = m_geom.size();
400  for (int lev = 0; lev < nlevs; lev++)
401  {
402  t_surf[lev]->setVal(surf_temp + surf_heating_rate*time);
403  amrex::Print() << "Surface temp at t=" << time
404  << ": "
405  << surf_temp + surf_heating_rate*time
406  << std::endl;
407  }
408  }
409  }
410 
411  void
412  update_mac_ptrs (const int& lev,
413  amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_old,
414  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Theta_prim,
415  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Qv_prim,
416  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Qr_prim)
417  { m_ma.update_field_ptrs(lev,vars_old,Theta_prim,Qv_prim,Qr_prim); }
418 
419  const amrex::MultiFab*
420  get_u_star (const int& lev) { return u_star[lev].get(); }
421 
422  const amrex::MultiFab*
423  get_w_star (const int& lev) { return w_star[lev].get(); }
424 
425  const amrex::MultiFab*
426  get_t_star (const int& lev) { return t_star[lev].get(); }
427 
428  const amrex::MultiFab*
429  get_q_star (const int& lev) { return q_star[lev].get(); }
430 
431  const amrex::MultiFab*
432  get_olen (const int& lev) { return olen[lev].get(); }
433 
434  const amrex::MultiFab*
435  get_pblh (const int& lev) { return pblh[lev].get(); }
436 
437  const amrex::MultiFab*
438  get_mac_avg (const int& lev, int comp) { return m_ma.get_average(lev,comp); }
439 
440  const amrex::MultiFab*
441  get_t_surf (const int& lev) { return t_surf[lev].get(); }
442 
443  amrex::Real
444  get_zref () { return m_ma.get_zref(); }
445 
446  const amrex::FArrayBox*
447  get_z0 (const int& lev) { return &z_0[lev]; }
448 
450 
451  const amrex::iMultiFab*
452  get_lmask(const int& lev) { return m_lmask_lev[lev][0]; }
453 
454  int
455  lmask_min_reduce (amrex::iMultiFab& lmask, const int& nghost)
456  {
457  int lmask_min = amrex::ReduceMin(lmask, nghost,
458  [=] AMREX_GPU_HOST_DEVICE (amrex::Box const& bx, amrex::Array4<int const> const& lm_arr) -> int
459  {
460  int locmin = std::numeric_limits<int>::max();
461  const auto lo = lbound(bx);
462  const auto hi = ubound(bx);
463  for (int j = lo.y; j <= hi.y; ++j) {
464  for (int i = lo.x; i <= hi.x; ++i) {
465  locmin = std::min(locmin, lm_arr(i,j,0));
466  }
467  }
468  return locmin;
469  });
470 
471  return lmask_min;
472  }
473 
474  enum struct FluxCalcType {
475  MOENG = 0, ///< Moeng functional form
476  DONELAN, ///< Donelan functional form
477  CUSTOM, ///< Custom constant flux functional form
478  ROTATE ///< Terrain rotation flux functional form
479  };
480 
481  enum struct ThetaCalcType {
482  ADIABATIC = 0,
483  HEAT_FLUX, ///< Heat-flux specified
484  SURFACE_TEMPERATURE ///< Surface temperature specified
485  };
486 
487  enum struct RoughCalcType {
488  CONSTANT = 0, ///< Constant z0
489  CHARNOCK,
491  DONELAN,
493  };
494 
495  enum struct PBLHeightCalcType {
496  None, MYNN25, YSU
497  };
498 
504 
505 private:
507  bool m_exp_most = false;
508  bool m_rotate = false;
509  bool m_include_wstar = false;
510  amrex::Real z0_const{0.1};
511  amrex::Real surf_temp;
512  amrex::Real surf_heating_rate{0};
513  amrex::Real surf_temp_flux{0};
514  amrex::Real custom_ustar{0};
515  amrex::Real custom_tstar{0};
516  amrex::Real custom_qstar{0};
517  amrex::Real cnk_a{0.0185};
518  bool cnk_visc{false};
519  amrex::Real depth{30.0};
520  amrex::Real m_start_bdy_time;
521  amrex::Real m_bdy_time_interval;
522  amrex::Vector<amrex::Geometry> m_geom;
523  amrex::Vector<amrex::FArrayBox> z_0;
524  bool m_var_z0{false};
525 
527  amrex::Vector<std::unique_ptr<amrex::MultiFab>> u_star;
528  amrex::Vector<std::unique_ptr<amrex::MultiFab>> w_star;
529  amrex::Vector<std::unique_ptr<amrex::MultiFab>> t_star;
530  amrex::Vector<std::unique_ptr<amrex::MultiFab>> q_star;
531  amrex::Vector<std::unique_ptr<amrex::MultiFab>> olen;
532  amrex::Vector<std::unique_ptr<amrex::MultiFab>> pblh;
533  amrex::Vector<std::unique_ptr<amrex::MultiFab>> t_surf;
534 
535  amrex::Vector<amrex::Vector<amrex::MultiFab*>> m_sst_lev;
536  amrex::Vector<amrex::Vector<amrex::iMultiFab*>> m_lmask_lev;
537  amrex::Vector<amrex::Vector<amrex::MultiFab*>> m_lsm_data_lev;
538  amrex::Vector<amrex::Vector<amrex::MultiFab*>> m_lsm_flux_lev;
539  amrex::Vector<amrex::MultiFab*> m_Hwave_lev;
540  amrex::Vector<amrex::MultiFab*> m_Lwave_lev;
541  amrex::Vector<amrex::MultiFab*> m_eddyDiffs_lev;
542 };
543 
544 #endif /* ABLMOST_H */
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:219
Definition: ERF_ABLMost.H:30
amrex::Real m_bdy_time_interval
Definition: ERF_ABLMost.H:521
amrex::Vector< amrex::MultiFab * > m_eddyDiffs_lev
Definition: ERF_ABLMost.H:541
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_star
Definition: ERF_ABLMost.H:530
const amrex::MultiFab * get_u_star(const int &lev)
Definition: ERF_ABLMost.H:420
amrex::Vector< std::unique_ptr< amrex::MultiFab > > olen
Definition: ERF_ABLMost.H:531
const amrex::MultiFab * get_t_star(const int &lev)
Definition: ERF_ABLMost.H:426
amrex::Vector< std::unique_ptr< amrex::MultiFab > > w_star
Definition: ERF_ABLMost.H:528
amrex::Vector< amrex::FArrayBox > z_0
Definition: ERF_ABLMost.H:523
amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_ABLMost.H:522
amrex::Real custom_ustar
Definition: ERF_ABLMost.H:514
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_surf
Definition: ERF_ABLMost.H:533
amrex::Real custom_qstar
Definition: ERF_ABLMost.H:516
amrex::Real depth
Definition: ERF_ABLMost.H:519
amrex::Real surf_heating_rate
Definition: ERF_ABLMost.H:512
void compute_pblh(const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const PBLHeightEstimator &est, const int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp)
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_flux_lev
Definition: ERF_ABLMost.H:538
MOSTAverage m_ma
Definition: ERF_ABLMost.H:526
ThetaCalcType
Definition: ERF_ABLMost.H:481
@ SURFACE_TEMPERATURE
Surface temperature specified.
@ HEAT_FLUX
Heat-flux specified.
const amrex::iMultiFab * get_lmask(const int &lev)
Definition: ERF_ABLMost.H:452
int lmask_min_reduce(amrex::iMultiFab &lmask, const int &nghost)
Definition: ERF_ABLMost.H:455
bool have_variable_sea_roughness()
Definition: ERF_ABLMost.H:449
bool cnk_visc
Definition: ERF_ABLMost.H:518
amrex::Real m_start_bdy_time
Definition: ERF_ABLMost.H:520
const amrex::MultiFab * get_q_star(const int &lev)
Definition: ERF_ABLMost.H:429
amrex::Vector< std::unique_ptr< amrex::MultiFab > > pblh
Definition: ERF_ABLMost.H:532
void update_mac_ptrs(const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qr_prim)
Definition: ERF_ABLMost.H:412
const amrex::FArrayBox * get_z0(const int &lev)
Definition: ERF_ABLMost.H:447
amrex::Vector< amrex::MultiFab * > m_Hwave_lev
Definition: ERF_ABLMost.H:539
void read_custom_roughness(const int &lev, const std::string &fname)
Definition: ERF_ABLMost.cpp:592
amrex::Real surf_temp
Definition: ERF_ABLMost.H:511
amrex::Real cnk_a
Definition: ERF_ABLMost.H:517
ThetaCalcType theta_type
Definition: ERF_ABLMost.H:500
void compute_fluxes(const int &lev, const int &max_iters, const FluxIter &most_flux, bool is_land)
Definition: ERF_ABLMost.cpp:145
PBLHeightCalcType pblh_type
Definition: ERF_ABLMost.H:503
FluxCalcType
Definition: ERF_ABLMost.H:474
@ MOENG
Moeng functional form.
@ CUSTOM
Custom constant flux functional form.
@ ROTATE
Terrain rotation flux functional form.
@ DONELAN
Donelan functional form.
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_data_lev
Definition: ERF_ABLMost.H:537
amrex::Real surf_temp_flux
Definition: ERF_ABLMost.H:513
void update_pblh(const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp)
Definition: ERF_ABLMost.cpp:561
const amrex::MultiFab * get_w_star(const int &lev)
Definition: ERF_ABLMost.H:423
void impose_most_bcs(const int &lev, const amrex::Vector< amrex::MultiFab * > &mfs, amrex::MultiFab *xxmom_flux, amrex::MultiFab *yymom_flux, amrex::MultiFab *zzmom_flux, amrex::MultiFab *xymom_flux, amrex::MultiFab *yxmom_flux, amrex::MultiFab *xzmom_flux, amrex::MultiFab *zxmom_flux, amrex::MultiFab *yzmom_flux, amrex::MultiFab *zymom_flux, amrex::MultiFab *xheat_flux, amrex::MultiFab *yheat_flux, amrex::MultiFab *zheat_flux, amrex::MultiFab *xqv_flux, amrex::MultiFab *yqv_flux, amrex::MultiFab *zqv_flux, amrex::MultiFab *z_phys)
Definition: ERF_ABLMost.cpp:211
amrex::Vector< amrex::Vector< amrex::iMultiFab * > > m_lmask_lev
Definition: ERF_ABLMost.H:536
const amrex::MultiFab * get_pblh(const int &lev)
Definition: ERF_ABLMost.H:435
void time_interp_sst(const int &lev, const amrex::Real &time)
Definition: ERF_ABLMost.cpp:493
RoughCalcType rough_type_land
Definition: ERF_ABLMost.H:501
bool m_rotate
Definition: ERF_ABLMost.H:508
amrex::Vector< amrex::MultiFab * > m_Lwave_lev
Definition: ERF_ABLMost.H:540
bool m_var_z0
Definition: ERF_ABLMost.H:524
RoughCalcType rough_type_sea
Definition: ERF_ABLMost.H:502
bool use_moisture
Definition: ERF_ABLMost.H:506
void compute_most_bcs(const int &lev, const amrex::Vector< amrex::MultiFab * > &mfs, amrex::MultiFab *xxmom_flux, amrex::MultiFab *yymom_flux, amrex::MultiFab *zzmom_flux, amrex::MultiFab *xymom_flux, amrex::MultiFab *yxmom_flux, amrex::MultiFab *xzmom_flux, amrex::MultiFab *zxmom_flux, amrex::MultiFab *yzmom_flux, amrex::MultiFab *zymom_flux, amrex::MultiFab *xheat_flux, amrex::MultiFab *yheat_flux, amrex::MultiFab *zheat_flux, amrex::MultiFab *xqv_flux, amrex::MultiFab *yqv_flux, amrex::MultiFab *zqv_flux, amrex::MultiFab *z_phys, const FluxCalc &flux_comp)
bool m_include_wstar
Definition: ERF_ABLMost.H:509
bool m_exp_most
Definition: ERF_ABLMost.H:507
ABLMost(const amrex::Vector< amrex::Geometry > &geom, bool &use_exp_most, bool &use_rot_most, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qr_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &z_phys_nd, amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab >>> &sst_lev, amrex::Vector< amrex::Vector< std::unique_ptr< amrex::iMultiFab >>> &lmask_lev, amrex::Vector< amrex::Vector< amrex::MultiFab * >> lsm_data, amrex::Vector< amrex::Vector< amrex::MultiFab * >> lsm_flux, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Hwave, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Lwave, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &eddyDiffs, amrex::Real start_bdy_time=0.0, amrex::Real bdy_time_interval=0.0)
Definition: ERF_ABLMost.H:35
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_star
Definition: ERF_ABLMost.H:529
const amrex::MultiFab * get_t_surf(const int &lev)
Definition: ERF_ABLMost.H:441
amrex::Vector< std::unique_ptr< amrex::MultiFab > > u_star
Definition: ERF_ABLMost.H:527
PBLHeightCalcType
Definition: ERF_ABLMost.H:495
amrex::Real z0_const
Definition: ERF_ABLMost.H:510
FluxCalcType flux_type
Definition: ERF_ABLMost.H:499
void update_surf_temp(const amrex::Real &time)
Definition: ERF_ABLMost.H:396
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_sst_lev
Definition: ERF_ABLMost.H:535
void get_lsm_tsurf(const int &lev)
Definition: ERF_ABLMost.cpp:528
amrex::Real get_zref()
Definition: ERF_ABLMost.H:444
RoughCalcType
Definition: ERF_ABLMost.H:487
@ CONSTANT
Constant z0.
const amrex::MultiFab * get_olen(const int &lev)
Definition: ERF_ABLMost.H:432
const amrex::MultiFab * get_mac_avg(const int &lev, int comp)
Definition: ERF_ABLMost.H:438
void update_fluxes(const int &lev, const amrex::Real &time, int max_iters=25)
Definition: ERF_ABLMost.cpp:12
amrex::Real custom_tstar
Definition: ERF_ABLMost.H:515
Definition: ERF_MOSTAverage.H:12
void update_field_ptrs(int lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qr_prim)
Definition: ERF_MOSTAverage.cpp:228
amrex::Real get_zref() const
Definition: ERF_MOSTAverage.H:95
const amrex::MultiFab * get_average(int lev, int comp) const
Definition: ERF_MOSTAverage.H:92
@ cons
Definition: ERF_IndexDefines.H:129