ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
MOSTAverage Class Reference

#include <ERF_MOSTAverage.H>

Collaboration diagram for MOSTAverage:

Public Member Functions

 MOSTAverage (amrex::Vector< amrex::Geometry > geom, const bool &has_zphys, std::string a_pp_prefix, const MeshType &m_mesh_type, const TerrainType &m_terrain_type)
 
 ~MOSTAverage ()
 
 MOSTAverage (MOSTAverage &&) noexcept=default
 
MOSTAverageoperator= (MOSTAverage &&other) noexcept=delete
 
 MOSTAverage (const MOSTAverage &other)=delete
 
MOSTAverageoperator= (const MOSTAverage &other)=delete
 
void make_MOSTAverage_at_level (const int &lev, const amrex::Vector< amrex::MultiFab * > &vars_old, std::unique_ptr< amrex::MultiFab > &Theta_prim, std::unique_ptr< amrex::MultiFab > &Qv_prim, std::unique_ptr< amrex::MultiFab > &Qr_prim, std::unique_ptr< amrex::MultiFab > &z_phys_nd)
 
void update_field_ptrs (const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qr_prim)
 
void set_rotated_fields (const int &lev)
 
void set_plane_normalization (const int &lev)
 
void set_region_normalization (const int &)
 
void set_k_indices_N (const int &lev)
 
void set_k_indices_EB (const int &lev)
 
void set_k_indices_T (const int &lev)
 
void set_norm_indices_T (const int &lev)
 
void set_z_positions_T (const int &lev)
 
void set_norm_positions_T (const int &lev)
 
void compute_averages (const int &lev)
 
void compute_plane_averages (const int &lev)
 
void compute_region_averages (const int &lev)
 
void write_k_indices (const int &lev)
 
void write_norm_indices (const int &lev)
 
void write_xz_positions (const int &lev, const int &j)
 
void write_averages (const int &lev)
 
const amrex::MultiFab * get_average (const int &lev, const int &comp) const
 
amrex::MultiFab * get_zref (const int &lev) const
 
const amrex::iMultiFab * get_k_indices (const int &lev) const
 

Static Public Member Functions

AMREX_GPU_HOST_DEVICE static AMREX_INLINE void trilinear_interp_T (const amrex::Real &xp, const amrex::Real &yp, const amrex::Real &zp, amrex::Real *interp_vals, amrex::Array4< amrex::Real const > const &interp_array, amrex::Array4< amrex::Real const > const &z_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &plo, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxi, const int interp_comp)
 

Protected Attributes

const amrex::Vector< amrex::Geometry > m_geom
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_fields
 
amrex::Vector< amrex::MultiFab * > m_z_phys_nd
 
std::string m_pp_prefix
 
MeshType m_mesh_type
 
TerrainType m_terrain_type
 
int m_nvar {6}
 
int m_navg {6}
 
int m_maxlev {0}
 
int m_policy {0}
 
bool m_rotate {false}
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_zref
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_x_pos
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_y_pos
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_z_pos
 
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_i_indx
 
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_j_indx
 
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_k_indx
 
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_averages
 
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_rot_fields
 
amrex::Vector< amrex::Vector< int > > m_ncell_plane
 
amrex::Vector< amrex::Vector< amrex::Real > > m_plane_average
 
int m_radius {0}
 
int m_ncell_region {1}
 
amrex::Vector< int > m_k_in
 
bool m_interp {false}
 
bool m_norm_vec {false}
 
bool m_t_avg {false}
 
amrex::Vector< int > m_t_init
 
amrex::Real m_time_window {amrex::Real(1.0e-16)}
 
amrex::Real m_fact_new
 
amrex::Real m_fact_old
 
bool include_subgrid_vel = false
 
amrex::Vector< amrex::Realm_Vsg
 
const amrex::Real zref_default = amrex::Real(10.0)
 

Constructor & Destructor Documentation

◆ MOSTAverage() [1/3]

MOSTAverage::MOSTAverage ( amrex::Vector< amrex::Geometry >  geom,
const bool &  has_zphys,
std::string  a_pp_prefix,
const MeshType &  m_mesh_type,
const TerrainType &  m_terrain_type 
)
explicit

◆ ~MOSTAverage()

MOSTAverage::~MOSTAverage ( )
inline
24  {}

◆ MOSTAverage() [2/3]

MOSTAverage::MOSTAverage ( MOSTAverage &&  )
defaultnoexcept

◆ MOSTAverage() [3/3]

MOSTAverage::MOSTAverage ( const MOSTAverage other)
delete

Member Function Documentation

◆ compute_averages()

void MOSTAverage::compute_averages ( const int &  lev)

Function to call the type of average computation.

Parameters
[in]levCurrent level
837 {
838  if (m_rotate) set_rotated_fields(lev);
839 
840  switch(m_policy) {
841  case 0: // Standard plane average
843  break;
844  case 1: // Local region/point
846  break;
847  default:
848  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false, "Unknown policy for MOSTAverage!");
849  }
850 
851  // We have initialized the averages
852  if (m_t_avg) m_t_init[lev] = 1;
853 }
bool m_t_avg
Definition: ERF_MOSTAverage.H:237
void compute_plane_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:862
int m_policy
Definition: ERF_MOSTAverage.H:207
void compute_region_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:1161
amrex::Vector< int > m_t_init
Definition: ERF_MOSTAverage.H:238
void set_rotated_fields(const int &lev)
Definition: ERF_MOSTAverage.cpp:288
bool m_rotate
Definition: ERF_MOSTAverage.H:208
Here is the call graph for this function:

◆ compute_plane_averages()

void MOSTAverage::compute_plane_averages ( const int &  lev)

Function to compute average over a plane.

