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

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

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 
)
156 {
157 
158  amrex::Real CRC_const1 = 11880.516; // m
159  amrex::Real CRC_const2 = 0.1902632;
160  amrex::Real CRC_const3 = 0.0065; // K km-1
161 
162  amrex::ignore_unused(i,j);
163 #ifndef AMREX_USE_GPU
164  bool debug = false;
165 #endif
166 
167  for (int new_k=0; new_k < new_n; new_k++) {
168 #ifndef AMREX_USE_GPU
169  if (debug) amrex::Print() << "new_k=" << new_k;
170 #endif
171  // Determine if interpolating or extrapolating.
172  // If interpolating, find bounding x values and store the indices.
173  bool extrapolating = true;
174  int kl, kr;
175  for (int ko=0; ko < orig_n-1; ko++) {
176  amrex::Real a = new_x_z[new_k]-orig_x_z[ko];
177  amrex::Real b = new_x_z[new_k]-orig_x_z[ko+1];
178  if (a*b <= 0.0) {
179  kl = ko;
180  kr = ko+1;
181  extrapolating = false;
182  break;
183  }
184  }
185 
186  if (extrapolating) {
187  if (var_type == 'T') {
188  // Assume a standard atmosphere -6.5 K km-1 lapse rate.
189  // Comparable to the WRF default, t_extrap_type=2.
190  amrex::Real depth_of_extrap_in_p = new_x_p[new_k]-orig_x_p[0];
191  amrex::Real avg_of_extrap_p = 0.5*(new_x_p[new_k]+orig_x_p[0]);
192  amrex::Real temp_extrap_starting_point = orig_y[0]*std::pow(orig_x_p[0]/100000.0, R_d/Cp_d);
193  amrex::Real dZdP = CRC_const1*CRC_const2*std::pow(avg_of_extrap_p/100.0, CRC_const2-1.0);
194  amrex::Real dZ = dZdP*(depth_of_extrap_in_p/100.0);
195  amrex::Real dT = dZ*CRC_const3;
196  new_y[new_k] = (temp_extrap_starting_point+dT)*std::pow(100000.0/new_x_p[new_k], R_d/Cp_d);
197  } else {
198  // Use a constant value below ground.
199  // Comparable to the WRF default, extrap_type=2.
200  new_y[new_k] = orig_y[0];
201  }
202  continue;
203  }
204 
205  if (order%2 != 0) {
206  if ((kl-((order+1)/2-1) >= 0) && (kr+((order+1)/2-1) <= orig_n-1)) {
207  // Odd order interpolation.
208  int ksta = kl-(((order+1)/2)-1);
209  int kend = ksta+order;
210 #ifndef AMREX_USE_GPU
211  int ksize = kend-ksta;
212  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;
213 #endif
214  amrex::Real new_x;
215  amrex::GpuArray<amrex::Real,9> orig_x_sub;
216  amrex::Real* orig_x_sub_p = orig_x_sub.data();
217  if (exp_interp) {
218  new_x = new_x_p[new_k];
219  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
220  } else {
221  new_x = new_x_z[new_k];
222  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
223  }
224  amrex::GpuArray<amrex::Real,9> orig_y_sub;
225  amrex::Real* orig_y_sub_p = orig_y_sub.data();
226  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
227 #ifndef AMREX_USE_GPU
228  if (debug) {
229  amrex::Print() << " orig_x_sub_p = [";
230  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_x_sub_p[k];
231  amrex::Print() << "]" << std::endl;
232  amrex::Print() << " orig_y_sub_p = [";
233  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_y_sub_p[k];
234  amrex::Print() << "]" << std::endl;
235  }
236 #endif
237  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
238  } else {
239  // Interpolation stencil is too big due to proximity to the surface or model top.
240  // We have bounding points so resort to a simple linear interpolation.
241  int ksta = kl;
242  int kend = kr;
243 #ifndef AMREX_USE_GPU
244  int ksize = kend-ksta+1;
245  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;
246 #endif
247  amrex::Real new_x;
248  amrex::GpuArray<amrex::Real,2> orig_x_sub;
249  amrex::GpuArray<amrex::Real,2> orig_y_sub;
250  amrex::Real* orig_x_sub_p = orig_x_sub.data();
251  amrex::Real* orig_y_sub_p = orig_y_sub.data();
252  if (exp_interp) {
253  new_x = new_x_p[new_k];
254  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
255  } else {
256  new_x = new_x_z[new_k];
257  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
258  }
259  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
260 #ifndef AMREX_USE_GPU
261  if (debug) {
262  amrex::Print() << " orig_x_sub_p = [";
263  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_x_sub_p[k];
264  amrex::Print() << "]" << std::endl;
265  amrex::Print() << " orig_y_sub_p = [";
266  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_y_sub_p[k];
267  amrex::Print() << "]" << std::endl;
268  }
269 #endif
270  lagrange_interp(1, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
271  }
272  } else if (order%2 == 0) {
273  if ((kl-(order/2) >= 0) && (kr+order/2 <= orig_n-1)) {
274  // Even order interpolation.
275  // Interpolate twice. Offset the range by 1 the 2nd time. Average the results.
276  amrex::Real new_y_l, new_y_r;
277  {
278  int ksta = kl-(order/2-1);
279  int kend = ksta+order;
280 #ifndef AMREX_USE_GPU
281  int ksize = kend-ksta;
282  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;
283 #endif
284  amrex::Real new_x;
285  amrex::GpuArray<amrex::Real,10> orig_x_sub;
286  amrex::GpuArray<amrex::Real,10> orig_y_sub;
287  amrex::Real* orig_x_sub_p = orig_x_sub.data();
288  amrex::Real* orig_y_sub_p = orig_y_sub.data();
289  if (exp_interp) {
290  new_x = new_x_p[new_k];
291  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
292  } else {
293  new_x = new_x_z[new_k];
294  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
295  }
296  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
297 #ifndef AMREX_USE_GPU
298  if (debug) {
299  amrex::Print() << " orig_x_sub_p = [";
300  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_x_sub_p[k];
301  amrex::Print() << "]" << std::endl;
302  amrex::Print() << " orig_y_sub_p = [";
303  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_y_sub_p[k];
304  amrex::Print() << "]" << std::endl;
305  }
306 #endif
307  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y_l);
308  }
309  {
310  int ksta = kl-order/2;
311  int kend = ksta+order;
312 #ifndef AMREX_USE_GPU
313  int ksize = kend-ksta;
314  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;
315 #endif
316  amrex::Real new_x;
317  amrex::GpuArray<amrex::Real,10> orig_x_sub;
318  amrex::GpuArray<amrex::Real,10> orig_y_sub;
319  amrex::Real* orig_x_sub_p = orig_x_sub.data();
320  amrex::Real* orig_y_sub_p = orig_y_sub.data();
321  if (exp_interp) {
322  new_x = new_x_p[new_k];
323  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
324  } else {
325  new_x = new_x_z[new_k];
326  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
327  }
328  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
329 #ifndef AMREX_USE_GPU
330  if (debug) {
331  amrex::Print() << " orig_x_sub_p = [";
332  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_x_sub_p[k];
333  amrex::Print() << "]" << std::endl;
334  amrex::Print() << " orig_y_sub_p = [";
335  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_y_sub_p[k];
336  amrex::Print() << "]" << std::endl;
337  }
338 #endif
339  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y_r);
340  }
341  new_y[new_k] = 0.5*(new_y_l+new_y_r);
342  } else if ((kl-(order/2-1) >= 0) && (kr+order/2 <= orig_n-1)) {
343  int ksta = kl-(order/2-1);
344  int kend = ksta+order;
345 #ifndef AMREX_USE_GPU
346  int ksize = kend-ksta;
347  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;
348 #endif
349  amrex::Real new_x;
350  amrex::GpuArray<amrex::Real,10> orig_x_sub;
351  amrex::GpuArray<amrex::Real,10> orig_y_sub;
352  amrex::Real* orig_x_sub_p = orig_x_sub.data();
353  amrex::Real* orig_y_sub_p = orig_y_sub.data();
354  if (exp_interp) {
355  new_x = new_x_p[new_k];
356  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
357  } else {
358  new_x = new_x_z[new_k];
359  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
360  }
361  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
362 #ifndef AMREX_USE_GPU
363  if (debug) {
364  amrex::Print() << " orig_x_sub_p = [";
365  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_x_sub_p[k];
366  amrex::Print() << "]" << std::endl;
367  amrex::Print() << " orig_y_sub_p = [";
368  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_y_sub_p[k];
369  amrex::Print() << "]" << std::endl;
370  }
371 #endif
372  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
373  } else if ((kl-order/2 >= 0) && (kr+order/2-1 <= orig_n-1)) {
374  int ksta = kl-order/2;
375  int kend = ksta+order;
376 #ifndef AMREX_USE_GPU
377  int ksize = kend-ksta;
378  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;
379 #endif
380  amrex::Real new_x;
381  amrex::GpuArray<amrex::Real,10> orig_x_sub;
382  amrex::GpuArray<amrex::Real,10> orig_y_sub;
383  amrex::Real* orig_x_sub_p = orig_x_sub.data();
384  amrex::Real* orig_y_sub_p = orig_y_sub.data();
385  if (exp_interp) {
386  new_x = new_x_p[new_k];
387  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
388  } else {
389  new_x = new_x_z[new_k];
390  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
391  }
392  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
393 #ifndef AMREX_USE_GPU
394  if (debug) {
395  amrex::Print() << " orig_x_sub_p = [";
396  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_x_sub_p[k];
397  amrex::Print() << "]" << std::endl;
398  amrex::Print() << " orig_y_sub_p = [";
399  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_y_sub_p[k];
400  amrex::Print() << "]" << std::endl;
401  }
402 #endif
403  lagrange_interp(order, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
404  } else {
405  // Linear interpolation.
406  int ksta = kl;
407  int kend = kr;
408 #ifndef AMREX_USE_GPU
409  int ksize = kend-ksta+1;
410  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;
411 #endif
412  amrex::Real new_x;
413  amrex::GpuArray<amrex::Real,2> orig_x_sub;
414  amrex::GpuArray<amrex::Real,2> orig_y_sub;
415  amrex::Real* orig_x_sub_p = orig_x_sub.data();
416  amrex::Real* orig_y_sub_p = orig_y_sub.data();
417  if (exp_interp) {
418  new_x = new_x_p[new_k];
419  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_p[k]; }
420  } else {
421  new_x = new_x_z[new_k];
422  for (int k=ksta; k <= kend; k++) { orig_x_sub_p[k-ksta] = orig_x_z[k]; }
423  }
424  for (int k=ksta; k <= kend; k++) { orig_y_sub_p[k-ksta] = orig_y[k]; }
425 #ifndef AMREX_USE_GPU
426  if (debug) {
427  amrex::Print() << " orig_x_sub_p = [";
428  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_x_sub_p[k];
429  amrex::Print() << "]" << std::endl;
430  amrex::Print() << " orig_y_sub_p = [";
431  for (int k=0; k < ksize; k++) amrex::Print() << " " << orig_y_sub_p[k];
432  amrex::Print() << "]" << std::endl;
433  }
434 #endif
435  lagrange_interp(1, orig_x_sub_p, orig_y_sub_p, new_x, new_y[new_k]);
436  }
437  }
438 #ifndef AMREX_USE_GPU
439  if (debug) amrex::Print() << " new_y[" << new_k << "]=" << new_y[new_k] << std::endl;
440 #endif
441  }
442 }
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:113

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