ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_MetgridUtils.H File Reference
#include <ERF.H>
#include <ERF_EOS.H>
#include <ERF_Utils.H>
#include <ERF_ProbCommon.H>
#include <ERF_HSEUtils.H>
Include dependency graph for ERF_MetgridUtils.H:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

void read_from_metgrid (int lev, 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)
 
void init_terrain_from_metgrid (amrex::FArrayBox &z_phys_nd_fab, const amrex::Vector< amrex::FArrayBox > &NC_hgt_fab)
 
void init_state_from_metgrid (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_cc_fab, const amrex::Vector< amrex::FArrayBox > &NC_hgt_fab, const amrex::Vector< amrex::FArrayBox > &NC_ght_fab, const amrex::Vector< amrex::FArrayBox > &NC_xvel_fab, const amrex::Vector< amrex::FArrayBox > &NC_yvel_fab, const amrex::Vector< amrex::FArrayBox > &NC_temp_fab, const amrex::Vector< amrex::FArrayBox > &NC_rhum_fab, const amrex::Vector< amrex::FArrayBox > &NC_pres_fab, amrex::FArrayBox &p_interp_fab, amrex::FArrayBox &t_interp_fab, amrex::FArrayBox &theta_fab, amrex::FArrayBox &mxrat_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)
 
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, const amrex::Vector< amrex::FArrayBox > &NC_MSFU_fab, const amrex::Vector< amrex::FArrayBox > &NC_MSFV_fab, const amrex::Vector< amrex::FArrayBox > &NC_MSFM_fab)
 
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 amrex::Vector< 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::Vector< amrex::FArrayBox > &NC_psfc_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)
 
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)
 
AMREX_FORCE_INLINE AMREX_GPU_DEVICE void calc_p_isothermal (const amrex::Real &z, amrex::Real &p)
 
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 &itime, 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, const bool &update_bc_data, const amrex::Array4< amrex::Real > &bc_data_xlo, const amrex::Array4< amrex::Real > &bc_data_xhi, const amrex::Array4< amrex::Real > &bc_data_ylo, const amrex::Array4< amrex::Real > &bc_data_yhi, const amrex::Box &bx_xlo, const amrex::Box &bx_xhi, const amrex::Box &bx_ylo, const amrex::Box &bx_yhi, const amrex::Array4< const int > &mask)
 
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)
 
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, const amrex::Array4< amrex::Real > &mxrat)
 

Function Documentation

◆ calc_p_isothermal()

AMREX_FORCE_INLINE AMREX_GPU_DEVICE void calc_p_isothermal ( const amrex::Real z,
amrex::Real p 
)
452 {
453  p = p_0*exp(-CONST_GRAV*z/(290.0*R_d));
454 }
constexpr amrex::Real p_0
Definition: ERF_Constants.H:18
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10

Referenced by interpolate_column_metgrid().

Here is the caller graph for this function:

◆ init_base_state_from_metgrid()

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 amrex::Vector< 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::Vector< amrex::FArrayBox > &  NC_psfc_fab 
)

◆ init_msfs_from_metgrid()

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,
const amrex::Vector< amrex::FArrayBox > &  NC_MSFU_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_MSFV_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_MSFM_fab 
)

◆ init_state_from_metgrid()

void init_state_from_metgrid ( 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_cc_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_hgt_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_ght_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_xvel_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_yvel_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_temp_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_rhum_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_pres_fab,
amrex::FArrayBox &  p_interp_fab,
amrex::FArrayBox &  t_interp_fab,
amrex::FArrayBox &  theta_fab,
amrex::FArrayBox &  mxrat_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 
)

◆ init_terrain_from_metgrid()

void init_terrain_from_metgrid ( amrex::FArrayBox &  z_phys_nd_fab,
const amrex::Vector< amrex::FArrayBox > &  NC_hgt_fab 
)

