ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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_nd_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 &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_nd_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 &  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 
)
486 {
487  // Here we closely follow WRF's vert_interp from
488  // dyn_em/module_initialize_real.F, although changes have been
489  // made to accommodate interpolation relative to height instead of
490  // pressure.
491  int imax_orig = amrex::ubound(amrex::Box(orig_data)).x;
492  int jmax_orig = amrex::ubound(amrex::Box(orig_data)).y;
493  int kmax_orig = amrex::ubound(amrex::Box(orig_data)).z;
494  int kmax_new = amrex::ubound(amrex::Box(new_z_full)).z;
495 
496  AMREX_ASSERT(kmax_orig < 256);
497  AMREX_ASSERT(kmax_new < 256);
498 
499  amrex::GpuArray<amrex::Real,256> new_z;
500  amrex::GpuArray<amrex::Real,256> new_p;
501  amrex::GpuArray<amrex::Real,256> new_data;
502  amrex::Real* new_z_p = new_z.data();
503  amrex::Real* new_p_p = new_p.data();
504  amrex::Real* new_data_p = new_data.data();
505  for (int k=0; k < kmax_new; k++) {
506  if (stag == 'X') {
507  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));
508  } else if (stag == 'Y') {
509  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));
510  } else if (stag == 'M') {
511  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 )+
512  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));
513  }
514  calc_p_isothermal(new_z_p[k], new_p_p[k]);
515  }
516 
517  amrex::GpuArray<amrex::Real,256> orig_z;
518  amrex::Real* orig_z_p = orig_z.data();
519  for (int k=0; k < kmax_orig; k++) {
520  if (stag == 'M') {
521  orig_z_p[k] = orig_z_full(i,j,k);
522  } else if (stag == 'X') {
523  if (i == 0) {
524  orig_z_p[k] = orig_z_full(i,j,k);
525  } else if (i == imax_orig) {
526  orig_z_p[k] = orig_z_full(imax_orig-1,j,k);
527  } else {
528  orig_z_p[k] = 0.5*(orig_z_full(i,j,k)+orig_z_full(i-1,j,k));
529  }
530  } else if (stag == 'Y') {
531  if (j == 0) {
532  orig_z_p[k] = orig_z_full(i,j,k);
533  } else if (j == jmax_orig) {
534  orig_z_p[k] = orig_z_full(i,jmax_orig-1,k);
535  } else {
536  orig_z_p[k] = 0.5*(orig_z_full(i,j,k)+orig_z_full(i,j-1,k));
537  }
538  }
539  }
540 
541  // Check if the data is top-down instead of bottom-up.
542  bool flip_data_required = false;
543  if (orig_z[1] > orig_z[kmax_orig-1]) flip_data_required = true;
544  if (flip_data_required) amrex::Abort("metgrid initialization flip_data_required. Not yet implemented.");
545 
546  // Search for the first level above the surface in the origin data.
547  // This is needed since the origin model topography will be
548  // different than the topography processed by WPS.
549  int k_above_sfc = 0;
550  for (int k=1; k < kmax_orig; k++) {
551  if (orig_z_p[k] > orig_z_p[0]) {
552  k_above_sfc = k;
553  break;
554  }
555  }
556 
557  int kend_order;
558  amrex::GpuArray<amrex::Real,256> ordered_z;
559  amrex::GpuArray<amrex::Real,256> ordered_data;
560  amrex::Real* ordered_z_p = ordered_z.data();
561  amrex::Real* ordered_data_p = ordered_data.data();
562  if (k_above_sfc > 1) {
563  // The levels are not monotonically increasing in height, so
564  // we sort and then make "artistic" quality control choices.
565  int count = 0;
566 
567  // Insert levels that are below the surface.
568  for (int k=1; k < k_above_sfc; k++) {
569  ordered_z_p[count] = orig_z_p[k];
570  ordered_data_p[count] = orig_data(i,j,k);
571  count++;
572  }
573 
574  // Check if the level that is nearest to and below the surface
575  // is "too close". If so, we'll ignore the upper level and keep
576  // the lower. Origin data is likely to be on pressure levels
577  // with higher spatial resolution near-surface, which supports
578  // the choice of eliminating levels that are "too close" in
579  // pressure-space. For simplicity, calculate delta P assuming a
580  // baroclinic atmosphere.
581  amrex::Real Pu, Pl;
582  calc_p_isothermal(orig_z_p[0], Pu);
583  calc_p_isothermal(ordered_z_p[count-1], Pl);
584  if (Pl-Pu < metgrid_proximity) {
585  count--;
586  }
587 
588  // Insert the surface level.
589  ordered_z_p[count] = orig_z_p[0];
590  ordered_data_p[count] = orig_data(i,j,0);
591  count++;
592 
593  // Quoting WRF's comments, the next level to use is at,
594  // "... ta da, the first level above the surface. I know, wow."
595  int knext = k_above_sfc;
596  // Conditionally more strongly use the surface data by removing
597  // levels between the surface and the height corresponding to a
598  // set number of ERF levels from the surface. This forces the
599  // interpolator to use the surface data up through a number of
600  // ERF levels from the surface.
601  if (metgrid_force_sfc_k > 0) {
602  for (int k=k_above_sfc; k < kmax_orig; k++) {
603  if (orig_z_p[k] > new_z_p[metgrid_force_sfc_k-1]) {
604  knext = k;
605  break;
606  }
607  }
608  }
609 
610  // Check if the level that is nearest to and above the surface
611  // is "too close". If so, we'll ignore that level.
612  calc_p_isothermal(orig_z_p[knext], Pu);
613  calc_p_isothermal(ordered_z_p[count-1], Pl);
614  if (Pl-Pu < metgrid_proximity) {
615  knext++;
616  }
617 
618  // Insert levels that are above the surface.
619  for (int k=knext; k < kmax_orig; k++) {
620  ordered_z_p[count] = orig_z_p[k];
621  ordered_data_p[count] = orig_data(i,j,k);
622  count++;
623  }
624 
625  kend_order = count;
626  } else {
627  // The surface is the lowest level in the origin data.
628 
629  // Insert the surface.
630  ordered_z_p[0] = orig_z[0];
631  ordered_data_p[0] = orig_data(i,j,0);
632 
633  // Similar to above, conditionally more strongly use the
634  // surface data.
635  int count = 1;
636  int knext = count;
637  if (metgrid_force_sfc_k > 0) {
638  for (int k=knext; k < kmax_orig; k++) {
639  if (orig_z_p[k] > new_z_p[metgrid_force_sfc_k]) {
640  knext = k;
641  break;
642  }
643  }
644  }
645 
646  // Insert the remaining levels, again ignoring levels that are
647  // "too close" to the prior valid level.
648  for (int k=knext; k < kmax_orig; k++) {
649  amrex::Real Pu, Pl;
650  calc_p_isothermal(orig_z_p[k], Pu);
651  calc_p_isothermal(ordered_z_p[count-1], Pl);
652  if (Pl-Pu < metgrid_proximity) {
653  continue;
654  }
655  ordered_z_p[count] = orig_z_p[k];
656  ordered_data_p[count] = orig_data(i,j,k);
657  count++;
658  }
659  kend_order = count;
660  }
661 
662  int ksta(0), kend(0);
663  if (metgrid_use_below_sfc && metgrid_use_sfc) {
664  // Use all levels.
665  ksta = 0;
666  kend = kend_order-1;
667  } else if (metgrid_use_below_sfc && !metgrid_use_sfc) {
668  // Use all levels except for the surface.
669  int ksfc = 0;
670  for (int k=0; k < kmax_orig; k++) {
671  if (ordered_z_p[k] == orig_z_p[0]) {
672  ksfc = k;
673  break;
674  }
675  }
676  for (int k=ksfc; k < kmax_orig-1; k++) {
677  ordered_z_p[k] = ordered_z_p[k+1];
678  ordered_data_p[k] = ordered_data_p[k+1];
679  }
680  ksta = 0;
681  kend = kend_order-2;
682  } else if (!metgrid_use_below_sfc && metgrid_use_sfc) {
683  // Use all levels above and including the surface.
684  int ksfc = 0;
685  for (int k=0; k < kend_order; k++) {
686  if (ordered_z_p[k] == orig_z_p[0]) {
687  ksfc = k;
688  break;
689  }
690  }
691  ksta = ksfc;
692  kend = kend_order-1;
693  } else {
694  // We shouldn't be in here!
695  amrex::Abort("metgrid initialization, !use_levels below_ground && !metgrid_use_sfc");
696  }
697 
698  // Insert the level of maximum winds.
699 // amrex::Real maxw_above_this_level = 30000.0;
700 // amrex::Real maxw_horiz_pres_diff = 5000.0;
701 // if ((flag_maxw == 1) && (use_maxw)) {
702 // amrex::Abort("metgrid initialization, use_maxw not yet implemented");
703 // }
704 
705  // Insert the level of the tropopause.
706 // amrex::Real trop_horiz_pres_diff = 5000.0;
707 // if ((flag_trop == 1) && (use_trop)) {
708 // amrex::Abort("metgrid initialization, use_trop not yet implemented");
709 // }
710 
711  amrex::GpuArray<amrex::Real,256> ordered_p;
712  amrex::Real* ordered_p_p = ordered_p.data();
713  for (int k=0; k < kend_order; k++) {
714  calc_p_isothermal(ordered_z_p[k], ordered_p_p[k]);
715  }
716 
717  int new_n = 0;
718  int zap_final = 0;
719  amrex::GpuArray<amrex::Real,256> final_z;
720  amrex::GpuArray<amrex::Real,256> final_p;
721  amrex::GpuArray<amrex::Real,256> final_data;
722  amrex::Real* final_z_p = final_z.data();
723  amrex::Real* final_p_p = final_p.data();
724  amrex::Real* final_data_p = final_data.data();
725  final_z_p[0] = ordered_z[ksta];
726  final_p_p[0] = ordered_p[ksta];
727  final_data_p[0] = ordered_data[ksta];
728  for (int k=ksta+1; k <= kend; k++) {
729  if ((final_p_p[new_n]-ordered_p_p[k]) < metgrid_proximity) {
730  zap_final++;
731  } else {
732  new_n++;
733  final_z_p[new_n] = ordered_z_p[k];
734  final_p_p[new_n] = ordered_p_p[k];
735  final_data_p[new_n] = ordered_data_p[k];
736  }
737  }
738  kend -= zap_final;
739 
740  // Call the interpolator.
741  lagrange_setup(var_type,
742  exp_interp,
743  kend-ksta,
744  kmax_new,
745  metgrid_order,
746  i,
747  j,
748  final_z_p,
749  final_p_p,
750  final_data_p,
751  new_z_p,
752  new_p_p,
753  new_data_p);
754 
755  // Optionally replace the lowest level of data with the surface
756  // field from the origin data.
757  if (metgrid_retain_sfc) new_data[0] = ordered_data[0];
758 
759  // Save the interpolated data.
760  for (int k=0; k < kmax_new; k++) {
761  if (mask(i,j,k) && update_bc_data && bx_xlo.contains(i,j,k)) bc_data_xlo(i,j,k,0) = new_data[k];
762  if (mask(i,j,k) && update_bc_data && bx_xhi.contains(i,j,k)) bc_data_xhi(i,j,k,0) = new_data[k];
763  if (mask(i,j,k) && update_bc_data && bx_ylo.contains(i,j,k)) bc_data_ylo(i,j,k,0) = new_data[k];
764  if (mask(i,j,k) && update_bc_data && bx_yhi.contains(i,j,k)) bc_data_yhi(i,j,k,0) = new_data[k];
765  if (itime == 0) new_data_full(i,j,k,src_comp) = new_data[k];
766  }
767 
768 }
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
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 
)
781 {
782  // This subroutine is a bit ham-handed and can be cleaned up later.
783  int imax_orig = amrex::ubound(amrex::Box(orig_data)).x;
784  int jmax_orig = amrex::ubound(amrex::Box(orig_data)).y;
785  int kmax_orig = amrex::ubound(amrex::Box(orig_data)).z;
786 
787  amrex::Real z;
788  if (stag == 'X') {
789  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));
790  }
791  else if (stag == 'Y') {
792  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));
793  }
794  else if (stag == 'M') {
795  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 )+
796  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));
797  }
798 
799  amrex::Real z0, z1;
800  int klow = -1;
801  int khi0 = -1;
802  amrex::Real dzlow = 1.0e12;
803  amrex::Real dzhi0 = -1.0e12;
804  for (int kk = 0; kk < kmax_orig; kk++) {
805  amrex::Real orig_z_stag = 0.0;
806  if (stag == 'M') {
807  orig_z_stag = orig_z(i,j,kk);
808  }
809  if (stag == 'X') {
810  if (i == 0) {
811  orig_z_stag = orig_z(i,j,kk);
812  }
813  else if (i == imax_orig) {
814  orig_z_stag = orig_z(imax_orig-1,j,kk);
815  }
816  else {
817  orig_z_stag = 0.5*(orig_z(i,j,kk)+orig_z(i-1,j,kk));
818  }
819  }
820  else if (stag == 'Y') {
821  if (j == 0) {
822  orig_z_stag = orig_z(i,j,kk);
823  }
824  else if (j == jmax_orig) {
825  orig_z_stag = orig_z(i,jmax_orig-1,kk);
826  }
827  else {
828  orig_z_stag = 0.5*(orig_z(i,j,kk)+orig_z(i,j-1,kk));
829  }
830  }
831 
832  amrex::Real dz = z - orig_z_stag;
833  if ((dz < 0.0) && (dz > dzhi0)) {
834  dzhi0 = dz;
835  khi0 = kk;
836  z1 = orig_z_stag;
837  }
838  if ((dz >= 0.0) && (dz < dzlow)) {
839  dzlow = dz;
840  klow = kk;
841  z0 = orig_z_stag;
842  }
843  } // kk
844 
845  // extrapolate below the bottom surface
846  if (klow == -1) {
847  int khi1 = -1;
848  amrex::Real dzhi1 = -1.0e12;
849  for (int kk = 0; kk < kmax_orig; kk++) {
850  amrex::Real orig_z_stag = 0.0;
851  if (stag == 'M') {
852  orig_z_stag = orig_z(i,j,kk);
853  }
854  else if (stag == 'X') {
855  if (i == 0) {
856  orig_z_stag = orig_z(i,j,kk);
857  }
858  else if (i == imax_orig) {
859  orig_z_stag = orig_z(imax_orig-1,j,kk);
860  }
861  else {
862  orig_z_stag = 0.5*(orig_z(i,j,kk)+orig_z(i-1,j,kk));
863  }
864  }
865  else if (stag == 'Y') {
866  if (j == 0) {
867  orig_z_stag = orig_z(i,j,kk);
868  }
869  else if (j == jmax_orig) {
870  orig_z_stag = orig_z(i,jmax_orig-1,kk);
871  }
872  else {
873  orig_z_stag = 0.5*(orig_z(i,j,kk)+orig_z(i,j-1,kk));
874  }
875  }
876  amrex::Real dz = z - orig_z_stag;
877  if ((dz < 0.0) && (dz > dzhi1) && (kk != khi0)) {
878  dzhi1 = dz;
879  khi1 = kk;
880  z1 = orig_z_stag;
881  }
882  }
883  amrex::Real y0 = orig_data(i,j,khi0,src_comp);
884  amrex::Real y1 = orig_data(i,j,khi1,src_comp);
885  return ( y0-(y1-y0)/(z1-z0)*(z0-z) );
886 
887  // Extrapolate above the top surface
888  } else if (khi0 == -1) {
889  khi0 = klow - 1;
890  int khi1 = klow;
891  if (stag == 'M') {
892  z0 = orig_z(i,j,khi0);
893  }
894  else if (stag == 'X') {
895  if (i == 0) {
896  z0 = orig_z(i,j,khi0);
897  }
898  else if (i == imax_orig) {
899  z0 = orig_z(imax_orig-1,j,khi0);
900  }
901  else {
902  z0 = 0.5*(orig_z(i,j,khi0)+orig_z(i-1,j,khi0));
903  }
904  }
905  else if (stag == 'Y') {
906  if (j == 0) {
907  z0 = orig_z(i,j,khi0);
908  }
909  else if (j == jmax_orig) {
910  z0 = orig_z(i,jmax_orig-1,khi0);
911  }
912  else {
913  z0 = 0.5*(orig_z(i,j,khi0)+orig_z(i,j-1,khi0));
914  }
915  }
916  amrex::Real y0 = orig_data(i,j,khi0,src_comp);
917  amrex::Real y1 = orig_data(i,j,khi1,src_comp);
918  return ( y0+(y1-y0)/(z1-z0)*(z-z0) );
919  } else {
920  // interpolate
921  amrex::Real y0 = orig_data(i,j,klow,src_comp);
922  amrex::Real y1 = orig_data(i,j,khi0,src_comp);
923  return ( y0+(y1-y0)/(z1-z0)*(z-z0) );
924 
925  }
926 }

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