Parameters
[in]levCurrent level
863 {
864  // Peel back the level
865  auto& fields = m_fields[lev];
866  auto& rot_fields = m_rot_fields[lev];
867  auto& averages = m_averages[lev];
868  const auto & geom = m_geom[lev];
869 
870  auto& z_phys = m_z_phys_nd[lev];
871  auto& x_pos = m_x_pos[lev];
872  auto& y_pos = m_y_pos[lev];
873  auto& z_pos = m_z_pos[lev];
874 
875  auto& i_indx = m_i_indx[lev];
876  auto& j_indx = m_j_indx[lev];
877  auto& k_indx = m_k_indx[lev];
878 
879  auto& ncell_plane = m_ncell_plane[lev];
880  auto& plane_average = m_plane_average[lev];
881 
882  // Set factors for time averaging
883  Real d_fact_new, d_fact_old;
884  if (m_t_avg && m_t_init[lev]) {
885  d_fact_new = m_fact_new;
886  d_fact_old = m_fact_old;
887  } else {
888  d_fact_new = one;
889  d_fact_old = zero;
890  }
891 
892  int klo = m_geom[lev].Domain().smallEnd(2);
893 
894  // GPU array to accumulate averages into
895  Gpu::DeviceVector<Real> pavg(plane_average.size(), zero);
896  Real* plane_avg = pavg.data();
897 
898  // Vectors for normalization and buffer storage
899  Vector<Real> denom(plane_average.size(),zero);
900  Vector<Real> val_old(plane_average.size(),zero);
901 
902  //
903  //----------------------------------------------------------
904  // Averages over all the fields
905  //----------------------------------------------------------
906  //
907  Box domain = geom.Domain();
908 
909  Array<int,AMREX_SPACEDIM> is_per = {0,0,0};
910  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
911  if (geom.isPeriodic(idim)) is_per[idim] = 1;
912  }
913 
914  // Averages for U,V,T,Qv (not Qc or W)
915  for (int imf(0); imf < 4; ++imf) {
916 
917  // Continue if no valid Qv pointer
918  if (!fields[imf]) continue;
919 
920  denom[imf] = one / (Real)ncell_plane[imf];
921  val_old[imf] = plane_average[imf]*d_fact_old;
922 
923 #ifdef _OPENMP
924 #pragma omp parallel if (Gpu::notInLaunchRegion())
925 #endif
926  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
927  Box vbx = mfi.validbox(); // This is the grid (not tile)
928  Box pbx = mfi.tilebox(); // This is the tile (not grid)
929 
930  if (pbx.smallEnd(2) != klo) { continue; }
931 
932  // Make planar since mfiter is over fields
933  pbx.makeSlab(2,klo);
934 
935  // Avoid double counting nodal data by changing the high end when we are
936  // at the high side of the grid (not just of the tile)
937  IndexType ixt = averages[imf]->boxArray().ixType();
938  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
939  if ( ixt.nodeCentered(idim) && (pbx.bigEnd(idim) == vbx.bigEnd(idim)) ) {
940  int dom_hi = domain.bigEnd(idim)+1;
941  if (pbx.bigEnd(idim) < dom_hi || is_per[idim]) {
942  pbx.growHi(idim,-1);
943  }
944  }
945  }
946 
947  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
948  fields[imf]->const_array(mfi);
949 
950  if (m_interp) {
951  const auto plo = geom.ProbLoArray();
952  const auto dxInv = geom.InvCellSizeArray();
953  const auto z_phys_arr = z_phys->const_array(mfi);
954  auto x_pos_arr = x_pos->array(mfi);
955  auto y_pos_arr = y_pos->array(mfi);
956  auto z_pos_arr = z_pos->array(mfi);
957  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
958  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
959  {
960  Real interp{0};
961  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
962  &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
963  Real val = interp;
964  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
965  });
966  } else {
967  auto k_arr = k_indx->const_array(mfi);
968  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
969  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
970  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
971  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
972  {
973  int mk = k_arr(i,j,0);
974  int mj = j_arr ? j_arr(i,j,0) : j;
975  int mi = i_arr ? i_arr(i,j,0) : i;
976  Real val = mf_arr(mi,mj,mk);
977  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
978  });
979  }
980  }
981  }
982 
983  //
984  //------------------------------------------------------------------------
985  // Averages for virtual potential temperature
986  // (This is cell-centered so we don't need to worry about double-counting)
987  //------------------------------------------------------------------------
988  //
989  if (fields[3]) // We have water vapor
990  {
991  int iavg = 4;
992  denom[iavg] = one / (Real)ncell_plane[iavg];
993  val_old[iavg] = plane_average[iavg]*d_fact_old;
994 
995 #ifdef _OPENMP
996 #pragma omp parallel if (Gpu::notInLaunchRegion())
997 #endif
998  for (MFIter mfi(*fields[3], TileNoZ()); mfi.isValid(); ++mfi)
999  {
1000  Box pbx = mfi.tilebox();
1001 
1002  if (pbx.smallEnd(2) != klo) { continue; }
1003 
1004  pbx.makeSlab(2,klo);
1005 
1006  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
1007  const Array4<Real const>& qv_mf_arr = fields[3]->const_array(mfi);
1008  const Array4<Real const>& qr_mf_arr = (fields[4]) ? fields[4]->const_array(mfi) :
1009  Array4<const Real> {};
1010 
1011  if (m_interp) {
1012  const auto plo = m_geom[lev].ProbLoArray();
1013  const auto dxInv = m_geom[lev].InvCellSizeArray();
1014  const auto z_phys_arr = z_phys->const_array(mfi);
1015  auto x_pos_arr = x_pos->array(mfi);
1016  auto y_pos_arr = y_pos->array(mfi);
1017  auto z_pos_arr = z_pos->array(mfi);
1018  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
1019  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
1020  {
1021  Real T_interp{0};
1022  Real qv_interp{0};
1023  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
1024  &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
1025  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
1026  &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
1027  Real vfac;
1028  if (qr_mf_arr) {
1029  // We also have liquid water
1030  Real qr_interp{0};
1031  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
1032  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
1033  vfac = one + Real(0.61)*qv_interp - qr_interp;
1034  } else {
1035  vfac = one + Real(0.61)*qv_interp;
1036  }
1037  const Real val = T_interp * vfac;
1038  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1039  });
1040  } else {
1041  auto k_arr = k_indx->const_array(mfi);
1042  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1043  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1044  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
1045  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
1046  {
1047  int mk = k_arr(i,j,0);
1048  int mj = j_arr ? j_arr(i,j,0) : j;
1049  int mi = i_arr ? i_arr(i,j,0) : i;
1050  Real vfac;
1051  if (qr_mf_arr) {
1052  // We also have liquid water
1053  vfac = one + Real(0.61)*qv_mf_arr(mi,mj,mk) - qr_mf_arr(mi,mj,mk);
1054  } else {
1055  vfac = one + Real(0.61)*qv_mf_arr(mi,mj,mk);
1056  }
1057  const Real val = T_mf_arr(mi,mj,mk) * vfac;
1058  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1059  });
1060  }
1061  }
1062  }
1063  else // copy temperature
1064  {
1065  int iavg = m_navg - 2;
1066  denom[iavg] = one / (Real)ncell_plane[iavg];
1067  // plane_avg[iavg] = plane_avg[2]
1068  Gpu::copy(Gpu::deviceToDevice, pavg.begin() + 2, pavg.begin() + 3,
1069  pavg.begin() + iavg);
1070  }
1071 
1072  //
1073  //------------------------------------------------------------------------
1074  // Averages for the tangential velocity magnitude
1075  // (This is cell-centered so we don't need to worry about double-counting)
1076  //------------------------------------------------------------------------
1077  //
1078  {
1079  int imf_cc = 2;
1080  int imf = 0;
1081  int iavg = m_navg - 1;
1082  denom[iavg] = one / (Real)ncell_plane[iavg];
1083  val_old[iavg] = plane_average[iavg]*d_fact_old;
1084 
1085  const Real Vsg = m_Vsg[lev];
1086 
1087 #ifdef _OPENMP
1088 #pragma omp parallel if (Gpu::notInLaunchRegion())
1089 #endif
1090  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi)
1091  {
1092  Box pbx = mfi.tilebox();
1093 
1094  if (pbx.smallEnd(2) != klo) { continue; }
1095 
1096  pbx.makeSlab(2,klo);
1097 
1098  // Last element is Umag and always cell centered
1099  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
1100  fields[imf ]->const_array(mfi);
1101  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
1102  fields[imf+1]->const_array(mfi);
1103 
1104  if (m_interp) {
1105  const auto plo = m_geom[lev].ProbLoArray();
1106  const auto dxInv = m_geom[lev].InvCellSizeArray();
1107  const auto z_phys_arr = z_phys->const_array(mfi);
1108  auto x_pos_arr = x_pos->array(mfi);
1109  auto y_pos_arr = y_pos->array(mfi);
1110  auto z_pos_arr = z_pos->array(mfi);
1111  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
1112  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
1113  {
1114  Real u_interp{0};
1115  Real v_interp{0};
1116  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
1117  &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
1118  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
1119  &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
1120  const Real val = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
1121  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1122  });
1123  } else {
1124  auto k_arr = k_indx->const_array(mfi);
1125  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1126  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1127  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
1128  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
1129  {
1130  int mk = k_arr(i,j,0);
1131  int mj = j_arr ? j_arr(i,j,0) : j;
1132  int mi = i_arr ? i_arr(i,j,0) : i;
1133  const Real u_val = myhalf * (u_mf_arr(mi,mj,mk) + u_mf_arr(mi+1,mj ,mk));
1134  const Real v_val = myhalf * (v_mf_arr(mi,mj,mk) + v_mf_arr(mi ,mj+1,mk));
1135  const Real val = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
1136  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1137  });
1138  }
1139  }
1140  }
1141 
1142  // Copy to host and sum across procs
1143  Gpu::copy(Gpu::deviceToHost, pavg.begin(), pavg.end(), plane_average.begin());
1144  ParallelDescriptor::ReduceRealSum(plane_average.data(), plane_average.size());
1145 
1146  // No spatial variation with plane averages
1147  for (int iavg(0); iavg < m_navg; ++iavg){
1148  plane_average[iavg] *= denom[iavg]*d_fact_new;
1149  plane_average[iavg] += val_old[iavg];
1150  averages[iavg]->setVal(plane_average[iavg]);
1151  }
1152 }
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
Real vfac
Definition: ERF_InitCustomPertVels_ABL.H:26
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+myhalf) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=one) { Real dT=T_pert *(std::cos(PI *L)+one)/two;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
int m_navg
Definition: ERF_MOSTAverage.H:205
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_averages
Definition: ERF_MOSTAverage.H:216
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_y_pos
Definition: ERF_MOSTAverage.H:211
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_i_indx
Definition: ERF_MOSTAverage.H:213
amrex::Vector< amrex::MultiFab * > m_z_phys_nd
Definition: ERF_MOSTAverage.H:197
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_x_pos
Definition: ERF_MOSTAverage.H:210
amrex::Vector< amrex::Real > m_Vsg
Definition: ERF_MOSTAverage.H:245
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_rot_fields
Definition: ERF_MOSTAverage.H:217
amrex::Vector< amrex::Vector< amrex::Real > > m_plane_average
Definition: ERF_MOSTAverage.H:222
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_z_pos
Definition: ERF_MOSTAverage.H:212
amrex::Vector< amrex::Vector< int > > m_ncell_plane
Definition: ERF_MOSTAverage.H:221
amrex::Real m_fact_new
Definition: ERF_MOSTAverage.H:240
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_j_indx
Definition: ERF_MOSTAverage.H:214
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_fields
Definition: ERF_MOSTAverage.H:196
amrex::Real m_fact_old
Definition: ERF_MOSTAverage.H:240
bool m_interp
Definition: ERF_MOSTAverage.H:232
const amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_MOSTAverage.H:195
AMREX_GPU_HOST_DEVICE static AMREX_INLINE void trilinear_interp_T(const amrex::Real &xp, const amrex::Real &yp, const amrex::Real &zp, amrex::Real *interp_vals, amrex::Array4< amrex::Real const > const &interp_array, amrex::Array4< amrex::Real const > const &z_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &plo, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxi, const int interp_comp)
Definition: ERF_MOSTAverage.H:126
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_k_indx
Definition: ERF_MOSTAverage.H:215

Referenced by compute_averages().

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

◆ compute_region_averages()

void MOSTAverage::compute_region_averages ( const int &  lev)

Function to compute average over local region.

Parameters
[in]levCurrent level
1162 {
1163  // Peel back the level
1164  auto& fields = m_fields[lev];
1165  auto& rot_fields = m_rot_fields[lev];
1166  auto& averages = m_averages[lev];
1167  const auto & geom = m_geom[lev];
1168 
1169  auto& z_phys = m_z_phys_nd[lev];
1170  auto& x_pos = m_x_pos[lev];
1171  auto& y_pos = m_y_pos[lev];
1172  auto& z_pos = m_z_pos[lev];
1173 
1174  auto& i_indx = m_i_indx[lev];
1175  auto& j_indx = m_j_indx[lev];
1176  auto& k_indx = m_k_indx[lev];
1177 
1178  int klo = m_geom[lev].Domain().smallEnd(2);
1179 
1180  // Set factors for time averaging
1181  Real d_fact_new, d_fact_old;
1182  if (m_t_avg && m_t_init[lev]) {
1183  d_fact_new = m_fact_new;
1184  d_fact_old = m_fact_old;
1185  } else {
1186  d_fact_new = one;
1187  d_fact_old = zero;
1188  }
1189 
1190  // Number of cells contained in the local average
1191  const Real denom = one / (Real) m_ncell_region;
1192 
1193  // Capture radius for device
1194  int d_radius = m_radius;
1195 
1196  //
1197  //----------------------------------------------------------
1198  // Averages for U,V,T,Qv
1199  //----------------------------------------------------------
1200  //
1201  for (int imf(0); imf < 4; ++imf) {
1202 
1203  // Continue if no valid Qv pointer
1204  if (!fields[imf]) continue;
1205 
1206 #ifdef _OPENMP
1207 #pragma omp parallel if (Gpu::notInLaunchRegion())
1208 #endif
1209  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
1210  Box pbx = mfi.tilebox();
1211 
1212  if (pbx.smallEnd(2) != klo) { continue; }
1213 
1214  // Make planar since mfiter is over fields
1215  pbx.makeSlab(2,klo);
1216 
1217  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
1218  fields[imf]->const_array(mfi);
1219  auto ma_arr = averages[imf]->array(mfi);
1220 
1221  if (m_interp) {
1222  const auto plo = geom.ProbLoArray();
1223  const auto dx = geom.CellSizeArray();
1224  const auto dxInv = geom.InvCellSizeArray();
1225  const auto z_phys_arr = z_phys->const_array(mfi);
1226  auto x_pos_arr = x_pos->array(mfi);
1227  auto y_pos_arr = y_pos->array(mfi);
1228  auto z_pos_arr = z_pos->array(mfi);
1229  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1230  {
1231  ma_arr(i,j,0) *= d_fact_old;
1232 
1233  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1234  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1235  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1236  for (int li(-d_radius); li <= (d_radius); ++li) {
1237  Real interp{0};
1238  Real xp = x_pos_arr(i+li,j+lj,0);
1239  Real yp = y_pos_arr(i+li,j+lj,0);
1240  Real zp = z_pos_arr(i+li,j+lj,0) + met_h_zeta*lk*dx[2];
1241  trilinear_interp_T(xp, yp, zp, &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
1242  Real val = denom * interp * d_fact_new;
1243  ma_arr(i,j,0) += val;
1244  }
1245  }
1246  }
1247  });
1248  } else {
1249  auto k_arr = k_indx->const_array(mfi);
1250  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1251  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1252  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int ) noexcept
1253  {
1254  ma_arr(i,j,0) *= d_fact_old;
1255 
1256  int mk = k_arr(i,j,0);
1257  int mj = j_arr ? j_arr(i,j,0) : j;
1258  int mi = i_arr ? i_arr(i,j,0) : i;
1259  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1260  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1261  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1262  Real val = denom * mf_arr(li, lj, lk) * d_fact_new;
1263  ma_arr(i,j,0) += val;
1264  }
1265  }
1266  }
1267  });
1268  }
1269  } // MFiter
1270 
1271  // Fill interior ghost cells and any ghost cells outside a periodic domain
1272  //***********************************************************************************
1273  averages[imf]->FillBoundary(geom.periodicity());
1274 
1275  } // imf
1276 
1277  //
1278  //----------------------------------------------------------
1279  // Averages for virtual potential temperature
1280  //----------------------------------------------------------
1281  //
1282  if (fields[3]) // We have water vapor
1283  {
1284  int iavg = 4;
1285 
1286 #ifdef _OPENMP
1287 #pragma omp parallel if (Gpu::notInLaunchRegion())
1288 #endif
1289  for (MFIter mfi(*fields[3], TileNoZ()); mfi.isValid(); ++mfi) {
1290  Box pbx = mfi.tilebox();
1291 
1292  if (pbx.smallEnd(2) != klo) { continue; }
1293 
1294  pbx.makeSlab(2,klo);
1295 
1296  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
1297  const Array4<Real const>& qv_mf_arr = fields[3]->const_array(mfi);
1298  const Array4<Real const>& qr_mf_arr = (fields[4]) ? fields[4]->const_array(mfi) :
1299  Array4<const Real> {};
1300  auto ma_arr = averages[iavg]->array(mfi);
1301 
1302  if (m_interp) {
1303  const auto plo = geom.ProbLoArray();
1304  const auto dx = geom.CellSizeArray();
1305  const auto dxInv = geom.InvCellSizeArray();
1306  const auto z_phys_arr = z_phys->const_array(mfi);
1307  auto x_pos_arr = x_pos->array(mfi);
1308  auto y_pos_arr = y_pos->array(mfi);
1309  auto z_pos_arr = z_pos->array(mfi);
1310  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1311  {
1312  ma_arr(i,j,0) *= d_fact_old;
1313 
1314  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1315  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1316  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1317  for (int li(-d_radius); li <= (d_radius); ++li) {
1318  Real T_interp{0};
1319  Real qv_interp{0};
1320  Real xp = x_pos_arr(i+li,j+lj,0);
1321  Real yp = y_pos_arr(i+li,j+lj,0);
1322  Real zp = z_pos_arr(i+li,j+lj,0) + met_h_zeta*lk*dx[2];
1323  trilinear_interp_T(xp, yp, zp, &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
1324  trilinear_interp_T(xp, yp, zp, &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
1325  Real vfac;
1326  if (qr_mf_arr) {
1327  // We also have liquid water
1328  Real qr_interp{0};
1329  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
1330  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
1331  vfac = one + Real(0.61)*qv_interp - qr_interp;
1332  } else {
1333  vfac = one + Real(0.61)*qv_interp;
1334  }
1335  const Real mag = T_interp * vfac;
1336  const Real val = denom * mag * d_fact_new;
1337  ma_arr(i,j,0) += val;
1338  }
1339  }
1340  }
1341  });
1342  } else {
1343  auto k_arr = k_indx->const_array(mfi);
1344  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1345  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1346  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int ) noexcept
1347  {
1348  ma_arr(i,j,0) *= d_fact_old;
1349 
1350  int mk = k_arr(i,j,0);
1351  int mj = j_arr ? j_arr(i,j,0) : j;
1352  int mi = i_arr ? i_arr(i,j,0) : i;
1353  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1354  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1355  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1356  Real vfac;
1357  if (qr_mf_arr) {
1358  // We also have liquid water
1359  vfac = one + Real(0.61)*qv_mf_arr(li,lj,lk) - qr_mf_arr(li,lj,lk);
1360  } else {
1361  vfac = one + Real(0.61)*qv_mf_arr(li,lj,lk);
1362  }
1363  const Real mag = T_mf_arr(li,lj,lk) * vfac;
1364  const Real val = denom * mag * d_fact_new;
1365  ma_arr(i,j,0) += val;
1366  }
1367  }
1368  }
1369  });
1370  }
1371  } // MFiter
1372 
1373  // Fill interior ghost cells and any ghost cells outside a periodic domain
1374  //***********************************************************************************
1375  averages[iavg]->FillBoundary(geom.periodicity());
1376 
1377  }
1378  else // copy temperature
1379  {
1380  int iavg = m_navg - 2;
1381  IntVect ng = averages[iavg]->nGrowVect();
1382  MultiFab::Copy(*(averages[iavg]),*(averages[2]),0,0,1,ng);
1383  }
1384 
1385  //
1386  //----------------------------------------------------------
1387  // Averages for the tangential velocity magnitude
1388  //----------------------------------------------------------
1389  //
1390  {
1391  int imf_cc = 2;
1392  int imf = 0;
1393  int iavg = m_navg - 1;
1394 
1395  const Real Vsg = m_Vsg[lev];
1396 
1397 #ifdef _OPENMP
1398 #pragma omp parallel if (Gpu::notInLaunchRegion())
1399 #endif
1400  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1401  Box pbx = mfi.tilebox();
1402 
1403  if (pbx.smallEnd(2) != klo) { continue; }
1404 
1405  pbx.makeSlab(2,klo);
1406 
1407  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
1408  fields[imf ]->const_array(mfi);
1409  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
1410  fields[imf+1]->const_array(mfi);
1411  auto ma_arr = averages[iavg]->array(mfi);
1412 
1413  if (m_interp) {
1414  const auto plo = geom.ProbLoArray();
1415  const auto dx = geom.CellSizeArray();
1416  const auto dxInv = geom.InvCellSizeArray();
1417  const auto z_phys_arr = z_phys->const_array(mfi);
1418  auto x_pos_arr = x_pos->array(mfi);
1419  auto y_pos_arr = y_pos->array(mfi);
1420  auto z_pos_arr = z_pos->array(mfi);
1421  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1422  {
1423  ma_arr(i,j,0) *= d_fact_old;
1424 
1425  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1426  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1427  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1428  for (int li(-d_radius); li <= (d_radius); ++li) {
1429  Real u_interp{0};
1430  Real v_interp{0};
1431  Real xp = x_pos_arr(i+li,j+lj,0);
1432  Real yp = y_pos_arr(i+li,j+lj,0);
1433  Real zp = z_pos_arr(i+li,j+lj,0) + met_h_zeta*lk*dx[2];
1434  trilinear_interp_T(xp, yp, zp, &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
1435  trilinear_interp_T(xp, yp, zp, &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
1436  const Real mag = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
1437  Real val = denom * mag * d_fact_new;
1438  ma_arr(i,j,0) += val;
1439  }
1440  }
1441  }
1442  });
1443  } else {
1444  auto k_arr = k_indx->const_array(mfi);
1445  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1446  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1447  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int ) noexcept
1448  {
1449  ma_arr(i,j,0) *= d_fact_old;
1450 
1451  int mk = k_arr(i,j,0);
1452  int mj = j_arr ? j_arr(i,j,0) : j;
1453  int mi = i_arr ? i_arr(i,j,0) : i;
1454  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1455  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1456  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1457  const Real u_val = myhalf * (u_mf_arr(li,lj,lk) + u_mf_arr(li+1,lj ,lk));
1458  const Real v_val = myhalf * (v_mf_arr(li,lj,lk) + v_mf_arr(li ,lj+1,lk));
1459  const Real mag = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
1460  Real val = denom * mag * d_fact_new;
1461  ma_arr(i,j,0) += val;
1462  }
1463  }
1464  }
1465  });
1466  }
1467  } // MFiter
1468 
1469  // Fill interior ghost cells and any ghost cells outside a periodic domain
1470  //***********************************************************************************
1471  averages[iavg]->FillBoundary(geom.periodicity());
1472 
1473  }
1474 
1475  // NOTE: Checking periodicity with the geom structure is not
1476  // sufficient at higher levels. The BA may be contained
1477  // within the domain and it's exterior ghost cells filled
1478  // from interpolation; yet the domain BCs are periodic.
1479 
1480  // Need to fill ghost cells outside the domain if not periodic
1481  bool not_per_x = !(geom.periodicity().isPeriodic(0));
1482  bool not_per_y = !(geom.periodicity().isPeriodic(1));
1483  Box cc_bnd_bx = (m_fields[lev][2]->boxArray()).minimalBox();
1484  Box domain = geom.Domain();
1485  if (domain.contains(cc_bnd_bx) || (not_per_x || not_per_y)) {
1486  for (int iavg(0); iavg < m_navg; ++iavg) {
1487  IntVect ng = averages[iavg]->nGrowVect(); ng[2]=0;
1488 
1489  // NOTE: Level 0 spans the whole domain, but finer
1490  // levels do not have such a restriction.
1491  // For now, use the bounding box of the boxArray.
1492 
1493  // NOTE2: The fields and averages have different indexing.
1494  // The averages are: U/V/T/Qv/Tv/Umag
1495  // The fields are: U/V/T/Qv/Qr/W
1496  // We clip iavg at 2 since all the remaining data is CC
1497 
1498  // Bounded box of CC data used for normalization
1499  int imf = min(iavg,2);
1500  Box bnd_bx = (fields[imf]->boxArray()).minimalBox();
1501 #ifdef _OPENMP
1502 #pragma omp parallel if (Gpu::notInLaunchRegion())
1503 #endif
1504  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
1505  Box gpbx = mfi.growntilebox(ng);
1506 
1507  if (gpbx.smallEnd(2) != klo) { continue; }
1508 
1509  gpbx.makeSlab(2,klo);
1510 
1511  if (bnd_bx.contains(gpbx)) continue;
1512 
1513  auto ma_arr = averages[iavg]->array(mfi);
1514 
1515  int i_lo = bnd_bx.smallEnd(0); int i_hi = bnd_bx.bigEnd(0);
1516  int j_lo = bnd_bx.smallEnd(1); int j_hi = bnd_bx.bigEnd(1);
1517  ParallelFor(gpbx, [=] AMREX_GPU_DEVICE(int i, int j, int ) noexcept
1518  {
1519  int li, lj;
1520  li = i < i_lo ? i_lo : i;
1521  li = li > i_hi ? i_hi : li;
1522  lj = j < j_lo ? j_lo : j;
1523  lj = lj > j_hi ? j_hi : lj;
1524 
1525  ma_arr(i,j,0) = ma_arr(li,lj,0);
1526  });
1527  } // MFiter
1528  } // iavg
1529  } // Not periodic
1530 }
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
AMREX_FORCE_INLINE AMREX_GPU_DEVICE amrex::Real Compute_h_zeta_AtCellCenter(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:55
int m_radius
Definition: ERF_MOSTAverage.H:226
int m_ncell_region
Definition: ERF_MOSTAverage.H:227
@ ng
Definition: ERF_Morrison.H:48

Referenced by compute_averages().

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

◆ get_average()

const amrex::MultiFab* MOSTAverage::get_average ( const int &  lev,
const int &  comp 
) const
inline
105 { return m_averages[lev][comp].get(); }

Referenced by SurfaceLayer::get_mac_avg().

Here is the caller graph for this function:

◆ get_k_indices()

const amrex::iMultiFab* MOSTAverage::get_k_indices ( const int &  lev) const
inline
111 { return m_k_indx[lev].get(); }

◆ get_zref()

amrex::MultiFab* MOSTAverage::get_zref ( const int &  lev) const
inline
108 { return m_zref[lev].get(); }
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_zref
Definition: ERF_MOSTAverage.H:209

Referenced by SurfaceLayer::get_zref().

Here is the caller graph for this function:

◆ make_MOSTAverage_at_level()

void MOSTAverage::make_MOSTAverage_at_level ( const int &  lev,
const amrex::Vector< amrex::MultiFab * > &  vars_old,
std::unique_ptr< amrex::MultiFab > &  Theta_prim,
std::unique_ptr< amrex::MultiFab > &  Qv_prim,
std::unique_ptr< amrex::MultiFab > &  Qr_prim,
std::unique_ptr< amrex::MultiFab > &  z_phys_nd 
)
91 {
92  m_fields[lev].resize(m_nvar);
93  m_rot_fields[lev].resize(m_nvar-1);
94  m_averages[lev].resize(m_navg);
95  m_z_phys_nd[lev] = z_phys_nd.get();
96 
97  bool use_terrain_fitted_coords = ( (m_terrain_type == TerrainType::StaticFittedMesh) ||
98  (m_terrain_type == TerrainType::MovingFittedMesh) );
99 
100  bool use_eb = (m_terrain_type == TerrainType::EB);
101 
102  { // Nodal in x
103  auto& mf = *vars_old[Vars::xvel];
104  // Create a 2D ba, dm, & ghost cells
105  const BoxArray& ba = mf.boxArray();
106  BoxList bl2d = ba.boxList();
107  for (auto& b : bl2d) { b.setRange(2,0); }
108  BoxArray ba2d(std::move(bl2d));
109  const DistributionMapping& dm = mf.DistributionMap();
110  const int ncomp = 1;
111  IntVect ng = mf.nGrowVect(); ng[2]=0;
112 
113  m_fields[lev][0] = vars_old[Vars::xvel];
114  m_averages[lev][0] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
115  m_averages[lev][0]->setVal(1.E34);
116  if (m_rotate) {
117  m_rot_fields[lev][0] = std::make_unique<MultiFab>(ba,dm,ncomp,ng);
118  MultiFab::Copy(*m_rot_fields[lev][0],mf,0,0,1,ng);
119  } else {
120  m_rot_fields[lev][0] = nullptr;
121  }
122  }
123  { // Nodal in y
124  auto& mf = *vars_old[Vars::yvel];
125  // Create a 2D ba, dm, & ghost cells
126  const BoxArray& ba = mf.boxArray();
127  BoxList bl2d = ba.boxList();
128  for (auto& b : bl2d) { b.setRange(2,0); }
129  BoxArray ba2d(std::move(bl2d));
130  const DistributionMapping& dm = mf.DistributionMap();
131  const int ncomp = 1;
132  IntVect ng = mf.nGrowVect(); ng[2]=0;
133 
134  m_fields[lev][1] = vars_old[Vars::yvel];
135  m_averages[lev][1] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
136  m_averages[lev][1]->setVal(1.E34);
137  if (m_rotate) {
138  m_rot_fields[lev][1] = std::make_unique<MultiFab>(ba,dm,ncomp,ng);
139  MultiFab::Copy(*m_rot_fields[lev][1],mf,0,0,1,ng);
140  } else {
141  m_rot_fields[lev][1] = nullptr;
142  }
143  }
144  { // CC vars
145  auto& mf = *Theta_prim;
146  // Create a 2D ba, dm, & ghost cells
147  const BoxArray& ba = mf.boxArray();
148  BoxList bl2d = ba.boxList();
149  for (auto& b : bl2d) { b.setRange(2,0); }
150  BoxArray ba2d(std::move(bl2d));
151  const DistributionMapping& dm = mf.DistributionMap();
152  const int ncomp = 1;
153  const int incomp = 1;
154  IntVect ng = mf.nGrowVect(); ng[2]=0;
155 
156  // Get field pointers
157  m_fields[lev][2] = Theta_prim.get();
158  m_fields[lev][3] = Qv_prim.get();
159  m_fields[lev][4] = Qr_prim.get();
160 
161  // Initialize remaining multifabs
162  for (int iavg(2); iavg < m_navg; ++iavg) {
163  m_averages[lev][iavg] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
164  m_averages[lev][iavg]->setVal(1.E34);
165  }
166 
167  // Default to dry
168  m_averages[lev][3]->setVal(0.0);
169 
170  if (m_rotate) {
171  m_rot_fields[lev][2] = std::make_unique<MultiFab>(ba,dm,ncomp,ng);
172  m_rot_fields[lev][3] = std::make_unique<MultiFab>(ba,dm,ncomp,ng);
173  MultiFab::Copy(*m_rot_fields[lev][2],*Theta_prim,0,0,1,ng);
174  if (Qv_prim) MultiFab::Copy(*m_rot_fields[lev][3],*Qv_prim,0,0,1,ng);
175  } else {
176  m_rot_fields[lev][2] = nullptr;
177  m_rot_fields[lev][3] = nullptr;
178  }
179 
180  // Default zref to 10 and fill will true values later
181  m_zref[lev] = std::make_unique<MultiFab>(ba2d,dm,1,ng);
182  m_zref[lev]->setVal(zref_default);
183 
184  if (use_terrain_fitted_coords && m_norm_vec && m_interp) {
185  m_x_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
186  m_y_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
187  m_z_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
188  } else if (use_terrain_fitted_coords && m_interp) {
189  m_x_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
190  m_y_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
191  m_z_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
192  } else if (use_terrain_fitted_coords && m_norm_vec) {
193  m_i_indx[lev] = std::make_unique<iMultiFab>(ba2d,dm,incomp,ng);
194  m_j_indx[lev] = std::make_unique<iMultiFab>(ba2d,dm,incomp,ng);
195  m_k_indx[lev] = std::make_unique<iMultiFab>(ba2d,dm,incomp,ng);
196  } else {
197  m_k_indx[lev] = std::make_unique<iMultiFab>(ba2d,dm,incomp,ng);
198  }
199  }
200  // Nodal in z (only used with terrain stress rotations)
201  m_fields[lev][5] = vars_old[Vars::zvel];
202 
203  // Setup auxiliary data for spatial configuration & policy
204  //--------------------------------------------------------
205  if (use_terrain_fitted_coords && m_norm_vec && m_interp) { // Terrain w/ norm & w/ interpolation
207  } else if (use_terrain_fitted_coords && m_interp) { // Terrain w/ interpolation
208  set_z_positions_T(lev);
209  } else if (use_terrain_fitted_coords && m_norm_vec) { // Terrain w/ norm & w/o interpolation
210  set_norm_indices_T(lev);
211  } else if (use_terrain_fitted_coords) { // Terrain
212  set_k_indices_T(lev);
213  } else if (use_eb) { // EB
214  set_k_indices_EB(lev);
215  } else { // No Terrain
216  set_k_indices_N(lev);
217  }
218 
219  // Setup normalization data for the chosen policy
220  //--------------------------------------------------------
221  switch(m_policy) {
222  case 0: // Plane average
224  break;
225  case 1: // Local region/point
227  break;
228  default:
229  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false, "Unknown policy for MOSTAverage!");
230  }
231 
232  // Set up the exponential time filtering
233  //--------------------------------------------------------
234  if (m_t_avg) {
235  // Exponential filter function
236  m_fact_old = std::exp(-one / m_time_window);
237 
238  // Enforce discrete normalization: (mfn*val_new + mfo*val_old)
240 
241  // None of the averages are initialized
242  m_t_init.resize(m_maxlev,0);
243  }
244 
245  // Corrections to the mean surface velocity
246  m_Vsg = Vector<Real>(m_maxlev, zero);
247  if (include_subgrid_vel) {
248  if (include_subgrid_vel) {
249  Print() << "Subgrid velocity scale correction at level : " << lev << ' ';
250  const auto dxArr = m_geom[lev].CellSizeArray();
251  Real dx = std::sqrt(dxArr[0]*dxArr[1]);
252  if (dx > Real(5000.)) {
253  m_Vsg[lev] = Real(0.32) * std::pow(dx/Real(5000.)-1, Real(0.33));
254  }
255  Print() << m_Vsg[lev] << std::endl;
256  }
257  }
258 }
void set_region_normalization(const int &)
Definition: ERF_MOSTAverage.H:61
void set_z_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:684
void set_norm_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:747
void set_k_indices_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:516
void set_plane_normalization(const int &lev)
Definition: ERF_MOSTAverage.cpp:347
TerrainType m_terrain_type
Definition: ERF_MOSTAverage.H:200
bool m_norm_vec
Definition: ERF_MOSTAverage.H:233
int m_nvar
Definition: ERF_MOSTAverage.H:204
void set_k_indices_EB(const int &lev)
Definition: ERF_MOSTAverage.cpp:461
void set_norm_indices_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:596
int m_maxlev
Definition: ERF_MOSTAverage.H:206
bool include_subgrid_vel
Definition: ERF_MOSTAverage.H:244
amrex::Real m_time_window
Definition: ERF_MOSTAverage.H:239
void set_k_indices_N(const int &lev)
Definition: ERF_MOSTAverage.cpp:405
const amrex::Real zref_default
Definition: ERF_MOSTAverage.H:249
@ xvel
Definition: ERF_IndexDefines.H:159
@ zvel
Definition: ERF_IndexDefines.H:161
@ yvel
Definition: ERF_IndexDefines.H:160