◆ interpolate_column_metgrid()

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 &  itime,
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,
const bool &  update_bc_data,
const amrex::Array4< amrex::Real > &  bc_data_xlo,
const amrex::Array4< amrex::Real > &  bc_data_xhi,
const amrex::Array4< amrex::Real > &  bc_data_ylo,
const amrex::Array4< amrex::Real > &  bc_data_yhi,
const amrex::Box &  bx_xlo,
const amrex::Box &  bx_xhi,
const amrex::Box &  bx_ylo,
const amrex::Box &  bx_yhi,
const amrex::Array4< const int > &  mask 
)
487 {
488  // Here we closely follow WRF's vert_interp from
489  // dyn_em/module_initialize_real.F, although changes have been
490  // made to accommodate interpolation relative to height instead of
491  // pressure.
492  int imax_orig = amrex::ubound(amrex::Box(orig_data)).x;
493  int jmax_orig = amrex::ubound(amrex::Box(orig_data)).y;
494  int kmax_orig = amrex::ubound(amrex::Box(orig_data)).z;
495  int kmax_new = kmax + 1;
496 
497  AMREX_ASSERT(kmax_orig < 256);
498  AMREX_ASSERT(kmax_new < 256);
499 
500  amrex::GpuArray<amrex::Real,256> new_z;
501  amrex::GpuArray<amrex::Real,256> new_p;
502  amrex::GpuArray<amrex::Real,256> new_data;
503  amrex::Real* new_z_p = new_z.data();
504  amrex::Real* new_p_p = new_p.data();
505  amrex::Real* new_data_p = new_data.data();
506  for (int k=0; k < kmax_new; k++) {
507  if (stag == 'X') {
508  new_z_p[k] = 0.25*(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));
509  } else if (stag == 'Y') {
510  new_z_p[k] = 0.25*(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));
511  } else if (stag == 'M') {
512  new_z_p[k] = 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 )+
513  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));
514  }
515  calc_p_isothermal(new_z_p[k], new_p_p[k]);
516  }
517 
518  amrex::GpuArray<amrex::Real,256> orig_z;
519  amrex::Real* orig_z_p = orig_z.data();
520  for (int k=0; k < kmax_orig; k++) {
521  if (stag == 'M') {
522  orig_z_p[k] = orig_z_full(i,j,k);
523  } else if (stag == 'X') {
524  if (i == 0) {
525  orig_z_p[k] = orig_z_full(i,j,k);
526  } else if (i == imax_orig) {
527  orig_z_p[k] = orig_z_full(imax_orig-1,j,k);
528  } else {
529  orig_z_p[k] = 0.5*(orig_z_full(i,j,k)+orig_z_full(i-1,j,k));
530  }
531  } else if (stag == 'Y') {
532  if (j == 0) {
533  orig_z_p[k] = orig_z_full(i,j,k);
534  } else if (j == jmax_orig) {
535  orig_z_p[k] = orig_z_full(i,jmax_orig-1,k);
536  } else {
537  orig_z_p[k] = 0.5*(orig_z_full(i,j,k)+orig_z_full(i,j-1,k));
538  }
539  }
540  }
541 
542  // Check if the data is top-down instead of bottom-up.
543  bool flip_data_required = false;
544  if (orig_z[1] > orig_z[kmax_orig-1]) flip_data_required = true;
545  if (flip_data_required) amrex::Abort("metgrid initialization flip_data_required. Not yet implemented.");
546 
547  // Search for the first level above the surface in the origin data.
548  // This is needed since the origin model topography will be
549  // different than the topography processed by WPS.
550  int k_above_sfc = 0;
551  for (int k=1; k < kmax_orig; k++) {
552  if (orig_z_p[k] > orig_z_p[0]) {
553  k_above_sfc = k;
554  break;
555  }
556  }
557 
558  int kend_order;
559  amrex::GpuArray<amrex::Real,256> ordered_z;
560  amrex::GpuArray<amrex::Real,256> ordered_data;
561  amrex::Real* ordered_z_p = ordered_z.data();
562  amrex::Real* ordered_data_p = ordered_data.data();
563  if (k_above_sfc > 1) {
564  // The levels are not monotonically increasing in height, so
565  // we sort and then make "artistic" quality control choices.
566  int count = 0;
567 
568  // Insert levels that are below the surface.
569  for (int k=1; k < k_above_sfc; k++) {
570  ordered_z_p[count] = orig_z_p[k];
571  ordered_data_p[count] = orig_data(i,j,k);
572  count++;
573  }
574 
575  // Check if the level that is nearest to and below the surface
576  // is "too close". If so, we'll ignore the upper level and keep
577  // the lower. Origin data is likely to be on pressure levels
578  // with higher spatial resolution near-surface, which supports
579  // the choice of eliminating levels that are "too close" in
580  // pressure-space. For simplicity, calculate delta P assuming a
581  // baroclinic atmosphere.
582  amrex::Real Pu, Pl;
583  calc_p_isothermal(orig_z_p[0], Pu);
584  calc_p_isothermal(ordered_z_p[count-1], Pl);
585  if (Pl-Pu < metgrid_proximity) {
586  count--;
587  }
588 
589  // Insert the surface level.
590  ordered_z_p[count] = orig_z_p[0];
591  ordered_data_p[count] = orig_data(i,j,0);
592  count++;
593 
594  // Quoting WRF's comments, the next level to use is at,
595  // "... ta da, the first level above the surface. I know, wow."
596  int knext = k_above_sfc;
597  // Conditionally more strongly use the surface data by removing
598  // levels between the surface and the height corresponding to a
599  // set number of ERF levels from the surface. This forces the
600  // interpolator to use the surface data up through a number of
601  // ERF levels from the surface.
602  if (metgrid_force_sfc_k > 0) {
603  for (int k=k_above_sfc; k < kmax_orig; k++) {
604  if (orig_z_p[k] > new_z_p[metgrid_force_sfc_k-1]) {
605  knext = k;
606  break;
607  }
608  }
609  }
610 
611  // Check if the level that is nearest to and above the surface
612  // is "too close". If so, we'll ignore that level.
613  calc_p_isothermal(orig_z_p[knext], Pu);
614  calc_p_isothermal(ordered_z_p[count-1], Pl);
615  if (Pl-Pu < metgrid_proximity) {
616  knext++;
617  }
618 
619  // Insert levels that are above the surface.
620  for (int k=knext; k < kmax_orig; k++) {
621  ordered_z_p[count] = orig_z_p[k];
622  ordered_data_p[count] = orig_data(i,j,k);
623  count++;
624  }
625 
626  kend_order = count;
627  } else {
628  // The surface is the lowest level in the origin data.
629 
630  // Insert the surface.
631  ordered_z_p[0] = orig_z[0];
632  ordered_data_p[0] = orig_data(i,j,0);
633 
634  // Similar to above, conditionally more strongly use the
635  // surface data.
636  int count = 1;
637  int knext = count;
638  if (metgrid_force_sfc_k > 0) {
639  for (int k=knext; k < kmax_orig; k++) {
640  if (orig_z_p[k] > new_z_p[metgrid_force_sfc_k]) {
641  knext = k;
642  break;
643  }
644  }
645  }
646 
647  // Insert the remaining levels, again ignoring levels that are
648  // "too close" to the prior valid level.
649  for (int k=knext; k < kmax_orig; k++) {
650  amrex::Real Pu, Pl;
651  calc_p_isothermal(orig_z_p[k], Pu);
652  calc_p_isothermal(ordered_z_p[count-1], Pl);
653  if (Pl-Pu < metgrid_proximity) {
654  continue;
655  }
656  ordered_z_p[count] = orig_z_p[k];
657  ordered_data_p[count] = orig_data(i,j,k);
658  count++;
659  }
660  kend_order = count;
661  }
662 
663  int ksta(0), kend(0);
664  if (metgrid_use_below_sfc && metgrid_use_sfc) {
665  // Use all levels.
666  ksta = 0;
667  kend = kend_order-1;
668  } else if (metgrid_use_below_sfc && !metgrid_use_sfc) {
669  // Use all levels except for the surface.
670  int ksfc = 0;
671  for (int k=0; k < kmax_orig; k++) {
672  if (ordered_z_p[k] == orig_z_p[0]) {
673  ksfc = k;
674  break;
675  }
676  }
677  for (int k=ksfc; k < kmax_orig-1; k++) {
678  ordered_z_p[k] = ordered_z_p[k+1];
679  ordered_data_p[k] = ordered_data_p[k+1];
680  }
681  ksta = 0;
682  kend = kend_order-2;
683  } else if (!metgrid_use_below_sfc && metgrid_use_sfc) {
684  // Use all levels above and including the surface.
685  int ksfc = 0;
686  for (int k=0; k < kend_order; k++) {
687  if (ordered_z_p[k] == orig_z_p[0]) {
688  ksfc = k;
689  break;
690  }
691  }
692  ksta = ksfc;
693  kend = kend_order-1;
694  } else {
695  // We shouldn't be in here!
696  amrex::Abort("metgrid initialization, !use_levels below_ground && !metgrid_use_sfc");
697  }
698 
699  // Insert the level of maximum winds.
700 // amrex::Real maxw_above_this_level = 30000.0;
701 // amrex::Real maxw_horiz_pres_diff = 5000.0;
702 // if ((flag_maxw == 1) && (use_maxw)) {
703 // amrex::Abort("metgrid initialization, use_maxw not yet implemented");
704 // }
705 
706  // Insert the level of the tropopause.
707 // amrex::Real trop_horiz_pres_diff = 5000.0;
708 // if ((flag_trop == 1) && (use_trop)) {
709 // amrex::Abort("metgrid initialization, use_trop not yet implemented");
710 // }
711 
712  amrex::GpuArray<amrex::Real,256> ordered_p;
713  amrex::Real* ordered_p_p = ordered_p.data();
714  for (int k=0; k < kend_order; k++) {
715  calc_p_isothermal(ordered_z_p[k], ordered_p_p[k]);
716  }
717 
718  int new_n = 0;
719  int zap_final = 0;
720  amrex::GpuArray<amrex::Real,256> final_z;
721  amrex::GpuArray<amrex::Real,256> final_p;
722  amrex::GpuArray<amrex::Real,256> final_data;
723  amrex::Real* final_z_p = final_z.data();
724  amrex::Real* final_p_p = final_p.data();
725  amrex::Real* final_data_p = final_data.data();
726  final_z_p[0] = ordered_z[ksta];
727  final_p_p[0] = ordered_p[ksta];
728  final_data_p[0] = ordered_data[ksta];
729  for (int k=ksta+1; k <= kend; k++) {
730  if ((final_p_p[new_n]-ordered_p_p[k]) < metgrid_proximity) {
731  zap_final++;
732  } else {
733  new_n++;
734  final_z_p[new_n] = ordered_z_p[k];
735  final_p_p[new_n] = ordered_p_p[k];
736  final_data_p[new_n] = ordered_data_p[k];
737  }
738  }
739  kend -= zap_final;
740 
741  // Call the interpolator.
742  lagrange_setup(var_type,
743  exp_interp,
744  kend-ksta,
745  kmax_new,
746  metgrid_order,
747  i,
748  j,
749  final_z_p,
750  final_p_p,
751  final_data_p,
752  new_z_p,
753  new_p_p,
754  new_data_p);
755 
756  // Optionally replace the lowest level of data with the surface
757  // field from the origin data.
758  if (metgrid_retain_sfc) new_data[0] = ordered_data[0];
759 
760  // Save the interpolated data.
761  for (int k=0; k < kmax_new; k++) {
762  if (mask(i,j,k) && update_bc_data && bx_xlo.contains(i,j,k)) bc_data_xlo(i,j,k,0) = new_data[k];
763  if (mask(i,j,k) && update_bc_data && bx_xhi.contains(i,j,k)) bc_data_xhi(i,j,k,0) = new_data[k];
764  if (mask(i,j,k) && update_bc_data && bx_ylo.contains(i,j,k)) bc_data_ylo(i,j,k,0) = new_data[k];
765  if (mask(i,j,k) && update_bc_data && bx_yhi.contains(i,j,k)) bc_data_yhi(i,j,k,0) = new_data[k];
766  if (itime == 0) new_data_full(i,j,k,src_comp) = new_data[k];
767  }
768 
769 }
AMREX_FORCE_INLINE AMREX_GPU_DEVICE void calc_p_isothermal(const amrex::Real &z, amrex::Real &p)
Definition: ERF_MetgridUtils.H:450
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::Real Real
Definition: ERF_ShocInterface.H:19
Here is the call graph for this function:

