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