Referenced by SurfaceLayer::make_SurfaceLayer_at_level().

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

◆ operator=() [1/2]

MOSTAverage& MOSTAverage::operator= ( const MOSTAverage other)
delete

◆ operator=() [2/2]

MOSTAverage& MOSTAverage::operator= ( MOSTAverage &&  other)
deletenoexcept

◆ set_k_indices_EB()

void MOSTAverage::set_k_indices_EB ( const int &  lev)

Function to set K indices for EB.

462 {
463  ParmParse pp(m_pp_prefix);
464  Real zref_tmp = zref_default;
465  auto read_z = pp.query("most.zref",zref_tmp);
466  auto read_k = pp.queryarr("most.k_arr_in",m_k_in);
467 
468  ParmParse pp_eb2("eb2");
469  std::string geometry;
470  pp_eb2.query("geometry", geometry);
471  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(geometry == "plane", "Only plane geometry supported for EB MOST");
472 
473  RealArray plane_point{zero, zero, zero};
474  pp_eb2.query("plane_point", plane_point);
475  Real z_eb = plane_point[2];
476 
477  // Default behavior is to use the first cell center
478  if (!read_z && !read_k) {
479  Real m_dz = m_geom[0].CellSize(2);
480  zref_tmp = myhalf * m_dz;
481  m_zref[lev]->setVal( zref_tmp );
482  Print() << "Reference height for MOST set to " << zref_tmp << std::endl;
483  read_z = true;
484  }
485 
486  // Specify z_ref & compute k_indx (z_ref takes precedence)
487  if (read_z) {
488  Real m_zlo = m_geom[lev].ProbLo(2);
489  Real m_zhi = m_geom[lev].ProbHi(2);
490  Real m_dz = m_geom[lev].CellSize(2);
491 
492  amrex::ignore_unused(m_zhi);
493 
494  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(zref_tmp >= m_zlo + myhalf * m_dz,
495  "Query point must be past first z-cell!");
496 
497  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(zref_tmp <= m_zhi - myhalf * m_dz,
498  "Query point must be below the last z-cell!");
499 
500  int lk = static_cast<int>(floor((zref_tmp + z_eb - m_zlo) / m_dz - myhalf));
501  int lk_phys = static_cast<int>(floor((zref_tmp - m_zlo) / m_dz - myhalf));
502 
503  m_zref[lev]->setVal( (lk_phys + 0.5) * m_dz + m_zlo );
504 
506 
507  m_k_indx[lev]->setVal(lk);
508  }
509 }
ParmParse pp("prob")
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
std::string m_pp_prefix
Definition: ERF_MOSTAverage.H:198
amrex::Vector< int > m_k_in
Definition: ERF_MOSTAverage.H:228

