1 #ifndef ERF_METGRIDUTIL_H_
2 #define ERF_METGRIDUTIL_H_
13 const amrex::Box& domain,
14 const std::string& fname,
15 std::string& NC_dateTime,
26 amrex::FArrayBox& NC_xvel_fab,
27 amrex::FArrayBox& NC_yvel_fab,
28 amrex::FArrayBox& NC_temp_fab,
29 amrex::FArrayBox& NC_rhum_fab,
30 amrex::FArrayBox& NC_pres_fab,
31 amrex::FArrayBox& NC_ght_fab,
32 amrex::FArrayBox& NC_hgt_fab,
33 amrex::FArrayBox& NC_psfc_fab,
34 amrex::FArrayBox& NC_msfu_fab,
35 amrex::FArrayBox& NC_msfv_fab,
36 amrex::FArrayBox& NC_msfm_fab,
37 amrex::FArrayBox& NC_sst_fab,
38 amrex::FArrayBox& NC_tsk_fab,
39 amrex::FArrayBox& NC_LAT_fab,
40 amrex::FArrayBox& NC_LON_fab,
41 amrex::IArrayBox& NC_lmask_iab,
42 amrex::Geometry& geom);
46 amrex::FArrayBox& NC_hgt_fab);
52 const bool interp_theta,
53 const bool metgrid_debug_quiescent,
54 const bool metgrid_debug_isothermal,
55 const bool metgrid_debug_dry,
56 const bool metgrid_basic_linear,
57 const bool metgrid_use_below_sfc,
58 const bool metgrid_use_sfc,
59 const bool metgrid_retain_sfc,
61 const int metgrid_order,
62 const int metgrid_metgrid_force_sfc_k,
67 amrex::FArrayBox& state_fab,
68 amrex::FArrayBox& x_vel_fab,
69 amrex::FArrayBox& y_vel_fab,
70 amrex::FArrayBox& z_vel_fab,
71 amrex::FArrayBox& z_phys_nd_fab,
72 const amrex::FArrayBox& NC_ght_fab,
73 const amrex::FArrayBox& NC_xvel_fab,
74 const amrex::FArrayBox& NC_yvel_fab,
75 const amrex::FArrayBox& NC_temp_fab,
76 const amrex::FArrayBox& NC_rhum_fab,
77 const amrex::FArrayBox& NC_pres_fab,
78 amrex::FArrayBox& tmp_src_fab,
79 amrex::FArrayBox& tmp_dst_fab,
80 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& fabs_for_bcs_xlo,
81 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& fabs_for_bcs_xhi,
82 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& fabs_for_bcs_ylo,
83 amrex::Vector<amrex::Vector<amrex::FArrayBox>>& fabs_for_bcs_yhi,
84 const amrex::Array4<const int>& mask_c_arr,
85 const amrex::Array4<const int>& mask_u_arr,
86 const amrex::Array4<const int>& mask_v_arr);
90 amrex::FArrayBox& msfu_fab,
91 amrex::FArrayBox& msfv_fab,
92 amrex::FArrayBox& msfm_fab,
94 amrex::FArrayBox& NC_MSFU_fab,
95 amrex::FArrayBox& NC_MSFV_fab,
96 amrex::FArrayBox& NC_MSFM_fab);
100 const bool metgrid_debug_psfc,
102 const amrex::Box& valid_bx,
103 const int& flag_psfc,
104 amrex::FArrayBox& state,
105 amrex::FArrayBox& r_hse_fab,
106 amrex::FArrayBox& p_hse_fab,
107 amrex::FArrayBox& pi_hse_fab,
108 amrex::FArrayBox& th_hse_fab,
109 amrex::FArrayBox& qv_hse_fab,
110 amrex::FArrayBox& z_phys_cc_fab,
111 const amrex::FArrayBox& NC_psfc_fab);
128 for (
int i=0; i <= order; i++) {
131 for (
int k=0; k <= order; k++) {
132 if (k == i)
continue;
147 const bool& exp_interp,
165 amrex::ignore_unused(i,j);
166 #ifndef AMREX_USE_GPU
170 for (
int new_k=0; new_k < new_n; new_k++) {
171 #ifndef AMREX_USE_GPU
172 if (debug) amrex::Print() <<
"new_k=" << new_k;
176 bool extrapolating =
true;
178 for (
int ko=0; ko < orig_n-1; ko++) {
184 extrapolating =
false;
190 if (var_type ==
'T') {
193 amrex::Real depth_of_extrap_in_p = new_x_p[new_k]-orig_x_p[0];
199 new_y[new_k] = (temp_extrap_starting_point+dT)*std::pow(
amrex::Real(100000.0)/new_x_p[new_k],
R_d/
Cp_d);
203 new_y[new_k] = orig_y[0];
209 if ((kl-((order+1)/2-1) >= 0) && (kr+((order+1)/2-1) <= orig_n-1)) {
211 int ksta = kl-(((order+1)/2)-1);
212 int kend = ksta+order;
213 #ifndef AMREX_USE_GPU
214 int ksize = kend-ksta;
215 if (debug) amrex::Print() <<
" (1a) order=" << order <<
" new_x_z=" << new_x_z[new_k] <<
" new_x_p=" << new_x_p[new_k] <<
" kl=" << kl <<
" kr=" << kr <<
" ksta=" << ksta <<
" kend=" << kend << std::endl;
218 amrex::GpuArray<amrex::Real,9> orig_x_sub;
221 new_x = new_x_p[new_k];
222 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
224 new_x = new_x_z[new_k];
225 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
227 amrex::GpuArray<amrex::Real,9> orig_y_sub;
229 for (
int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
230 #ifndef AMREX_USE_GPU
232 amrex::Print() <<
" orig_x_sub_p = [";
233 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_x_sub_p[k];
234 amrex::Print() <<
"]" << std::endl;
235 amrex::Print() <<
" orig_y_sub_p = [";
236 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_y_sub_p[k];
237 amrex::Print() <<
"]" << std::endl;
240 lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
246 #ifndef AMREX_USE_GPU
247 int ksize = kend-ksta+1;
248 if (debug) amrex::Print() <<
" (1b) order=" << order <<
" new_x=" << new_x_z[new_k] <<
" new_x_p=" << new_x_p[new_k] <<
" kl=" << kl <<
" kr=" << kr <<
" ksta=" << ksta <<
" kend=" << kend << std::endl;
251 amrex::GpuArray<amrex::Real,2> orig_x_sub;
252 amrex::GpuArray<amrex::Real,2> orig_y_sub;
256 new_x = new_x_p[new_k];
257 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
259 new_x = new_x_z[new_k];
260 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
262 for (
int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
263 #ifndef AMREX_USE_GPU
265 amrex::Print() <<
" orig_x_sub_p = [";
266 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_x_sub_p[k];
267 amrex::Print() <<
"]" << std::endl;
268 amrex::Print() <<
" orig_y_sub_p = [";
269 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_y_sub_p[k];
270 amrex::Print() <<
"]" << std::endl;
275 }
else if (order%2 == 0) {
276 if ((kl-(order/2) >= 0) && (kr+order/2 <= orig_n-1)) {
281 int ksta = kl-(order/2-1);
282 int kend = ksta+order;
283 #ifndef AMREX_USE_GPU
284 int ksize = kend-ksta;
285 if (debug) amrex::Print() <<
" (2a) order=" << order <<
" new_x_z=" << new_x_z[new_k] <<
" new_x_p=" << new_x_p[new_k] <<
" kl=" << kl <<
" kr=" << kr <<
" ksta=" << ksta <<
" kend=" << kend << std::endl;
288 amrex::GpuArray<amrex::Real,10> orig_x_sub;
289 amrex::GpuArray<amrex::Real,10> orig_y_sub;
293 new_x = new_x_p[new_k];
294 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
296 new_x = new_x_z[new_k];
297 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
299 for (
int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
300 #ifndef AMREX_USE_GPU
302 amrex::Print() <<
" orig_x_sub_p = [";
303 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_x_sub_p[k];
304 amrex::Print() <<
"]" << std::endl;
305 amrex::Print() <<
" orig_y_sub_p = [";
306 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_y_sub_p[k];
307 amrex::Print() <<
"]" << std::endl;
313 int ksta = kl-order/2;
314 int kend = ksta+order;
315 #ifndef AMREX_USE_GPU
316 int ksize = kend-ksta;
317 if (debug) amrex::Print() <<
"new_k=" << new_k <<
" (2b) order=" << order <<
" new_x_z=" << new_x_z[new_k] <<
" new_x_p=" << new_x_p[new_k] <<
" kl=" << kl <<
" kr=" << kr <<
" ksta=" << ksta <<
" kend=" << kend << std::endl;
320 amrex::GpuArray<amrex::Real,10> orig_x_sub;
321 amrex::GpuArray<amrex::Real,10> orig_y_sub;
325 new_x = new_x_p[new_k];
326 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
328 new_x = new_x_z[new_k];
329 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
331 for (
int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
332 #ifndef AMREX_USE_GPU
334 amrex::Print() <<
" orig_x_sub_p = [";
335 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_x_sub_p[k];
336 amrex::Print() <<
"]" << std::endl;
337 amrex::Print() <<
" orig_y_sub_p = [";
338 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_y_sub_p[k];
339 amrex::Print() <<
"]" << std::endl;
344 new_y[new_k] =
myhalf*(new_y_l+new_y_r);
345 }
else if ((kl-(order/2-1) >= 0) && (kr+order/2 <= orig_n-1)) {
346 int ksta = kl-(order/2-1);
347 int kend = ksta+order;
348 #ifndef AMREX_USE_GPU
349 int ksize = kend-ksta;
350 if (debug) amrex::Print() <<
" (3) order=" << order <<
" new_x_z=" << new_x_z[new_k] <<
" new_x_p=" << new_x_p[new_k] <<
" kl=" << kl <<
" kr=" << kr <<
" ksta=" << ksta <<
" kend=" << kend << std::endl;
353 amrex::GpuArray<amrex::Real,10> orig_x_sub;
354 amrex::GpuArray<amrex::Real,10> orig_y_sub;
358 new_x = new_x_p[new_k];
359 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
361 new_x = new_x_z[new_k];
362 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
364 for (
int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
365 #ifndef AMREX_USE_GPU
367 amrex::Print() <<
" orig_x_sub_p = [";
368 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_x_sub_p[k];
369 amrex::Print() <<
"]" << std::endl;
370 amrex::Print() <<
" orig_y_sub_p = [";
371 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_y_sub_p[k];
372 amrex::Print() <<
"]" << std::endl;
375 lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
376 }
else if ((kl-order/2 >= 0) && (kr+order/2-1 <= orig_n-1)) {
377 int ksta = kl-order/2;
378 int kend = ksta+order;
379 #ifndef AMREX_USE_GPU
380 int ksize = kend-ksta;
381 if (debug) amrex::Print() <<
" (4) order=" << order <<
" new_x_z=" << new_x_z[new_k] <<
" new_x_p=" << new_x_p[new_k] <<
" kl=" << kl <<
" kr=" << kr <<
" ksta=" << ksta <<
" kend=" << kend << std::endl;
384 amrex::GpuArray<amrex::Real,10> orig_x_sub;
385 amrex::GpuArray<amrex::Real,10> orig_y_sub;
389 new_x = new_x_p[new_k];
390 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
392 new_x = new_x_z[new_k];
393 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
395 for (
int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
396 #ifndef AMREX_USE_GPU
398 amrex::Print() <<
" orig_x_sub_p = [";
399 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_x_sub_p[k];
400 amrex::Print() <<
"]" << std::endl;
401 amrex::Print() <<
" orig_y_sub_p = [";
402 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_y_sub_p[k];
403 amrex::Print() <<
"]" << std::endl;
406 lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
411 #ifndef AMREX_USE_GPU
412 int ksize = kend-ksta+1;
413 if (debug) amrex::Print() <<
" (5) order=" << order <<
" new_x=" << new_x_z[new_k] <<
" new_x_p=" << new_x_p[new_k] <<
" kl=" << kl <<
" kr=" << kr <<
" ksta=" << ksta <<
" kend=" << kend << std::endl;
416 amrex::GpuArray<amrex::Real,2> orig_x_sub;
417 amrex::GpuArray<amrex::Real,2> orig_y_sub;
421 new_x = new_x_p[new_k];
422 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
424 new_x = new_x_z[new_k];
425 for (
int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
427 for (
int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
428 #ifndef AMREX_USE_GPU
430 amrex::Print() <<
" orig_x_sub_p = [";
431 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_x_sub_p[k];
432 amrex::Print() <<
"]" << std::endl;
433 amrex::Print() <<
" orig_y_sub_p = [";
434 for (
int k=0; k < ksize; k++) amrex::Print() <<
" " << orig_y_sub_p[k];
435 amrex::Print() <<
"]" << std::endl;
441 #ifndef AMREX_USE_GPU
442 if (debug) amrex::Print() <<
" new_y[" << new_k <<
"]=" << new_y[new_k] << std::endl;
460 const bool& metgrid_use_sfc,
461 const bool& exp_interp,
462 const bool& metgrid_retain_sfc,
464 const int& metgrid_order,
465 const int& metgrid_force_sfc_k,
473 const amrex::Array4<amrex::Real const>& orig_z_full,
474 const amrex::Array4<amrex::Real const>& orig_data,
475 const amrex::Array4<amrex::Real const>& new_z_full,
476 const amrex::Array4<amrex::Real>& new_data_full)
482 int imin_orig = amrex::lbound(amrex::Box(orig_data)).x;
483 int imax_orig = amrex::ubound(amrex::Box(orig_data)).x;
484 int jmin_orig = amrex::lbound(amrex::Box(orig_data)).y;
485 int jmax_orig = amrex::ubound(amrex::Box(orig_data)).y;
486 int kmax_orig = amrex::ubound(amrex::Box(orig_data)).z;
487 int kmax_new = kmax + 1;
489 AMREX_ASSERT(kmax_orig < 256);
490 AMREX_ASSERT(kmax_new < 256);
492 amrex::GpuArray<amrex::Real,256> new_z;
493 amrex::GpuArray<amrex::Real,256> new_p;
494 amrex::GpuArray<amrex::Real,256> new_data;
498 for (
int k=0; k < kmax_new; k++) {
500 new_z_p[k] =
fourth*(new_z_full(i,j,k)+new_z_full(i,j+1,k)+new_z_full(i,j,k+1)+new_z_full(i,j+1,k+1));
501 }
else if (stag ==
'Y') {
502 new_z_p[k] =
fourth*(new_z_full(i,j,k)+new_z_full(i+1,j,k)+new_z_full(i,j,k+1)+new_z_full(i+1,j,k+1));
503 }
else if (stag ==
'M') {
504 new_z_p[k] =
amrex::Real(0.125)*(new_z_full(i,j,k )+new_z_full(i,j+1,k )+new_z_full(i+1,j,k )+new_z_full(i+1,j+1,k )+
505 new_z_full(i,j,k+1)+new_z_full(i,j+1,k+1)+new_z_full(i+1,j,k+1)+new_z_full(i+1,j+1,k+1));
510 amrex::GpuArray<amrex::Real,256> orig_z;
512 for (
int k=0; k < kmax_orig; k++) {
514 orig_z_p[k] = orig_z_full(i,j,k);
515 }
else if (stag ==
'X') {
516 if (i <= imin_orig) {
517 orig_z_p[k] = orig_z_full(i,j,k);
518 }
else if (i >= imax_orig) {
519 orig_z_p[k] = orig_z_full(imax_orig-1,j,k);
521 orig_z_p[k] =
myhalf*(orig_z_full(i,j,k)+orig_z_full(i-1,j,k));
523 }
else if (stag ==
'Y') {
524 if (j <= jmin_orig) {
525 orig_z_p[k] = orig_z_full(i,j,k);
526 }
else if (j >= jmax_orig) {
527 orig_z_p[k] = orig_z_full(i,jmax_orig-1,k);
529 orig_z_p[k] =
myhalf*(orig_z_full(i,j,k)+orig_z_full(i,j-1,k));
535 bool flip_data_required =
false;
536 if (orig_z[1] > orig_z[kmax_orig-1]) flip_data_required =
true;
537 if (flip_data_required) amrex::Abort(
"metgrid initialization flip_data_required. Not yet implemented.");
543 for (
int k=1; k < kmax_orig; k++) {
544 if (orig_z_p[k] > orig_z_p[0]) {
551 amrex::GpuArray<amrex::Real,256> ordered_z;
552 amrex::GpuArray<amrex::Real,256> ordered_data;
555 if (k_above_sfc > 1) {
561 for (
int k=1; k < k_above_sfc; k++) {
562 ordered_z_p[count] = orig_z_p[k];
563 ordered_data_p[count] = orig_data(i,j,k,src_comp);
577 if (Pl-Pu < metgrid_proximity) {
582 ordered_z_p[count] = orig_z_p[0];
583 ordered_data_p[count] = orig_data(i,j,0,src_comp);
588 int knext = k_above_sfc;
594 if (metgrid_force_sfc_k > 0) {
595 for (
int k=k_above_sfc; k < kmax_orig; k++) {
596 if (orig_z_p[k] > new_z_p[metgrid_force_sfc_k-1]) {
607 if (Pl-Pu < metgrid_proximity) {
612 for (
int k=knext; k < kmax_orig; k++) {
613 ordered_z_p[count] = orig_z_p[k];
614 ordered_data_p[count] = orig_data(i,j,k,src_comp);
623 ordered_z_p[0] = orig_z[0];
624 ordered_data_p[0] = orig_data(i,j,0,src_comp);
630 if (metgrid_force_sfc_k > 0) {
631 for (
int k=knext; k < kmax_orig; k++) {
632 if (orig_z_p[k] > new_z_p[metgrid_force_sfc_k]) {
641 for (
int k=knext; k < kmax_orig; k++) {
645 if (Pl-Pu < metgrid_proximity) {
648 ordered_z_p[count] = orig_z_p[k];
649 ordered_data_p[count] = orig_data(i,j,k,src_comp);
655 int ksta(0), kend(0);
656 if (metgrid_use_below_sfc && metgrid_use_sfc) {
660 }
else if (metgrid_use_below_sfc && !metgrid_use_sfc) {
663 for (
int k=0; k < kmax_orig; k++) {
664 if (ordered_z_p[k] == orig_z_p[0]) {
669 for (
int k=ksfc; k < kmax_orig-1; k++) {
670 ordered_z_p[k] = ordered_z_p[k+1];
671 ordered_data_p[k] = ordered_data_p[k+1];
675 }
else if (!metgrid_use_below_sfc && metgrid_use_sfc) {
678 for (
int k=0; k < kend_order; k++) {
679 if (ordered_z_p[k] == orig_z_p[0]) {
688 amrex::Abort(
"metgrid initialization, !use_levels below_ground && !metgrid_use_sfc");
704 amrex::GpuArray<amrex::Real,256> ordered_p;
706 for (
int k=0; k < kend_order; k++) {
712 amrex::GpuArray<amrex::Real,256> final_z;
713 amrex::GpuArray<amrex::Real,256> final_p;
714 amrex::GpuArray<amrex::Real,256> final_data;
718 final_z_p[0] = ordered_z[ksta];
719 final_p_p[0] = ordered_p[ksta];
720 final_data_p[0] = ordered_data[ksta];
721 for (
int k=ksta+1; k <= kend; k++) {
722 if ((final_p_p[new_n]-ordered_p_p[k]) < metgrid_proximity) {
726 final_z_p[new_n] = ordered_z_p[k];
727 final_p_p[new_n] = ordered_p_p[k];
728 final_data_p[new_n] = ordered_data_p[k];
750 if (metgrid_retain_sfc) new_data[0] = ordered_data[0];
753 for (
int k=0; k < kmax_new; k++) {
754 new_data_full(i,j,k,dst_comp) = new_data[k];
767 const amrex::Array4<amrex::Real const>& orig_z,
768 const amrex::Array4<amrex::Real const>& orig_data,
769 const amrex::Array4<amrex::Real const>& new_z)
772 int imax_orig = amrex::ubound(amrex::Box(orig_data)).x;
773 int jmax_orig = amrex::ubound(amrex::Box(orig_data)).y;
774 int kmax_orig = amrex::ubound(amrex::Box(orig_data)).z;
778 z =
fourth*(new_z(i,j,k)+new_z(i,j+1,k)+new_z(i,j,k+1)+new_z(i,j+1,k+1));
780 else if (stag ==
'Y') {
781 z =
fourth*(new_z(i,j,k)+new_z(i+1,j,k)+new_z(i,j,k+1)+new_z(i+1,j,k+1));
783 else if (stag ==
'M') {
784 z =
amrex::Real(0.125)*(new_z(i,j,k )+new_z(i,j+1,k )+new_z(i+1,j,k )+new_z(i+1,j+1,k )+
785 new_z(i,j,k+1)+new_z(i,j+1,k+1)+new_z(i+1,j,k+1)+new_z(i+1,j+1,k+1));
793 for (
int kk = 0; kk < kmax_orig; kk++) {
796 orig_z_stag = orig_z(i,j,kk);
800 orig_z_stag = orig_z(i,j,kk);
802 else if (i == imax_orig) {
803 orig_z_stag = orig_z(imax_orig-1,j,kk);
806 orig_z_stag =
myhalf*(orig_z(i,j,kk)+orig_z(i-1,j,kk));
809 else if (stag ==
'Y') {
811 orig_z_stag = orig_z(i,j,kk);
813 else if (j == jmax_orig) {
814 orig_z_stag = orig_z(i,jmax_orig-1,kk);
817 orig_z_stag =
myhalf*(orig_z(i,j,kk)+orig_z(i,j-1,kk));
838 for (
int kk = 0; kk < kmax_orig; kk++) {
841 orig_z_stag = orig_z(i,j,kk);
843 else if (stag ==
'X') {
845 orig_z_stag = orig_z(i,j,kk);
847 else if (i == imax_orig) {
848 orig_z_stag = orig_z(imax_orig-1,j,kk);
851 orig_z_stag =
myhalf*(orig_z(i,j,kk)+orig_z(i-1,j,kk));
854 else if (stag ==
'Y') {
856 orig_z_stag = orig_z(i,j,kk);
858 else if (j == jmax_orig) {
859 orig_z_stag = orig_z(i,jmax_orig-1,kk);
862 orig_z_stag =
myhalf*(orig_z(i,j,kk)+orig_z(i,j-1,kk));
866 if ((
dz <
zero) && (
dz > dzhi1) && (kk != khi0)) {
874 return ( y0-(y1-y0)/(z1-
z0)*(
z0-
z) );
877 }
else if (khi0 == -1) {
881 z0 = orig_z(i,j,khi0);
883 else if (stag ==
'X') {
885 z0 = orig_z(i,j,khi0);
887 else if (i == imax_orig) {
888 z0 = orig_z(imax_orig-1,j,khi0);
891 z0 =
myhalf*(orig_z(i,j,khi0)+orig_z(i-1,j,khi0));
894 else if (stag ==
'Y') {
896 z0 = orig_z(i,j,khi0);
898 else if (j == jmax_orig) {
899 z0 = orig_z(i,jmax_orig-1,khi0);
902 z0 =
myhalf*(orig_z(i,j,khi0)+orig_z(i,j-1,khi0));
907 return ( y0+(y1-y0)/(z1-
z0)*(
z-
z0) );
912 return ( y0+(y1-y0)/(z1-
z0)*(
z-
z0) );
923 const amrex::Array4<amrex::Real const>& rhum,
924 const amrex::Array4<amrex::Real const>& temp,
925 const amrex::Array4<amrex::Real const>&
pres,
927 const amrex::Array4<amrex::Real>& mxrat)
941 if (temp(i,j,k) !=
zero) {
945 mxrat(i,j,k,src_indx) = std::pow(
amrex::Real(10.0), -6);
962 if (
pres(i,j,k) < qv_max_p_safe) {
963 if (mxrat(i,j,k,src_indx) > qv_max_flag) {
964 mxrat(i,j,k,src_indx) = qv_max_value;
967 if (
pres(i,j,k) < qv_min_p_safe) {
968 if (mxrat(i,j,k,src_indx) < qv_min_flag) {
969 mxrat(i,j,k,src_indx) = qv_min_value;
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:22
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real p_0
Definition: ERF_Constants.H:28
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:31
constexpr amrex::Real R_d
Definition: ERF_Constants.H:20
Real z0
Definition: ERF_InitCustomPertVels_ScalarAdvDiff.H:8
const bool use_moisture
Definition: ERF_InitCustomPert_Bomex.H:14
const Real dz
Definition: ERF_InitCustomPert_Bubble.H:25
Real * p
Definition: ERF_InitCustomPert_SquallLine.H:61
void init_terrain_from_metgrid(amrex::FArrayBox &z_phys_nd_fab, amrex::FArrayBox &NC_hgt_fab)
AMREX_FORCE_INLINE AMREX_GPU_DEVICE void lagrange_interp(const int &order, amrex::Real *x, amrex::Real *y, amrex::Real &new_x, amrex::Real &new_y)
Definition: ERF_MetgridUtils.H:116
AMREX_FORCE_INLINE AMREX_GPU_DEVICE amrex::Real interpolate_column_metgrid_linear(const int &i, const int &j, const int &k, char stag, int src_comp, const amrex::Array4< amrex::Real const > &orig_z, const amrex::Array4< amrex::Real const > &orig_data, const amrex::Array4< amrex::Real const > &new_z)
Definition: ERF_MetgridUtils.H:762
AMREX_FORCE_INLINE AMREX_GPU_DEVICE void interpolate_column_metgrid(const bool &metgrid_use_below_sfc, const bool &metgrid_use_sfc, const bool &exp_interp, const bool &metgrid_retain_sfc, const amrex::Real &metgrid_proximity, const int &metgrid_order, const int &metgrid_force_sfc_k, const int &i, const int &j, const int &kmax, const int &src_comp, const int &dst_comp, char var_type, char stag, const amrex::Array4< amrex::Real const > &orig_z_full, const amrex::Array4< amrex::Real const > &orig_data, const amrex::Array4< amrex::Real const > &new_z_full, const amrex::Array4< amrex::Real > &new_data_full)
Definition: ERF_MetgridUtils.H:459
void read_from_metgrid(int lev, int itime, const amrex::Box &domain, const std::string &fname, std::string &NC_dateTime, amrex::Real &NC_epochTime, int &flag_psfc, int &flag_msf, int &flag_sst, int &flag_tsk, int &flag_lmask, int &NC_nx, int &NC_ny, amrex::Real &NC_dx, amrex::Real &NC_dy, amrex::FArrayBox &NC_xvel_fab, amrex::FArrayBox &NC_yvel_fab, amrex::FArrayBox &NC_temp_fab, amrex::FArrayBox &NC_rhum_fab, amrex::FArrayBox &NC_pres_fab, amrex::FArrayBox &NC_ght_fab, amrex::FArrayBox &NC_hgt_fab, amrex::FArrayBox &NC_psfc_fab, amrex::FArrayBox &NC_msfu_fab, amrex::FArrayBox &NC_msfv_fab, amrex::FArrayBox &NC_msfm_fab, amrex::FArrayBox &NC_sst_fab, amrex::FArrayBox &NC_tsk_fab, amrex::FArrayBox &NC_LAT_fab, amrex::FArrayBox &NC_LON_fab, amrex::IArrayBox &NC_lmask_iab, amrex::Geometry &geom)
AMREX_FORCE_INLINE AMREX_GPU_DEVICE void calc_p_isothermal(const amrex::Real &z, amrex::Real &p)
Definition: ERF_MetgridUtils.H:450
void init_base_state_from_metgrid(const bool use_moisture, const bool metgrid_debug_psfc, const amrex::Real l_rdOcp, const amrex::Box &valid_bx, const int &flag_psfc, amrex::FArrayBox &state, amrex::FArrayBox &r_hse_fab, amrex::FArrayBox &p_hse_fab, amrex::FArrayBox &pi_hse_fab, amrex::FArrayBox &th_hse_fab, amrex::FArrayBox &qv_hse_fab, amrex::FArrayBox &z_phys_cc_fab, const amrex::FArrayBox &NC_psfc_fab)
void init_msfs_from_metgrid(const bool metgrid_debug_msf, amrex::FArrayBox &msfu_fab, amrex::FArrayBox &msfv_fab, amrex::FArrayBox &msfm_fab, const int &flag_msf, amrex::FArrayBox &NC_MSFU_fab, amrex::FArrayBox &NC_MSFV_fab, amrex::FArrayBox &NC_MSFM_fab)
void init_state_from_metgrid(const int lev, const int itime, const bool use_moisture, const bool interp_theta, const bool metgrid_debug_quiescent, const bool metgrid_debug_isothermal, const bool metgrid_debug_dry, const bool metgrid_basic_linear, const bool metgrid_use_below_sfc, const bool metgrid_use_sfc, const bool metgrid_retain_sfc, const amrex::Real metgrid_proximity, const int metgrid_order, const int metgrid_metgrid_force_sfc_k, const amrex::Real l_rdOcp, amrex::Box &tbxc, amrex::Box &tbxu, amrex::Box &tbxv, amrex::FArrayBox &state_fab, amrex::FArrayBox &x_vel_fab, amrex::FArrayBox &y_vel_fab, amrex::FArrayBox &z_vel_fab, amrex::FArrayBox &z_phys_nd_fab, const amrex::FArrayBox &NC_ght_fab, const amrex::FArrayBox &NC_xvel_fab, const amrex::FArrayBox &NC_yvel_fab, const amrex::FArrayBox &NC_temp_fab, const amrex::FArrayBox &NC_rhum_fab, const amrex::FArrayBox &NC_pres_fab, amrex::FArrayBox &tmp_src_fab, amrex::FArrayBox &tmp_dst_fab, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &fabs_for_bcs_xlo, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &fabs_for_bcs_xhi, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &fabs_for_bcs_ylo, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &fabs_for_bcs_yhi, const amrex::Array4< const int > &mask_c_arr, const amrex::Array4< const int > &mask_u_arr, const amrex::Array4< const int > &mask_v_arr)
AMREX_FORCE_INLINE AMREX_GPU_DEVICE void lagrange_setup(char var_type, const bool &exp_interp, const int &orig_n, const int &new_n, const int &order, const int &i, const int &j, amrex::Real *orig_x_z, amrex::Real *orig_x_p, amrex::Real *orig_y, amrex::Real *new_x_z, amrex::Real *new_x_p, amrex::Real *new_y)
Definition: ERF_MetgridUtils.H:146
AMREX_FORCE_INLINE AMREX_GPU_DEVICE void rh_to_mxrat(int i, int j, int k, const amrex::Array4< amrex::Real const > &rhum, const amrex::Array4< amrex::Real const > &temp, const amrex::Array4< amrex::Real const > &pres, int src_indx, const amrex::Array4< amrex::Real > &mxrat)
Definition: ERF_MetgridUtils.H:920
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ pres
Definition: ERF_Kessler.H:25
real(c_double), parameter svp1
Definition: ERF_module_model_constants.F90:78
real(c_double), parameter svp3
Definition: ERF_module_model_constants.F90:80
real(c_double), parameter svp2
Definition: ERF_module_model_constants.F90:79
real(c_double), parameter svpt0
Definition: ERF_module_model_constants.F90:81