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_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
 

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 {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 = 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
779 {
780  if (m_rotate) set_rotated_fields(lev);
781 
782  switch(m_policy) {
783  case 0: // Standard plane average
785  break;
786  case 1: // Local region/point
788  break;
789  default:
790  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false, "Unknown policy for MOSTAverage!");
791  }
792 
793  // We have initialized the averages
794  if (m_t_avg) m_t_init[lev] = 1;
795 }
bool m_t_avg
Definition: ERF_MOSTAverage.H:231
void compute_plane_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:804
int m_policy
Definition: ERF_MOSTAverage.H:201
void compute_region_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:1103
amrex::Vector< int > m_t_init
Definition: ERF_MOSTAverage.H:232
void set_rotated_fields(const int &lev)
Definition: ERF_MOSTAverage.cpp:284
bool m_rotate
Definition: ERF_MOSTAverage.H:202
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
805 {
806  // Peel back the level
807  auto& fields = m_fields[lev];
808  auto& rot_fields = m_rot_fields[lev];
809  auto& averages = m_averages[lev];
810  const auto & geom = m_geom[lev];
811 
812  auto& z_phys = m_z_phys_nd[lev];
813  auto& x_pos = m_x_pos[lev];
814  auto& y_pos = m_y_pos[lev];
815  auto& z_pos = m_z_pos[lev];
816 
817  auto& i_indx = m_i_indx[lev];
818  auto& j_indx = m_j_indx[lev];
819  auto& k_indx = m_k_indx[lev];
820 
821  auto& ncell_plane = m_ncell_plane[lev];
822  auto& plane_average = m_plane_average[lev];
823 
824  // Set factors for time averaging
825  Real d_fact_new, d_fact_old;
826  if (m_t_avg && m_t_init[lev]) {
827  d_fact_new = m_fact_new;
828  d_fact_old = m_fact_old;
829  } else {
830  d_fact_new = 1.0;
831  d_fact_old = 0.0;
832  }
833 
834  int klo = m_geom[lev].Domain().smallEnd(2);
835 
836  // GPU array to accumulate averages into
837  Gpu::DeviceVector<Real> pavg(plane_average.size(), 0.0);
838  Real* plane_avg = pavg.data();
839 
840  // Vectors for normalization and buffer storage
841  Vector<Real> denom(plane_average.size(),0.0);
842  Vector<Real> val_old(plane_average.size(),0.0);
843 
844  //
845  //----------------------------------------------------------
846  // Averages over all the fields
847  //----------------------------------------------------------
848  //
849  Box domain = geom.Domain();
850 
851  Array<int,AMREX_SPACEDIM> is_per = {0,0,0};
852  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
853  if (geom.isPeriodic(idim)) is_per[idim] = 1;
854  }
855 
856  // Averages for U,V,T,Qv (not Qc or W)
857  for (int imf(0); imf < 4; ++imf) {
858 
859  // Continue if no valid Qv pointer
860  if (!fields[imf]) continue;
861 
862  denom[imf] = 1.0 / (Real)ncell_plane[imf];
863  val_old[imf] = plane_average[imf]*d_fact_old;
864 
865 #ifdef _OPENMP
866 #pragma omp parallel if (Gpu::notInLaunchRegion())
867 #endif
868  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
869  Box vbx = mfi.validbox(); // This is the grid (not tile)
870  Box pbx = mfi.tilebox(); // This is the tile (not grid)
871 
872  if (pbx.smallEnd(2) != klo) { continue; }
873 
874  // Make planar since mfiter is over fields
875  pbx.makeSlab(2,klo);
876 
877  // Avoid double counting nodal data by changing the high end when we are
878  // at the high side of the grid (not just of the tile)
879  IndexType ixt = averages[imf]->boxArray().ixType();
880  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
881  if ( ixt.nodeCentered(idim) && (pbx.bigEnd(idim) == vbx.bigEnd(idim)) ) {
882  int dom_hi = domain.bigEnd(idim)+1;
883  if (pbx.bigEnd(idim) < dom_hi || is_per[idim]) {
884  pbx.growHi(idim,-1);
885  }
886  }
887  }
888 
889  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
890  fields[imf]->const_array(mfi);
891 
892  if (m_interp) {
893  const auto plo = geom.ProbLoArray();
894  const auto dxInv = geom.InvCellSizeArray();
895  const auto z_phys_arr = z_phys->const_array(mfi);
896  auto x_pos_arr = x_pos->array(mfi);
897  auto y_pos_arr = y_pos->array(mfi);
898  auto z_pos_arr = z_pos->array(mfi);
899  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
900  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
901  {
902  Real interp{0};
903  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
904  &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
905  Real val = interp;
906  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
907  });
908  } else {
909  auto k_arr = k_indx->const_array(mfi);
910  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
911  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
912  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
913  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
914  {
915  int mk = k_arr(i,j,0);
916  int mj = j_arr ? j_arr(i,j,0) : j;
917  int mi = i_arr ? i_arr(i,j,0) : i;
918  Real val = mf_arr(mi,mj,mk);
919  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
920  });
921  }
922  }
923  }
924 
925  //
926  //------------------------------------------------------------------------
927  // Averages for virtual potential temperature
928  // (This is cell-centered so we don't need to worry about double-counting)
929  //------------------------------------------------------------------------
930  //
931  if (fields[3]) // We have water vapor
932  {
933  int iavg = 4;
934  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
935  val_old[iavg] = plane_average[iavg]*d_fact_old;
936 
937 #ifdef _OPENMP
938 #pragma omp parallel if (Gpu::notInLaunchRegion())
939 #endif
940  for (MFIter mfi(*fields[3], TileNoZ()); mfi.isValid(); ++mfi)
941  {
942  Box pbx = mfi.tilebox();
943 
944  if (pbx.smallEnd(2) != klo) { continue; }
945 
946  pbx.makeSlab(2,klo);
947 
948  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
949  const Array4<Real const>& qv_mf_arr = fields[3]->const_array(mfi);
950  const Array4<Real const>& qr_mf_arr = (fields[4]) ? fields[4]->const_array(mfi) :
951  Array4<const Real> {};
952 
953  if (m_interp) {
954  const auto plo = m_geom[lev].ProbLoArray();
955  const auto dxInv = m_geom[lev].InvCellSizeArray();
956  const auto z_phys_arr = z_phys->const_array(mfi);
957  auto x_pos_arr = x_pos->array(mfi);
958  auto y_pos_arr = y_pos->array(mfi);
959  auto z_pos_arr = z_pos->array(mfi);
960  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
961  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
962  {
963  Real T_interp{0};
964  Real qv_interp{0};
965  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
966  &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
967  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
968  &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
969  Real vfac;
970  if (qr_mf_arr) {
971  // We also have liquid water
972  Real qr_interp{0};
973  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
974  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
975  vfac = 1.0 + 0.61*qv_interp - qr_interp;
976  } else {
977  vfac = 1.0 + 0.61*qv_interp;
978  }
979  const Real val = T_interp * vfac;
980  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
981  });
982  } else {
983  auto k_arr = k_indx->const_array(mfi);
984  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
985  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
986  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
987  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
988  {
989  int mk = k_arr(i,j,0);
990  int mj = j_arr ? j_arr(i,j,0) : j;
991  int mi = i_arr ? i_arr(i,j,0) : i;
992  Real vfac;
993  if (qr_mf_arr) {
994  // We also have liquid water
995  vfac = 1.0 + 0.61*qv_mf_arr(mi,mj,mk) - qr_mf_arr(mi,mj,mk);
996  } else {
997  vfac = 1.0 + 0.61*qv_mf_arr(mi,mj,mk);
998  }
999  const Real val = T_mf_arr(mi,mj,mk) * vfac;
1000  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1001  });
1002  }
1003  }
1004  }
1005  else // copy temperature
1006  {
1007  int iavg = m_navg - 2;
1008  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
1009  // plane_avg[iavg] = plane_avg[2]
1010  Gpu::copy(Gpu::deviceToDevice, pavg.begin() + 2, pavg.begin() + 3,
1011  pavg.begin() + iavg);
1012  }
1013 
1014  //
1015  //------------------------------------------------------------------------
1016  // Averages for the tangential velocity magnitude
1017  // (This is cell-centered so we don't need to worry about double-counting)
1018  //------------------------------------------------------------------------
1019  //
1020  {
1021  int imf_cc = 2;
1022  int imf = 0;
1023  int iavg = m_navg - 1;
1024  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
1025  val_old[iavg] = plane_average[iavg]*d_fact_old;
1026 
1027  const Real Vsg = m_Vsg[lev];
1028 
1029 #ifdef _OPENMP
1030 #pragma omp parallel if (Gpu::notInLaunchRegion())
1031 #endif
1032  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi)
1033  {
1034  Box pbx = mfi.tilebox();
1035 
1036  if (pbx.smallEnd(2) != klo) { continue; }
1037 
1038  pbx.makeSlab(2,klo);
1039 
1040  // Last element is Umag and always cell centered
1041  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
1042  fields[imf ]->const_array(mfi);
1043  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
1044  fields[imf+1]->const_array(mfi);
1045 
1046  if (m_interp) {
1047  const auto plo = m_geom[lev].ProbLoArray();
1048  const auto dxInv = m_geom[lev].InvCellSizeArray();
1049  const auto z_phys_arr = z_phys->const_array(mfi);
1050  auto x_pos_arr = x_pos->array(mfi);
1051  auto y_pos_arr = y_pos->array(mfi);
1052  auto z_pos_arr = z_pos->array(mfi);
1053  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
1054  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
1055  {
1056  Real u_interp{0};
1057  Real v_interp{0};
1058  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
1059  &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
1060  trilinear_interp_T(x_pos_arr(i,j,0), y_pos_arr(i,j,0), z_pos_arr(i,j,0),
1061  &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
1062  const Real val = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
1063  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1064  });
1065  } else {
1066  auto k_arr = k_indx->const_array(mfi);
1067  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1068  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1069  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
1070  AMREX_GPU_DEVICE(int i, int j, int , Gpu::Handler const& handler) noexcept
1071  {
1072  int mk = k_arr(i,j,0);
1073  int mj = j_arr ? j_arr(i,j,0) : j;
1074  int mi = i_arr ? i_arr(i,j,0) : i;
1075  const Real u_val = 0.5 * (u_mf_arr(mi,mj,mk) + u_mf_arr(mi+1,mj ,mk));
1076  const Real v_val = 0.5 * (v_mf_arr(mi,mj,mk) + v_mf_arr(mi ,mj+1,mk));
1077  const Real val = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
1078  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1079  });
1080  }
1081  }
1082  }
1083 
1084  // Copy to host and sum across procs
1085  Gpu::copy(Gpu::deviceToHost, pavg.begin(), pavg.end(), plane_average.begin());
1086  ParallelDescriptor::ReduceRealSum(plane_average.data(), plane_average.size());
1087 
1088  // No spatial variation with plane averages
1089  for (int iavg(0); iavg < m_navg; ++iavg){
1090  plane_average[iavg] *= denom[iavg]*d_fact_new;
1091  plane_average[iavg] += val_old[iavg];
1092  averages[iavg]->setVal(plane_average[iavg]);
1093  }
1094 }
Real vfac
Definition: ERF_InitCustomPertVels_ABL.H:26
const int klo
Definition: ERF_InitCustomPertVels_ParticleTests.H:4
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:16
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+0.5) *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<=1.0) { Real dT=T_pert *(std::cos(PI *L)+1.0)/2.0;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:199
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_averages
Definition: ERF_MOSTAverage.H:210
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_y_pos
Definition: ERF_MOSTAverage.H:205
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_i_indx
Definition: ERF_MOSTAverage.H:207
amrex::Vector< amrex::MultiFab * > m_z_phys_nd
Definition: ERF_MOSTAverage.H:191
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_x_pos
Definition: ERF_MOSTAverage.H:204
amrex::Vector< amrex::Real > m_Vsg
Definition: ERF_MOSTAverage.H:239
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_rot_fields
Definition: ERF_MOSTAverage.H:211
amrex::Vector< amrex::Vector< amrex::Real > > m_plane_average
Definition: ERF_MOSTAverage.H:216
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_z_pos
Definition: ERF_MOSTAverage.H:206
amrex::Vector< amrex::Vector< int > > m_ncell_plane
Definition: ERF_MOSTAverage.H:215
amrex::Real m_fact_new
Definition: ERF_MOSTAverage.H:234
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_j_indx
Definition: ERF_MOSTAverage.H:208
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_fields
Definition: ERF_MOSTAverage.H:190
amrex::Real m_fact_old
Definition: ERF_MOSTAverage.H:234
bool m_interp
Definition: ERF_MOSTAverage.H:226
const amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_MOSTAverage.H:189
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:120
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_k_indx
Definition: ERF_MOSTAverage.H:209

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

Referenced by SurfaceLayer::get_mac_avg().

Here is the caller graph for this function:

◆ get_zref()

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

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

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_N()

void MOSTAverage::set_k_indices_N ( const int &  lev)

Function to set K indices without terrain.

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

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).