Referenced by make_MOSTAverage_at_level().

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

◆ set_k_indices_N()

void MOSTAverage::set_k_indices_N ( const int &  lev)

Function to set K indices without terrain.

406 {
407  ParmParse pp(m_pp_prefix);
408  Real zref_tmp = zref_default;
409  auto read_z = pp.query("most.zref",zref_tmp);
410  auto read_k = pp.queryarr("most.k_arr_in",m_k_in);
411 
412  // Default behavior is to use the first cell center
413  if (!read_z && !read_k) {
414  Real m_zlo = m_geom[0].ProbLo(2);
415  Real m_dz = m_geom[0].CellSize(2);
416  zref_tmp = m_zlo + myhalf * m_dz;
417  m_zref[lev]->setVal( zref_tmp );
418  Print() << "Reference height for MOST set to " << zref_tmp << std::endl;
419  read_z = true;
420  }
421 
422  // Specify z_ref & compute k_indx (z_ref takes precedence)
423  if (read_z) {
424  Real m_zlo = m_geom[lev].ProbLo(2);
425  Real m_zhi = m_geom[lev].ProbHi(2);
426  Real m_dz = m_geom[lev].CellSize(2);
427 
428  amrex::ignore_unused(m_zhi);
429 
430  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(zref_tmp >= m_zlo + myhalf * m_dz,
431  "Query point must be past first z-cell!");
432 
433  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(zref_tmp <= m_zhi - myhalf * m_dz,
434  "Query point must be below the last z-cell!");
435 
436  int lk = static_cast<int>(floor((zref_tmp - m_zlo) / m_dz - myhalf));
437 
438  m_zref[lev]->setVal( (lk + 0.5) * m_dz + m_zlo );
439 
441 
442  m_k_indx[lev]->setVal(lk);
443  // Specified k_indx & compute z_ref
444  } else if (read_k) {
445  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(m_k_in[lev] >= m_radius,
446  "K index must be larger than averaging radius!");
447  m_k_indx[lev]->setVal(m_k_in[lev]);
448 
449  // TODO: check that z_ref is constant across levels
450  Real m_zlo = m_geom[0].ProbLo(2);
451  Real m_dz = m_geom[0].CellSize(2);
452  m_zref[lev]->setVal( ((Real)m_k_in[0] + 0.5) * m_dz + m_zlo );
453  }
454 }