◆ interpolate_column_metgrid_linear()

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 
)
782 {
783  // This subroutine is a bit ham-handed and can be cleaned up later.
784  int imax_orig = amrex::ubound(amrex::Box(orig_data)).x;
785  int jmax_orig = amrex::ubound(amrex::Box(orig_data)).y;
786  int kmax_orig = amrex::ubound(amrex::Box(orig_data)).z;
787 
788  amrex::Real z;
789  if (stag == 'X') {
790  z = 0.25*(new_z(i,j,k)+new_z(i,j+1,k)+new_z(i,j,k+1)+new_z(i,j+1,k+1));
791  }
792  else if (stag == 'Y') {
793  z = 0.25*(new_z(i,j,k)+new_z(i+1,j,k)+new_z(i,j,k+1)+new_z(i+1,j,k+1));
794  }
795  else if (stag == 'M') {
796  z = 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 )+
797  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));
798  }
799 
800  amrex::Real z0, z1;
801  int klow = -1;
802  int khi0 = -1;
803  amrex::Real dzlow = 1.0e12;
804  amrex::Real dzhi0 = -1.0e12;
805  for (int kk = 0; kk < kmax_orig; kk++) {
806  amrex::Real orig_z_stag = 0.0;
807  if (stag == 'M') {
808  orig_z_stag = orig_z(i,j,kk);
809  }
810  if (stag == 'X') {
811  if (i == 0) {
812  orig_z_stag = orig_z(i,j,kk);
813  }
814  else if (i == imax_orig) {
815  orig_z_stag = orig_z(imax_orig-1,j,kk);
816  }
817  else {
818  orig_z_stag = 0.5*(orig_z(i,j,kk)+orig_z(i-1,j,kk));
819  }
820  }
821  else if (stag == 'Y') {
822  if (j == 0) {
823  orig_z_stag = orig_z(i,j,kk);
824  }
825  else if (j == jmax_orig) {
826  orig_z_stag = orig_z(i,jmax_orig-1,kk);
827  }
828  else {
829  orig_z_stag = 0.5*(orig_z(i,j,kk)+orig_z(i,j-1,kk));
830  }
831  }
832 
833  amrex::Real dz = z - orig_z_stag;
834  if ((dz < 0.0) && (dz > dzhi0)) {
835  dzhi0 = dz;
836  khi0 = kk;
837  z1 = orig_z_stag;
838  }
839  if ((dz >= 0.0) && (dz < dzlow)) {
840  dzlow = dz;
841  klow = kk;
842  z0 = orig_z_stag;
843  }
844  } // kk
845 
846  // extrapolate below the bottom surface
847  if (klow == -1) {
848  int khi1 = -1;
849  amrex::Real dzhi1 = -1.0e12;
850  for (int kk = 0; kk < kmax_orig; kk++) {
851  amrex::Real orig_z_stag = 0.0;
852  if (stag == 'M') {
853  orig_z_stag = orig_z(i,j,kk);
854  }
855  else if (stag == 'X') {
856  if (i == 0) {
857  orig_z_stag = orig_z(i,j,kk);
858  }
859  else if (i == imax_orig) {
860  orig_z_stag = orig_z(imax_orig-1,j,kk);
861  }
862  else {
863  orig_z_stag = 0.5*(orig_z(i,j,kk)+orig_z(i-1,j,kk));
864  }
865  }
866  else if (stag == 'Y') {
867  if (j == 0) {
868  orig_z_stag = orig_z(i,j,kk);
869  }
870  else if (j == jmax_orig) {
871  orig_z_stag = orig_z(i,jmax_orig-1,kk);
872  }
873  else {
874  orig_z_stag = 0.5*(orig_z(i,j,kk)+orig_z(i,j-1,kk));
875  }
876  }
877  amrex::Real dz = z - orig_z_stag;
878  if ((dz < 0.0) && (dz > dzhi1) && (kk != khi0)) {
879  dzhi1 = dz;
880  khi1 = kk;
881  z1 = orig_z_stag;
882  }
883  }
884  amrex::Real y0 = orig_data(i,j,khi0,src_comp);
885  amrex::Real y1 = orig_data(i,j,khi1,src_comp);
886  return ( y0-(y1-y0)/(z1-z0)*(z0-z) );
887 
888  // Extrapolate above the top surface
889  } else if (khi0 == -1) {
890  khi0 = klow - 1;
891  int khi1 = klow;
892  if (stag == 'M') {
893  z0 = orig_z(i,j,khi0);
894  }
895  else if (stag == 'X') {
896  if (i == 0) {
897  z0 = orig_z(i,j,khi0);
898  }
899  else if (i == imax_orig) {
900  z0 = orig_z(imax_orig-1,j,khi0);
901  }
902  else {
903  z0 = 0.5*(orig_z(i,j,khi0)+orig_z(i-1,j,khi0));
904  }
905  }
906  else if (stag == 'Y') {
907  if (j == 0) {
908  z0 = orig_z(i,j,khi0);
909  }
910  else if (j == jmax_orig) {
911  z0 = orig_z(i,jmax_orig-1,khi0);
912  }
913  else {
914  z0 = 0.5*(orig_z(i,j,khi0)+orig_z(i,j-1,khi0));
915  }
916  }
917  amrex::Real y0 = orig_data(i,j,khi0,src_comp);
918  amrex::Real y1 = orig_data(i,j,khi1,src_comp);
919  return ( y0+(y1-y0)/(z1-z0)*(z-z0) );
920  } else {
921  // interpolate
922  amrex::Real y0 = orig_data(i,j,klow,src_comp);
923  amrex::Real y1 = orig_data(i,j,khi0,src_comp);
924  return ( y0+(y1-y0)/(z1-z0)*(z-z0) );
925 
926  }
927 }
Real z0
Definition: ERF_InitCustomPertVels_ScalarAdvDiff.H:8