459 {
460  // Peel back the level
461  auto& fields = m_fields[lev];
462 
463  // MFIter over CC data
464  int imf_cc = 2;
465 
466  ParmParse pp(m_pp_prefix);
467  Real zref_tmp = zref_default;
468  auto read_z = pp.query("most.zref",zref_tmp);
469  auto read_k = pp.queryarr("most.k_arr_in",m_k_in);
470  int klo = m_geom[lev].Domain().smallEnd(2);
471 
472  // Allow default zref
473  if (!read_z) {
474  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
475  read_z = true;
476  }
477 
478  // No default behavior with terrain (we can't tell the difference between
479  // vertical grid stretching and true terrain)
480  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(read_z != read_k,
481  "Need to specify zref or k_arr_in for MOST");
482 
483  // Capture for device
484  Real d_zref = zref_tmp;
485  Real d_radius = m_radius;
486  amrex::ignore_unused(d_radius);
487 
488  // Specify z_ref & compute k_indx (z_ref takes precedence)
489  if (read_z) {
490  int kmax = m_geom[lev].Domain().bigEnd(2);
491  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
492  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
493 
494  if (npbx.smallEnd(2) != klo) { continue; }
495 
496  npbx.makeSlab(2,klo);
497 
498  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
499  auto k_arr = m_k_indx[lev]->array(mfi);
500  auto zref_arr = m_zref[lev]->array(mfi);
501  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
502  {
503  k_arr(i,j,k) = klo;
504  bool found = false;
505  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
506  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
507  Real z_target = z_bot_face + d_zref;
508  for (int lk(klo); lk<=kmax; ++lk) {
509  Real z_lo = 0.25 * ( z_phys_arr(i,j ,lk ) + z_phys_arr(i+1,j ,lk )
510  + z_phys_arr(i,j+1,lk ) + z_phys_arr(i+1,j+1,lk ) );
511  Real z_hi = 0.25 * ( z_phys_arr(i,j ,lk+1) + z_phys_arr(i+1,j ,lk+1)
512  + z_phys_arr(i,j+1,lk+1) + z_phys_arr(i+1,j+1,lk+1) );
513  if (z_target > z_lo && z_target < z_hi){
514  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(lk >= d_radius,
515  "K index must be larger than averaging radius!");
516  k_arr(i,j,0) = lk;
517  zref_arr(i,j,0) = 0.5 * (z_hi + z_lo) - z_bot_face;
518  found = true;
519  break;
520  }
521  }
522  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(found,
523  "zref not found with terrain!");
524  });
525  }
526  // Specified k_indx & compute z_ref
527  } else if (read_k) {
528  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false, "Specified k-indx with terrain not implemented!");
529  }
530 }

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).