Referenced by make_MOSTAverage_at_level().

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

◆ set_k_indices_T()

void MOSTAverage::set_k_indices_T ( const int &  lev)

Function to set K indices with terrain (w/o terrain normals or interpolation).

517 {
518  // Peel back the level
519  auto& fields = m_fields[lev];
520 
521  // MFIter over CC data
522  int imf_cc = 2;
523 
524  ParmParse pp(m_pp_prefix);
525  Real zref_tmp = zref_default;
526  auto read_z = pp.query("most.zref",zref_tmp);
527  auto read_k = pp.queryarr("most.k_arr_in",m_k_in);
528  int klo = m_geom[lev].Domain().smallEnd(2);
529 
530  // Allow default zref
531  if (!read_z) {
532  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
533  read_z = true;
534  }
535 
536  // No default behavior with terrain (we can't tell the difference between
537  // vertical grid stretching and true terrain)
538  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(read_z != read_k,
539  "Need to specify zref or k_arr_in for MOST");
540 
541  // Capture for device
542  Real d_zref = zref_tmp;
543  Real d_radius = m_radius;
544  amrex::ignore_unused(d_radius);
545 
546  // Specify z_ref & compute k_indx (z_ref takes precedence)
547  if (read_z) {
548  int kmax = m_geom[lev].Domain().bigEnd(2);
549  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
550  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
551 
552  if (npbx.smallEnd(2) != klo) { continue; }
553 
554  npbx.makeSlab(2,klo);
555 
556  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
557  auto k_arr = m_k_indx[lev]->array(mfi);
558  auto zref_arr = m_zref[lev]->array(mfi);
559  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
560  {
561  k_arr(i,j,k) = klo;
562  bool found = false;
563  Real z_bot_face = fourth * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
564  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
565  Real z_target = z_bot_face + d_zref;
566  for (int lk(klo); lk<=kmax; ++lk) {
567  Real z_lo = fourth * ( z_phys_arr(i,j ,lk ) + z_phys_arr(i+1,j ,lk )
568  + z_phys_arr(i,j+1,lk ) + z_phys_arr(i+1,j+1,lk ) );
569  Real z_hi = fourth * ( z_phys_arr(i,j ,lk+1) + z_phys_arr(i+1,j ,lk+1)
570  + z_phys_arr(i,j+1,lk+1) + z_phys_arr(i+1,j+1,lk+1) );
571  if (z_target > z_lo && z_target < z_hi){
572  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(lk >= d_radius,
573  "K index must be larger than averaging radius!");
574  k_arr(i,j,0) = lk;
575  zref_arr(i,j,0) = myhalf * (z_hi + z_lo) - z_bot_face;
576  found = true;
577  break;
578  }
579  }
580  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(found,
581  "zref not found with terrain!");
582  });
583  }
584  // Specified k_indx & compute z_ref
585  } else if (read_k) {
586  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false, "Specified k-indx with terrain not implemented!");
587  }
588 }
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12

Referenced by make_MOSTAverage_at_level().

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

◆ set_norm_indices_T()

void MOSTAverage::set_norm_indices_T ( const int &  lev)

Function to set I,J,K indices with terrain normals (w/o interpolation).

597 {
598  // Peel back the level
599  auto& fields = m_fields[lev];
600 
601  // MFIter over CC data
602  int imf_cc = 2;
603 
604  ParmParse pp(m_pp_prefix);
605  Real zref_tmp = zref_default;
606  auto read_zref = pp.query("most.zref",zref_tmp);
607  if (!read_zref) {
608  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
609  }
610  int klo = m_geom[lev].Domain().smallEnd(2);
611 
612  // Capture for device
613  Real d_zref = zref_tmp;
614  Real d_radius = m_radius;
615 
616  const auto dxInv = m_geom[lev].InvCellSizeArray();
617  IntVect ng = m_k_indx[lev]->nGrowVect(); ng[2]=0;
618  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
619  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
620 
621  if (npbx.smallEnd(2) != klo) { continue; }
622 
623  int kmax = npbx.bigEnd(2);
624 
625  npbx.makeSlab(2,klo);
626 
627  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
628  auto i_arr = m_i_indx[lev]->array(mfi);
629  auto j_arr = m_j_indx[lev]->array(mfi);
630  auto k_arr = m_k_indx[lev]->array(mfi);
631  auto zref_arr = m_zref[lev]->array(mfi);
632  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
633  {
634  // Elements of normal vector
635  Real met_h_xi = Compute_h_xi_AtCellCenter (i,j,k,dxInv,z_phys_arr);
636  Real met_h_eta = Compute_h_eta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
637  Real mag = std::sqrt(met_h_xi*met_h_xi + met_h_eta*met_h_eta + one);
638 
639  // Unit-normal vector scaled by z_ref
640  Real delta_x = -met_h_xi/mag * d_zref;
641  Real delta_y = -met_h_eta/mag * d_zref;
642  Real delta_z = one/mag * d_zref;
643 
644  // Compute i & j as displacements (no grid stretching)
645  int delta_i = static_cast<int>(std::round(delta_x*dxInv[0]));
646  int delta_j = static_cast<int>(std::round(delta_y*dxInv[1]));
647  int i_new = i + delta_i;
648  int j_new = j + delta_j;
649  i_arr(i,j,0) = i_new;
650  j_arr(i,j,0) = j_new;
651 
652  // Search for k (grid is stretched in z)
653  Real z_bot_face = fourth * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
654  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
655  Real z_target = z_bot_face + delta_z;
656  k_arr(i,j,0) = klo;
657  zref_arr(i,j,0) = myhalf * z_bot_face +
658  Real(0.125) * ( z_phys_arr(i ,j ,k+1) + z_phys_arr(i+1,j ,k+1)
659  + z_phys_arr(i ,j+1,k+1) + z_phys_arr(i+1,j+1,k+1) );
660  for (int lk(klo); lk<=kmax; ++lk) {
661  Real z_lo = fourth * ( z_phys_arr(i_new,j_new ,lk ) + z_phys_arr(i_new+1,j_new ,lk )
662  + z_phys_arr(i_new,j_new+1,lk ) + z_phys_arr(i_new+1,j_new+1,lk ) );
663  Real z_hi = fourth * ( z_phys_arr(i_new,j_new ,lk+1) + z_phys_arr(i_new+1,j_new ,lk+1)
664  + z_phys_arr(i_new,j_new+1,lk+1) + z_phys_arr(i_new+1,j_new+1,lk+1) );
665  if (z_target > z_lo && z_target < z_hi){
666  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(lk >= d_radius,
667  "K index must be larger than averaging radius!");
668  amrex::ignore_unused(d_radius);
669  k_arr(i,j,0) = lk;
670  zref_arr(i,j,0) = myhalf * (z_hi + z_lo) - z_bot_face;
671  break;
672  }
673  }
674  });
675  }
676 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtCellCenter(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:85
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtCellCenter(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:70

Referenced by make_MOSTAverage_at_level().

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

◆ set_norm_positions_T()

void MOSTAverage::set_norm_positions_T ( const int &  lev)

Function to set positions with terrain and normal vector (with interpolation).

748 {
749  // Peel back the level
750  auto& fields = m_fields[lev];
751 
752  // MFIter over CC data
753  int imf_cc = 2;
754 
755  ParmParse pp(m_pp_prefix);
756  Real zref_tmp = zref_default;
757  auto read_zref = pp.query("most.zref",zref_tmp);
758  if (!read_zref) {
759  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
760  }
761  int klo = m_geom[lev].Domain().smallEnd(2);
762 
763  // Capture for device
764  Real d_zref = zref_tmp;
765  const auto plo = m_geom[lev].ProbLoArray();
766 
767  RealVect base;
768  const auto dx = m_geom[lev].CellSizeArray();
769  const auto dxInv = m_geom[lev].InvCellSizeArray();
770  IntVect ng = m_x_pos[lev]->nGrowVect(); ng[2]=0;
771  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
772  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
773  Box gtbx = mfi.growntilebox(ng);
774  RealBox grb{gtbx,dx.data(),base.dataPtr()};
775 
776  if (npbx.smallEnd(2) != klo) { continue; }
777 
778  npbx.makeSlab(2,klo);
779 
780  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
781  auto x_pos_arr = m_x_pos[lev]->array(mfi);
782  auto y_pos_arr = m_y_pos[lev]->array(mfi);
783  auto z_pos_arr = m_z_pos[lev]->array(mfi);
784  auto zref_arr = m_zref[lev]->array(mfi);
785  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
786  {
787  // Elements of normal vector
788  Real met_h_xi = Compute_h_xi_AtCellCenter (i,j,k,dxInv,z_phys_arr);
789  Real met_h_eta = Compute_h_eta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
790  Real imag = one / std::sqrt(met_h_xi*met_h_xi + met_h_eta*met_h_eta + one);
791 
792  // Unit-normal vector scaled by z_ref
793  Real delta_x = -met_h_xi * imag * d_zref;
794  Real delta_y = -met_h_eta * imag * d_zref;
795  Real delta_z = imag * d_zref;
796 
797  // Position of the current node (indx:0,0,1)
798  Real x0 = plo[0] + ((Real) i + myhalf) * dx[0];
799  Real y0 = plo[1] + ((Real) j + myhalf) * dx[1];
800 
801  // Final position at end of vector
802  x_pos_arr(i,j,0) = x0 + delta_x;
803  y_pos_arr(i,j,0) = y0 + delta_y;
804  Real z_bot_face = fourth * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
805  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
806  z_pos_arr(i,j,0) = z_bot_face + delta_z;
807 
808  // NOTE: Normal vector end point can be below the surface for concave regions.
809  // Here we protect against that by augmenting the normal if needed.
810  int i_new = (int) ((x_pos_arr(i,j,0) - plo[0]) / dx[0] - myhalf);
811  int j_new = (int) ((y_pos_arr(i,j,0) - plo[1]) / dx[1] - myhalf);
812  Real z_new_bot_face = fourth * ( z_phys_arr(i_new,j_new ,k) + z_phys_arr(i_new+1,j_new ,k)
813  + z_phys_arr(i_new,j_new+1,k) + z_phys_arr(i_new+1,j_new+1,k) );
814  if (z_pos_arr(i,j,0) < z_new_bot_face) {
815  z_pos_arr(i,j,0) = z_new_bot_face + delta_z;
816  }
817 
818  zref_arr(i,j,0) = delta_z;
819 
820  // Destination position must be contained on the current process!
821  Real pos[] = {x_pos_arr(i,j,0)-plo[0],y_pos_arr(i,j,0)-plo[1],myhalf*dx[2]};
822  amrex::ignore_unused(pos);
823  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(grb.contains(&pos[0]),
824  "Query point outside of proc domain!");
825  });
826  }
827 }

Referenced by make_MOSTAverage_at_level().

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

◆ set_plane_normalization()

void MOSTAverage::set_plane_normalization ( const int &  lev)

Function to compute normalization for plane average.

348 {
349  // Cells per plane and temp avg storage
350  m_ncell_plane.resize(m_maxlev);
351  m_plane_average.resize(m_maxlev);
352 
353  // True domain not used for normalization
354  Box domain = m_geom[lev].Domain();
355 
356  // NOTE: Level 0 spans the whole domain, but finer
357  // levels do not have such a restriction.
358  // For now, use the bounding box of the boxArray
359  // for normalization, consistent with avg routine.
360 
361  // Bounded box of CC data used for normalization
362  Box bnd_bx = (m_fields[lev][2]->boxArray()).minimalBox();
363 
364  // NOTE: Bounding box must lie on the periodic boundaries
365  // in order to trip the is_per flag
366 
367  // Num components, plane avg, cells per plane
368  Array<int,AMREX_SPACEDIM> is_per = {0,0,0};
369  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
370  if ( m_geom[lev].isPeriodic(idim) &&
371  bnd_bx.bigEnd(idim)==domain.bigEnd(idim) &&
372  bnd_bx.smallEnd(idim)==domain.smallEnd(idim) ) { is_per[idim] = 1; }
373  }
374 
375  m_ncell_plane[lev].resize(m_navg);
376  m_plane_average[lev].resize(m_navg);
377  for (int iavg(0); iavg < m_navg; ++iavg) {
378  // Convert bnd_bx to current index type
379  IndexType ixt = m_averages[lev][iavg]->boxArray().ixType();
380  bnd_bx.convert(ixt);
381  IntVect bnd_bx_lo(bnd_bx.loVect());
382  IntVect bnd_bx_hi(bnd_bx.hiVect());
383 
384  m_plane_average[lev][iavg] = zero;
385 
386  m_ncell_plane[lev][iavg] = 1;
387  for (int idim(0); idim < AMREX_SPACEDIM; ++idim) {
388  if (idim != 2) {
389  if (ixt.nodeCentered(idim) && is_per[idim]) {
390  m_ncell_plane[lev][iavg] *= (bnd_bx_hi[idim] - bnd_bx_lo[idim]);
391  } else {
392  m_ncell_plane[lev][iavg] *= (bnd_bx_hi[idim] - bnd_bx_lo[idim] + 1);
393  }
394  }
395  } // idim
396  } // iavg
397 }

Referenced by make_MOSTAverage_at_level().

Here is the caller graph for this function:

◆ set_region_normalization()

void MOSTAverage::set_region_normalization ( const int &  )
inline
62  {m_ncell_region = (2 * m_radius + 1) * (2 * m_radius + 1) * (2 * m_radius + 1);}

Referenced by make_MOSTAverage_at_level().

Here is the caller graph for this function:

◆ set_rotated_fields()

void MOSTAverage::set_rotated_fields ( const int &  lev)

Function to set the rotated velocities.

289 {
290  // Peel back the level
291  auto& fields = m_fields[lev];
292  auto& rot_fields = m_rot_fields[lev];
293  auto z_phys_nd = m_z_phys_nd[lev];
294 
295  // Inverse grid size
296  const auto dxInv = m_geom[lev].InvCellSizeArray();
297 
298  // Single MFIter over CC data
299  int imf_cc = 2;
300 
301  // Populate rotated U & V for terrain
302 #ifdef _OPENMP
303 #pragma omp parallel if (Gpu::notInLaunchRegion())
304 #endif
305  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
306  Box ubx = mfi.tilebox(IntVect(1,0,0));
307  Box vbx = mfi.tilebox(IntVect(0,1,0));
308 
309  const Array4<const Real>& z_phys_arr = z_phys_nd->const_array(mfi);
310 
311  const Array4<const Real>& u_arr = fields[0]->const_array(mfi);
312  const Array4<const Real>& v_arr = fields[1]->const_array(mfi);
313  const Array4<const Real>& w_arr = fields[5]->const_array(mfi);
314 
315  const Array4<Real>& u_rot_arr = rot_fields[0]->array(mfi);
316  const Array4<Real>& v_rot_arr = rot_fields[1]->array(mfi);
317 
318  // U rotated magnitude
319  ParallelFor(ubx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
320  {
321  // Elements of first tangent vector
322  Real met_h_xi = Compute_h_xi_AtIface(i,j,k,dxInv,z_phys_arr);
323  u_rot_arr(i,j,k) = (u_arr(i,j,k) + met_h_xi*w_arr(i,j,k))
324  / std::sqrt(met_h_xi*met_h_xi + one);
325  });
326 
327  // V rotated magnitude
328  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
329  {
330  // Elements of second tangent vector
331  Real met_h_eta = Compute_h_eta_AtJface(i,j,k,dxInv,z_phys_arr);
332  v_rot_arr(i,j,k) = (v_arr(i,j,k) + met_h_eta*w_arr(i,j,k))
333  / std::sqrt(met_h_eta*met_h_eta + one);
334  });
335  }
336 
337  // Direct copy of other scalar variables
338  MultiFab::Copy(*rot_fields[2],*fields[2],0,0,1,rot_fields[2]->nGrowVect());
339  if (fields[3]) MultiFab::Copy(*rot_fields[3],*fields[3],0,0,1,rot_fields[3]->nGrowVect());
340 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:117
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:170

Referenced by compute_averages().

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

◆ set_z_positions_T()

void MOSTAverage::set_z_positions_T ( const int &  lev)

Function to set positions with terrain and e_z vector (with interpolation but no terrain normals)

685 {
686  // Peel back the level
687  auto& fields = m_fields[lev];
688 
689  // MFIter over CC data
690  int imf_cc = 2;
691 
692  ParmParse pp(m_pp_prefix);
693  Real zref_tmp = zref_default;
694  auto read_zref = pp.query("most.zref",zref_tmp);
695  if (!read_zref) {
696  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
697  } else {
698  m_zref[lev]->setVal(zref_tmp);
699  }
700  int klo = m_geom[lev].Domain().smallEnd(2);
701 
702  // Capture for device
703  Real d_zref = zref_tmp;
704  const auto plo = m_geom[lev].ProbLoArray();
705 
706  RealVect base;
707  const auto dx = m_geom[lev].CellSizeArray();
708  IntVect ng = m_x_pos[lev]->nGrowVect(); ng[2]=0;
709  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
710  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
711  Box gtbx = mfi.growntilebox(ng);
712 
713  if (npbx.smallEnd(2) != klo) { continue; }
714 
715  npbx.makeSlab(2,klo);
716 
717  RealBox grb{gtbx,dx.data(),base.dataPtr()};
718 
719  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
720  auto x_pos_arr = m_x_pos[lev]->array(mfi);
721  auto y_pos_arr = m_y_pos[lev]->array(mfi);
722  auto z_pos_arr = m_z_pos[lev]->array(mfi);
723  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
724  {
725  // Final position at end of vector
726  x_pos_arr(i,j,0) = plo[0] + ((Real) i + myhalf) * dx[0];
727  y_pos_arr(i,j,0) = plo[1] + ((Real) j + myhalf) * dx[1];
728  Real z_bot_face = fourth * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
729  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
730  z_pos_arr(i,j,0) = z_bot_face + d_zref;
731 
732  // Destination position must be contained on the current process!
733  Real pos[] = {x_pos_arr(i,j,0)-plo[0],y_pos_arr(i,j,0)-plo[1],myhalf*dx[2]};
734  amrex::ignore_unused(pos);
735  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(grb.contains(&pos[0]),
736  "Query point outside of proc domain!");
737  });
738  }
739 }

Referenced by make_MOSTAverage_at_level().

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

◆ trilinear_interp_T()

AMREX_GPU_HOST_DEVICE static AMREX_INLINE void MOSTAverage::trilinear_interp_T ( const amrex::Real xp,
const amrex::Real yp,
const amrex::Real zp,
amrex::Real interp_vals,
amrex::Array4< amrex::Real const > const &  interp_array,
amrex::Array4< amrex::Real const > const &  z_arr,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  plo,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  dxi,
const int  interp_comp 
)
inlinestatic

Function to compute trilinear interpolation with terrain.

Parameters
[in]xpX-position
[in]ypY-position
[out]interp_valsValues interpolated
[in]interp_arrayArray to interpolate on
[in]z_arrPhysical heights
[in]ploProblem lower bounds
[in]dxiInverse cell size array
[in]interp_compNumber of components to interpolate
135  {
136  // Search to get z/k
137  bool found = false;
138  int kmax = ubound(z_arr).z;
139  amrex::Real zval = zero;
140  amrex::Real z_target = zp;
141 
142  // Map position to i,j (must be same mapping in cpp file)
143  amrex::Real ireal = (xp - plo[0]) * dxi[0];
144  amrex::Real jreal = (yp - plo[1]) * dxi[1];
145  int i_new = (int) (ireal - myhalf);
146  int j_new = (int) (jreal - myhalf);
147 
148  for (int lk(0); lk<kmax; ++lk) {
149  amrex::Real z_lo = fourth * ( z_arr(i_new,j_new ,lk ) + z_arr(i_new+1,j_new ,lk )
150  + z_arr(i_new,j_new+1,lk ) + z_arr(i_new+1,j_new+1,lk ) );
151  amrex::Real z_hi = fourth * ( z_arr(i_new,j_new ,lk+1) + z_arr(i_new+1,j_new ,lk+1)
152  + z_arr(i_new,j_new+1,lk+1) + z_arr(i_new+1,j_new+1,lk+1) );
153  if (z_target > z_lo && z_target < z_hi){
154  found = true;
155  zval = (amrex::Real) lk + ((z_target - z_lo) / (z_hi - z_lo)) + myhalf;
156  break;
157  }
158  }
159 
160  amrex::ignore_unused(found);
161  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(found, "MOSTAverage: Height above terrain not found, try increasing z_ref!");
162 
163  // NOTE: This is the point ahead of the current i,j (e.g. i/j_new + 1)
164  const amrex::RealVect lx(ireal + myhalf, jreal + myhalf, zval);
165 
166  const amrex::IntVect ijk = lx.floor();
167 
168  int i = ijk[0]; int j = ijk[1]; int k = ijk[2];
169 
170  // Convert ijk (IntVect) to a RealVect explicitly
171  amrex::RealVect ijk_r(static_cast<amrex::Real>(ijk[0]),
172  static_cast<amrex::Real>(ijk[1]),
173  static_cast<amrex::Real>(ijk[2]));
174 
175  // Weights
176  const amrex::RealVect sx_hi = lx - ijk_r;
177  const amrex::RealVect sx_lo = one - sx_hi;
178 
179  for (int n = 0; n < interp_comp; n++) {
180  interp_vals[n] = sx_lo[0]*sx_lo[1]*sx_lo[2]*interp_array(i-1, j-1, k-1,n) +
181  sx_lo[0]*sx_lo[1]*sx_hi[2]*interp_array(i-1, j-1, k ,n) +
182  sx_lo[0]*sx_hi[1]*sx_lo[2]*interp_array(i-1, j , k-1,n) +
183  sx_lo[0]*sx_hi[1]*sx_hi[2]*interp_array(i-1, j , k ,n) +
184  sx_hi[0]*sx_lo[1]*sx_lo[2]*interp_array(i , j-1, k-1,n) +
185  sx_hi[0]*sx_lo[1]*sx_hi[2]*interp_array(i , j-1, k ,n) +
186  sx_hi[0]*sx_hi[1]*sx_lo[2]*interp_array(i , j , k-1,n) +
187  sx_hi[0]*sx_hi[1]*sx_hi[2]*interp_array(i , j , k ,n);
188  }
189  }

Referenced by compute_plane_averages(), and compute_region_averages().

Here is the caller graph for this function:

◆ update_field_ptrs()

void MOSTAverage::update_field_ptrs ( const int &  lev,
amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars_old,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Theta_prim,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Qv_prim,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Qr_prim 
)

Function to reset the pointers to field variables.

