ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_MOSTStress.H
Go to the documentation of this file.
1 #ifndef ERF_MOSTStress_H
2 #define ERF_MOSTStress_H
3 
4 #include <ERF_Constants.H>
5 #include <ERF_IndexDefines.H>
6 #include <ERF_MOSTRoughness.H>
7 #include <ERF_Wstar.H>
8 
9 /**
10  * Structure of plain old data relevant to MOST BCs
11  */
12 struct most_data
13 {
14 public:
15  amrex::Real z0_const{amrex::Real(0.1)}; ///< Roughness height -- default constant value(m)
16  amrex::Real kappa{KAPPA}; ///< von Karman constant
17  amrex::Real gravity{CONST_GRAV}; ///< Acceleration due to gravity (m/s^2)
18  amrex::Real surf_temp_flux{zero}; ///< Heat flux TODO: decide whether this is <θ'w'> or <θv'w'> under moist conditions
19  amrex::Real surf_moist_flux{zero}; ///< Moisture flux
20 
21  amrex::Real Cnk_a{amrex::Real(0.0185)}; ///< Standard Charnock constant https://doi.org/amrex::Real(10.1175)/JAMC-D-17-amrex::Real(0137.1)
22  amrex::Real Cnk_b1{one/amrex::Real(30.0)}; ///< Modified Charnock Eq (4) https://doi.org/amrex::Real(10.1175)/JAMC-D-17-amrex::Real(0137.1)
23  amrex::Real Cnk_b2{amrex::Real(1260.0)}; ///< Modified Charnock Eq (4) https://doi.org/amrex::Real(10.1175)/JAMC-D-17-amrex::Real(0137.1)
24  amrex::Real Cnk_d{amrex::Real(30.0)}; ///< Modified Charnock Eq (4) https://doi.org/amrex::Real(10.1175)/JAMC-D-17-amrex::Real(0137.1)
26  bool visc{false}; ///< Use viscous Charnock formulation
27 
31 
32  const amrex::Real Bjr_beta = amrex::Real(1.2); // Empirical parameter from Beljaars 1995 QJRMS
33 };
34 
35 
36 /**
37  * Structure of similarity functions for Moeng formulation
38  */
40 {
41  //
42  // Similarity functions from Jimenez (2012)
43  //
44  AMREX_GPU_HOST_DEVICE
45  AMREX_FORCE_INLINE
47  calc_psi_m2 (amrex::Real zeta) const
48  {
49  if (zeta > 0) {
50  amrex::Real x = std::pow(one + std::pow(zeta, amrex::Real(2.5)), one/amrex::Real(2.5));
51  return ( -amrex::Real(6.1)*std::log(zeta + x) );
52  } else {
53  amrex::Real x = std::pow(one - amrex::Real(16.0)*zeta, fourth);
54  amrex::Real psi_k_m = two * std::log(myhalf * (one + x)) + std::log(myhalf * (one + x * x)) -
55  two * std::atan(x) + PIoTwo;
56  amrex::Real y = std::pow(one - amrex::Real(10.0)*zeta, one/three);
57  amrex::Real psi_c_m = (three/two)*std::log((y*y + y + one)/three)
58  - std::sqrt(three)*std::atan((two*y + one)/std::sqrt(three))
59  + PI/std::sqrt(three);
60  return ( (psi_k_m + zeta*zeta*psi_c_m) / (one + zeta*zeta) );
61  }
62  }
63 
64  AMREX_GPU_HOST_DEVICE
65  AMREX_FORCE_INLINE
67  calc_psi_h2 (amrex::Real zeta) const
68  {
69  if (zeta > 0) {
70  amrex::Real x = std::pow(one + std::pow(zeta, amrex::Real(1.1)), one/amrex::Real(1.1));
71  return ( -amrex::Real(6.1)*std::log(zeta + x) );
72  } else {
73  amrex::Real x = std::sqrt(one - amrex::Real(16.0)*zeta);
74  amrex::Real psi_k_h = two * std::log(myhalf * (one + x));
75  amrex::Real y = std::pow(one - amrex::Real(34.0)*zeta, one/three);
76  amrex::Real psi_c_h = (three/two)*std::log((y*y + y + one)/three)
77  - std::sqrt(three)*std::atan((two*y + one)/std::sqrt(three))
78  + PI/std::sqrt(three);
79  return ( (psi_k_h + zeta*zeta*psi_c_h) / (one + zeta*zeta) );
80  }
81  }
82 
83 
84  //
85  // Similarity functions from Businger & Dyer (1966)
86  //
87  AMREX_GPU_HOST_DEVICE
88  AMREX_FORCE_INLINE
90  calc_psi_m (amrex::Real zeta) const
91  {
92  if (zeta > 0) {
93  return -beta_m * zeta;
94  } else {
95  amrex::Real x = std::sqrt(std::sqrt(one - gamma_m * zeta));
96  return two * std::log(myhalf * (one + x)) + std::log(myhalf * (one + x * x)) -
97  two * std::atan(x) + PIoTwo;
98  }
99  }
100 
101  AMREX_GPU_HOST_DEVICE
102  AMREX_FORCE_INLINE
104  calc_psi_h (amrex::Real zeta) const
105  {
106  if (zeta > 0) {
107  return -beta_h * zeta;
108  } else {
109  amrex::Real x = std::sqrt(one - gamma_h * zeta);
110  return two * std::log(myhalf * (one + x));
111  }
112  }
113 
114 private:
115  amrex::Real beta_m{amrex::Real(5.0)}; ///< Constants from Dyer, BLM, 1974
116  amrex::Real beta_h{amrex::Real(5.0)}; ///< https://doi.org/amrex::Real(10.1007)/BF00240838
119 };
120 
121 
122 /**
123  * Empirical kinematic viscosity [m2/s] formula from Andreas (1989) CRREL Rep.
124  * 89-11, valid between -173 and 277 deg C.
125  */
126 AMREX_GPU_DEVICE
127 AMREX_FORCE_INLINE
130 {
131  amrex::Real TC = T_degK - amrex::Real(273.15);
132  return amrex::Real(1.326e-5)*(one + amrex::Real(6.542e-3)*TC + amrex::Real(8.301e-6)*TC*TC - amrex::Real(4.84e-9)*TC*TC*TC);
133 }
134 
135 
136 /**
137  * Adiabatic with constant roughness
138  */
139 struct adiabatic
140 {
142  amrex::Real Qvflux)
143  {
144  mdata.surf_temp_flux = Tflux;
145  mdata.surf_moist_flux = Qvflux;
146  }
147 
148  AMREX_GPU_DEVICE
149  AMREX_FORCE_INLINE
150  void
151  iterate_flux (const int& i,
152  const int& j,
153  const int& k,
154  const int& /*max_iters*/,
155  const amrex::Array4<const amrex::Real>& zref_arr,
156  const amrex::Array4<const amrex::Real>& z0_arr,
157  const amrex::Array4<const amrex::Real>& umm_arr,
158  const amrex::Array4<const amrex::Real>& /*tm_arr*/,
159  const amrex::Array4<const amrex::Real>& /*tvm_arr*/,
160  const amrex::Array4<const amrex::Real>& /*qvm_arr*/,
161  const amrex::Array4<amrex::Real>& u_star_arr,
162  const amrex::Array4<amrex::Real>& /*w_star_arr*/,
163  const amrex::Array4<amrex::Real>& t_star_arr,
164  const amrex::Array4<amrex::Real>& q_star_arr,
165  const amrex::Array4<amrex::Real>& /*t_surf_arr*/,
166  const amrex::Array4<amrex::Real>& /*q_surf_arr*/,
167  const amrex::Array4<amrex::Real>& olen_arr,
168  const amrex::Array4<amrex::Real>& /*pblh_arr*/,
169  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
170  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
171  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
172  {
173  olen_arr(i,j,k) = amrex::Real(1.0e16);
174  u_star_arr(i,j,k) = mdata.kappa * umm_arr(i,j,k) / std::log(zref_arr(i,j,k) / z0_arr(i,j,k));
175  t_star_arr(i,j,k) = zero;
176  q_star_arr(i,j,k) = zero;
177  }
178 
179 private:
182 };
183 
184 
185 /**
186  * Adiabatic with charnock roughness
187  */
189 {
191  amrex::Real Qvflux,
192  amrex::Real cnk_a,
193  bool cnk_visc)
194  {
195  mdata.surf_temp_flux = Tflux;
196  mdata.surf_moist_flux = Qvflux;
197  mdata.Cnk_a = cnk_a;
198  mdata.visc = cnk_visc;
199  }
200 
201  AMREX_GPU_DEVICE
202  AMREX_FORCE_INLINE
203  void
204  iterate_flux (const int& i,
205  const int& j,
206  const int& k,
207  const int& max_iters,
208  const amrex::Array4<const amrex::Real>& zref_arr,
209  const amrex::Array4<amrex::Real>& z0_arr,
210  const amrex::Array4<const amrex::Real>& umm_arr,
211  const amrex::Array4<const amrex::Real>& tm_arr,
212  const amrex::Array4<const amrex::Real>& /*tvm_arr*/,
213  const amrex::Array4<const amrex::Real>& /*qvm_arr*/,
214  const amrex::Array4<amrex::Real>& u_star_arr,
215  const amrex::Array4<amrex::Real>& /*w_star_arr*/,
216  const amrex::Array4<amrex::Real>& t_star_arr,
217  const amrex::Array4<amrex::Real>& q_star_arr,
218  const amrex::Array4<amrex::Real>& /*t_surf_arr*/,
219  const amrex::Array4<amrex::Real>& /*q_surf_arr*/,
220  const amrex::Array4<amrex::Real>& olen_arr,
221  const amrex::Array4<amrex::Real>& /*pblh_arr*/,
222  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
223  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
224  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
225  {
226  amrex::Real ustar = zero;
227  amrex::Real zref = zref_arr(i,j,k);
228  amrex::Real z0 = z0_arr(i,j,k);
229  amrex::Real z0_old = z0;
230  amrex::Real psi_m = zero;
231  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
232  amrex::Real C = std::log(zref / z0);
233 
234  // Fixed point iteration on roughness
235  int iter_z = 0;
236  do {
237  // Transfer curr to old
238  z0_old = z0;
239 
240  // Update
241  C = std::log(zref / z0_old);
242  ustar = mdata.kappa * umm / (C - psi_m);
243  if (mdata.Cnk_a > 0) {
244  z0 = (mdata.Cnk_a / mdata.gravity) * ustar * ustar;
245  if (mdata.visc) {
246  z0 += air_viscosity(tm_arr(i,j,k)) / std::max(ustar, amrex::Real(0.05));
247  }
248  } else {
249  z0 = COARE3_roughness(zref, umm, ustar);
250  }
251 
252  ++iter_z;
253  } while ( (std::abs(z0 - z0_old) > tol_z) && (iter_z <= max_iters) );
254  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter_z < max_iters,
255  "Maximum number of MOST roughness iterations reached.");
256  C = std::log(zref / z0);
257 
258  // Populate the stored MOST arrays
259  z0_arr(i,j,k) = z0;
260  olen_arr(i,j,k) = amrex::Real(1.0e16);
261  u_star_arr(i,j,k) = mdata.kappa * umm / (C - psi_m);
262  t_star_arr(i,j,k) = zero;
263  q_star_arr(i,j,k) = zero;
264  }
265 
266 private:
269  const amrex::Real tol_z = amrex::Real(1.0e-10);
270  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
271 };
272 
273 
274 /**
275  * Adiabatic with modified charnock roughness
276  */
278 {
280  amrex::Real Qvflux,
281  amrex::Real depth)
282  {
283  mdata.surf_temp_flux = Tflux;
284  mdata.surf_moist_flux = Qvflux;
285  mdata.Cnk_d = depth;
286  mdata.Cnk_b = mdata.Cnk_b1 * std::log(mdata.Cnk_b2 / mdata.Cnk_d);
287  }
288 
289  AMREX_GPU_DEVICE
290  AMREX_FORCE_INLINE
291  void
292  iterate_flux (const int& i,
293  const int& j,
294  const int& k,
295  const int& max_iters,
296  const amrex::Array4<const amrex::Real>& zref_arr,
297  const amrex::Array4<amrex::Real>& z0_arr,
298  const amrex::Array4<const amrex::Real>& umm_arr,
299  const amrex::Array4<const amrex::Real>& /*tm_arr*/,
300  const amrex::Array4<const amrex::Real>& /*tvm_arr*/,
301  const amrex::Array4<const amrex::Real>& /*qvm_arr*/,
302  const amrex::Array4<amrex::Real>& u_star_arr,
303  const amrex::Array4<amrex::Real>& /*w_star_arr*/,
304  const amrex::Array4<amrex::Real>& t_star_arr,
305  const amrex::Array4<amrex::Real>& q_star_arr,
306  const amrex::Array4<amrex::Real>& /*t_surf_arr*/,
307  const amrex::Array4<amrex::Real>& /*q_surf_arr*/,
308  const amrex::Array4<amrex::Real>& olen_arr,
309  const amrex::Array4<amrex::Real>& /*pblh_arr*/,
310  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
311  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
312  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
313  {
314  amrex::Real ustar = zero;
315  amrex::Real zref = zref_arr(i,j,k);
316  amrex::Real z0 = z0_arr(i,j,k);
317  amrex::Real z0_old = z0;
318  amrex::Real psi_m = zero;
319  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
320  amrex::Real C = std::log(zref / z0);
321 
322  // Fixed point iteration on roughness
323  int iter_z = 0;
324  do {
325  // Transfer curr to old
326  z0_old = z0;
327 
328  // Update
329  C = std::log(zref / z0_old);
330  ustar = mdata.kappa * umm / (C - psi_m);
331  z0 = std::exp( (amrex::Real(2.7)*ustar - amrex::Real(1.8)/mdata.Cnk_b) / (ustar + amrex::Real(0.17)/mdata.Cnk_b) );
332 
333  ++iter_z;
334  } while ( (std::abs(z0 - z0_old) > tol_z) && (iter_z <= max_iters) );
335  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter_z < max_iters,
336  "Maximum number of MOST roughness iterations reached.");
337  C = std::log(zref / z0);
338 
339  // Populate the stored MOST arrays
340  z0_arr(i,j,k) = z0;
341  olen_arr(i,j,k) = amrex::Real(1.0e16);
342  u_star_arr(i,j,k) = mdata.kappa * umm / (C - psi_m);
343  t_star_arr(i,j,k) = zero;
344  q_star_arr(i,j,k) = zero;
345  }
346 
347 private:
350  const amrex::Real tol_z = amrex::Real(1.0e-10);
351  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
352 };
353 
354 
355 /**
356  * Adiabatic with Donelan roughness
357  */
359 {
361  amrex::Real Qvflux)
362  {
363  mdata.surf_temp_flux = Tflux;
364  mdata.surf_moist_flux = Qvflux;
365  }
366 
367  AMREX_GPU_DEVICE
368  AMREX_FORCE_INLINE
369  void
370  iterate_flux (const int& i,
371  const int& j,
372  const int& k,
373  const int& max_iters,
374  const amrex::Array4<const amrex::Real>& zref_arr,
375  const amrex::Array4<amrex::Real>& z0_arr,
376  const amrex::Array4<const amrex::Real>& umm_arr,
377  const amrex::Array4<const amrex::Real>& /*tm_arr*/,
378  const amrex::Array4<const amrex::Real>& /*tvm_arr*/,
379  const amrex::Array4<const amrex::Real>& /*qvm_arr*/,
380  const amrex::Array4<amrex::Real>& u_star_arr,
381  const amrex::Array4<amrex::Real>& /*w_star_arr*/,
382  const amrex::Array4<amrex::Real>& t_star_arr,
383  const amrex::Array4<amrex::Real>& q_star_arr,
384  const amrex::Array4<amrex::Real>& /*t_surf_arr*/,
385  const amrex::Array4<amrex::Real>& /*q_surf_arr*/,
386  const amrex::Array4<amrex::Real>& olen_arr,
387  const amrex::Array4<amrex::Real>& /*pblh_arr*/,
388  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
389  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
390  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
391  {
392  int iter = 0;
393  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
394  amrex::Real ustar = zero;
395  amrex::Real zref = zref_arr(i,j,k);
396  amrex::Real z0 = zero;
397  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
398  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) {
399  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
400  }
401  do {
402  ustar = u_star_arr(i,j,k);
403  z0 = Donelan_roughness(ustar);
404  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0);
405  ++iter;
406  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
407  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
408  "Maximum number of MOST iterations reached.");
409 
410  // Populate the stored MOST arrays
411  z0_arr(i,j,k) = z0;
412  olen_arr(i,j,k) = amrex::Real(1.0e16);
413  t_star_arr(i,j,k) = zero;
414  q_star_arr(i,j,k) = zero;
415  }
416 
417 private:
420  const amrex::Real tol = amrex::Real(1.0e-5);
421  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
422 };
423 
424 
425 /**
426  * Adiabatic with wave-coupled roughness
427  */
429 {
431  amrex::Real Qvflux)
432  {
433  mdata.surf_temp_flux = Tflux;
434  mdata.surf_moist_flux = Qvflux;
435  }
436 
437  AMREX_GPU_DEVICE
438  AMREX_FORCE_INLINE
439  void
440  iterate_flux (const int& i,
441  const int& j,
442  const int& k,
443  const int& max_iters,
444  const amrex::Array4<const amrex::Real>& zref_arr,
445  const amrex::Array4<amrex::Real>& z0_arr,
446  const amrex::Array4<const amrex::Real>& umm_arr,
447  const amrex::Array4<const amrex::Real>& /*tm_arr*/,
448  const amrex::Array4<const amrex::Real>& /*tvm_arr*/,
449  const amrex::Array4<const amrex::Real>& /*qvm_arr*/,
450  const amrex::Array4<amrex::Real>& u_star_arr,
451  const amrex::Array4<amrex::Real>& /*w_star_arr*/,
452  const amrex::Array4<amrex::Real>& t_star_arr,
453  const amrex::Array4<amrex::Real>& q_star_arr,
454  const amrex::Array4<amrex::Real>& /*t_surf_arr*/,
455  const amrex::Array4<amrex::Real>& /*q_surf_arr*/,
456  const amrex::Array4<amrex::Real>& olen_arr,
457  const amrex::Array4<amrex::Real>& /*pblh_arr*/,
458  const amrex::Array4<amrex::Real>& Hwave_arr,
459  const amrex::Array4<amrex::Real>& Lwave_arr,
460  const amrex::Array4<amrex::Real>& eta_arr) const
461  {
462  int iter = 0;
463  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
464  amrex::Real ustar = zero;
465  amrex::Real zref = zref_arr(i,j,k);
466  amrex::Real z0 = zero;
467  int ie, je;
468  ie = i < lbound(eta_arr).x ? lbound(eta_arr).x : i;
469  je = j < lbound(eta_arr).y ? lbound(eta_arr).y : j;
470  ie = ie > ubound(eta_arr).x ? ubound(eta_arr).x : ie;
471  je = je > ubound(eta_arr).y ? ubound(eta_arr).y : je;
472  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
473  do {
474  ustar = u_star_arr(i,j,k);
475  z0 = std::min( std::max(amrex::Real(1200.0) * Hwave_arr(i,j,k) * std::pow( Hwave_arr(i,j,k)/(Lwave_arr(i,j,k)+eps), amrex::Real(4.5) )
476  + amrex::Real(0.11) * eta_arr(ie,je,k,EddyDiff::Mom_v) / ustar, z0_eps), z0_max );
477  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0);
478  ++iter;
479  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
480  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
481  "Maximum number of MOST iterations reached.");
482 
483  // Populate the stored MOST arrays
484  z0_arr(i,j,k) = z0;
485  olen_arr(i,j,k) = amrex::Real(1.0e16);
486  t_star_arr(i,j,k) = zero;
487  q_star_arr(i,j,k) = zero;
488  }
489 
490 private:
493  const amrex::Real tol = amrex::Real(1.0e-5);
494 #ifdef AMREX_USE_FLOAT
495  const amrex::Real eps = amrex::Real(1e-8);
496 #else
497  const amrex::Real eps = amrex::Real(1e-15);
498 #endif
499  const amrex::Real z0_eps = amrex::Real(1.0e-6);
501  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
502 };
503 
504 
505 /**
506  * Surface flux with constant roughness
507  */
509 {
511  amrex::Real Qvflux,
512  bool cons_qflux)
513  {
514  mdata.surf_temp_flux = Tflux;
515  mdata.surf_moist_flux = Qvflux;
516  spec_qflux = cons_qflux;
517  }
518 
519  AMREX_GPU_DEVICE
520  AMREX_FORCE_INLINE
521  void
522  iterate_flux (const int& i,
523  const int& j,
524  const int& k,
525  const int& max_iters,
526  const amrex::Array4<const amrex::Real>& zref_arr,
527  const amrex::Array4<const amrex::Real>& z0_arr,
528  const amrex::Array4<const amrex::Real>& umm_arr,
529  const amrex::Array4<const amrex::Real>& tm_arr,
530  const amrex::Array4<const amrex::Real>& tvm_arr,
531  const amrex::Array4<const amrex::Real>& qvm_arr,
532  const amrex::Array4<amrex::Real>& u_star_arr,
533  const amrex::Array4<amrex::Real>& w_star_arr,
534  const amrex::Array4<amrex::Real>& t_star_arr,
535  const amrex::Array4<amrex::Real>& q_star_arr,
536  const amrex::Array4<amrex::Real>& t_surf_arr,
537  const amrex::Array4<amrex::Real>& q_surf_arr,
538  const amrex::Array4<amrex::Real>& olen_arr,
539  const amrex::Array4<amrex::Real>& pblh_arr,
540  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
541  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
542  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
543  {
544  int iter = 0;
545  amrex::Real ustar = zero;
546  amrex::Real wstar = zero;
547  amrex::Real tflux = zero;
548  amrex::Real qflux = zero;
549  amrex::Real zeta = zero;
550  amrex::Real psi_m = zero;
551  amrex::Real psi_h = zero;
552  amrex::Real Olen = zero;
553  amrex::Real zref = zref_arr(i,j,k);
554  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
555  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) {
556  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
557  } else {
558  Olen = olen_arr(i,j,k);
559  zeta = zref / Olen;
560  psi_m = sfuns.calc_psi_m(zeta);
561  psi_h = sfuns.calc_psi_h(zeta);
562  }
563  do {
564  ustar = u_star_arr(i,j,k);
565  qflux = (spec_qflux) ? mdata.surf_moist_flux :
566  -(qvm_arr(i,j,k) - q_surf_arr(i,j,k)) * ustar * mdata.kappa /
567  (std::log(zref / z0_arr(i,j,k)) - psi_h); // <w'Qv'>
568  tflux = mdata.surf_temp_flux*(1 + amrex::Real(0.61)*qvm_arr(i,j,k)) + qflux*amrex::Real(0.61)*tm_arr(i,j,k);
569  if (w_star_arr) {
570  // update w* and Umagmean
571  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
572  wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
573  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
574  umm = std::max(umm, WSMIN);
575  }
576  Olen = -ustar * ustar * ustar * tvm_arr(i,j,k) / (mdata.kappa * mdata.gravity * tflux);
577  zeta = zref / Olen;
578  psi_m = sfuns.calc_psi_m(zeta);
579  psi_h = sfuns.calc_psi_h(zeta);
580  u_star_arr(i,j,k) = mdata.kappa * umm / (std::log(zref / z0_arr(i,j,k)) - psi_m);
581  ++iter;
582  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
583  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
584  "Maximum number of MOST iterations reached.");
585 
586  // Populate the stored MOST arrays
587  olen_arr(i,j,k) = Olen;
588  t_surf_arr(i,j,k) = mdata.surf_temp_flux * (std::log(zref / z0_arr(i,j,k)) - psi_h) /
589  (u_star_arr(i,j,k) * mdata.kappa) + tm_arr(i,j,k);
590  t_star_arr(i,j,k) = -mdata.surf_temp_flux / u_star_arr(i,j,k);
591  if (spec_qflux) {
592  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (std::log(zref / z0_arr(i,j,k)) - psi_h) /
593  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
594  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
595  } else {
596  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) /
597  (std::log(zref / z0_arr(i,j,k)) - psi_h);
598  }
599  }
600 
601 private:
605  const amrex::Real tol = amrex::Real(1.0e-5);
606  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
607 };
608 
609 
610 /**
611  * Surface flux with charnock roughness
612  */
614 {
616  amrex::Real Qvflux,
617  amrex::Real cnk_a,
618  bool cnk_visc,
619  bool cons_qflux)
620  {
621  mdata.surf_temp_flux = Tflux;
622  mdata.surf_moist_flux = Qvflux;
623  mdata.Cnk_a = cnk_a;
624  mdata.visc = cnk_visc;
625  spec_qflux = cons_qflux;
626  }
627 
628  AMREX_GPU_DEVICE
629  AMREX_FORCE_INLINE
630  void
631  iterate_flux (const int& i,
632  const int& j,
633  const int& k,
634  const int& max_iters,
635  const amrex::Array4<const amrex::Real>& zref_arr,
636  const amrex::Array4<amrex::Real>& z0_arr,
637  const amrex::Array4<const amrex::Real>& umm_arr,
638  const amrex::Array4<const amrex::Real>& tm_arr,
639  const amrex::Array4<const amrex::Real>& tvm_arr,
640  const amrex::Array4<const amrex::Real>& qvm_arr,
641  const amrex::Array4<amrex::Real>& u_star_arr,
642  const amrex::Array4<amrex::Real>& w_star_arr,
643  const amrex::Array4<amrex::Real>& t_star_arr,
644  const amrex::Array4<amrex::Real>& q_star_arr,
645  const amrex::Array4<amrex::Real>& t_surf_arr,
646  const amrex::Array4<amrex::Real>& q_surf_arr,
647  const amrex::Array4<amrex::Real>& olen_arr,
648  const amrex::Array4<amrex::Real>& pblh_arr,
649  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
650  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
651  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
652  {
653  int iter = 0;
654  amrex::Real ustar = zero;
655  amrex::Real wstar = zero;
656  amrex::Real tflux = zero;
657  amrex::Real qflux = zero;
658  amrex::Real z0 = zero;
659  amrex::Real zeta = zero;
660  amrex::Real psi_m = zero;
661  amrex::Real psi_h = zero;
662  amrex::Real Olen = zero;
663  amrex::Real zref = zref_arr(i,j,k);
664  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
665  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) {
666  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
667  } else {
668  Olen = olen_arr(i,j,k);
669  zeta = zref / Olen;
670  psi_m = sfuns.calc_psi_m(zeta);
671  psi_h = sfuns.calc_psi_h(zeta);
672  }
673  do {
674  ustar = u_star_arr(i,j,k);
675  if (mdata.Cnk_a > 0) {
676  z0 = (mdata.Cnk_a / mdata.gravity) * ustar * ustar;
677  if (mdata.visc) {
678  z0 += air_viscosity(tm_arr(i,j,k)) / std::max(ustar, amrex::Real(0.05));
679  }
680  } else {
681  z0 = COARE3_roughness(zref, umm, ustar);
682  }
683  qflux = (spec_qflux) ? mdata.surf_moist_flux :
684  -(qvm_arr(i,j,k) - q_surf_arr(i,j,k)) * ustar * mdata.kappa /
685  (std::log(zref / z0_arr(i,j,k)) - psi_h); // <w'Qv'>
686  tflux = mdata.surf_temp_flux*(1 + amrex::Real(0.61)*qvm_arr(i,j,k)) + qflux*amrex::Real(0.61)*tm_arr(i,j,k);
687  if (w_star_arr) {
688  // update w* and Umagmean
689  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
690  wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
691  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
692  umm = std::max(umm, WSMIN);
693  }
694  Olen = -ustar * ustar * ustar * tvm_arr(i,j,k) / (mdata.kappa * mdata.gravity * tflux);
695  zeta = zref / Olen;
696  psi_m = sfuns.calc_psi_m(zeta);
697  psi_h = sfuns.calc_psi_h(zeta);
698  u_star_arr(i,j,k) = mdata.kappa * umm / (std::log(zref / z0) - psi_m);
699  ++iter;
700  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
701  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
702  "Maximum number of MOST iterations reached.");
703 
704  // Populate the stored MOST arrays
705  z0_arr(i,j,k) = z0;
706  olen_arr(i,j,k) = Olen;
707  t_surf_arr(i,j,k) = mdata.surf_temp_flux * (std::log(zref / z0) - psi_h) /
708  (u_star_arr(i,j,k) * mdata.kappa) + tm_arr(i,j,k);
709  t_star_arr(i,j,k) = -mdata.surf_temp_flux / u_star_arr(i,j,k);
710  if (spec_qflux) {
711  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (std::log(zref / z0_arr(i,j,k)) - psi_h) /
712  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
713  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
714  } else {
715  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) /
716  (std::log(zref / z0_arr(i,j,k)) - psi_h);
717  }
718  }
719 
720 private:
724  const amrex::Real tol = amrex::Real(1.0e-5);
725  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
726 };
727 
728 
729 /**
730  * Surface flux with modified charnock roughness
731  */
733 {
735  amrex::Real Qvflux,
736  amrex::Real depth,
737  bool cons_qflux)
738  {
739  mdata.surf_temp_flux = Tflux;
740  mdata.surf_moist_flux = Qvflux;
741  mdata.Cnk_d = depth;
742  mdata.Cnk_b = mdata.Cnk_b1 * std::log(mdata.Cnk_b2 / mdata.Cnk_d);
743  spec_qflux = cons_qflux;
744  }
745 
746  AMREX_GPU_DEVICE
747  AMREX_FORCE_INLINE
748  void
749  iterate_flux (const int& i,
750  const int& j,
751  const int& k,
752  const int& max_iters,
753  const amrex::Array4<const amrex::Real>& zref_arr,
754  const amrex::Array4<amrex::Real>& z0_arr,
755  const amrex::Array4<const amrex::Real>& umm_arr,
756  const amrex::Array4<const amrex::Real>& tm_arr,
757  const amrex::Array4<const amrex::Real>& tvm_arr,
758  const amrex::Array4<const amrex::Real>& qvm_arr,
759  const amrex::Array4<amrex::Real>& u_star_arr,
760  const amrex::Array4<amrex::Real>& w_star_arr,
761  const amrex::Array4<amrex::Real>& t_star_arr,
762  const amrex::Array4<amrex::Real>& q_star_arr,
763  const amrex::Array4<amrex::Real>& t_surf_arr,
764  const amrex::Array4<amrex::Real>& q_surf_arr,
765  const amrex::Array4<amrex::Real>& olen_arr,
766  const amrex::Array4<amrex::Real>& pblh_arr,
767  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
768  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
769  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
770  {
771  int iter = 0;
772  amrex::Real ustar = zero;
773  amrex::Real wstar = zero;
774  amrex::Real tflux = zero;
775  amrex::Real qflux = zero;
776  amrex::Real z0 = zero;
777  amrex::Real zeta = zero;
778  amrex::Real psi_m = zero;
779  amrex::Real psi_h = zero;
780  amrex::Real Olen = zero;
781  amrex::Real zref = zref_arr(i,j,k);
782  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
783  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) {
784  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
785  } else {
786  Olen = olen_arr(i,j,k);
787  zeta = zref / Olen;
788  psi_m = sfuns.calc_psi_m(zeta);
789  psi_h = sfuns.calc_psi_h(zeta);
790  }
791  do {
792  ustar = u_star_arr(i,j,k);
793  z0 = std::exp( (amrex::Real(2.7)*ustar - amrex::Real(1.8)/mdata.Cnk_b) / (ustar + amrex::Real(0.17)/mdata.Cnk_b) );
794  qflux = (spec_qflux) ? mdata.surf_moist_flux :
795  -(qvm_arr(i,j,k) - q_surf_arr(i,j,k)) * ustar * mdata.kappa /
796  (std::log(zref / z0_arr(i,j,k)) - psi_h); // <w'Qv'>
797  tflux = mdata.surf_temp_flux*(1 + amrex::Real(0.61)*qvm_arr(i,j,k)) + qflux*amrex::Real(0.61)*tm_arr(i,j,k);
798  if (w_star_arr) {
799  // update w* and Umagmean
800  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
801  wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
802  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
803  umm = std::max(umm, WSMIN);
804  }
805  Olen = -ustar * ustar * ustar * tvm_arr(i,j,k) / (mdata.kappa * mdata.gravity * tflux);
806  zeta = zref / Olen;
807  psi_m = sfuns.calc_psi_m(zeta);
808  psi_h = sfuns.calc_psi_h(zeta);
809  u_star_arr(i,j,k) = mdata.kappa * umm / (std::log(zref / z0) - psi_m);
810  ++iter;
811  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
812  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
813  "Maximum number of MOST iterations reached.");
814 
815  // Populate the stored MOST arrays
816  z0_arr(i,j,k) = z0;
817  olen_arr(i,j,k) = Olen;
818  t_surf_arr(i,j,k) = mdata.surf_temp_flux * (std::log(zref / z0) - psi_h) /
819  (u_star_arr(i,j,k) * mdata.kappa) + tm_arr(i,j,k);
820  t_star_arr(i,j,k) = -mdata.surf_temp_flux / u_star_arr(i,j,k);
821  if (spec_qflux) {
822  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (std::log(zref / z0_arr(i,j,k)) - psi_h) /
823  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
824  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
825  } else {
826  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) /
827  (std::log(zref / z0_arr(i,j,k)) - psi_h);
828  }
829  }
830 
831 private:
835  const amrex::Real tol = amrex::Real(1.0e-5);
836  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
837 };
838 
839 
840 /**
841  * Surface flux with donelan roughness
842  */
844 {
846  amrex::Real Qvflux,
847  bool cons_qflux)
848  {
849  mdata.surf_temp_flux = Tflux;
850  mdata.surf_moist_flux = Qvflux;
851  spec_qflux = cons_qflux;
852  }
853 
854  AMREX_GPU_DEVICE
855  AMREX_FORCE_INLINE
856  void
857  iterate_flux (const int& i,
858  const int& j,
859  const int& k,
860  const int& max_iters,
861  const amrex::Array4<const amrex::Real>& zref_arr,
862  const amrex::Array4<amrex::Real>& z0_arr,
863  const amrex::Array4<const amrex::Real>& umm_arr,
864  const amrex::Array4<const amrex::Real>& tm_arr,
865  const amrex::Array4<const amrex::Real>& tvm_arr,
866  const amrex::Array4<const amrex::Real>& qvm_arr,
867  const amrex::Array4<amrex::Real>& u_star_arr,
868  const amrex::Array4<amrex::Real>& w_star_arr,
869  const amrex::Array4<amrex::Real>& t_star_arr,
870  const amrex::Array4<amrex::Real>& q_star_arr,
871  const amrex::Array4<amrex::Real>& t_surf_arr,
872  const amrex::Array4<amrex::Real>& q_surf_arr,
873  const amrex::Array4<amrex::Real>& olen_arr,
874  const amrex::Array4<amrex::Real>& pblh_arr,
875  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
876  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
877  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
878  {
879  int iter = 0;
880  amrex::Real ustar = zero;
881  amrex::Real wstar = zero;
882  amrex::Real tflux = zero;
883  amrex::Real qflux = zero;
884  amrex::Real z0 = zero;
885  amrex::Real zeta = zero;
886  amrex::Real psi_m = zero;
887  amrex::Real psi_h = zero;
888  amrex::Real Olen = zero;
889  amrex::Real zref = zref_arr(i,j,k);
890  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
891  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) {
892  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
893  } else {
894  Olen = olen_arr(i,j,k);
895  zeta = zref / Olen;
896  psi_m = sfuns.calc_psi_m(zeta);
897  psi_h = sfuns.calc_psi_h(zeta);
898  }
899  do {
900  ustar = u_star_arr(i,j,k);
901  z0 = Donelan_roughness(ustar);
902  qflux = (spec_qflux) ? mdata.surf_moist_flux :
903  -(qvm_arr(i,j,k) - q_surf_arr(i,j,k)) * ustar * mdata.kappa /
904  (std::log(zref / z0_arr(i,j,k)) - psi_h); // <w'Qv'>
905  tflux = mdata.surf_temp_flux*(1 + amrex::Real(0.61)*qvm_arr(i,j,k)) + qflux*amrex::Real(0.61)*tm_arr(i,j,k);
906  if (w_star_arr) {
907  // update w* and Umagmean
908  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
909  wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
910  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
911  umm = std::max(umm, WSMIN);
912  }
913  Olen = -ustar * ustar * ustar * tvm_arr(i,j,k) / (mdata.kappa * mdata.gravity * tflux);
914  zeta = zref / Olen;
915  psi_m = sfuns.calc_psi_m(zeta);
916  psi_h = sfuns.calc_psi_h(zeta);
917  u_star_arr(i,j,k) = mdata.kappa * umm / (std::log(zref / z0) - psi_m);
918  ++iter;
919  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
920  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
921  "Maximum number of MOST iterations reached.");
922 
923  // Populate the stored MOST arrays
924  z0_arr(i,j,k) = z0;
925  olen_arr(i,j,k) = Olen;
926  t_surf_arr(i,j,k) = mdata.surf_temp_flux * (std::log(zref / z0) - psi_h) /
927  (u_star_arr(i,j,k) * mdata.kappa) + tm_arr(i,j,k);
928  t_star_arr(i,j,k) = -mdata.surf_temp_flux / u_star_arr(i,j,k);
929  if (spec_qflux) {
930  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (std::log(zref / z0_arr(i,j,k)) - psi_h) /
931  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
932  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
933  } else {
934  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) /
935  (std::log(zref / z0_arr(i,j,k)) - psi_h);
936  }
937  }
938 
939 private:
943  const amrex::Real tol = amrex::Real(1.0e-5);
944  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
945 };
946 
947 
948 /**
949  * Surface flux with wave-coupled roughness
950  */
952 {
954  amrex::Real Qvflux,
955  bool cons_qflux)
956  {
957  mdata.surf_temp_flux = Tflux;
958  mdata.surf_moist_flux = Qvflux;
959  spec_qflux = cons_qflux;
960  }
961 
962  AMREX_GPU_DEVICE
963  AMREX_FORCE_INLINE
964  void
965  iterate_flux (const int& i,
966  const int& j,
967  const int& k,
968  const int& max_iters,
969  const amrex::Array4<const amrex::Real>& zref_arr,
970  const amrex::Array4<amrex::Real>& z0_arr,
971  const amrex::Array4<const amrex::Real>& umm_arr,
972  const amrex::Array4<const amrex::Real>& tm_arr,
973  const amrex::Array4<const amrex::Real>& tvm_arr,
974  const amrex::Array4<const amrex::Real>& qvm_arr,
975  const amrex::Array4<amrex::Real>& u_star_arr,
976  const amrex::Array4<amrex::Real>& w_star_arr,
977  const amrex::Array4<amrex::Real>& t_star_arr,
978  const amrex::Array4<amrex::Real>& q_star_arr,
979  const amrex::Array4<amrex::Real>& t_surf_arr,
980  const amrex::Array4<amrex::Real>& q_surf_arr,
981  const amrex::Array4<amrex::Real>& olen_arr,
982  const amrex::Array4<amrex::Real>& pblh_arr,
983  const amrex::Array4<amrex::Real>& Hwave_arr,
984  const amrex::Array4<amrex::Real>& Lwave_arr,
985  const amrex::Array4<amrex::Real>& eta_arr) const
986  {
987  int iter = 0;
988  amrex::Real ustar = zero;
989  amrex::Real wstar = zero;
990  amrex::Real tflux = zero;
991  amrex::Real qflux = zero;
992  amrex::Real z0 = zero;
993  amrex::Real zeta = zero;
994  amrex::Real psi_m = zero;
995  amrex::Real psi_h = zero;
996  amrex::Real Olen = zero;
997  amrex::Real zref = zref_arr(i,j,k);
998  int ie, je;
999  ie = i < lbound(eta_arr).x ? lbound(eta_arr).x : i;
1000  je = j < lbound(eta_arr).y ? lbound(eta_arr).y : j;
1001  ie = ie > ubound(eta_arr).x ? ubound(eta_arr).x : ie;
1002  je = je > ubound(eta_arr).y ? ubound(eta_arr).y : je;
1003  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
1004  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) {
1005  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
1006  } else {
1007  Olen = olen_arr(i,j,k);
1008  zeta = zref / Olen;
1009  psi_m = sfuns.calc_psi_m(zeta);
1010  psi_h = sfuns.calc_psi_h(zeta);
1011  }
1012  do {
1013  ustar = u_star_arr(i,j,k);
1014  z0 = std::min( std::max(amrex::Real(1200.0) * Hwave_arr(i,j,k) * std::pow( Hwave_arr(i,j,k)/(Lwave_arr(i,j,k)+eps), amrex::Real(4.5) )
1015  + amrex::Real(0.11) * eta_arr(ie,je,k,EddyDiff::Mom_v) / ustar, z0_eps), z0_max );
1016  qflux = (spec_qflux) ? mdata.surf_moist_flux :
1017  -(qvm_arr(i,j,k) - q_surf_arr(i,j,k)) * ustar * mdata.kappa /
1018  (std::log(zref / z0_arr(i,j,k)) - psi_h); // <w'Qv'>
1019  tflux = mdata.surf_temp_flux*(1 + amrex::Real(0.61)*qvm_arr(i,j,k)) + qflux*amrex::Real(0.61)*tm_arr(i,j,k);
1020  if (w_star_arr) {
1021  // update w* and Umagmean
1022  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
1023  wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
1024  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
1025  umm = std::max(umm, WSMIN);
1026  }
1027  Olen = -ustar * ustar * ustar * tvm_arr(i,j,k) / (mdata.kappa * mdata.gravity * tflux);
1028  zeta = zref / Olen;
1029  psi_m = sfuns.calc_psi_m(zeta);
1030  psi_h = sfuns.calc_psi_h(zeta);
1031  u_star_arr(i,j,k) = mdata.kappa * umm / (std::log(zref / z0) - psi_m);
1032  ++iter;
1033  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
1034  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
1035  "Maximum number of MOST iterations reached.");
1036 
1037  // Populate the stored MOST arrays
1038  z0_arr(i,j,k) = z0;
1039  olen_arr(i,j,k) = Olen;
1040  t_surf_arr(i,j,k) = mdata.surf_temp_flux * (std::log(zref / z0) - psi_h) /
1041  (u_star_arr(i,j,k) * mdata.kappa) + tm_arr(i,j,k);
1042  t_star_arr(i,j,k) = -mdata.surf_temp_flux / u_star_arr(i,j,k);
1043  if (spec_qflux) {
1044  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (std::log(zref / z0_arr(i,j,k)) - psi_h) /
1045  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
1046  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
1047  } else {
1048  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) /
1049  (std::log(zref / z0_arr(i,j,k)) - psi_h);
1050  }
1051  }
1052 
1053 private:
1057  const amrex::Real tol = amrex::Real(1.0e-5);
1058 #ifdef AMREX_USE_FLOAT
1059  const amrex::Real eps = amrex::Real(1e-8);
1060 #else
1061  const amrex::Real eps = amrex::Real(1e-15);
1062 #endif
1063  const amrex::Real z0_eps = amrex::Real(1.0e-6);
1065  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
1066 };
1067 
1068 
1069 /**
1070  * Surface temperature with constant roughness
1071  */
1073 {
1075  amrex::Real Qvflux,
1076  bool cons_qflux)
1077  {
1078  mdata.surf_temp_flux = Tflux;
1079  mdata.surf_moist_flux = Qvflux;
1080  spec_qflux = cons_qflux;
1081  }
1082 
1083  AMREX_GPU_DEVICE
1084  AMREX_FORCE_INLINE
1085  void
1086  iterate_flux (const int& i,
1087  const int& j,
1088  const int& k,
1089  const int& max_iters,
1090  const amrex::Array4<const amrex::Real>& zref_arr,
1091  const amrex::Array4<const amrex::Real>& z0_arr,
1092  const amrex::Array4<const amrex::Real>& umm_arr,
1093  const amrex::Array4<const amrex::Real>& tm_arr,
1094  const amrex::Array4<const amrex::Real>& tvm_arr,
1095  const amrex::Array4<const amrex::Real>& qvm_arr,
1096  const amrex::Array4<amrex::Real>& u_star_arr,
1097  const amrex::Array4<amrex::Real>& w_star_arr,
1098  const amrex::Array4<amrex::Real>& t_star_arr,
1099  const amrex::Array4<amrex::Real>& q_star_arr,
1100  const amrex::Array4<amrex::Real>& t_surf_arr,
1101  const amrex::Array4<amrex::Real>& q_surf_arr,
1102  const amrex::Array4<amrex::Real>& olen_arr,
1103  const amrex::Array4<amrex::Real>& pblh_arr,
1104  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
1105  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
1106  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
1107  {
1108  amrex::Real Rib = zero;
1109  amrex::Real zeta = zero;
1110  amrex::Real zeta_old = zero;
1111  amrex::Real psi_m = zero;
1112  amrex::Real psi_h = zero;
1113  amrex::Real num = zero;
1114  amrex::Real den = zero;
1115  amrex::Real zref = zref_arr(i,j,k);
1116  amrex::Real z0 = z0_arr(i,j,k);
1117  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
1118  amrex::Real C = std::log(zref / z0);
1119 
1120  // First iteration we assume neutral (L -> inf)
1121  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) { olen_arr(i,j,k) = amrex::Real(1.0e3); }
1122  zeta = zref / olen_arr(i,j,k);
1123 
1124  // Water vapor in atmos and surface
1125  amrex::Real qv_s, qv_a;
1126  if (q_surf_arr(i,j,k) > zero) {
1127  qv_s = q_surf_arr(i,j,k);
1128  } else {
1129  // First iteration and no qv_surf was specified
1130  // Use mean since there will be no flux
1131  qv_s = qvm_arr(i,j,k);
1132  }
1133  qv_a = qvm_arr(i,j,k);
1134 
1135  // update w* and Umagmean from Beljaars (1995)
1136  if (w_star_arr) {
1137  // NOTE: Thv flux is lagged, similar to WRF
1138  psi_m = sfuns.calc_psi_m2(zeta);
1139  psi_h = sfuns.calc_psi_h2(zeta);
1140  amrex::Real ustar = mdata.kappa * umm / (C - psi_m);
1141  amrex::Real tstar = mdata.kappa * (tm_arr(i,j,k) - t_surf_arr(i,j,k)) / (C - psi_h);
1143  -ustar * mdata.kappa * (qv_a - qv_s) / (C - psi_h);
1144  amrex::Real tflux = -ustar*tstar*(1 + amrex::Real(0.61)*qvm_arr(i,j,k)) + amrex::Real(0.61)*tm_arr(i,j,k)*qflux;
1145  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
1146  amrex::Real wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
1147  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
1148  umm = std::max(umm, WSMIN);
1149  }
1150 
1151  // Bulk Richardson number w/ moisture
1152  amrex::Real thv_s = t_surf_arr(i,j,k) * (one + amrex::Real(0.61)*qv_s);
1153  amrex::Real thv_a = tm_arr(i,j,k) * (one + amrex::Real(0.61)*qv_a);
1154  Rib = ( (mdata.gravity * zref) / tm_arr(i,j,k) ) *
1155  ( (thv_a - thv_s) / (umm * umm) );
1156  Rib = std::min(std::max(Rib,-amrex::Real(4.0)),amrex::Real(4.0));
1157 
1158  // Fixed point iteration on zeta
1159  int iter = 0;
1160  do {
1161  // Transfer curr to old
1162  zeta_old = zeta;
1163 
1164  // Stability functions
1165  psi_m = sfuns.calc_psi_m2(zeta_old);
1166  psi_h = sfuns.calc_psi_h2(zeta_old);
1167 
1168  // Limiting
1169  num = std::max(C - psi_m, amrex::Real(1.0));
1170  den = std::max(C - psi_h, amrex::Real(1.0));
1171 
1172  // Update with under relaxation
1173  zeta = (one - alpha) * zeta_old + alpha * Rib * num * num / den;
1174 
1175  ++iter;
1176  } while ( (std::abs(zeta - zeta_old) > tol) && (iter <= max_iters) );
1177  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
1178  "Maximum number of MOST iterations reached.");
1179 
1180  // Populate the stored MOST arrays
1181  olen_arr(i,j,k) = zref / zeta;
1182  u_star_arr(i,j,k) = mdata.kappa * umm / (C - psi_m);
1183  t_star_arr(i,j,k) = mdata.kappa * (tm_arr(i,j,k) - t_surf_arr(i,j,k)) / (C - psi_h);
1184  if (spec_qflux) {
1185  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (C - psi_h) /
1186  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
1187  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
1188  } else {
1189  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) / (C - psi_h);
1190  }
1191  }
1192 
1193 private:
1197  const amrex::Real tol = amrex::Real(1.0e-3);
1199  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
1200 };
1201 
1202 
1203 /**
1204  * Surface temperature with charnock roughness
1205  */
1207 {
1209  amrex::Real Qvflux,
1210  amrex::Real cnk_a,
1211  bool cnk_visc,
1212  bool cons_qflux)
1213  {
1214  mdata.surf_temp_flux = Tflux;
1215  mdata.surf_moist_flux = Qvflux;
1216  mdata.Cnk_a = cnk_a;
1217  mdata.visc = cnk_visc;
1218  spec_qflux = cons_qflux;
1219  }
1220 
1221  AMREX_GPU_DEVICE
1222  AMREX_FORCE_INLINE
1223  void
1224  iterate_flux (const int& i,
1225  const int& j,
1226  const int& k,
1227  const int& max_iters,
1228  const amrex::Array4<const amrex::Real>& zref_arr,
1229  const amrex::Array4<amrex::Real>& z0_arr,
1230  const amrex::Array4<const amrex::Real>& umm_arr,
1231  const amrex::Array4<const amrex::Real>& tm_arr,
1232  const amrex::Array4<const amrex::Real>& tvm_arr,
1233  const amrex::Array4<const amrex::Real>& qvm_arr,
1234  const amrex::Array4<amrex::Real>& u_star_arr,
1235  const amrex::Array4<amrex::Real>& w_star_arr,
1236  const amrex::Array4<amrex::Real>& t_star_arr,
1237  const amrex::Array4<amrex::Real>& q_star_arr,
1238  const amrex::Array4<amrex::Real>& t_surf_arr,
1239  const amrex::Array4<amrex::Real>& q_surf_arr,
1240  const amrex::Array4<amrex::Real>& olen_arr,
1241  const amrex::Array4<amrex::Real>& pblh_arr,
1242  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
1243  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
1244  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
1245  {
1246  amrex::Real Rib = zero;
1247  amrex::Real zeta = zero;
1248  amrex::Real zeta_old = zero;
1249  amrex::Real psi_m = zero;
1250  amrex::Real psi_h = zero;
1251  amrex::Real num = zero;
1252  amrex::Real den = zero;
1253  amrex::Real z0 = z0_arr(i,j,k);
1254  amrex::Real z0_old = z0;
1255  amrex::Real ustar = zero;
1256  amrex::Real zref = zref_arr(i,j,k);
1257  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
1258  amrex::Real C = std::log(zref / z0);
1259 
1260  // First iteration we assume neutral (L -> inf)
1261  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) { olen_arr(i,j,k) = amrex::Real(1.0e3); }
1262  zeta = zref / olen_arr(i,j,k);
1263 
1264  // Water vapor in atmos and surface
1265  amrex::Real qv_s, qv_a;
1266  if (q_surf_arr(i,j,k) > zero) {
1267  qv_s = q_surf_arr(i,j,k);
1268  } else {
1269  // First iteration and no qv_surf was specified
1270  // Use mean since there will be no flux
1271  qv_s = qvm_arr(i,j,k);
1272  }
1273  qv_a = qvm_arr(i,j,k);
1274 
1275  // update w* and Umagmean from Beljaars (1995)
1276  if (w_star_arr) {
1277  // NOTE: Thv flux is lagged, similar to WRF
1278  psi_m = sfuns.calc_psi_m2(zeta);
1279  psi_h = sfuns.calc_psi_h2(zeta);
1280  ustar = mdata.kappa * umm / (C - psi_m);
1281  amrex::Real tstar = mdata.kappa * (tm_arr(i,j,k) - t_surf_arr(i,j,k)) / (C - psi_h);
1283  -ustar * mdata.kappa * (qv_a - qv_s) / (C - psi_h);
1284  amrex::Real tflux = -ustar*tstar*(1 + amrex::Real(0.61)*qvm_arr(i,j,k)) + amrex::Real(0.61)*tm_arr(i,j,k)*qflux;
1285  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
1286  amrex::Real wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
1287  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
1288  umm = std::max(umm, WSMIN);
1289  }
1290 
1291  // Bulk Richardson number w/ moisture
1292  amrex::Real thv_s = t_surf_arr(i,j,k) * (one + amrex::Real(0.61)*qv_s);
1293  amrex::Real thv_a = tm_arr(i,j,k) * (one + amrex::Real(0.61)*qv_a);
1294  Rib = ( (mdata.gravity * zref) / tm_arr(i,j,k) ) *
1295  ( (thv_a - thv_s) / (umm * umm) );
1296  Rib = std::min(std::max(Rib,-amrex::Real(4.0)),amrex::Real(4.0));
1297 
1298  // Fixed point iteration on zeta
1299  int iter = 0;
1300  do {
1301  // Transfer curr to old
1302  zeta_old = zeta;
1303 
1304  // Stability functions
1305  psi_m = sfuns.calc_psi_m2(zeta_old);
1306  psi_h = sfuns.calc_psi_h2(zeta_old);
1307 
1308  // Fixed point iteration on roughness
1309  int iter_z = 0;
1310  do {
1311  // Transfer curr to old
1312  z0_old = z0;
1313 
1314  // Update
1315  C = std::log(zref / z0_old);
1316  ustar = mdata.kappa * umm / (C - psi_m);
1317  if (mdata.Cnk_a > 0) {
1318  z0 = (mdata.Cnk_a / mdata.gravity) * ustar * ustar;
1319  if (mdata.visc) {
1320  z0 += air_viscosity(tm_arr(i,j,k)) / std::max(ustar, amrex::Real(0.05));
1321  }
1322  } else {
1323  z0 = COARE3_roughness(zref, umm, ustar);
1324  }
1325 
1326  ++iter_z;
1327  } while ( (std::abs(z0 - z0_old) > tol_z) && (iter_z <= max_iters) );
1328  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter_z < max_iters,
1329  "Maximum number of MOST roughness iterations reached.");
1330  C = std::log(zref / z0);
1331 
1332  // Limiting
1333  num = std::max(C - psi_m, amrex::Real(1.0));
1334  den = std::max(C - psi_h, amrex::Real(1.0));
1335 
1336  // Update with under relaxation
1337  zeta = (amrex::Real(1.0) - alpha) * zeta_old + alpha * Rib * num * num / den;
1338 
1339  ++iter;
1340  } while ( (std::abs(zeta - zeta_old) > tol) && (iter <= max_iters) );
1341  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
1342  "Maximum number of MOST iterations reached.");
1343 
1344  // Populate the stored MOST arrays
1345  z0_arr(i,j,k) = z0;
1346  olen_arr(i,j,k) = zref / zeta;
1347  u_star_arr(i,j,k) = mdata.kappa * umm / (C - psi_m);
1348  t_star_arr(i,j,k) = mdata.kappa * (tm_arr(i,j,k) - t_surf_arr(i,j,k)) / (C - psi_h);
1349  if (spec_qflux) {
1350  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (C - psi_h) /
1351  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
1352  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
1353  } else {
1354  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) / (C - psi_h);
1355  }
1356  }
1357 
1358 private:
1362  const amrex::Real tol = amrex::Real(1.0e-3);
1363  const amrex::Real tol_z = amrex::Real(1.0e-10);
1365  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
1366 };
1367 
1368 
1369 /**
1370  * Surface temperature with modified charnock roughness
1371  */
1373 {
1375  amrex::Real Qvflux,
1376  amrex::Real depth,
1377  bool cons_qflux)
1378  {
1379  mdata.surf_temp_flux = Tflux;
1380  mdata.surf_moist_flux = Qvflux;
1381  mdata.Cnk_d = depth;
1382  mdata.Cnk_b = mdata.Cnk_b1 * std::log(mdata.Cnk_b2 / mdata.Cnk_d);
1383  spec_qflux = cons_qflux;
1384  }
1385 
1386  AMREX_GPU_DEVICE
1387  AMREX_FORCE_INLINE
1388  void
1389  iterate_flux (const int& i,
1390  const int& j,
1391  const int& k,
1392  const int& max_iters,
1393  const amrex::Array4<const amrex::Real>& zref_arr,
1394  const amrex::Array4<amrex::Real>& z0_arr,
1395  const amrex::Array4<const amrex::Real>& umm_arr,
1396  const amrex::Array4<const amrex::Real>& tm_arr,
1397  const amrex::Array4<const amrex::Real>& tvm_arr,
1398  const amrex::Array4<const amrex::Real>& qvm_arr,
1399  const amrex::Array4<amrex::Real>& u_star_arr,
1400  const amrex::Array4<amrex::Real>& w_star_arr,
1401  const amrex::Array4<amrex::Real>& t_star_arr,
1402  const amrex::Array4<amrex::Real>& q_star_arr,
1403  const amrex::Array4<amrex::Real>& t_surf_arr,
1404  const amrex::Array4<amrex::Real>& q_surf_arr,
1405  const amrex::Array4<amrex::Real>& olen_arr,
1406  const amrex::Array4<amrex::Real>& pblh_arr,
1407  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
1408  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
1409  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
1410  {
1411  amrex::Real Rib = zero;
1412  amrex::Real zeta = zero;
1413  amrex::Real zeta_old = zero;
1414  amrex::Real psi_m = zero;
1415  amrex::Real psi_h = zero;
1416  amrex::Real num = zero;
1417  amrex::Real den = zero;
1418  amrex::Real z0 = z0_arr(i,j,k);
1419  amrex::Real z0_old = z0;
1420  amrex::Real ustar = zero;
1421  amrex::Real zref = zref_arr(i,j,k);
1422  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
1423  amrex::Real C = std::log(zref / z0);
1424 
1425  // First iteration we assume neutral (L -> inf)
1426  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) { olen_arr(i,j,k) = amrex::Real(1.0e3); }
1427  zeta = zref / olen_arr(i,j,k);
1428 
1429  // Water vapor in atmos and surface
1430  amrex::Real qv_s, qv_a;
1431  if (q_surf_arr(i,j,k) > zero) {
1432  qv_s = q_surf_arr(i,j,k);
1433  } else {
1434  // First iteration and no qv_surf was specified
1435  // Use mean since there will be no flux
1436  qv_s = qvm_arr(i,j,k);
1437  }
1438  qv_a = qvm_arr(i,j,k);
1439 
1440  // update w* and Umagmean from Beljaars (1995)
1441  if (w_star_arr) {
1442  // NOTE: Thv flux is lagged, similar to WRF
1443  psi_m = sfuns.calc_psi_m2(zeta);
1444  psi_h = sfuns.calc_psi_h2(zeta);
1445  ustar = mdata.kappa * umm / (C - psi_m);
1446  amrex::Real tstar = mdata.kappa * (tm_arr(i,j,k) - t_surf_arr(i,j,k)) / (C - psi_h);
1448  -ustar * mdata.kappa * (qv_a - qv_s) / (C - psi_h);
1449  amrex::Real tflux = -ustar*tstar*(1 + amrex::Real(0.61)*qvm_arr(i,j,k)) + amrex::Real(0.61)*tm_arr(i,j,k)*qflux;
1450  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
1451  amrex::Real wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
1452  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
1453  umm = std::max(umm, WSMIN);
1454  }
1455 
1456  // Bulk Richardson number w/ moisture
1457  amrex::Real thv_s = t_surf_arr(i,j,k) * (amrex::Real(1.0) + amrex::Real(0.61)*qv_s);
1458  amrex::Real thv_a = tm_arr(i,j,k) * (amrex::Real(1.0) + amrex::Real(0.61)*qv_a);
1459  Rib = ( (mdata.gravity * zref) / tm_arr(i,j,k) ) *
1460  ( (thv_a - thv_s) / (umm * umm) );
1461  Rib = std::min(std::max(Rib,-amrex::Real(4.0)),amrex::Real(4.0));
1462 
1463  // Fixed point iteration on zeta
1464  int iter = 0;
1465  do {
1466  // Transfer curr to old
1467  zeta_old = zeta;
1468 
1469  // Stability functions
1470  psi_m = sfuns.calc_psi_m2(zeta_old);
1471  psi_h = sfuns.calc_psi_h2(zeta_old);
1472 
1473  // Fixed point iteration on roughness
1474  int iter_z = 0;
1475  do {
1476  // Transfer curr to old
1477  z0_old = z0;
1478 
1479  // Update
1480  C = std::log(zref / z0_old);
1481  ustar = mdata.kappa * umm / (C - psi_m);
1482  z0 = std::exp( (amrex::Real(2.7)*ustar - amrex::Real(1.8)/mdata.Cnk_b) / (ustar + amrex::Real(0.17)/mdata.Cnk_b) );
1483 
1484  ++iter_z;
1485  } while ( (std::abs(z0 - z0_old) > tol_z) && (iter_z <= max_iters) );
1486  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter_z < max_iters,
1487  "Maximum number of MOST roughness iterations reached.");
1488  C = std::log(zref / z0);
1489 
1490  // Limiting
1491  num = std::max(C - psi_m, amrex::Real(1.0));
1492  den = std::max(C - psi_h, amrex::Real(1.0));
1493 
1494  // Update with under relaxation
1495  zeta = (amrex::Real(1.0) - alpha) * zeta_old + alpha * Rib * num * num / den;
1496 
1497  ++iter;
1498  } while ( (std::abs(zeta - zeta_old) > tol) && (iter <= max_iters) );
1499  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
1500  "Maximum number of MOST iterations reached.");
1501 
1502  // Populate the stored MOST arrays
1503  z0_arr(i,j,k) = z0;
1504  olen_arr(i,j,k) = zref / zeta;
1505  u_star_arr(i,j,k) = mdata.kappa * umm / (C - psi_m);
1506  t_star_arr(i,j,k) = mdata.kappa * (tm_arr(i,j,k) - t_surf_arr(i,j,k)) / (C - psi_h);
1507  if (spec_qflux) {
1508  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (C - psi_h) /
1509  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
1510  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
1511  } else {
1512  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) / (C - psi_h);
1513  }
1514  }
1515 
1516 private:
1520  const amrex::Real tol = amrex::Real(1.0e-3);
1521  const amrex::Real tol_z = amrex::Real(1.0e-10);
1523  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
1524 };
1525 
1526 
1527 /**
1528  * Surface temperature with donelan roughness
1529  */
1531 {
1533  amrex::Real Qvflux,
1534  bool cons_qflux)
1535  {
1536  mdata.surf_temp_flux = Tflux;
1537  mdata.surf_moist_flux = Qvflux;
1538  spec_qflux = cons_qflux;
1539  }
1540 
1541  AMREX_GPU_DEVICE
1542  AMREX_FORCE_INLINE
1543  void
1544  iterate_flux (const int& i,
1545  const int& j,
1546  const int& k,
1547  const int& max_iters,
1548  const amrex::Array4<const amrex::Real>& zref_arr,
1549  const amrex::Array4<amrex::Real>& z0_arr,
1550  const amrex::Array4<const amrex::Real>& umm_arr,
1551  const amrex::Array4<const amrex::Real>& tm_arr,
1552  const amrex::Array4<const amrex::Real>& tvm_arr,
1553  const amrex::Array4<const amrex::Real>& qvm_arr,
1554  const amrex::Array4<amrex::Real>& u_star_arr,
1555  const amrex::Array4<amrex::Real>& w_star_arr,
1556  const amrex::Array4<amrex::Real>& t_star_arr,
1557  const amrex::Array4<amrex::Real>& q_star_arr,
1558  const amrex::Array4<amrex::Real>& t_surf_arr,
1559  const amrex::Array4<amrex::Real>& q_surf_arr,
1560  const amrex::Array4<amrex::Real>& olen_arr,
1561  const amrex::Array4<amrex::Real>& pblh_arr,
1562  const amrex::Array4<amrex::Real>& /*Hwave_arr*/,
1563  const amrex::Array4<amrex::Real>& /*Lwave_arr*/,
1564  const amrex::Array4<amrex::Real>& /*eta_arr*/) const
1565  {
1566  int iter = 0;
1567  amrex::Real ustar = zero;
1568  amrex::Real wstar = zero;
1569  amrex::Real z0 = zero;
1570  amrex::Real tflux = zero;
1571  amrex::Real qflux = zero;
1572  amrex::Real zeta = zero;
1573  amrex::Real psi_m = zero;
1574  amrex::Real psi_h = zero;
1575  amrex::Real Olen = zero;
1576  amrex::Real Oleno = zero;
1577  amrex::Real zref = zref_arr(i,j,k);
1578  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
1579  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) {
1580  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
1581  } else {
1582  Olen = olen_arr(i,j,k);
1583  Oleno = Olen;
1584  zeta = zref / Olen;
1585  psi_m = sfuns.calc_psi_m(zeta);
1586  psi_h = sfuns.calc_psi_h(zeta);
1587  }
1588  do {
1589  ustar = u_star_arr(i,j,k);
1590  z0 = Donelan_roughness(ustar);
1591  tflux = -(tm_arr(i,j,k) - t_surf_arr(i,j,k)) * ustar * mdata.kappa /
1592  (std::log(zref / z0) - psi_h); // <w'T'>
1593  tflux *= (1 + amrex::Real(0.61)*qvm_arr(i,j,k));
1594  qflux = (spec_qflux) ? mdata.surf_moist_flux :
1595  -(qvm_arr(i,j,k) - q_surf_arr(i,j,k)) * ustar * mdata.kappa /
1596  (std::log(zref / z0_arr(i,j,k)) - psi_h); // <w'Qv'>
1597  tflux += amrex::Real(0.61)*tm_arr(i,j,k) * qflux; // ~= <w'Tv'>
1598  if (w_star_arr) {
1599  // update w* and Umagmean
1600  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
1601  wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
1602  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
1603  umm = std::max(umm, WSMIN);
1604  }
1605  Olen = -ustar * ustar * ustar * tvm_arr(i,j,k) / (mdata.kappa * mdata.gravity * tflux);
1606  if ( (((Olen >= zero) && (Oleno <= zero)) ||
1607  ((Olen <= zero) && (Oleno >= zero))) &&
1608  std::fabs(Olen) + std::fabs(Oleno) < amrex::Real(1.0)) {
1609  Olen = myhalf * (Olen + Oleno);
1610  }
1611  Oleno = Olen;
1612  zeta = zref / Olen;
1613  psi_m = sfuns.calc_psi_m(zeta);
1614  psi_h = sfuns.calc_psi_h(zeta);
1615  u_star_arr(i,j,k) = mdata.kappa * umm / (std::log(zref / z0) - psi_m);
1616  ++iter;
1617  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
1618  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
1619  "Maximum number of MOST iterations reached.");
1620 
1621  // Populate the stored MOST arrays
1622  z0_arr(i,j,k) = z0;
1623  olen_arr(i,j,k) = Olen;
1624  t_star_arr(i,j,k) = mdata.kappa * (tm_arr(i,j,k) - t_surf_arr(i,j,k)) /
1625  (std::log(zref / z0) - psi_h);
1626  if (spec_qflux) {
1627  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (std::log(zref / z0_arr(i,j,k)) - psi_h) /
1628  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
1629  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
1630  } else {
1631  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) /
1632  (std::log(zref / z0_arr(i,j,k)) - psi_h);
1633  }
1634  }
1635 
1636 private:
1640  const amrex::Real tol = amrex::Real(1.0e-5);
1641  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
1642 };
1643 
1644 
1645 /**
1646  * Surface temperature with wave-coupled roughness
1647  */
1649 {
1651  amrex::Real Qvflux,
1652  bool cons_qflux)
1653  {
1654  mdata.surf_temp_flux = Tflux;
1655  mdata.surf_moist_flux = Qvflux;
1656  spec_qflux = cons_qflux;
1657  }
1658 
1659  AMREX_GPU_DEVICE
1660  AMREX_FORCE_INLINE
1661  void
1662  iterate_flux (const int& i,
1663  const int& j,
1664  const int& k,
1665  const int& max_iters,
1666  const amrex::Array4<const amrex::Real>& zref_arr,
1667  const amrex::Array4<amrex::Real>& z0_arr,
1668  const amrex::Array4<const amrex::Real>& umm_arr,
1669  const amrex::Array4<const amrex::Real>& tm_arr,
1670  const amrex::Array4<const amrex::Real>& tvm_arr,
1671  const amrex::Array4<const amrex::Real>& qvm_arr,
1672  const amrex::Array4<amrex::Real>& u_star_arr,
1673  const amrex::Array4<amrex::Real>& w_star_arr,
1674  const amrex::Array4<amrex::Real>& t_star_arr,
1675  const amrex::Array4<amrex::Real>& q_star_arr,
1676  const amrex::Array4<amrex::Real>& t_surf_arr,
1677  const amrex::Array4<amrex::Real>& q_surf_arr,
1678  const amrex::Array4<amrex::Real>& olen_arr,
1679  const amrex::Array4<amrex::Real>& pblh_arr,
1680  const amrex::Array4<amrex::Real>& Hwave_arr,
1681  const amrex::Array4<amrex::Real>& Lwave_arr,
1682  const amrex::Array4<amrex::Real>& eta_arr) const
1683  {
1684  int iter = 0;
1685  amrex::Real ustar = zero;
1686  amrex::Real wstar = zero;
1687  amrex::Real z0 = zero;
1688  amrex::Real tflux = zero;
1689  amrex::Real qflux = zero;
1690  amrex::Real zeta = zero;
1691  amrex::Real psi_m = zero;
1692  amrex::Real psi_h = zero;
1693  amrex::Real Olen = zero;
1694  amrex::Real Oleno = zero;
1695  amrex::Real zref = zref_arr(i,j,k);
1696  int ie, je;
1697  ie = i < lbound(eta_arr).x ? lbound(eta_arr).x : i;
1698  je = j < lbound(eta_arr).y ? lbound(eta_arr).y : j;
1699  ie = ie > ubound(eta_arr).x ? ubound(eta_arr).x : ie;
1700  je = je > ubound(eta_arr).y ? ubound(eta_arr).y : je;
1701  amrex::Real umm = std::max(umm_arr(i,j,k), WSMIN);
1702  if (u_star_arr(i,j,k) == amrex::Real(1.e34)) {
1703  u_star_arr(i,j,k) = mdata.kappa * umm / std::log(zref / z0_arr(i,j,k));
1704  } else {
1705  Olen = olen_arr(i,j,k);
1706  Oleno = Olen;
1707  zeta = zref / Olen;
1708  psi_m = sfuns.calc_psi_m(zeta);
1709  psi_h = sfuns.calc_psi_h(zeta);
1710  }
1711  do {
1712  ustar = u_star_arr(i,j,k);
1713  z0 = std::min( std::max(amrex::Real(1200.0) * Hwave_arr(i,j,k) * std::pow( Hwave_arr(i,j,k)/(Lwave_arr(i,j,k)+eps), amrex::Real(4.5) )
1714  + amrex::Real(0.11) * eta_arr(ie,je,k,EddyDiff::Mom_v) / ustar, z0_eps), z0_max );
1715  tflux = -(tm_arr(i,j,k) - t_surf_arr(i,j,k)) * ustar * mdata.kappa /
1716  (std::log(zref / z0) - psi_h); // <w'T'>
1717  tflux *= (1 + amrex::Real(0.61)*qvm_arr(i,j,k));
1718  qflux = (spec_qflux) ? mdata.surf_moist_flux :
1719  -(qvm_arr(i,j,k) - q_surf_arr(i,j,k)) * ustar * mdata.kappa /
1720  (std::log(zref / z0_arr(i,j,k)) - psi_h); // <w'Qv'>
1721  tflux += amrex::Real(0.61)*tm_arr(i,j,k) * qflux; // ~= <w'Tv'>
1722  if (w_star_arr) {
1723  // update w* and Umagmean
1724  w_star_arr(i,j,k) = calc_wstar(tflux, pblh_arr(i,j,k), tvm_arr(i,j,k));
1725  wstar = mdata.Bjr_beta * w_star_arr(i,j,k);
1726  umm = std::sqrt(umm_arr(i,j,k)*umm_arr(i,j,k) + wstar*wstar);
1727  umm = std::max(umm, WSMIN);
1728  }
1729  Olen = -ustar * ustar * ustar * tvm_arr(i,j,k) / (mdata.kappa * mdata.gravity * tflux);
1730  if ( (((Olen >= zero) && (Oleno <= zero)) ||
1731  ((Olen <= zero) && (Oleno >= zero))) &&
1732  std::fabs(Olen) + std::fabs(Oleno) < amrex::Real(1.0)) {
1733  Olen = myhalf * (Olen + Oleno);
1734  }
1735  Oleno = Olen;
1736  zeta = zref / Olen;
1737  psi_m = sfuns.calc_psi_m(zeta);
1738  psi_h = sfuns.calc_psi_h(zeta);
1739  u_star_arr(i,j,k) = mdata.kappa * umm / (std::log(zref / z0) - psi_m);
1740  ++iter;
1741  } while ((std::abs(u_star_arr(i,j,k) - ustar) > tol) && iter <= max_iters);
1742  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < max_iters,
1743  "Maximum number of MOST iterations reached.");
1744 
1745  // Populate the stored MOST arrays
1746  z0_arr(i,j,k) = z0;
1747  olen_arr(i,j,k) = Olen;
1748  t_star_arr(i,j,k) = mdata.kappa * (tm_arr(i,j,k) - t_surf_arr(i,j,k)) /
1749  (std::log(zref / z0) - psi_h);
1750  if (spec_qflux) {
1751  q_surf_arr(i,j,k) = mdata.surf_moist_flux * (std::log(zref / z0_arr(i,j,k)) - psi_h) /
1752  (u_star_arr(i,j,k) * mdata.kappa) + qvm_arr(i,j,k);
1753  q_star_arr(i,j,k) = -mdata.surf_moist_flux / u_star_arr(i,j,k);
1754  } else {
1755  q_star_arr(i,j,k) = mdata.kappa * (qvm_arr(i,j,k) - q_surf_arr(i,j,k)) /
1756  (std::log(zref / z0_arr(i,j,k)) - psi_h);
1757  }
1758  }
1759 
1760 private:
1764  const amrex::Real tol = amrex::Real(1.0e-5);
1765 #ifdef AMREX_USE_FLOAT
1766  const amrex::Real eps = amrex::Real(1e-8);
1767 #else
1768  const amrex::Real eps = amrex::Real(1e-15);
1769 #endif
1770  const amrex::Real z0_eps = amrex::Real(1.0e-6);
1772  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
1773 };
1774 
1775 
1776 /**
1777  * Moeng flux formulation
1778  */
1780 {
1782 
1783  AMREX_GPU_DEVICE
1784  AMREX_FORCE_INLINE
1785  amrex::Real
1786  compute_q_flux (const int& i,
1787  const int& j,
1788  const int& k,
1789  const amrex::Array4<const amrex::Real>& cons_arr,
1790  const amrex::Array4<const amrex::Real>& velx_arr,
1791  const amrex::Array4<const amrex::Real>& vely_arr,
1792  const amrex::Array4<const amrex::Real>& umm_arr,
1793  const amrex::Array4<const amrex::Real>& qvm_arr,
1794  const amrex::Array4<const amrex::Real>& u_star_arr,
1795  const amrex::Array4<const amrex::Real>& q_star_arr,
1796  const amrex::Array4<const amrex::Real>& q_surf_arr) const
1797  {
1798  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
1799  amrex::Real qv = cons_arr(i,j,k,RhoQ1_comp) / rho;
1800  amrex::Real velx = myhalf * ( velx_arr(i,j,k) + velx_arr(i+1,j ,k) );
1801  amrex::Real vely = myhalf * ( vely_arr(i,j,k) + vely_arr(i ,j+1,k) );
1802 
1803  amrex::Real qv_mean = qvm_arr(i,j,0);
1804  amrex::Real ustar = u_star_arr(i,j,0);
1805  amrex::Real qstar = q_star_arr(i,j,0);
1806  amrex::Real qv_surf = q_surf_arr(i,j,0);
1807  amrex::Real wsp_mean = umm_arr(i,j,0);
1808  wsp_mean = std::max(wsp_mean, WSMIN);
1809 
1810  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
1811  amrex::Real num1 = wsp * (qv_mean-qv_surf);
1812  amrex::Real num2 = wsp_mean * (qv-qv_mean);
1813 
1814  // NOTE: this is rho*<Qv'w'> = -K dQvdz
1815  amrex::Real moflux = (std::abs(qstar) > eps) ?
1816  -rho*qstar*ustar*(num1+num2)/((qv_mean-qv_surf)*wsp_mean) : zero;
1817 
1818  return moflux;
1819  }
1820 
1821  AMREX_GPU_DEVICE
1822  AMREX_FORCE_INLINE
1823  amrex::Real
1824  compute_t_flux (const int& i,
1825  const int& j,
1826  const int& k,
1827  const amrex::Array4<const amrex::Real>& cons_arr,
1828  const amrex::Array4<const amrex::Real>& velx_arr,
1829  const amrex::Array4<const amrex::Real>& vely_arr,
1830  const amrex::Array4<const amrex::Real>& umm_arr,
1831  const amrex::Array4<const amrex::Real>& tm_arr,
1832  const amrex::Array4<const amrex::Real>& u_star_arr,
1833  const amrex::Array4<const amrex::Real>& t_star_arr,
1834  const amrex::Array4<const amrex::Real>& t_surf_arr) const
1835  {
1836  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
1837  amrex::Real theta = cons_arr(i,j,k,RhoTheta_comp) / rho;
1838  amrex::Real velx = myhalf * ( velx_arr(i,j,k) + velx_arr(i+1,j ,k) );
1839  amrex::Real vely = myhalf * ( vely_arr(i,j,k) + vely_arr(i ,j+1,k) );
1840 
1841  amrex::Real theta_mean = tm_arr(i,j,0);
1842  amrex::Real ustar = u_star_arr(i,j,0);
1843  amrex::Real tstar = t_star_arr(i,j,0);
1844  amrex::Real theta_surf = t_surf_arr(i,j,0);
1845  amrex::Real wsp_mean = umm_arr(i,j,0);
1846  wsp_mean = std::max(wsp_mean, WSMIN);
1847 
1848  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
1849  amrex::Real num1 = wsp * (theta_mean-theta_surf);
1850  amrex::Real num2 = wsp_mean * (theta-theta_mean);
1851 
1852  // NOTE: this is rho*<T'w'> = -K dTdz
1853  amrex::Real moflux = (std::abs(tstar) > eps) ?
1854  -rho*tstar*ustar*(num1+num2)/((theta_mean-theta_surf)*wsp_mean) : zero;
1855 
1856  return moflux;
1857  }
1858 
1859  AMREX_GPU_DEVICE
1860  AMREX_FORCE_INLINE
1861  amrex::Real
1862  compute_u_flux (const int& i,
1863  const int& j,
1864  const int& k,
1865  const amrex::Array4<const amrex::Real>& cons_arr,
1866  const amrex::Array4<const amrex::Real>& velx_arr,
1867  const amrex::Array4<const amrex::Real>& vely_arr,
1868  const amrex::Array4<const amrex::Real>& umm_arr,
1869  const amrex::Array4<const amrex::Real>& um_arr,
1870  const amrex::Array4<const amrex::Real>& u_star_arr) const
1871  {
1872  amrex::Real velx = velx_arr(i,j,k);
1873  amrex::Real vely = fourth * ( vely_arr(i ,j,k) + vely_arr(i ,j+1,k)
1874  + vely_arr(i-1,j,k) + vely_arr(i-1,j+1,k) );
1875  amrex::Real rho = myhalf * ( cons_arr(i-1,j,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
1876 
1877  amrex::Real umean = um_arr(i,j,0);
1878  amrex::Real ustar = myhalf * ( u_star_arr(i-1,j,0) + u_star_arr(i,j,0) );
1879  amrex::Real wsp_mean = myhalf * ( umm_arr(i-1,j,0) + umm_arr(i,j,0) );
1880  wsp_mean = std::max(wsp_mean, WSMIN);
1881 
1882  // Note: The surface mean shear stress is decomposed into tau_xz by
1883  // multiplying the modeled shear stress (rho*ustar^2) with
1884  // a factor of umean/wsp_mean for directionality; this factor
1885  // modifies the denominator from what is in Moeng amrex::Real(1984.)
1886  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
1887  amrex::Real num1 = wsp * umean;
1888  amrex::Real num2 = wsp_mean * (velx-umean);
1889 
1890  // NOTE: this is rho*<u'w'> = -K dudz
1891  amrex::Real stressx = -rho*ustar*ustar * (num1+num2)/(wsp_mean*wsp_mean);
1892 
1893  return stressx;
1894  }
1895 
1896  AMREX_GPU_DEVICE
1897  AMREX_FORCE_INLINE
1898  amrex::Real
1899  compute_v_flux (const int& i,
1900  const int& j,
1901  const int& k,
1902  const amrex::Array4<const amrex::Real>& cons_arr,
1903  const amrex::Array4<const amrex::Real>& velx_arr,
1904  const amrex::Array4<const amrex::Real>& vely_arr,
1905  const amrex::Array4<const amrex::Real>& umm_arr,
1906  const amrex::Array4<const amrex::Real>& vm_arr,
1907  const amrex::Array4<const amrex::Real>& u_star_arr) const
1908  {
1909  amrex::Real velx = fourth * ( velx_arr(i,j ,k) + velx_arr(i+1,j ,k)
1910  + velx_arr(i,j-1,k) + velx_arr(i+1,j-1,k) );
1911  amrex::Real vely = vely_arr(i,j,k);
1912  amrex::Real rho = myhalf * ( cons_arr(i,j-1,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
1913 
1914  amrex::Real vmean = vm_arr(i,j,0);
1915  amrex::Real ustar = myhalf * ( u_star_arr(i,j-1,0) + u_star_arr(i,j,0) );
1916  amrex::Real wsp_mean = myhalf * ( umm_arr(i,j-1,0) + umm_arr(i,j,0) );
1917  wsp_mean = std::max(wsp_mean, WSMIN);
1918 
1919  // Note: The surface mean shear stress is decomposed into tau_yz by
1920  // multiplying the modeled shear stress (rho*ustar^2) with
1921  // a factor of vmean/wsp_mean for directionality; this factor
1922  // modifies the denominator from what is in Moeng amrex::Real(1984.)
1923  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
1924  amrex::Real num1 = wsp * vmean;
1925  amrex::Real num2 = wsp_mean * (vely-vmean);
1926 
1927  // NOTE: this is rho*<v'w'> = -K dvdz
1928  amrex::Real stressy = -rho*ustar*ustar * (num1+num2)/(wsp_mean*wsp_mean);
1929 
1930  return stressy;
1931  }
1932 
1933 private:
1934 #ifdef AMREX_USE_FLOAT
1935  const amrex::Real eps = amrex::Real(1e-8);
1936 #else
1937  const amrex::Real eps = amrex::Real(1e-15);
1938 #endif
1939  const amrex::Real WSMIN = amrex::Real(0.1); // minimum wind speed
1940 };
1941 
1942 
1943 /**
1944  * Donelan flux formulation
1945  */
1947 {
1949 
1950  AMREX_GPU_DEVICE
1951  AMREX_FORCE_INLINE
1952  amrex::Real
1953  compute_q_flux (const int& /*i*/,
1954  const int& /*j*/,
1955  const int& /*k*/,
1956  const amrex::Array4<const amrex::Real>& /*cons_arr*/,
1957  const amrex::Array4<const amrex::Real>& /*velx_arr*/,
1958  const amrex::Array4<const amrex::Real>& /*vely_arr*/,
1959  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
1960  const amrex::Array4<const amrex::Real>& /*qm_arr*/,
1961  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
1962  const amrex::Array4<const amrex::Real>& /*q_star_arr*/,
1963  const amrex::Array4<const amrex::Real>& /*q_surf_arr*/) const
1964  {
1965  // NOTE: this is rho*<q'w'> = -K dqdz
1966  amrex::Real moflux = zero;
1967 
1968  return moflux;
1969  }
1970 
1971  AMREX_GPU_DEVICE
1972  AMREX_FORCE_INLINE
1973  amrex::Real
1974  compute_t_flux (const int& i,
1975  const int& j,
1976  const int& k,
1977  const amrex::Array4<const amrex::Real>& cons_arr,
1978  const amrex::Array4<const amrex::Real>& /*velx_arr*/,
1979  const amrex::Array4<const amrex::Real>& /*vely_arr*/,
1980  const amrex::Array4<const amrex::Real>& umm_arr,
1981  const amrex::Array4<const amrex::Real>& tm_arr,
1982  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
1983  const amrex::Array4<const amrex::Real>& /*t_star_arr*/,
1984  const amrex::Array4<const amrex::Real>& t_surf_arr) const
1985  {
1986  amrex::Real Ch = amrex::Real(0.0012);
1987  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
1988 
1989  amrex::Real theta_surf = t_surf_arr(i,j,0);
1990  amrex::Real theta_mean = tm_arr(i,j,0);
1991  amrex::Real wsp_mean = umm_arr(i,j,0);
1992 
1993  // NOTE: this is rho*<T'w'> = -K dTdz
1994  amrex::Real moflux = -rho * Ch * wsp_mean * (theta_mean - theta_surf);
1995 
1996  return moflux;
1997  }
1998 
1999  AMREX_GPU_DEVICE
2000  AMREX_FORCE_INLINE
2001  amrex::Real
2002  compute_u_flux (const int& i,
2003  const int& j,
2004  const int& k,
2005  const amrex::Array4<const amrex::Real>& cons_arr,
2006  const amrex::Array4<const amrex::Real>& velx_arr,
2007  const amrex::Array4<const amrex::Real>& vely_arr,
2008  const amrex::Array4<const amrex::Real>& umm_arr,
2009  const amrex::Array4<const amrex::Real>& /*um_arr*/,
2010  const amrex::Array4<const amrex::Real>& /*u_star_arr*/) const
2011  {
2012  amrex::Real velx = velx_arr(i,j,k);
2013  amrex::Real vely = fourth * ( vely_arr(i ,j,k) + vely_arr(i ,j+1,k)
2014  + vely_arr(i-1,j,k) + vely_arr(i-1,j+1,k) );
2015  amrex::Real rho = myhalf * ( cons_arr(i-1,j,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
2016 
2017  amrex::Real Cd = amrex::Real(0.001);
2018  const amrex::Real c = amrex::Real(7e-5);
2019  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2020  amrex::Real wsp_mean = myhalf * ( umm_arr(i-1,j,0) + umm_arr(i,j,0) );
2021  if (wsp_mean <= amrex::Real(5.0)) {
2022  Cd = amrex::Real(0.001);
2023  } else if (wsp_mean < amrex::Real(25.0) && wsp_mean > amrex::Real(5.0)) {
2024  Cd = amrex::Real(0.001) + c * (wsp_mean - amrex::Real(5.0));
2025  } else {
2026  Cd = amrex::Real(0.0024);
2027  }
2028 
2029  // NOTE: this is rho*<u'w'> = -K dudz
2030  amrex::Real stressx = -rho * Cd * velx * wsp;
2031 
2032  return stressx;
2033  }
2034 
2035  AMREX_GPU_DEVICE
2036  AMREX_FORCE_INLINE
2037  amrex::Real
2038  compute_v_flux (const int& i,
2039  const int& j,
2040  const int& k,
2041  const amrex::Array4<const amrex::Real>& cons_arr,
2042  const amrex::Array4<const amrex::Real>& velx_arr,
2043  const amrex::Array4<const amrex::Real>& vely_arr,
2044  const amrex::Array4<const amrex::Real>& umm_arr,
2045  const amrex::Array4<const amrex::Real>& /*vm_arr*/,
2046  const amrex::Array4<const amrex::Real>& /*u_star_arr*/) const
2047  {
2048  amrex::Real velx = fourth * ( velx_arr(i,j ,k) + velx_arr(i+1,j ,k)
2049  + velx_arr(i,j-1,k) + velx_arr(i+1,j-1,k) );
2050  amrex::Real vely = vely_arr(i,j,k);
2051  amrex::Real rho = myhalf * ( cons_arr(i,j-1,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
2052 
2053  amrex::Real Cd = amrex::Real(0.001);
2054  const amrex::Real c = amrex::Real(7e-5);
2055  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2056  amrex::Real wsp_mean = myhalf * ( umm_arr(i,j-1,0) + umm_arr(i,j,0) );
2057  if (wsp_mean <= amrex::Real(5.0)) {
2058  Cd = amrex::Real(0.001);
2059  } else if (wsp_mean < amrex::Real(25.0) && wsp_mean > amrex::Real(5.0)) {
2060  Cd = amrex::Real(0.001) + c * (wsp_mean - amrex::Real(5.0));
2061  } else {
2062  Cd = amrex::Real(0.0024);
2063  }
2064 
2065  // NOTE: this is rho*<v'w'> = -K dvdz
2066  amrex::Real stressy = -rho * Cd * vely * wsp;
2067 
2068  return stressy;
2069  }
2070 
2071 private:
2072 
2073 };
2074 
2075 
2076 /**
2077  * Custom flux formulation
2078  */
2080 {
2081  custom_flux (bool specified_rho_surf)
2082  : fluxes_include_rho(specified_rho_surf)
2083  {}
2084 
2085  AMREX_GPU_DEVICE
2086  AMREX_FORCE_INLINE
2087  amrex::Real
2088  compute_q_flux (const int& i,
2089  const int& j,
2090  const int& k,
2091  const amrex::Array4<const amrex::Real>& cons_arr,
2092  const amrex::Array4<const amrex::Real>& /*velx_arr*/,
2093  const amrex::Array4<const amrex::Real>& /*vely_arr*/,
2094  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2095  const amrex::Array4<const amrex::Real>& /*qm_arr*/,
2096  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
2097  const amrex::Array4<const amrex::Real>& q_star_arr,
2098  const amrex::Array4<const amrex::Real>& /*q_surf_arr*/) const
2099  {
2100  amrex::Real rho = (fluxes_include_rho) ? one : cons_arr(i,j,k,Rho_comp);
2101  amrex::Real qstar = q_star_arr(i,j,0);
2102 
2103  // NOTE: this is rho*<q'w'> = -K dqdz
2104  amrex::Real moflux = (std::abs(qstar) > eps) ? rho * qstar : zero;
2105 
2106  return moflux;
2107  }
2108 
2109  AMREX_GPU_DEVICE
2110  AMREX_FORCE_INLINE
2111  amrex::Real
2112  compute_t_flux (const int& i,
2113  const int& j,
2114  const int& k,
2115  const amrex::Array4<const amrex::Real>& cons_arr,
2116  const amrex::Array4<const amrex::Real>& /*velx_arr*/,
2117  const amrex::Array4<const amrex::Real>& /*vely_arr*/,
2118  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2119  const amrex::Array4<const amrex::Real>& /*tm_arr*/,
2120  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
2121  const amrex::Array4<const amrex::Real>& t_star_arr,
2122  const amrex::Array4<const amrex::Real>& /*t_surf_arr*/) const
2123  {
2124  amrex::Real rho = (fluxes_include_rho) ? one : cons_arr(i,j,k,Rho_comp);
2125  amrex::Real tstar = t_star_arr(i,j,0);
2126 
2127  // NOTE: this is rho*<T'w'> = -K dTdz
2128  amrex::Real moflux = (std::abs(tstar) > eps) ? rho * tstar : zero;
2129 
2130  return moflux;
2131  }
2132 
2133  AMREX_GPU_DEVICE
2134  AMREX_FORCE_INLINE
2135  amrex::Real
2136  compute_u_flux (const int& i,
2137  const int& j,
2138  const int& k,
2139  const amrex::Array4<const amrex::Real>& cons_arr,
2140  const amrex::Array4<const amrex::Real>& velx_arr,
2141  const amrex::Array4<const amrex::Real>& vely_arr,
2142  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2143  const amrex::Array4<const amrex::Real>& /*um_arr*/,
2144  const amrex::Array4<const amrex::Real>& u_star_arr) const
2145  {
2146  amrex::Real velx = velx_arr(i,j,k);
2147  amrex::Real vely = fourth * ( vely_arr(i ,j,k) + vely_arr(i ,j+1,k)
2148  + vely_arr(i-1,j,k) + vely_arr(i-1,j+1,k) );
2149  amrex::Real rho = (fluxes_include_rho) ? one : myhalf * ( cons_arr(i-1,j,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
2150 
2151  amrex::Real ustar = myhalf * ( u_star_arr(i-1,j,0) + u_star_arr(i,j,0) );
2152  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2153 
2154  // NOTE: this is rho*<u'w'> = -K dudz
2155  amrex::Real stressx = -rho * ustar * ustar * velx / wsp;
2156 
2157  return stressx;
2158  }
2159 
2160  AMREX_GPU_DEVICE
2161  AMREX_FORCE_INLINE
2162  amrex::Real
2163  compute_v_flux (const int& i,
2164  const int& j,
2165  const int& k,
2166  const amrex::Array4<const amrex::Real>& cons_arr,
2167  const amrex::Array4<const amrex::Real>& velx_arr,
2168  const amrex::Array4<const amrex::Real>& vely_arr,
2169  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2170  const amrex::Array4<const amrex::Real>& /*vm_arr*/,
2171  const amrex::Array4<const amrex::Real>& u_star_arr) const
2172  {
2173  amrex::Real velx = fourth * ( velx_arr(i,j ,k) + velx_arr(i+1,j ,k)
2174  + velx_arr(i,j-1,k) + velx_arr(i+1,j-1,k) );
2175  amrex::Real vely = vely_arr(i,j,k);
2176  amrex::Real rho = (fluxes_include_rho) ? one : myhalf * ( cons_arr(i,j-1,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
2177 
2178  amrex::Real ustar = myhalf * ( u_star_arr(i,j-1,0) + u_star_arr(i,j,0) );
2179  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2180 
2181  // NOTE: this is rho*<v'w'> = -K dvdz
2182  amrex::Real stressy = -rho * ustar * ustar * vely / wsp;
2183 
2184  return stressy;
2185  }
2186 
2187 private:
2188 #ifdef AMREX_USE_FLOAT
2189  const amrex::Real eps = amrex::Real(1e-8);
2190 #else
2191  const amrex::Real eps = amrex::Real(1e-15);
2192 #endif
2193  const bool fluxes_include_rho{false};
2194 };
2195 
2196 
2197 /**
2198  * Bulk coefficient flux formulation
2199  */
2201 {
2203  amrex::Real m_Ch,
2204  amrex::Real m_Cq)
2205  {
2206  mdata.Cd = m_Cd;
2207  mdata.Ch = m_Ch;
2208  mdata.Cq = m_Cq;
2209  }
2210 
2211  AMREX_GPU_DEVICE
2212  AMREX_FORCE_INLINE
2213  amrex::Real
2214  compute_q_flux (const int& i,
2215  const int& j,
2216  const int& k,
2217  const amrex::Array4<const amrex::Real>& cons_arr,
2218  const amrex::Array4<const amrex::Real>& velx_arr,
2219  const amrex::Array4<const amrex::Real>& vely_arr,
2220  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2221  const amrex::Array4<const amrex::Real>& /*qm_arr*/,
2222  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
2223  const amrex::Array4<const amrex::Real>& /*q_star_arr*/,
2224  const amrex::Array4<const amrex::Real>& q_surf_arr) const
2225  {
2226  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
2227  amrex::Real qv = cons_arr(i,j,k,RhoQ1_comp)/cons_arr(i,j,k,Rho_comp);
2228  amrex::Real qvsurf = q_surf_arr(i,j,0);
2229  amrex::Real velx = myhalf * ( velx_arr(i,j,k) + velx_arr(i+1,j ,k) );
2230  amrex::Real vely = myhalf * ( vely_arr(i,j,k) + vely_arr(i ,j+1,k) );
2231  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2232 
2233  // NOTE: this is rho*<q'w'> = -K dqdz
2234  amrex::Real moflux = rho * mdata.Cq * wsp * (qvsurf - qv);
2235 
2236  return moflux;
2237  }
2238 
2239  AMREX_GPU_DEVICE
2240  AMREX_FORCE_INLINE
2241  amrex::Real
2242  compute_t_flux (const int& i,
2243  const int& j,
2244  const int& k,
2245  const amrex::Array4<const amrex::Real>& cons_arr,
2246  const amrex::Array4<const amrex::Real>& velx_arr,
2247  const amrex::Array4<const amrex::Real>& vely_arr,
2248  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2249  const amrex::Array4<const amrex::Real>& /*tm_arr*/,
2250  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
2251  const amrex::Array4<const amrex::Real>& /*t_star_arr*/,
2252  const amrex::Array4<const amrex::Real>& t_surf_arr) const
2253  {
2254  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
2255  amrex::Real th = cons_arr(i,j,k,RhoTheta_comp)/cons_arr(i,j,k,Rho_comp);
2256  amrex::Real thsurf = t_surf_arr(i,j,0);
2257  amrex::Real velx = myhalf * ( velx_arr(i,j,k) + velx_arr(i+1,j ,k) );
2258  amrex::Real vely = myhalf * ( vely_arr(i,j,k) + vely_arr(i ,j+1,k) );
2259  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2260 
2261  // NOTE: this is rho*<T'w'> = -K dTdz
2262  amrex::Real moflux = rho * mdata.Ch * wsp * (thsurf - th);
2263 
2264  return moflux;
2265  }
2266 
2267  AMREX_GPU_DEVICE
2268  AMREX_FORCE_INLINE
2269  amrex::Real
2270  compute_u_flux (const int& i,
2271  const int& j,
2272  const int& k,
2273  const amrex::Array4<const amrex::Real>& cons_arr,
2274  const amrex::Array4<const amrex::Real>& velx_arr,
2275  const amrex::Array4<const amrex::Real>& vely_arr,
2276  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2277  const amrex::Array4<const amrex::Real>& /*um_arr*/,
2278  const amrex::Array4<const amrex::Real>& /*u_star_arr*/) const
2279  {
2280  amrex::Real rho = myhalf * ( cons_arr(i-1,j,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
2281  amrex::Real velx = velx_arr(i,j,k);
2282  amrex::Real vely = fourth * ( vely_arr(i ,j,k) + vely_arr(i ,j+1,k)
2283  + vely_arr(i-1,j,k) + vely_arr(i-1,j+1,k) );
2284  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2285 
2286  // NOTE: this is rho*<u'w'> = -K dudz
2287  // NOTE: tau_tot = rho * Cd * wsp^2, multiply by u/wsp to get tau_13
2288  amrex::Real stressx = -rho * mdata.Cd * wsp * velx;
2289 
2290  return stressx;
2291  }
2292 
2293  AMREX_GPU_DEVICE
2294  AMREX_FORCE_INLINE
2295  amrex::Real
2296  compute_v_flux (const int& i,
2297  const int& j,
2298  const int& k,
2299  const amrex::Array4<const amrex::Real>& cons_arr,
2300  const amrex::Array4<const amrex::Real>& velx_arr,
2301  const amrex::Array4<const amrex::Real>& vely_arr,
2302  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2303  const amrex::Array4<const amrex::Real>& /*vm_arr*/,
2304  const amrex::Array4<const amrex::Real>& /*u_star_arr*/) const
2305  {
2306  amrex::Real rho = myhalf * ( cons_arr(i,j-1,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
2307  amrex::Real velx = fourth * ( velx_arr(i,j ,k) + velx_arr(i+1,j ,k)
2308  + velx_arr(i,j-1,k) + velx_arr(i+1,j-1,k) );
2309  amrex::Real vely = vely_arr(i,j,k);
2310  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2311 
2312  // NOTE: this is rho*<v'w'> = -K dvdz
2313  // NOTE: tau_tot = rho * Cd * wsp^2, multiply by v/wsp to get tau_13
2314  amrex::Real stressy = -rho * mdata.Cd * wsp * vely;
2315 
2316  return stressy;
2317  }
2318 
2319 private:
2321 };
2322 
2323 
2324 /**
2325  * RICO flux formulation
2326  */
2328 {
2329  rico_flux (amrex::Real l_theta_z0, amrex::Real l_qsat_z0)
2330  : theta_z0{l_theta_z0}, qsat_z0{l_qsat_z0} {}
2331 
2332 
2333  AMREX_GPU_DEVICE
2334  AMREX_FORCE_INLINE
2335  amrex::Real
2336  compute_q_flux (const int& i,
2337  const int& j,
2338  const int& k,
2339  const amrex::Array4<const amrex::Real>& cons_arr,
2340  const amrex::Array4<const amrex::Real>& velx_arr,
2341  const amrex::Array4<const amrex::Real>& vely_arr,
2342  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2343  const amrex::Array4<const amrex::Real>& /*qm_arr*/,
2344  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
2345  const amrex::Array4<const amrex::Real>& q_star_arr,
2346  const amrex::Array4<const amrex::Real>& /*t_surf_arr*/) const
2347  {
2348  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
2349  amrex::Real qstar = q_star_arr(i,j,0);
2350  amrex::Real q = cons_arr(i, j, k, RhoQ1_comp) / rho;
2351 
2352  amrex::Real velx = myhalf * (velx_arr(i,j,k) + velx_arr(i+1,j,k));
2353  amrex::Real vely = myhalf * (vely_arr(i,j,k) + vely_arr(i,j+1,k));
2354  amrex::Real wsp = std::sqrt(velx*velx + vely*vely);
2355 
2356  // NOTE: this is rho*<q'w'> = -K dqdz
2357  amrex::Real moflux = (std::abs(qstar) > eps) ? - rho * qstar * wsp * (q - qsat_z0) : zero;
2358 
2359  return moflux;
2360  }
2361 
2362  AMREX_GPU_DEVICE
2363  AMREX_FORCE_INLINE
2364  amrex::Real
2365  compute_t_flux (const int& i,
2366  const int& j,
2367  const int& k,
2368  const amrex::Array4<const amrex::Real>& cons_arr,
2369  const amrex::Array4<const amrex::Real>& velx_arr,
2370  const amrex::Array4<const amrex::Real>& vely_arr,
2371  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2372  const amrex::Array4<const amrex::Real>& /*tm_arr*/,
2373  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
2374  const amrex::Array4<const amrex::Real>& t_star_arr,
2375  const amrex::Array4<const amrex::Real>& /*t_surf_arr*/) const
2376  {
2377  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
2378  amrex::Real tstar = t_star_arr(i, j, 0);
2379  amrex::Real theta = cons_arr(i, j, k, RhoTheta_comp) / rho;
2380 
2381  amrex::Real velx = myhalf * (velx_arr(i,j,k) + velx_arr(i+1,j,k));
2382  amrex::Real vely = myhalf * (vely_arr(i,j,k) + vely_arr(i,j+1,k));
2383  amrex::Real wsp = std::sqrt(velx*velx + vely*vely);
2384 
2385  // NOTE: this is rho*<T'w'> = -K dTdz
2386  amrex::Real moflux = (std::abs(tstar) > eps) ? - rho * tstar * wsp * (theta - theta_z0) : zero;
2387 
2388  return moflux;
2389  }
2390 
2391  AMREX_GPU_DEVICE
2392  AMREX_FORCE_INLINE
2393  amrex::Real
2394  compute_u_flux (const int& i,
2395  const int& j,
2396  const int& k,
2397  const amrex::Array4<const amrex::Real>& cons_arr,
2398  const amrex::Array4<const amrex::Real>& velx_arr,
2399  const amrex::Array4<const amrex::Real>& vely_arr,
2400  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2401  const amrex::Array4<const amrex::Real>& /*um_arr*/,
2402  const amrex::Array4<const amrex::Real>& u_star_arr) const
2403  {
2404  amrex::Real velx = velx_arr(i,j,k);
2405  amrex::Real vely = fourth * ( vely_arr(i ,j,k) + vely_arr(i ,j+1,k)
2406  + vely_arr(i-1,j,k) + vely_arr(i-1,j+1,k) );
2407  amrex::Real rho = myhalf * ( cons_arr(i-1,j,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
2408 
2409  amrex::Real ustar = myhalf * ( u_star_arr(i-1,j,0) + u_star_arr(i,j,0) );
2410  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2411 
2412  // NOTE: this is rho*<u'w'> = -K dudz
2413  amrex::Real stressx = -rho * ustar * wsp * velx;
2414 
2415  return stressx;
2416  }
2417 
2418  AMREX_GPU_DEVICE
2419  AMREX_FORCE_INLINE
2420  amrex::Real
2421  compute_v_flux (const int& i,
2422  const int& j,
2423  const int& k,
2424  const amrex::Array4<const amrex::Real>& cons_arr,
2425  const amrex::Array4<const amrex::Real>& velx_arr,
2426  const amrex::Array4<const amrex::Real>& vely_arr,
2427  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2428  const amrex::Array4<const amrex::Real>& /*vm_arr*/,
2429  const amrex::Array4<const amrex::Real>& u_star_arr) const
2430  {
2431  amrex::Real velx = fourth * ( velx_arr(i,j ,k) + velx_arr(i+1,j ,k)
2432  + velx_arr(i,j-1,k) + velx_arr(i+1,j-1,k) );
2433  amrex::Real vely = vely_arr(i,j,k);
2434  amrex::Real rho = myhalf * ( cons_arr(i,j-1,k,Rho_comp) + cons_arr(i,j,k,Rho_comp) );
2435 
2436  amrex::Real ustar = myhalf * ( u_star_arr(i,j-1,0) + u_star_arr(i,j,0) );
2437  amrex::Real wsp = std::sqrt(velx*velx+vely*vely);
2438 
2439  // NOTE: this is rho*<v'w'> = -K dvdz
2440  amrex::Real stressy = -rho * ustar * wsp * vely;
2441 
2442  return stressy;
2443  }
2444 
2445 private:
2446 #ifdef AMREX_USE_FLOAT
2447  const amrex::Real eps = amrex::Real(1e-8);
2448 #else
2449  const amrex::Real eps = amrex::Real(1e-15);
2450 #endif
2453 };
2454 
2455 
2456 /**
2457  * Rotate flux formulation
2458  */
2460 {
2462 
2463  AMREX_GPU_DEVICE
2464  AMREX_FORCE_INLINE
2465  amrex::Real
2466  compute_q_flux (const int& i,
2467  const int& j,
2468  const int& k,
2469  const amrex::Array4<const amrex::Real>& cons_arr,
2470  const amrex::Array4<const amrex::Real>& /*velx_arr*/,
2471  const amrex::Array4<const amrex::Real>& /*vely_arr*/,
2472  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2473  const amrex::Array4<const amrex::Real>& qvm_arr,
2474  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
2475  const amrex::Array4<const amrex::Real>& q_star_arr,
2476  const amrex::Array4<const amrex::Real>& q_surf_arr) const
2477  {
2478  // NOTE: this is the total stress
2479  amrex::Real qv_mean = qvm_arr(i,j,0);
2480  amrex::Real qv_surf = q_surf_arr(i,j,0);
2481  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
2482  amrex::Real qstar = q_star_arr(i,j,0);
2483  amrex::Real moflux = (std::abs(qstar) > eps) ? -rho * qstar * (qv_mean-qv_surf): zero;
2484 
2485  return moflux;
2486  }
2487 
2488  AMREX_GPU_DEVICE
2489  AMREX_FORCE_INLINE
2490  amrex::Real
2491  compute_t_flux (const int& i,
2492  const int& j,
2493  const int& k,
2494  const amrex::Array4<const amrex::Real>& cons_arr,
2495  const amrex::Array4<const amrex::Real>& /*velx_arr*/,
2496  const amrex::Array4<const amrex::Real>& /*vely_arr*/,
2497  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2498  const amrex::Array4<const amrex::Real>& tm_arr,
2499  const amrex::Array4<const amrex::Real>& /*u_star_arr*/,
2500  const amrex::Array4<const amrex::Real>& t_star_arr,
2501  const amrex::Array4<const amrex::Real>& t_surf_arr) const
2502  {
2503  // NOTE: this is the total stress
2504  amrex::Real theta_mean = tm_arr(i,j,0);
2505  amrex::Real theta_surf = t_surf_arr(i,j,0);
2506  amrex::Real rho = cons_arr(i,j,k,Rho_comp);
2507  amrex::Real tstar = t_star_arr(i,j,0);
2508  amrex::Real moflux = (std::abs(tstar) > eps) ? -rho * tstar * (theta_mean-theta_surf) : zero;
2509 
2510  return moflux;
2511  }
2512 
2513  AMREX_GPU_DEVICE
2514  AMREX_FORCE_INLINE
2515  amrex::Real
2516  compute_u_flux (const int& i,
2517  const int& j,
2518  const int& k,
2519  const amrex::Array4<const amrex::Real>& cons_arr,
2520  const amrex::Array4<const amrex::Real>& /*velx_arr*/,
2521  const amrex::Array4<const amrex::Real>& /*vely_arr*/,
2522  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2523  const amrex::Array4<const amrex::Real>& /*um_arr*/,
2524  const amrex::Array4<const amrex::Real>& u_star_arr) const
2525  {
2526  // NOTE: this is the total stress
2527  amrex::Real rho = myhalf *( cons_arr(i-1,j,k,Rho_comp) + cons_arr(i ,j,k,Rho_comp) );
2528  amrex::Real ustar = myhalf * ( u_star_arr(i-1,j,0) + u_star_arr(i,j,0) );
2529  amrex::Real stressx = -rho * ustar * ustar;
2530 
2531  return stressx;
2532  }
2533 
2534  AMREX_GPU_DEVICE
2535  AMREX_FORCE_INLINE
2536  amrex::Real
2537  compute_v_flux (const int& /*i*/,
2538  const int& /*j*/,
2539  const int& /*k*/,
2540  const amrex::Array4<const amrex::Real>& /*cons_arr*/,
2541  const amrex::Array4<const amrex::Real>& /*velx_arr*/,
2542  const amrex::Array4<const amrex::Real>& /*vely_arr*/,
2543  const amrex::Array4<const amrex::Real>& /*umm_arr*/,
2544  const amrex::Array4<const amrex::Real>& /*vm_arr*/,
2545  const amrex::Array4<const amrex::Real>& /*u_star_arr*/) const
2546  {
2547  // NOTE: this is the total stress
2548  amrex::Real stressy = zero;
2549 
2550  return stressy;
2551  }
2552 
2553 private:
2554 #ifdef AMREX_USE_FLOAT
2555  const amrex::Real eps = amrex::Real(1e-8);
2556 #else
2557  const amrex::Real eps = amrex::Real(1e-15);
2558 #endif
2559 };
2560 #endif
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:31
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real PI
Definition: ERF_Constants.H:16
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:32
constexpr amrex::Real PIoTwo
Definition: ERF_Constants.H:17
@ num
Definition: ERF_DataStruct.H:24
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
Real z0
Definition: ERF_InitCustomPertVels_ScalarAdvDiff.H:8
rho
Definition: ERF_InitCustomPert_Bubble.H:106
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Donelan_roughness(amrex::Real ustar)
Definition: ERF_MOSTRoughness.H:29
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real COARE3_roughness(amrex::Real zref, amrex::Real umm, amrex::Real ustar)
Definition: ERF_MOSTRoughness.H:9
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real air_viscosity(amrex::Real T_degK)
Definition: ERF_MOSTStress.H:129
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real calc_wstar(const amrex::Real &ust, const amrex::Real &tst, const amrex::Real &qst, const amrex::Real &pblh, const amrex::Real &th, const amrex::Real &thv, const amrex::Real &qv=amrex::Real(0))
Definition: ERF_Wstar.H:13
@ Mom_v
Definition: ERF_IndexDefines.H:193
@ theta
Definition: ERF_MM5.H:20
@ qv
Definition: ERF_Kessler.H:28
Definition: ERF_MOSTStress.H:189
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:270
most_data mdata
Definition: ERF_MOSTStress.H:267
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:204
const amrex::Real tol_z
Definition: ERF_MOSTStress.H:269
similarity_funs sfuns
Definition: ERF_MOSTStress.H:268
adiabatic_charnock(amrex::Real Tflux, amrex::Real Qvflux, amrex::Real cnk_a, bool cnk_visc)
Definition: ERF_MOSTStress.H:190
Definition: ERF_MOSTStress.H:359
similarity_funs sfuns
Definition: ERF_MOSTStress.H:419
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:421
most_data mdata
Definition: ERF_MOSTStress.H:418
adiabatic_donelan(amrex::Real Tflux, amrex::Real Qvflux)
Definition: ERF_MOSTStress.H:360
const amrex::Real tol
Definition: ERF_MOSTStress.H:420
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:370
Definition: ERF_MOSTStress.H:278
similarity_funs sfuns
Definition: ERF_MOSTStress.H:349
const amrex::Real tol_z
Definition: ERF_MOSTStress.H:350
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:351
adiabatic_mod_charnock(amrex::Real Tflux, amrex::Real Qvflux, amrex::Real depth)
Definition: ERF_MOSTStress.H:279
most_data mdata
Definition: ERF_MOSTStress.H:348
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:292
Definition: ERF_MOSTStress.H:429
const amrex::Real eps
Definition: ERF_MOSTStress.H:497
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &Hwave_arr, const amrex::Array4< amrex::Real > &Lwave_arr, const amrex::Array4< amrex::Real > &eta_arr) const
Definition: ERF_MOSTStress.H:440
const amrex::Real z0_eps
Definition: ERF_MOSTStress.H:499
most_data mdata
Definition: ERF_MOSTStress.H:491
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:501
similarity_funs sfuns
Definition: ERF_MOSTStress.H:492
const amrex::Real tol
Definition: ERF_MOSTStress.H:493
const amrex::Real z0_max
Definition: ERF_MOSTStress.H:500
adiabatic_wave_coupled(amrex::Real Tflux, amrex::Real Qvflux)
Definition: ERF_MOSTStress.H:430
Definition: ERF_MOSTStress.H:140
adiabatic(amrex::Real Tflux, amrex::Real Qvflux)
Definition: ERF_MOSTStress.H:141
similarity_funs sfuns
Definition: ERF_MOSTStress.H:181
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< const amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:151
most_data mdata
Definition: ERF_MOSTStress.H:180
Definition: ERF_MOSTStress.H:2201
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_u_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2270
bulk_coeff_flux(amrex::Real m_Cd, amrex::Real m_Ch, amrex::Real m_Cq)
Definition: ERF_MOSTStress.H:2202
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_q_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &q_surf_arr) const
Definition: ERF_MOSTStress.H:2214
most_data mdata
Definition: ERF_MOSTStress.H:2320
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_t_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &t_surf_arr) const
Definition: ERF_MOSTStress.H:2242
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_v_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2296
Definition: ERF_MOSTStress.H:2080
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_u_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &u_star_arr) const
Definition: ERF_MOSTStress.H:2136
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_q_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &q_star_arr, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2088
custom_flux(bool specified_rho_surf)
Definition: ERF_MOSTStress.H:2081
const bool fluxes_include_rho
Definition: ERF_MOSTStress.H:2193
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_t_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &t_star_arr, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2112
const amrex::Real eps
Definition: ERF_MOSTStress.H:2191
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_v_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &u_star_arr) const
Definition: ERF_MOSTStress.H:2163
Definition: ERF_MOSTStress.H:1947
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_u_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2002
donelan_flux()
Definition: ERF_MOSTStress.H:1948
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_v_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2038
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_t_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &t_surf_arr) const
Definition: ERF_MOSTStress.H:1974
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_q_flux(const int &, const int &, const int &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:1953
Definition: ERF_MOSTStress.H:1780
moeng_flux()
Definition: ERF_MOSTStress.H:1781
const amrex::Real eps
Definition: ERF_MOSTStress.H:1937
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_u_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &um_arr, const amrex::Array4< const amrex::Real > &u_star_arr) const
Definition: ERF_MOSTStress.H:1862
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:1939
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_t_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &u_star_arr, const amrex::Array4< const amrex::Real > &t_star_arr, const amrex::Array4< const amrex::Real > &t_surf_arr) const
Definition: ERF_MOSTStress.H:1824
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_q_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< const amrex::Real > &u_star_arr, const amrex::Array4< const amrex::Real > &q_star_arr, const amrex::Array4< const amrex::Real > &q_surf_arr) const
Definition: ERF_MOSTStress.H:1786
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_v_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &vm_arr, const amrex::Array4< const amrex::Real > &u_star_arr) const
Definition: ERF_MOSTStress.H:1899
Definition: ERF_MOSTStress.H:13
amrex::Real surf_moist_flux
Moisture flux.
Definition: ERF_MOSTStress.H:19
amrex::Real Cnk_b2
Modified Charnock Eq (4) https://doi.org/amrex::Real(10.1175)/JAMC-D-17-amrex::Real(0137....
Definition: ERF_MOSTStress.H:23
amrex::Real Cnk_b
Definition: ERF_MOSTStress.H:25
amrex::Real Ch
Definition: ERF_MOSTStress.H:29
amrex::Real Cnk_d
Modified Charnock Eq (4) https://doi.org/amrex::Real(10.1175)/JAMC-D-17-amrex::Real(0137....
Definition: ERF_MOSTStress.H:24
amrex::Real kappa
von Karman constant
Definition: ERF_MOSTStress.H:16
amrex::Real gravity
Acceleration due to gravity (m/s^2)
Definition: ERF_MOSTStress.H:17
amrex::Real Cnk_a
Standard Charnock constant https://doi.org/amrex::Real(10.1175)/JAMC-D-17-amrex::Real(0137....
Definition: ERF_MOSTStress.H:21
amrex::Real Cd
Definition: ERF_MOSTStress.H:28
amrex::Real Cq
Definition: ERF_MOSTStress.H:30
const amrex::Real Bjr_beta
Definition: ERF_MOSTStress.H:32
amrex::Real Cnk_b1
Modified Charnock Eq (4) https://doi.org/amrex::Real(10.1175)/JAMC-D-17-amrex::Real(0137....
Definition: ERF_MOSTStress.H:22
amrex::Real z0_const
Roughness height – default constant value(m)
Definition: ERF_MOSTStress.H:15
bool visc
Use viscous Charnock formulation.
Definition: ERF_MOSTStress.H:26
amrex::Real surf_temp_flux
Heat flux TODO: decide whether this is <θ'w'> or <θv'w'> under moist conditions.
Definition: ERF_MOSTStress.H:18
Definition: ERF_MOSTStress.H:2328
rico_flux(amrex::Real l_theta_z0, amrex::Real l_qsat_z0)
Definition: ERF_MOSTStress.H:2329
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_q_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &q_star_arr, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2336
amrex::Real qsat_z0
Definition: ERF_MOSTStress.H:2452
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_t_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &t_star_arr, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2365
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_u_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &u_star_arr) const
Definition: ERF_MOSTStress.H:2394
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_v_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &velx_arr, const amrex::Array4< const amrex::Real > &vely_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &u_star_arr) const
Definition: ERF_MOSTStress.H:2421
const amrex::Real eps
Definition: ERF_MOSTStress.H:2449
amrex::Real theta_z0
Definition: ERF_MOSTStress.H:2451
Definition: ERF_MOSTStress.H:2460
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_v_flux(const int &, const int &, const int &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &) const
Definition: ERF_MOSTStress.H:2537
rotate_flux()
Definition: ERF_MOSTStress.H:2461
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_u_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &u_star_arr) const
Definition: ERF_MOSTStress.H:2516
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_t_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &t_star_arr, const amrex::Array4< const amrex::Real > &t_surf_arr) const
Definition: ERF_MOSTStress.H:2491
const amrex::Real eps
Definition: ERF_MOSTStress.H:2557
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real compute_q_flux(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &cons_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &q_star_arr, const amrex::Array4< const amrex::Real > &q_surf_arr) const
Definition: ERF_MOSTStress.H:2466
Definition: ERF_MOSTStress.H:40
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_m(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:90
amrex::Real beta_m
Constants from Dyer, BLM, 1974.
Definition: ERF_MOSTStress.H:115
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_h2(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:67
amrex::Real beta_h
https://doi.org/amrex::Real(10.1007)/BF00240838
Definition: ERF_MOSTStress.H:116
amrex::Real gamma_h
Definition: ERF_MOSTStress.H:118
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_h(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:104
amrex::Real gamma_m
Definition: ERF_MOSTStress.H:117
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_m2(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:47
Definition: ERF_MOSTStress.H:614
const amrex::Real tol
Definition: ERF_MOSTStress.H:724
surface_flux_charnock(amrex::Real Tflux, amrex::Real Qvflux, amrex::Real cnk_a, bool cnk_visc, bool cons_qflux)
Definition: ERF_MOSTStress.H:615
similarity_funs sfuns
Definition: ERF_MOSTStress.H:723
most_data mdata
Definition: ERF_MOSTStress.H:721
bool spec_qflux
Definition: ERF_MOSTStress.H:722
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:631
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:725
Definition: ERF_MOSTStress.H:844
bool spec_qflux
Definition: ERF_MOSTStress.H:941
const amrex::Real tol
Definition: ERF_MOSTStress.H:943
surface_flux_donelan(amrex::Real Tflux, amrex::Real Qvflux, bool cons_qflux)
Definition: ERF_MOSTStress.H:845
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:944
most_data mdata
Definition: ERF_MOSTStress.H:940
similarity_funs sfuns
Definition: ERF_MOSTStress.H:942
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:857
Definition: ERF_MOSTStress.H:733
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:749
surface_flux_mod_charnock(amrex::Real Tflux, amrex::Real Qvflux, amrex::Real depth, bool cons_qflux)
Definition: ERF_MOSTStress.H:734
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:836
most_data mdata
Definition: ERF_MOSTStress.H:832
bool spec_qflux
Definition: ERF_MOSTStress.H:833
const amrex::Real tol
Definition: ERF_MOSTStress.H:835
similarity_funs sfuns
Definition: ERF_MOSTStress.H:834
Definition: ERF_MOSTStress.H:952
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:1065
surface_flux_wave_coupled(amrex::Real Tflux, amrex::Real Qvflux, bool cons_qflux)
Definition: ERF_MOSTStress.H:953
const amrex::Real z0_max
Definition: ERF_MOSTStress.H:1064
const amrex::Real z0_eps
Definition: ERF_MOSTStress.H:1063
most_data mdata
Definition: ERF_MOSTStress.H:1054
similarity_funs sfuns
Definition: ERF_MOSTStress.H:1056
const amrex::Real tol
Definition: ERF_MOSTStress.H:1057
bool spec_qflux
Definition: ERF_MOSTStress.H:1055
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &Hwave_arr, const amrex::Array4< amrex::Real > &Lwave_arr, const amrex::Array4< amrex::Real > &eta_arr) const
Definition: ERF_MOSTStress.H:965
const amrex::Real eps
Definition: ERF_MOSTStress.H:1061
Definition: ERF_MOSTStress.H:509
similarity_funs sfuns
Definition: ERF_MOSTStress.H:604
surface_flux(amrex::Real Tflux, amrex::Real Qvflux, bool cons_qflux)
Definition: ERF_MOSTStress.H:510
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< const amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:522
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:606
most_data mdata
Definition: ERF_MOSTStress.H:602
const amrex::Real tol
Definition: ERF_MOSTStress.H:605
bool spec_qflux
Definition: ERF_MOSTStress.H:603
Definition: ERF_MOSTStress.H:1207
most_data mdata
Definition: ERF_MOSTStress.H:1359
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:1224
bool spec_qflux
Definition: ERF_MOSTStress.H:1360
surface_temp_charnock(amrex::Real Tflux, amrex::Real Qvflux, amrex::Real cnk_a, bool cnk_visc, bool cons_qflux)
Definition: ERF_MOSTStress.H:1208
const amrex::Real tol_z
Definition: ERF_MOSTStress.H:1363
const amrex::Real alpha
Definition: ERF_MOSTStress.H:1364
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:1365
const amrex::Real tol
Definition: ERF_MOSTStress.H:1362
similarity_funs sfuns
Definition: ERF_MOSTStress.H:1361
Definition: ERF_MOSTStress.H:1531
similarity_funs sfuns
Definition: ERF_MOSTStress.H:1639
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:1544
most_data mdata
Definition: ERF_MOSTStress.H:1637
const amrex::Real tol
Definition: ERF_MOSTStress.H:1640
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:1641
bool spec_qflux
Definition: ERF_MOSTStress.H:1638
surface_temp_donelan(amrex::Real Tflux, amrex::Real Qvflux, bool cons_qflux)
Definition: ERF_MOSTStress.H:1532
Definition: ERF_MOSTStress.H:1373
surface_temp_mod_charnock(amrex::Real Tflux, amrex::Real Qvflux, amrex::Real depth, bool cons_qflux)
Definition: ERF_MOSTStress.H:1374
similarity_funs sfuns
Definition: ERF_MOSTStress.H:1519
const amrex::Real tol
Definition: ERF_MOSTStress.H:1520
const amrex::Real tol_z
Definition: ERF_MOSTStress.H:1521
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:1523
most_data mdata
Definition: ERF_MOSTStress.H:1517
bool spec_qflux
Definition: ERF_MOSTStress.H:1518
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:1389
const amrex::Real alpha
Definition: ERF_MOSTStress.H:1522
Definition: ERF_MOSTStress.H:1649
const amrex::Real eps
Definition: ERF_MOSTStress.H:1768
surface_temp_wave_coupled(amrex::Real Tflux, amrex::Real Qvflux, bool cons_qflux)
Definition: ERF_MOSTStress.H:1650
const amrex::Real tol
Definition: ERF_MOSTStress.H:1764
most_data mdata
Definition: ERF_MOSTStress.H:1761
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &Hwave_arr, const amrex::Array4< amrex::Real > &Lwave_arr, const amrex::Array4< amrex::Real > &eta_arr) const
Definition: ERF_MOSTStress.H:1662
const amrex::Real z0_eps
Definition: ERF_MOSTStress.H:1770
similarity_funs sfuns
Definition: ERF_MOSTStress.H:1763
bool spec_qflux
Definition: ERF_MOSTStress.H:1762
const amrex::Real z0_max
Definition: ERF_MOSTStress.H:1771
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:1772
Definition: ERF_MOSTStress.H:1073
surface_temp(amrex::Real Tflux, amrex::Real Qvflux, bool cons_qflux)
Definition: ERF_MOSTStress.H:1074
similarity_funs sfuns
Definition: ERF_MOSTStress.H:1196
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void iterate_flux(const int &i, const int &j, const int &k, const int &max_iters, const amrex::Array4< const amrex::Real > &zref_arr, const amrex::Array4< const amrex::Real > &z0_arr, const amrex::Array4< const amrex::Real > &umm_arr, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::Array4< const amrex::Real > &tvm_arr, const amrex::Array4< const amrex::Real > &qvm_arr, const amrex::Array4< amrex::Real > &u_star_arr, const amrex::Array4< amrex::Real > &w_star_arr, const amrex::Array4< amrex::Real > &t_star_arr, const amrex::Array4< amrex::Real > &q_star_arr, const amrex::Array4< amrex::Real > &t_surf_arr, const amrex::Array4< amrex::Real > &q_surf_arr, const amrex::Array4< amrex::Real > &olen_arr, const amrex::Array4< amrex::Real > &pblh_arr, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &, const amrex::Array4< amrex::Real > &) const
Definition: ERF_MOSTStress.H:1086
const amrex::Real tol
Definition: ERF_MOSTStress.H:1197
bool spec_qflux
Definition: ERF_MOSTStress.H:1195
const amrex::Real alpha
Definition: ERF_MOSTStress.H:1198
const amrex::Real WSMIN
Definition: ERF_MOSTStress.H:1199
most_data mdata
Definition: ERF_MOSTStress.H:1194