539 {
540  // Peel back the level
541  auto& fields = m_fields[lev];
542 
543  // MFIter over CC data
544  int imf_cc = 2;
545 
546  ParmParse pp(m_pp_prefix);
547  Real zref_tmp = zref_default;
548  auto read_zref = pp.query("most.zref",zref_tmp);
549  if (!read_zref) {
550  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
551  }
552  int klo = m_geom[lev].Domain().smallEnd(2);
553 
554  // Capture for device
555  Real d_zref = zref_tmp;
556  Real d_radius = m_radius;
557 
558  const auto dxInv = m_geom[lev].InvCellSizeArray();
559  IntVect ng = m_k_indx[lev]->nGrowVect(); ng[2]=0;
560  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
561  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
562 
563  if (npbx.smallEnd(2) != klo) { continue; }
564 
565  int kmax = npbx.bigEnd(2);
566 
567  npbx.makeSlab(2,klo);
568 
569  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
570  auto i_arr = m_i_indx[lev]->array(mfi);
571  auto j_arr = m_j_indx[lev]->array(mfi);
572  auto k_arr = m_k_indx[lev]->array(mfi);
573  auto zref_arr = m_zref[lev]->array(mfi);
574  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
575  {
576  // Elements of normal vector
577  Real met_h_xi = Compute_h_xi_AtCellCenter (i,j,k,dxInv,z_phys_arr);
578  Real met_h_eta = Compute_h_eta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
579  Real mag = std::sqrt(met_h_xi*met_h_xi + met_h_eta*met_h_eta + 1.0);
580 
581  // Unit-normal vector scaled by z_ref
582  Real delta_x = -met_h_xi/mag * d_zref;
583  Real delta_y = -met_h_eta/mag * d_zref;
584  Real delta_z = 1.0/mag * d_zref;
585 
586  // Compute i & j as displacements (no grid stretching)
587  int delta_i = static_cast<int>(std::round(delta_x*dxInv[0]));
588  int delta_j = static_cast<int>(std::round(delta_y*dxInv[1]));
589  int i_new = i + delta_i;
590  int j_new = j + delta_j;
591  i_arr(i,j,0) = i_new;
592  j_arr(i,j,0) = j_new;
593 
594  // Search for k (grid is stretched in z)
595  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
596  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
597  Real z_target = z_bot_face + delta_z;
598  k_arr(i,j,0) = klo;
599  zref_arr(i,j,0) = 0.5 * z_bot_face +
600  0.125 * ( z_phys_arr(i ,j ,k+1) + z_phys_arr(i+1,j ,k+1)
601  + z_phys_arr(i ,j+1,k+1) + z_phys_arr(i+1,j+1,k+1) );
602  for (int lk(klo); lk<=kmax; ++lk) {
603  Real z_lo = 0.25 * ( z_phys_arr(i_new,j_new ,lk ) + z_phys_arr(i_new+1,j_new ,lk )
604  + z_phys_arr(i_new,j_new+1,lk ) + z_phys_arr(i_new+1,j_new+1,lk ) );
605  Real z_hi = 0.25 * ( z_phys_arr(i_new,j_new ,lk+1) + z_phys_arr(i_new+1,j_new ,lk+1)
606  + z_phys_arr(i_new,j_new+1,lk+1) + z_phys_arr(i_new+1,j_new+1,lk+1) );
607  if (z_target > z_lo && z_target < z_hi){
608  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(lk >= d_radius,
609  "K index must be larger than averaging radius!");
610  amrex::ignore_unused(d_radius);
611  k_arr(i,j,0) = lk;
612  zref_arr(i,j,0) = 0.5 * (z_hi + z_lo) - z_bot_face;
613  break;
614  }
615  }
616  });
617  }
618 }
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:83
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:68

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).