Parameters
[in]levCurrent level
[in]vars_oldConserved variables at each level
[in]Theta_primPrimitive theta component at each level
274 {
275  m_fields[lev][0] = &vars_old[lev][Vars::xvel];
276  m_fields[lev][1] = &vars_old[lev][Vars::yvel];
277  m_fields[lev][2] = Theta_prim[lev].get();
278  m_fields[lev][3] = Qv_prim[lev].get();
279  m_fields[lev][4] = Qr_prim[lev].get();
280  m_fields[lev][5] = &vars_old[lev][Vars::zvel];
281 }

Referenced by SurfaceLayer::update_mac_ptrs().

Here is the caller graph for this function:

◆ write_averages()

void MOSTAverage::write_averages ( const int &  lev)

Function to write averages to text file.

Parameters
[in]levCurrent level
1686 {
1687  // Peel back the level
1688  auto& fields = m_fields[lev];
1689  auto& averages = m_averages[lev];
1690 
1691  // MFIter on CC
1692  int imf_cc = 2;
1693 
1694  int klo = m_geom[lev].Domain().smallEnd(2);
1695 
1696  int navg = m_navg - 1;
1697 
1698  std::ofstream ofile;
1699  ofile.open ("MOST_averages.txt");
1700  ofile << "Averages computed via MOSTAverages class:\n";
1701 
1702  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1703  Box pbx = mfi.tilebox();
1704 
1705  if(pbx.smallEnd(2) != klo) { continue; }
1706 
1707  pbx.makeSlab(2,klo);
1708 
1709  int il = pbx.smallEnd(0); int iu = pbx.bigEnd(0);
1710  int jl = pbx.smallEnd(1); int ju = pbx.bigEnd(1);
1711 
1712  for (int j(jl); j <= ju; ++j) {
1713  for (int i(il); i <= iu; ++i) {
1714  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1715  int k = 0;
1716  for (int iavg(0); iavg <= navg; ++iavg) {
1717  auto mf_arr = averages[iavg]->array(mfi);
1718  ofile << "iavg val: "
1719  << iavg << ' '
1720  << mf_arr(i,j,k) << "\n";
1721  }
1722  ofile << "\n";
1723  }
1724  }
1725  }
1726  ofile.close();
1727 }
Here is the call graph for this function:

◆ write_k_indices()

void MOSTAverage::write_k_indices ( const int &  lev)

Function to write the K indices to text file.

Parameters
[in]levCurrent level
1540 {
1541  // Peel back the level
1542  auto& fields = m_fields[lev];
1543  auto& k_indx = m_k_indx[lev];
1544 
1545  // MFIter on CC
1546  int imf_cc = 2;
1547 
1548  int klo = m_geom[lev].Domain().smallEnd(2);
1549 
1550  std::ofstream ofile;
1551  ofile.open ("MOST_k_indices.txt");
1552  ofile << "K indices used to compute averages via MOSTAverages class:\n";
1553 
1554  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1555  Box pbx = mfi.tilebox();
1556 
1557  if(pbx.smallEnd(2) != klo) { continue; }
1558 
1559  pbx.makeSlab(2,klo);
1560 
1561  int il = pbx.smallEnd(0); int iu = pbx.bigEnd(0);
1562  int jl = pbx.smallEnd(1); int ju = pbx.bigEnd(1);
1563 
1564  auto k_arr = k_indx->array(mfi);
1565 
1566  for (int j(jl); j <= ju; ++j) {
1567  for (int i(il); i <= iu; ++i) {
1568  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1569  int k = 0;
1570  ofile << "K_ind: "
1571  << k_arr(i,j,k) << "\n";
1572  ofile << "\n";
1573  }
1574  }
1575  }
1576  ofile.close();
1577 }
Here is the call graph for this function:

◆ write_norm_indices()

void MOSTAverage::write_norm_indices ( const int &  lev)

Function to write I,J,K indices to text file.

Parameters
[in]levCurrent level
1587 {
1588  // Peel back the level
1589  auto& fields = m_fields[lev];
1590  auto& k_indx = m_k_indx[lev];
1591  auto& j_indx = m_j_indx[lev];
1592  auto& i_indx = m_i_indx[lev];
1593 
1594  // MFIter on CC
1595  int imf_cc = 2;
1596 
1597  int klo = m_geom[lev].Domain().smallEnd(2);
1598 
1599  std::ofstream ofile;
1600  ofile.open ("MOST_ijk_indices.txt");
1601  ofile << "IJK indices used to compute averages via MOSTAverages class:\n";
1602 
1603  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1604  Box pbx = mfi.tilebox();
1605 
1606  if(pbx.smallEnd(2) != klo) { continue; }
1607 
1608  pbx.makeSlab(2,klo);
1609 
1610  int il = pbx.smallEnd(0); int iu = pbx.bigEnd(0);
1611  int jl = pbx.smallEnd(1); int ju = pbx.bigEnd(1);
1612 
1613  auto k_arr = k_indx->array(mfi);
1614  auto j_arr = j_indx ? j_indx->array(mfi) : Array4<int> {};
1615  auto i_arr = i_indx ? i_indx->array(mfi) : Array4<int> {};
1616 
1617  for (int j(jl); j <= ju; ++j) {
1618  for (int i(il); i <= iu; ++i) {
1619  ofile << "(I1,J1,K1): " << "(" << i << "," << j << "," << 0 << ")" << "\n";
1620 
1621  int k = 0;
1622  int km = k_arr(i,j,k);
1623  int jm = j_arr ? j_arr(i,j,k) : j;
1624  int im = i_arr ? i_arr(i,j,k) : i;
1625 
1626  ofile << "(I2,J2,K2): "
1627  << "(" << im << "," << jm << "," << km << ")" << "\n";
1628  ofile << "\n";
1629  }
1630  }
1631  }
1632  ofile.close();
1633 }
Here is the call graph for this function:

◆ write_xz_positions()

void MOSTAverage::write_xz_positions ( const int &  lev,
const int &  j 
)

Function to write X & Z positions to text file.

Parameters
[in]levCurrent level
[in]jIndex in y-dir
1645 {
1646  // Peel back the level
1647  auto& fields = m_fields[lev];
1648  auto& x_pos_mf = m_x_pos[lev];
1649  auto& z_pos_mf = m_z_pos[lev];
1650 
1651  // MFIter on CC
1652  int imf_cc = 2;
1653 
1654  int klo = m_geom[lev].Domain().smallEnd(2);
1655 
1656  std::ofstream ofile;
1657  ofile.open ("MOST_xz_positions.txt");
1658 
1659  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1660  Box pbx = mfi.tilebox();
1661 
1662  if(pbx.smallEnd(2) != klo) { continue; }
1663 
1664  pbx.makeSlab(2,klo);
1665 
1666  int il = pbx.smallEnd(0); int iu = pbx.bigEnd(0);
1667 
1668  auto x_pos_arr = x_pos_mf->array(mfi);
1669  auto z_pos_arr = z_pos_mf->array(mfi);
1670 
1671  int k = 0;
1672  for (int i(il); i <= iu; ++i)
1673  ofile << x_pos_arr(i,j,k) << ' ' << z_pos_arr(i,j,k) << "\n";
1674  }
1675  ofile.close();
1676 }
Here is the call graph for this function:

Member Data Documentation

◆ include_subgrid_vel

bool MOSTAverage::include_subgrid_vel = false
protected

◆ m_averages

amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab> > > MOSTAverage::m_averages
protected

◆ m_fact_new

amrex::Real MOSTAverage::m_fact_new
protected

◆ m_fact_old

amrex::Real MOSTAverage::m_fact_old
protected

◆ m_fields

◆ m_geom

◆ m_i_indx

amrex::Vector<std::unique_ptr<amrex::iMultiFab> > MOSTAverage::m_i_indx
protected

◆ m_interp

bool MOSTAverage::m_interp {false}
protected

◆ m_j_indx

amrex::Vector<std::unique_ptr<amrex::iMultiFab> > MOSTAverage::m_j_indx
protected

◆ m_k_in

amrex::Vector<int> MOSTAverage::m_k_in
protected

◆ m_k_indx

◆ m_maxlev

int MOSTAverage::m_maxlev {0}
protected

◆ m_mesh_type

MeshType MOSTAverage::m_mesh_type
protected

◆ m_navg

◆ m_ncell_plane

amrex::Vector<amrex::Vector<int> > MOSTAverage::m_ncell_plane
protected

◆ m_ncell_region

int MOSTAverage::m_ncell_region {1}
protected

◆ m_norm_vec

bool MOSTAverage::m_norm_vec {false}
protected

◆ m_nvar

int MOSTAverage::m_nvar {6}
protected

◆ m_plane_average

amrex::Vector<amrex::Vector<amrex::Real> > MOSTAverage::m_plane_average
protected

◆ m_policy

int MOSTAverage::m_policy {0}
protected

◆ m_pp_prefix

std::string MOSTAverage::m_pp_prefix
protected

◆ m_radius

◆ m_rot_fields

amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab> > > MOSTAverage::m_rot_fields
protected

◆ m_rotate

bool MOSTAverage::m_rotate {false}
protected

◆ m_t_avg

bool MOSTAverage::m_t_avg {false}
protected

◆ m_t_init

amrex::Vector<int> MOSTAverage::m_t_init
protected

◆ m_terrain_type

TerrainType MOSTAverage::m_terrain_type
protected

◆ m_time_window

amrex::Real MOSTAverage::m_time_window {amrex::Real(1.0e-16)}
protected

◆ m_Vsg

amrex::Vector<amrex::Real> MOSTAverage::m_Vsg
protected

◆ m_x_pos

amrex::Vector<std::unique_ptr<amrex::MultiFab> > MOSTAverage::m_x_pos
protected

◆ m_y_pos

amrex::Vector<std::unique_ptr<amrex::MultiFab> > MOSTAverage::m_y_pos
protected

◆ m_z_phys_nd

◆ m_z_pos

amrex::Vector<std::unique_ptr<amrex::MultiFab> > MOSTAverage::m_z_pos
protected

◆ m_zref

amrex::Vector<std::unique_ptr<amrex::MultiFab> > MOSTAverage::m_zref
protected

◆ zref_default


The documentation for this class was generated from the following files: