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