690 {
691  // Peel back the level
692  auto& fields = m_fields[lev];
693 
694  // MFIter over CC data
695  int imf_cc = 2;
696 
697  ParmParse pp(m_pp_prefix);
698  Real zref_tmp = zref_default;
699  auto read_zref = pp.query("most.zref",zref_tmp);
700  if (!read_zref) {
701  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
702  }
703  int klo = m_geom[lev].Domain().smallEnd(2);
704 
705  // Capture for device
706  Real d_zref = zref_tmp;
707  const auto plo = m_geom[lev].ProbLoArray();
708 
709  RealVect base;
710  const auto dx = m_geom[lev].CellSizeArray();
711  const auto dxInv = m_geom[lev].InvCellSizeArray();
712  IntVect ng = m_x_pos[lev]->nGrowVect(); ng[2]=0;
713  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
714  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
715  Box gtbx = mfi.growntilebox(ng);
716  RealBox grb{gtbx,dx.data(),base.dataPtr()};
717 
718  if (npbx.smallEnd(2) != klo) { continue; }
719 
720  npbx.makeSlab(2,klo);
721 
722  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
723  auto x_pos_arr = m_x_pos[lev]->array(mfi);
724  auto y_pos_arr = m_y_pos[lev]->array(mfi);
725  auto z_pos_arr = m_z_pos[lev]->array(mfi);
726  auto zref_arr = m_zref[lev]->array(mfi);
727  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
728  {
729  // Elements of normal vector
730  Real met_h_xi = Compute_h_xi_AtCellCenter (i,j,k,dxInv,z_phys_arr);
731  Real met_h_eta = Compute_h_eta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
732  Real imag = 1.0 / std::sqrt(met_h_xi*met_h_xi + met_h_eta*met_h_eta + 1.0);
733 
734  // Unit-normal vector scaled by z_ref
735  Real delta_x = -met_h_xi * imag * d_zref;
736  Real delta_y = -met_h_eta * imag * d_zref;
737  Real delta_z = imag * d_zref;
738 
739  // Position of the current node (indx:0,0,1)
740  Real x0 = plo[0] + ((Real) i + 0.5) * dx[0];
741  Real y0 = plo[1] + ((Real) j + 0.5) * dx[1];
742 
743  // Final position at end of vector
744  x_pos_arr(i,j,0) = x0 + delta_x;
745  y_pos_arr(i,j,0) = y0 + delta_y;
746  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
747  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
748  z_pos_arr(i,j,0) = z_bot_face + delta_z;
749 
750  // NOTE: Normal vector end point can be below the surface for concave regions.
751  // Here we protect against that by augmenting the normal if needed.
752  int i_new = (int) ((x_pos_arr(i,j,0) - plo[0]) / dx[0] - 0.5);
753  int j_new = (int) ((y_pos_arr(i,j,0) - plo[1]) / dx[1] - 0.5);
754  Real z_new_bot_face = 0.25 * ( z_phys_arr(i_new,j_new ,k) + z_phys_arr(i_new+1,j_new ,k)
755  + z_phys_arr(i_new,j_new+1,k) + z_phys_arr(i_new+1,j_new+1,k) );
756  if (z_pos_arr(i,j,0) < z_new_bot_face) {
757  z_pos_arr(i,j,0) = z_new_bot_face + delta_z;
758  }
759 
760  zref_arr(i,j,0) = delta_z;
761 
762  // Destination position must be contained on the current process!
763  Real pos[] = {x_pos_arr(i,j,0)-plo[0],y_pos_arr(i,j,0)-plo[1],0.5*dx[2]};
764  amrex::ignore_unused(pos);
765  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(grb.contains(&pos[0]),
766  "Query point outside of proc domain!");
767  });
768  }
769 }

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.

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

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.

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

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)