◆ lagrange_interp()

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 
)
121 {
122  // Interpolation using Lagrange polynomials.
123  // P(x) = f(x0)Ln0(x) + ... + f(xn)Lnn(x)
124  // where Lnk(x) = (x -x0)(x -x1)...(x -xk-1)(x -xk+1)...(x -xn)
125  // ---------------------------------------------
126  // (xk-x0)(xk-x1)...(xk-xk-1)(xk-xk+1)...(xk-xn)
127  amrex::Real Px = 0.0;
128  for (int i=0; i <= order; i++) {
129  amrex::Real n = 1.0;
130  amrex::Real d = 1.0;
131  for (int k=0; k <= order; k++) {
132  if (k == i) continue;
133  n *= new_x-x[k];
134  d *= x[i]-x[k];
135  }
136  if (d != 0.0) {
137  Px += y[i]*n/d;
138  }
139  }
140  new_y = Px;
141 }

Referenced by lagrange_setup().

Here is the caller graph for this function:

◆ lagrange_setup()

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 
)
159 {
160 
161  amrex::Real CRC_const1 = 11880.516; // m
162  amrex::Real CRC_const2 = 0.1902632;
163  amrex::Real CRC_const3 = 0.0065; // K km-1
164 
165  amrex::ignore_unused(i,j);
166 #ifndef AMREX_USE_GPU
167  bool debug = false;
168 #endif
169 
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;
173 #endif
174  // Determine if interpolating or extrapolating.
175  // If interpolating, find bounding x values and store the indices.
176  bool extrapolating = true;
177  int kl, kr;
178  for (int ko=0; ko < orig_n-1; ko++) {
179  amrex::Real a = new_x_z[new_k]-orig_x_z[ko];
180  amrex::Real b = new_x_z[new_k]-orig_x_z[ko+1];
181  if (a*b <= 0.0) {
182  kl = ko;
183  kr = ko+1;
184  extrapolating = false;
185  break;
186  }
187  }
188 
189  if (extrapolating) {
190  if (var_type == 'T') {
191  // Assume a standard atmosphere -6.5 K km-1 lapse rate.
192  // Comparable to the WRF default, t_extrap_type=2.
193  amrex::Real depth_of_extrap_in_p = new_x_p[new_k]-orig_x_p[0];
194  amrex::Real avg_of_extrap_p = 0.5*(new_x_p[new_k]+orig_x_p[0]);
195  amrex::Real temp_extrap_starting_point = orig_y[0]*std::pow(orig_x_p[0]/100000.0, R_d/Cp_d);
196  amrex::Real dZdP = CRC_const1*CRC_const2*std::pow(avg_of_extrap_p/100.0, CRC_const2-1.0);
197  amrex::Real dZ = dZdP*(depth_of_extrap_in_p/100.0);
198  amrex::Real dT = dZ*CRC_const3;
199  new_y[new_k] = (temp_extrap_starting_point+dT)*std::pow(100000.0/new_x_p[new_k], R_d/Cp_d);
200  } else {
201  // Use a constant value below ground.
202  // Comparable to the WRF default, extrap_type=2.
203  new_y[new_k] = orig_y[0];
204  }
205  continue;
206  }
207 
208  if (order%2 != 0) {
209  if ((kl-((order+1)/2-1) >= 0) && (kr+((order+1)/2-1) <= orig_n-1)) {
210  // Odd order interpolation.
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;
216 #endif
217  amrex::Real new_x;
218  amrex::GpuArray<amrex::Real,9> orig_x_sub;
219  amrex::Real* orig_x_sub_p = orig_x_sub.data();
220  if (exp_interp) {
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]; }
223  } else {
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]; }
226  }
227  amrex::GpuArray<amrex::Real,9> orig_y_sub;
228  amrex::Real* orig_y_sub_p = orig_y_sub.data();
229  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
230 #ifndef AMREX_USE_GPU
231  if (debug) {
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;
238  }
239 #endif
240  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
241  } else {
242  // Interpolation stencil is too big due to proximity to the surface or model top.
243  // We have bounding points so resort to a simple linear interpolation.
244  int ksta = kl;
245  int kend = kr;
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;
249 #endif
250  amrex::Real new_x;
251  amrex::GpuArray<amrex::Real,2> orig_x_sub;
252  amrex::GpuArray<amrex::Real,2> orig_y_sub;
253  amrex::Real* orig_x_sub_p = orig_x_sub.data();
254  amrex::Real* orig_y_sub_p = orig_y_sub.data();
255  if (exp_interp) {
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]; }
258  } else {
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]; }
261  }
262  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
263 #ifndef AMREX_USE_GPU
264  if (debug) {
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;
271  }
272 #endif
273  lagrange_interp(1, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
274  }
275  } else if (order%2 == 0) {
276  if ((kl-(order/2) >= 0) && (kr+order/2 <= orig_n-1)) {
277  // Even order interpolation.
278  // Interpolate twice. Offset the range by 1 the 2nd time. Average the results.
279  amrex::Real new_y_l, new_y_r;
280  {
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;
286 #endif
287  amrex::Real new_x;
288  amrex::GpuArray<amrex::Real,10> orig_x_sub;
289  amrex::GpuArray<amrex::Real,10> orig_y_sub;
290  amrex::Real* orig_x_sub_p = orig_x_sub.data();
291  amrex::Real* orig_y_sub_p = orig_y_sub.data();
292  if (exp_interp) {
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]; }
295  } else {
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]; }
298  }
299  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
300 #ifndef AMREX_USE_GPU
301  if (debug) {
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;
308  }
309 #endif
310  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y_l);
311  }
312  {
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;
318 #endif
319  amrex::Real new_x;
320  amrex::GpuArray<amrex::Real,10> orig_x_sub;
321  amrex::GpuArray<amrex::Real,10> orig_y_sub;
322  amrex::Real* orig_x_sub_p = orig_x_sub.data();
323  amrex::Real* orig_y_sub_p = orig_y_sub.data();
324  if (exp_interp) {
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]; }
327  } else {
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]; }
330  }
331  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
332 #ifndef AMREX_USE_GPU
333  if (debug) {
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;
340  }
341 #endif
342  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y_r);
343  }
344  new_y[new_k] = 0.5*(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;
351 #endif
352  amrex::Real new_x;
353  amrex::GpuArray<amrex::Real,10> orig_x_sub;
354  amrex::GpuArray<amrex::Real,10> orig_y_sub;
355  amrex::Real* orig_x_sub_p = orig_x_sub.data();
356  amrex::Real* orig_y_sub_p = orig_y_sub.data();
357  if (exp_interp) {
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]; }
360  } else {
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]; }
363  }
364  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
365 #ifndef AMREX_USE_GPU
366  if (debug) {
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;
373  }
374 #endif
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;
382 #endif
383  amrex::Real new_x;
384  amrex::GpuArray<amrex::Real,10> orig_x_sub;
385  amrex::GpuArray<amrex::Real,10> orig_y_sub;
386  amrex::Real* orig_x_sub_p = orig_x_sub.data();
387  amrex::Real* orig_y_sub_p = orig_y_sub.data();
388  if (exp_interp) {
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]; }
391  } else {
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]; }
394  }
395  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
396 #ifndef AMREX_USE_GPU
397  if (debug) {
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;
404  }
405 #endif
406  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
407  } else {
408  // Linear interpolation.
409  int ksta = kl;
410  int kend = kr;
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;
414 #endif
415  amrex::Real new_x;
416  amrex::GpuArray<amrex::Real,2> orig_x_sub;
417  amrex::GpuArray<amrex::Real,2> orig_y_sub;
418  amrex::Real* orig_x_sub_p = orig_x_sub.data();
419  amrex::Real* orig_y_sub_p = orig_y_sub.data();
420  if (exp_interp) {
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]; }
423  } else {
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]; }
426  }
427  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
428 #ifndef AMREX_USE_GPU
429  if (debug) {
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;
436  }
437 #endif
438  lagrange_interp(1, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
439  }
440  }
441 #ifndef AMREX_USE_GPU
442  if (debug) amrex::Print() << " new_y[" << new_k << "]=" << new_y[new_k] << std::endl;
443 #endif
444  }
445 }
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:12
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

Referenced by interpolate_column_metgrid().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_from_metgrid()

void read_from_metgrid ( int  lev,
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 
)

◆ rh_to_mxrat()

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,
const amrex::Array4< amrex::Real > &  mxrat 
)
939 {
940  amrex::Real qv_max_p_safe = 10000.0; // WRF default value
941  amrex::Real qv_max_flag = 1.0e-5; // WRF default value
942  amrex::Real qv_max_value = 3.0e-6; // WRF default value
943  amrex::Real qv_min_p_safe = 110000.0; // WRF default value
944  amrex::Real qv_min_flag = 1.0e-6; // WRF default value
945  amrex::Real qv_min_value = 1.0e-6; // WRF default value
946  amrex::Real eps = 0.622;
947  amrex::Real svp1 = 0.6112;
948  amrex::Real svp2 = 17.67;
949  amrex::Real svp3 = 29.65;
950  amrex::Real svpt0 = 273.15;
951  // WRF's method when model_config_rec%rh2qv_wrt_liquid=.true. (default behavior)
952  if (temp(i,j,k) != 0.0) {
953  amrex::Real es=0.01*rhum(i,j,k)*svp1*10.0*exp(svp2*(temp(i,j,k)-svpt0)/(temp(i,j,k)-svp3));
954  if (es >= pres(i,j,k)/100.0) {
955  // vapor pressure exceeds total pressure
956  mxrat(i,j,k) = std::pow(10.0, -6);
957  }
958  else {
959  mxrat(i,j,k) = amrex::max(eps*es/(pres(i,j,k)/100.0-es), 1.0e-6);
960  }
961  }
962  else {
963  // I don't know why there's a fringe case handled in WRF where T is absolute zero...
964  // Let's just deal with it here in case we also end up needing it.
965  mxrat(i,j,k) = 1.0e-6;
966  }
967  // See the below comment from WRF dyn_em/module_initialize_real.F rh_to_mxrat1.
968  // For pressures above a defined level, reasonable Qv values should be
969  // a certain value or smaller. If they are larger than this, the input data
970  // probably had "missing" RH, and we filled in some values. This is an
971  // attempt to catch those. Also, set the minimum value for the entire
972  // domain that is above the selected pressure level.
973  if (pres(i,j,k) < qv_max_p_safe) {
974  if (mxrat(i,j,k) > qv_max_flag) {
975  mxrat(i,j,k) = qv_max_value;
976  }
977  }
978  if (pres(i,j,k) < qv_min_p_safe) {
979  if (mxrat(i,j,k) < qv_min_flag) {
980  mxrat(i,j,k) = qv_min_value;
981  }
982  }
983 }
@ 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