627 {
628  // Peel back the level
629  auto& fields = m_fields[lev];
630 
631  // MFIter over CC data
632  int imf_cc = 2;
633 
634  ParmParse pp(m_pp_prefix);
635  Real zref_tmp = zref_default;
636  auto read_zref = pp.query("most.zref",zref_tmp);
637  if (!read_zref) {
638  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
639  } else {
640  m_zref[lev]->setVal(zref_tmp);
641  }
642  int klo = m_geom[lev].Domain().smallEnd(2);
643 
644  // Capture for device
645  Real d_zref = zref_tmp;
646  const auto plo = m_geom[lev].ProbLoArray();
647 
648  RealVect base;
649  const auto dx = m_geom[lev].CellSizeArray();
650  IntVect ng = m_x_pos[lev]->nGrowVect(); ng[2]=0;
651  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
652  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
653  Box gtbx = mfi.growntilebox(ng);
654 
655  if (npbx.smallEnd(2) != klo) { continue; }
656 
657  npbx.makeSlab(2,klo);
658 
659  RealBox grb{gtbx,dx.data(),base.dataPtr()};
660 
661  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
662  auto x_pos_arr = m_x_pos[lev]->array(mfi);
663  auto y_pos_arr = m_y_pos[lev]->array(mfi);
664  auto z_pos_arr = m_z_pos[lev]->array(mfi);
665  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
666  {
667  // Final position at end of vector
668  x_pos_arr(i,j,0) = plo[0] + ((Real) i + 0.5) * dx[0];
669  y_pos_arr(i,j,0) = plo[1] + ((Real) j + 0.5) * dx[1];
670  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
671  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
672  z_pos_arr(i,j,0) = z_bot_face + d_zref;
673 
674  // Destination position must be contained on the current process!
675  Real pos[] = {x_pos_arr(i,j,0)-plo[0],y_pos_arr(i,j,0)-plo[1],0.5*dx[2]};
676  amrex::ignore_unused(pos);
677  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(grb.contains(&pos[0]),
678  "Query point outside of proc domain!");
679  });
680  }
681 }

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
129  {
130  // Search to get z/k
131  bool found = false;
132  int kmax = ubound(z_arr).z;
133  amrex::Real zval = 0.0;
134  amrex::Real z_target = zp;
135 
136  // Map position to i,j (must be same mapping in cpp file)
137  amrex::Real ireal = (xp - plo[0]) * dxi[0];
138  amrex::Real jreal = (yp - plo[1]) * dxi[1];
139  int i_new = (int) (ireal - 0.5);
140  int j_new = (int) (jreal - 0.5);
141 
142  for (int lk(0); lk<kmax; ++lk) {
143  amrex::Real z_lo = 0.25 * ( z_arr(i_new,j_new ,lk ) + z_arr(i_new+1,j_new ,lk )
144  + z_arr(i_new,j_new+1,lk ) + z_arr(i_new+1,j_new+1,lk ) );
145  amrex::Real z_hi = 0.25 * ( z_arr(i_new,j_new ,lk+1) + z_arr(i_new+1,j_new ,lk+1)
146  + z_arr(i_new,j_new+1,lk+1) + z_arr(i_new+1,j_new+1,lk+1) );
147  if (z_target > z_lo && z_target < z_hi){
148  found = true;
149  zval = (amrex::Real) lk + ((z_target - z_lo) / (z_hi - z_lo)) + 0.5;
150  break;
151  }
152  }
153 
154  amrex::ignore_unused(found);
155  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(found, "MOSTAverage: Height above terrain not found, try increasing z_ref!");
156 
157  // NOTE: This is the point ahead of the current i,j (e.g. i/j_new + 1)
158  const amrex::RealVect lx(ireal + 0.5, jreal + 0.5, zval);
159 
160  const amrex::IntVect ijk = lx.floor();
161 
162  int i = ijk[0]; int j = ijk[1]; int k = ijk[2];
163 
164  // Convert ijk (IntVect) to a RealVect explicitly
165  amrex::RealVect ijk_r(static_cast<amrex::Real>(ijk[0]),
166  static_cast<amrex::Real>(ijk[1]),
167  static_cast<amrex::Real>(ijk[2]));
168 
169  // Weights
170  const amrex::RealVect sx_hi = lx - ijk_r;
171  const amrex::RealVect sx_lo = 1.0 - sx_hi;
172 
173  for (int n = 0; n < interp_comp; n++) {
174  interp_vals[n] = sx_lo[0]*sx_lo[1]*sx_lo[2]*interp_array(i-1, j-1, k-1,n) +
175  sx_lo[0]*sx_lo[1]*sx_hi[2]*interp_array(i-1, j-1, k ,n) +
176  sx_lo[0]*sx_hi[1]*sx_lo[2]*interp_array(i-1, j , k-1,n) +
177  sx_lo[0]*sx_hi[1]*sx_hi[2]*interp_array(i-1, j , k ,n) +
178  sx_hi[0]*sx_lo[1]*sx_lo[2]*interp_array(i , j-1, k-1,n) +
179  sx_hi[0]*sx_lo[1]*sx_hi[2]*interp_array(i , j-1, k ,n) +
180  sx_hi[0]*sx_hi[1]*sx_lo[2]*interp_array(i , j , k-1,n) +
181  sx_hi[0]*sx_hi[1]*sx_hi[2]*interp_array(i , j , k ,n);
182  }
183  }

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
270 {
271  m_fields[lev][0] = &vars_old[lev][Vars::xvel];
272  m_fields[lev][1] = &vars_old[lev][Vars::yvel];
273  m_fields[lev][2] = Theta_prim[lev].get();
274  m_fields[lev][3] = Qv_prim[lev].get();
275  m_fields[lev][4] = Qr_prim[lev].get();
276  m_fields[lev][5] = &vars_old[lev][Vars::zvel];
277 }

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
1628 {
1629  // Peel back the level
1630  auto& fields = m_fields[lev];
1631  auto& averages = m_averages[lev];
1632 
1633  // MFIter on CC
1634  int imf_cc = 2;
1635 
1636  int klo = m_geom[lev].Domain().smallEnd(2);
1637 
1638  int navg = m_navg - 1;
1639 
1640  std::ofstream ofile;
1641  ofile.open ("MOST_averages.txt");
1642  ofile << "Averages computed via MOSTAverages class:\n";
1643 
1644  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1645  Box pbx = mfi.tilebox();
1646 
1647  if(pbx.smallEnd(2) != klo) { continue; }
1648 
1649  pbx.makeSlab(2,klo);
1650 
1651  int il = pbx.smallEnd(0); int iu = pbx.bigEnd(0);
1652  int jl = pbx.smallEnd(1); int ju = pbx.bigEnd(1);
1653 
1654  for (int j(jl); j <= ju; ++j) {
1655  for (int i(il); i <= iu; ++i) {
1656  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1657  int k = 0;
1658  for (int iavg(0); iavg <= navg; ++iavg) {
1659  auto mf_arr = averages[iavg]->array(mfi);
1660  ofile << "iavg val: "
1661  << iavg << ' '
1662  << mf_arr(i,j,k) << "\n";
1663  }
1664  ofile << "\n";
1665  }
1666  }
1667  }
1668  ofile.close();
1669 }
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
1482 {
1483  // Peel back the level
1484  auto& fields = m_fields[lev];
1485  auto& k_indx = m_k_indx[lev];
1486 
1487  // MFIter on CC
1488  int imf_cc = 2;
1489 
1490  int klo = m_geom[lev].Domain().smallEnd(2);
1491 
1492  std::ofstream ofile;
1493  ofile.open ("MOST_k_indices.txt");
1494  ofile << "K indices used to compute averages via MOSTAverages class:\n";
1495 
1496  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1497  Box pbx = mfi.tilebox();
1498 
1499  if(pbx.smallEnd(2) != klo) { continue; }
1500 
1501  pbx.makeSlab(2,klo);
1502 
1503  int il = pbx.smallEnd(0); int iu = pbx.bigEnd(0);
1504  int jl = pbx.smallEnd(1); int ju = pbx.bigEnd(1);
1505 
1506  auto k_arr = k_indx->array(mfi);
1507 
1508  for (int j(jl); j <= ju; ++j) {
1509  for (int i(il); i <= iu; ++i) {
1510  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1511  int k = 0;
1512  ofile << "K_ind: "
1513  << k_arr(i,j,k) << "\n";
1514  ofile << "\n";
1515  }
1516  }
1517  }
1518  ofile.close();
1519 }
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
1529 {
1530  // Peel back the level
1531  auto& fields = m_fields[lev];
1532  auto& k_indx = m_k_indx[lev];
1533  auto& j_indx = m_j_indx[lev];
1534  auto& i_indx = m_i_indx[lev];
1535 
1536  // MFIter on CC
1537  int imf_cc = 2;
1538 
1539  int klo = m_geom[lev].Domain().smallEnd(2);
1540 
1541  std::ofstream ofile;
1542  ofile.open ("MOST_ijk_indices.txt");
1543  ofile << "IJK indices used to compute averages via MOSTAverages class:\n";
1544 
1545  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1546  Box pbx = mfi.tilebox();
1547 
1548  if(pbx.smallEnd(2) != klo) { continue; }
1549 
1550  pbx.makeSlab(2,klo);
1551 
1552  int il = pbx.smallEnd(0); int iu = pbx.bigEnd(0);
1553  int jl = pbx.smallEnd(1); int ju = pbx.bigEnd(1);
1554 
1555  auto k_arr = k_indx->array(mfi);
1556  auto j_arr = j_indx ? j_indx->array(mfi) : Array4<int> {};
1557  auto i_arr = i_indx ? i_indx->array(mfi) : Array4<int> {};
1558 
1559  for (int j(jl); j <= ju; ++j) {
1560  for (int i(il); i <= iu; ++i) {
1561  ofile << "(I1,J1,K1): " << "(" << i << "," << j << "," << 0 << ")" << "\n";
1562 
1563  int k = 0;
1564  int km = k_arr(i,j,k);
1565  int jm = j_arr ? j_arr(i,j,k) : j;
1566  int im = i_arr ? i_arr(i,j,k) : i;
1567 
1568  ofile << "(I2,J2,K2): "
1569  << "(" << im << "," << jm << "," << km << ")" << "\n";
1570  ofile << "\n";
1571  }
1572  }
1573  }
1574  ofile.close();
1575 }
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
1587 {
1588  // Peel back the level
1589  auto& fields = m_fields[lev];
1590  auto& x_pos_mf = m_x_pos[lev];
1591  auto& z_pos_mf = m_z_pos[lev];
1592 
1593  // MFIter on CC
1594  int imf_cc = 2;
1595 
1596  int klo = m_geom[lev].Domain().smallEnd(2);
1597 
1598  std::ofstream ofile;
1599  ofile.open ("MOST_xz_positions.txt");
1600 
1601  for (MFIter mfi(*fields[imf_cc], TileNoZ()); mfi.isValid(); ++mfi) {
1602  Box pbx = mfi.tilebox();
1603 
1604  if(pbx.smallEnd(2) != klo) { continue; }
1605 
1606  pbx.makeSlab(2,klo);
1607 
1608  int il = pbx.smallEnd(0); int iu = pbx.bigEnd(0);
1609 
1610  auto x_pos_arr = x_pos_mf->array(mfi);
1611  auto z_pos_arr = z_pos_mf->array(mfi);
1612 
1613  int k = 0;
1614  for (int i(il); i <= iu; ++i)
1615  ofile << x_pos_arr(i,j,k) << ' ' << z_pos_arr(i,j,k) << "\n";
1616  }
1617  ofile.close();
1618 }
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

Referenced by set_k_indices_N(), and set_k_indices_T().

◆ m_k_indx

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

◆ 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 {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: