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
733 {
734  if (m_rotate) set_rotated_fields(lev);
735 
736  switch(m_policy) {
737  case 0: // Standard plane average
739  break;
740  case 1: // Local region/point
742  break;
743  default:
744  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(false, "Unknown policy for MOSTAverage!");
745  }
746 
747  // We have initialized the averages
748  if (m_t_avg) m_t_init[lev] = 1;
749 }
bool m_t_avg
Definition: ERF_MOSTAverage.H:231
void compute_plane_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:758
int m_policy
Definition: ERF_MOSTAverage.H:201
void compute_region_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:1043
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
759 {
760  // Peel back the level
761  auto& fields = m_fields[lev];
762  auto& rot_fields = m_rot_fields[lev];
763  auto& averages = m_averages[lev];
764  const auto & geom = m_geom[lev];
765 
766  auto& z_phys = m_z_phys_nd[lev];
767  auto& x_pos = m_x_pos[lev];
768  auto& y_pos = m_y_pos[lev];
769  auto& z_pos = m_z_pos[lev];
770 
771  auto& i_indx = m_i_indx[lev];
772  auto& j_indx = m_j_indx[lev];
773  auto& k_indx = m_k_indx[lev];
774 
775  auto& ncell_plane = m_ncell_plane[lev];
776  auto& plane_average = m_plane_average[lev];
777 
778  // Set factors for time averaging
779  Real d_fact_new, d_fact_old;
780  if (m_t_avg && m_t_init[lev]) {
781  d_fact_new = m_fact_new;
782  d_fact_old = m_fact_old;
783  } else {
784  d_fact_new = 1.0;
785  d_fact_old = 0.0;
786  }
787 
788  // GPU array to accumulate averages into
789  Gpu::DeviceVector<Real> pavg(plane_average.size(), 0.0);
790  Real* plane_avg = pavg.data();
791 
792  // Vectors for normalization and buffer storage
793  Vector<Real> denom(plane_average.size(),0.0);
794  Vector<Real> val_old(plane_average.size(),0.0);
795 
796  //
797  //----------------------------------------------------------
798  // Averages over all the fields
799  //----------------------------------------------------------
800  //
801  Box domain = geom.Domain();
802 
803  Array<int,AMREX_SPACEDIM> is_per = {0,0,0};
804  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
805  if (geom.isPeriodic(idim)) is_per[idim] = 1;
806  }
807 
808  // Averages for U,V,T,Qv (not Qr or W)
809  for (int imf(0); imf < 4; ++imf) {
810 
811  // Continue if no valid Qv pointer
812  if (!fields[imf]) continue;
813 
814  denom[imf] = 1.0 / (Real)ncell_plane[imf];
815  val_old[imf] = plane_average[imf]*d_fact_old;
816 
817 #ifdef _OPENMP
818 #pragma omp parallel if (Gpu::notInLaunchRegion())
819 #endif
820  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
821  Box vbx = mfi.validbox(); // This is the grid (not tile)
822  Box pbx = mfi.tilebox(); // This is the tile (not grid)
823  pbx.setSmall(2,0); pbx.setBig(2,0);
824 
825  // Avoid double counting nodal data by changing the high end when we are
826  // at the high side of the grid (not just of the tile)
827  IndexType ixt = averages[imf]->boxArray().ixType();
828  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
829  if ( ixt.nodeCentered(idim) && (pbx.bigEnd(idim) == vbx.bigEnd(idim)) ) {
830  int dom_hi = domain.bigEnd(idim)+1;
831  if (pbx.bigEnd(idim) < dom_hi || is_per[idim]) {
832  pbx.growHi(idim,-1);
833  }
834  }
835  }
836 
837  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
838  fields[imf]->const_array(mfi);
839 
840  if (m_interp) {
841  const auto plo = geom.ProbLoArray();
842  const auto dxInv = geom.InvCellSizeArray();
843  const auto z_phys_arr = z_phys->const_array(mfi);
844  auto x_pos_arr = x_pos->array(mfi);
845  auto y_pos_arr = y_pos->array(mfi);
846  auto z_pos_arr = z_pos->array(mfi);
847  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
848  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
849  {
850  Real interp{0};
851  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
852  &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
853  Real val = interp;
854  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
855  });
856  } else {
857  auto k_arr = k_indx->const_array(mfi);
858  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
859  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
860  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
861  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
862  {
863  int mk = k_arr(i,j,k);
864  int mj = j_arr ? j_arr(i,j,k) : j;
865  int mi = i_arr ? i_arr(i,j,k) : i;
866  Real val = mf_arr(mi,mj,mk);
867  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
868  });
869  }
870  }
871  }
872 
873  //
874  //------------------------------------------------------------------------
875  // Averages for virtual potential temperature
876  // (This is cell-centered so we don't need to worry about double-counting)
877  //------------------------------------------------------------------------
878  //
879  if (fields[3]) // We have water vapor
880  {
881  int iavg = 4;
882  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
883  val_old[iavg] = plane_average[iavg]*d_fact_old;
884 
885 #ifdef _OPENMP
886 #pragma omp parallel if (Gpu::notInLaunchRegion())
887 #endif
888  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi)
889  {
890  Box pbx = mfi.tilebox();
891  pbx.setSmall(2,0); pbx.setBig(2,0);
892 
893  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
894  const Array4<Real const>& qv_mf_arr = (fields[3])? fields[3]->const_array(mfi) : Array4<const Real>{};
895  const Array4<Real const>& qr_mf_arr = (fields[4])? fields[4]->const_array(mfi) : Array4<const Real>{};
896 
897  if (m_interp) {
898  const auto plo = m_geom[lev].ProbLoArray();
899  const auto dxInv = m_geom[lev].InvCellSizeArray();
900  const auto z_phys_arr = z_phys->const_array(mfi);
901  auto x_pos_arr = x_pos->array(mfi);
902  auto y_pos_arr = y_pos->array(mfi);
903  auto z_pos_arr = z_pos->array(mfi);
904  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
905  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
906  {
907  Real T_interp{0};
908  Real qv_interp{0};
909  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
910  &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
911  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
912  &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
913  Real vfac;
914  if (qr_mf_arr) {
915  // We also have liquid water
916  Real qr_interp{0};
917  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
918  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
919  vfac = 1.0 + 0.61*qv_interp - qr_interp;
920  } else {
921  vfac = 1.0 + 0.61*qv_interp;
922  }
923  const Real val = T_interp * vfac;
924  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
925  });
926  } else {
927  auto k_arr = k_indx->const_array(mfi);
928  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
929  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
930  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
931  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
932  {
933  int mk = k_arr(i,j,k);
934  int mj = j_arr ? j_arr(i,j,k) : j;
935  int mi = i_arr ? i_arr(i,j,k) : i;
936  Real vfac;
937  if (qr_mf_arr) {
938  // We also have liquid water
939  vfac = 1.0 + 0.61*qv_mf_arr(mi,mj,mk) - qr_mf_arr(mi,mj,mk);
940  } else {
941  vfac = 1.0 + 0.61*qv_mf_arr(mi,mj,mk);
942  }
943  const Real val = T_mf_arr(mi,mj,mk) * vfac;
944  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
945  });
946  }
947  }
948  }
949  else // copy temperature
950  {
951  int iavg = m_navg - 2;
952  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
953  // plane_avg[iavg] = plane_avg[2]
954  Gpu::copy(Gpu::deviceToDevice, pavg.begin() + 2, pavg.begin() + 3,
955  pavg.begin() + iavg);
956  }
957 
958  //
959  //------------------------------------------------------------------------
960  // Averages for the tangential velocity magnitude
961  // (This is cell-centered so we don't need to worry about double-counting)
962  //------------------------------------------------------------------------
963  //
964  {
965  int imf = 0;
966  int iavg = m_navg - 1;
967  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
968  val_old[iavg] = plane_average[iavg]*d_fact_old;
969 
970  const Real Vsg = m_Vsg[lev];
971 
972 #ifdef _OPENMP
973 #pragma omp parallel if (Gpu::notInLaunchRegion())
974 #endif
975  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi)
976  {
977  Box pbx = mfi.tilebox();
978  pbx.setSmall(2,0); pbx.setBig(2,0);
979 
980  // Last element is Umag and always cell centered
981  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
982  fields[imf ]->const_array(mfi);
983  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
984  fields[imf+1]->const_array(mfi);
985 
986  if (m_interp) {
987  const auto plo = m_geom[lev].ProbLoArray();
988  const auto dxInv = m_geom[lev].InvCellSizeArray();
989  const auto z_phys_arr = z_phys->const_array(mfi);
990  auto x_pos_arr = x_pos->array(mfi);
991  auto y_pos_arr = y_pos->array(mfi);
992  auto z_pos_arr = z_pos->array(mfi);
993  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
994  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
995  {
996  Real u_interp{0};
997  Real v_interp{0};
998  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
999  &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
1000  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
1001  &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
1002  const Real val = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
1003  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1004  });
1005  } else {
1006  auto k_arr = k_indx->const_array(mfi);
1007  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1008  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1009  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
1010  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
1011  {
1012  int mk = k_arr(i,j,k);
1013  int mj = j_arr ? j_arr(i,j,k) : j;
1014  int mi = i_arr ? i_arr(i,j,k) : i;
1015  const Real u_val = 0.5 * (u_mf_arr(mi,mj,mk) + u_mf_arr(mi+1,mj ,mk));
1016  const Real v_val = 0.5 * (v_mf_arr(mi,mj,mk) + v_mf_arr(mi ,mj+1,mk));
1017  const Real val = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
1018  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
1019  });
1020  }
1021  }
1022  }
1023 
1024  // Copy to host and sum across procs
1025  Gpu::copy(Gpu::deviceToHost, pavg.begin(), pavg.end(), plane_average.begin());
1026  ParallelDescriptor::ReduceRealSum(plane_average.data(), plane_average.size());
1027 
1028  // No spatial variation with plane averages
1029  for (int iavg(0); iavg < m_navg; ++iavg){
1030  plane_average[iavg] *= denom[iavg]*d_fact_new;
1031  plane_average[iavg] += val_old[iavg];
1032  averages[iavg]->setVal(plane_average[iavg]);
1033  }
1034 }
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
1044 {
1045  // Peel back the level
1046  auto& fields = m_fields[lev];
1047  auto& rot_fields = m_rot_fields[lev];
1048  auto& averages = m_averages[lev];
1049  const auto & geom = m_geom[lev];
1050 
1051  auto& z_phys = m_z_phys_nd[lev];
1052  auto& x_pos = m_x_pos[lev];
1053  auto& y_pos = m_y_pos[lev];
1054  auto& z_pos = m_z_pos[lev];
1055 
1056  auto& i_indx = m_i_indx[lev];
1057  auto& j_indx = m_j_indx[lev];
1058  auto& k_indx = m_k_indx[lev];
1059 
1060  // Set factors for time averaging
1061  Real d_fact_new, d_fact_old;
1062  if (m_t_avg && m_t_init[lev]) {
1063  d_fact_new = m_fact_new;
1064  d_fact_old = m_fact_old;
1065  } else {
1066  d_fact_new = 1.0;
1067  d_fact_old = 0.0;
1068  }
1069 
1070  // Number of cells contained in the local average
1071  const Real denom = 1.0 / (Real) m_ncell_region;
1072 
1073  // Capture radius for device
1074  int d_radius = m_radius;
1075 
1076  //
1077  //----------------------------------------------------------
1078  // Averages for U,V,T,Qv
1079  //----------------------------------------------------------
1080  //
1081  for (int imf(0); imf < 4; ++imf) {
1082 
1083  // Continue if no valid Qv pointer
1084  if (!fields[imf]) continue;
1085 
1086 #ifdef _OPENMP
1087 #pragma omp parallel if (Gpu::notInLaunchRegion())
1088 #endif
1089  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
1090  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1091 
1092  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
1093  fields[imf]->const_array(mfi);
1094  auto ma_arr = averages[imf]->array(mfi);
1095 
1096  if (m_interp) {
1097  const auto plo = geom.ProbLoArray();
1098  const auto dx = geom.CellSizeArray();
1099  const auto dxInv = geom.InvCellSizeArray();
1100  const auto z_phys_arr = z_phys->const_array(mfi);
1101  auto x_pos_arr = x_pos->array(mfi);
1102  auto y_pos_arr = y_pos->array(mfi);
1103  auto z_pos_arr = z_pos->array(mfi);
1104  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1105  {
1106  ma_arr(i,j,k) *= d_fact_old;
1107 
1108  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1109  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1110  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1111  for (int li(-d_radius); li <= (d_radius); ++li) {
1112  Real interp{0};
1113  Real xp = x_pos_arr(i+li,j+lj,k);
1114  Real yp = y_pos_arr(i+li,j+lj,k);
1115  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1116  trilinear_interp_T(xp, yp, zp, &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
1117  Real val = denom * interp * d_fact_new;
1118  ma_arr(i,j,k) += val;
1119  }
1120  }
1121  }
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(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1128  {
1129  ma_arr(i,j,k) *= d_fact_old;
1130 
1131  int mk = k_arr(i,j,k);
1132  int mj = j_arr ? j_arr(i,j,k) : j;
1133  int mi = i_arr ? i_arr(i,j,k) : i;
1134  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1135  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1136  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1137  Real val = denom * mf_arr(li, lj, lk) * d_fact_new;
1138  ma_arr(i,j,k) += val;
1139  }
1140  }
1141  }
1142  });
1143  }
1144  } // MFiter
1145 
1146  // Fill interior ghost cells and any ghost cells outside a periodic domain
1147  //***********************************************************************************
1148  averages[imf]->FillBoundary(geom.periodicity());
1149 
1150  } // imf
1151 
1152  //
1153  //----------------------------------------------------------
1154  // Averages for virtual potential temperature
1155  //----------------------------------------------------------
1156  //
1157  if (fields[3]) // We have water vapor
1158  {
1159  int iavg = 4;
1160 
1161 #ifdef _OPENMP
1162 #pragma omp parallel if (Gpu::notInLaunchRegion())
1163 #endif
1164  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1165  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1166 
1167  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
1168  const Array4<Real const>& qv_mf_arr = (fields[3])? fields[3]->const_array(mfi) : Array4<const Real>{};
1169  const Array4<Real const>& qr_mf_arr = (fields[4])? fields[4]->const_array(mfi) : Array4<const Real>{};
1170  auto ma_arr = averages[iavg]->array(mfi);
1171 
1172  if (m_interp) {
1173  const auto plo = geom.ProbLoArray();
1174  const auto dx = geom.CellSizeArray();
1175  const auto dxInv = geom.InvCellSizeArray();
1176  const auto z_phys_arr = z_phys->const_array(mfi);
1177  auto x_pos_arr = x_pos->array(mfi);
1178  auto y_pos_arr = y_pos->array(mfi);
1179  auto z_pos_arr = z_pos->array(mfi);
1180  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1181  {
1182  ma_arr(i,j,k) *= d_fact_old;
1183 
1184  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1185  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1186  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1187  for (int li(-d_radius); li <= (d_radius); ++li) {
1188  Real T_interp{0};
1189  Real qv_interp{0};
1190  Real xp = x_pos_arr(i+li,j+lj,k);
1191  Real yp = y_pos_arr(i+li,j+lj,k);
1192  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1193  trilinear_interp_T(xp, yp, zp, &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
1194  trilinear_interp_T(xp, yp, zp, &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
1195  Real vfac;
1196  if (qr_mf_arr) {
1197  // We also have liquid water
1198  Real qr_interp{0};
1199  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
1200  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
1201  vfac = 1.0 + 0.61*qv_interp - qr_interp;
1202  } else {
1203  vfac = 1.0 + 0.61*qv_interp;
1204  }
1205  const Real mag = T_interp * vfac;
1206  const Real val = denom * mag * d_fact_new;
1207  ma_arr(i,j,k) += val;
1208  }
1209  }
1210  }
1211  });
1212  } else {
1213  auto k_arr = k_indx->const_array(mfi);
1214  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1215  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1216  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1217  {
1218  ma_arr(i,j,k) *= d_fact_old;
1219 
1220  int mk = k_arr(i,j,k);
1221  int mj = j_arr ? j_arr(i,j,k) : j;
1222  int mi = i_arr ? i_arr(i,j,k) : i;
1223  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1224  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1225  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1226  Real vfac;
1227  if (qr_mf_arr) {
1228  // We also have liquid water
1229  vfac = 1.0 + 0.61*qv_mf_arr(li,lj,lk) - qr_mf_arr(li,lj,lk);
1230  } else {
1231  vfac = 1.0 + 0.61*qv_mf_arr(li,lj,lk);
1232  }
1233  const Real mag = T_mf_arr(li,lj,lk) * vfac;
1234  const Real val = denom * mag * d_fact_new;
1235  ma_arr(i,j,k) += val;
1236  }
1237  }
1238  }
1239  });
1240  }
1241  } // MFiter
1242 
1243  // Fill interior ghost cells and any ghost cells outside a periodic domain
1244  //***********************************************************************************
1245  averages[iavg]->FillBoundary(geom.periodicity());
1246 
1247  }
1248  else // copy temperature
1249  {
1250  int iavg = m_navg - 2;
1251  IntVect ng = averages[iavg]->nGrowVect();
1252  MultiFab::Copy(*(averages[iavg]),*(averages[2]),0,0,1,ng);
1253  }
1254 
1255  //
1256  //----------------------------------------------------------
1257  // Averages for the tangential velocity magnitude
1258  //----------------------------------------------------------
1259  //
1260  {
1261  int imf = 0;
1262  int iavg = m_navg - 1;
1263 
1264  const Real Vsg = m_Vsg[lev];
1265 
1266 #ifdef _OPENMP
1267 #pragma omp parallel if (Gpu::notInLaunchRegion())
1268 #endif
1269  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1270  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1271 
1272  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
1273  fields[imf ]->const_array(mfi);
1274  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
1275  fields[imf+1]->const_array(mfi);
1276  auto ma_arr = averages[iavg]->array(mfi);
1277 
1278  if (m_interp) {
1279  const auto plo = geom.ProbLoArray();
1280  const auto dx = geom.CellSizeArray();
1281  const auto dxInv = geom.InvCellSizeArray();
1282  const auto z_phys_arr = z_phys->const_array(mfi);
1283  auto x_pos_arr = x_pos->array(mfi);
1284  auto y_pos_arr = y_pos->array(mfi);
1285  auto z_pos_arr = z_pos->array(mfi);
1286  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1287  {
1288  ma_arr(i,j,k) *= d_fact_old;
1289 
1290  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1291  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1292  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1293  for (int li(-d_radius); li <= (d_radius); ++li) {
1294  Real u_interp{0};
1295  Real v_interp{0};
1296  Real xp = x_pos_arr(i+li,j+lj,k);
1297  Real yp = y_pos_arr(i+li,j+lj,k);
1298  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1299  trilinear_interp_T(xp, yp, zp, &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
1300  trilinear_interp_T(xp, yp, zp, &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
1301  const Real mag = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
1302  Real val = denom * mag * d_fact_new;
1303  ma_arr(i,j,k) += val;
1304  }
1305  }
1306  }
1307  });
1308  } else {
1309  auto k_arr = k_indx->const_array(mfi);
1310  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1311  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1312  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1313  {
1314  ma_arr(i,j,k) *= d_fact_old;
1315 
1316  int mk = k_arr(i,j,k);
1317  int mj = j_arr ? j_arr(i,j,k) : j;
1318  int mi = i_arr ? i_arr(i,j,k) : i;
1319  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1320  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1321  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1322  const Real u_val = 0.5 * (u_mf_arr(li,lj,lk) + u_mf_arr(li+1,lj ,lk));
1323  const Real v_val = 0.5 * (v_mf_arr(li,lj,lk) + v_mf_arr(li ,lj+1,lk));
1324  const Real mag = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
1325  Real val = denom * mag * d_fact_new;
1326  ma_arr(i,j,k) += val;
1327  }
1328  }
1329  }
1330  });
1331  }
1332  } // MFiter
1333 
1334  // Fill interior ghost cells and any ghost cells outside a periodic domain
1335  //***********************************************************************************
1336  averages[iavg]->FillBoundary(geom.periodicity());
1337 
1338  }
1339 
1340  // NOTE: Checking periodicity with the geom structure is not
1341  // sufficient at higher levels. The BA may be contained
1342  // within the domain and it's exterior ghost cells filled
1343  // from interpolation; yet the domain BCs are periodic.
1344 
1345  // Need to fill ghost cells outside the domain if not periodic
1346  bool not_per_x = !(geom.periodicity().isPeriodic(0));
1347  bool not_per_y = !(geom.periodicity().isPeriodic(1));
1348  Box cc_bnd_bx = (m_fields[lev][2]->boxArray()).minimalBox();
1349  Box domain = geom.Domain();
1350  if (domain.contains(cc_bnd_bx) || (not_per_x || not_per_y)) {
1351  for (int iavg(0); iavg < m_navg; ++iavg) {
1352  IntVect ng = averages[iavg]->nGrowVect(); ng[2]=0;
1353 
1354  // NOTE: Level 0 spans the whole domain, but finer
1355  // levels do not have such a restriction.
1356  // For now, use the bounding box of the boxArray.
1357 
1358  // NOTE2: The fields and averages have different indexing.
1359  // The averages are: U/V/T/Qv/Tv/Umag
1360  // The fields are: U/V/T/Qv/Qr/W
1361  // We clip iavg at 2 since all the remaining data is CC
1362 
1363  // Bounded box of CC data used for normalization
1364  int imf = min(iavg,2);
1365  Box bnd_bx = (m_fields[lev][imf]->boxArray()).minimalBox();
1366 #ifdef _OPENMP
1367 #pragma omp parallel if (Gpu::notInLaunchRegion())
1368 #endif
1369  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1370  Box gpbx = mfi.growntilebox(ng); gpbx.setSmall(2,0); gpbx.setBig(2,0);
1371 
1372  if (bnd_bx.contains(gpbx)) continue;
1373 
1374  auto ma_arr = averages[iavg]->array(mfi);
1375 
1376  int i_lo = bnd_bx.smallEnd(0); int i_hi = bnd_bx.bigEnd(0);
1377  int j_lo = bnd_bx.smallEnd(1); int j_hi = bnd_bx.bigEnd(1);
1378  ParallelFor(gpbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1379  {
1380  int li, lj;
1381  li = i < i_lo ? i_lo : i;
1382  li = li > i_hi ? i_hi : li;
1383  lj = j < j_lo ? j_lo : j;
1384  lj = lj > j_hi ? j_hi : lj;
1385 
1386  ma_arr(i,j,k) = ma_arr(li,lj,k);
1387  });
1388  } // MFiter
1389  } // iavg
1390  } // Not periodic
1391 }
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>(ba,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][4] = 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:603
void set_norm_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:654
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:526
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:141
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142

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

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

527 {
528  ParmParse pp(m_pp_prefix);
529  Real zref_tmp = zref_default;
530  auto read_zref = pp.query("most.zref",zref_tmp);
531  if (!read_zref) {
532  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
533  }
534 
535  // Capture for device
536  Real d_zref = zref_tmp;
537  Real d_radius = m_radius;
538 
539  int kmax = m_geom[lev].Domain().bigEnd(2);
540  const auto dxInv = m_geom[lev].InvCellSizeArray();
541  IntVect ng = m_k_indx[lev]->nGrowVect(); ng[2]=0;
542  for (MFIter mfi(*m_k_indx[lev], TileNoZ()); mfi.isValid(); ++mfi) {
543  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
544  Box gpbx = mfi.growntilebox(ng);
545  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
546  auto i_arr = m_i_indx[lev]->array(mfi);
547  auto j_arr = m_j_indx[lev]->array(mfi);
548  auto k_arr = m_k_indx[lev]->array(mfi);
549  auto zref_arr = m_zref[lev]->array(mfi);
550  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
551  {
552  // Elements of normal vector
553  Real met_h_xi = Compute_h_xi_AtCellCenter (i,j,k,dxInv,z_phys_arr);
554  Real met_h_eta = Compute_h_eta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
555  Real mag = std::sqrt(met_h_xi*met_h_xi + met_h_eta*met_h_eta + 1.0);
556 
557  // Unit-normal vector scaled by z_ref
558  Real delta_x = -met_h_xi/mag * d_zref;
559  Real delta_y = -met_h_eta/mag * d_zref;
560  Real delta_z = 1.0/mag * d_zref;
561 
562  // Compute i & j as displacements (no grid stretching)
563  int delta_i = static_cast<int>(std::round(delta_x*dxInv[0]));
564  int delta_j = static_cast<int>(std::round(delta_y*dxInv[1]));
565  int i_new = i + delta_i;
566  int j_new = j + delta_j;
567  i_arr(i,j,k) = i_new;
568  j_arr(i,j,k) = j_new;
569 
570  // Search for k (grid is stretched in z)
571  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
572  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
573  Real z_target = z_bot_face + delta_z;
574  for (int lk(0); lk<=kmax; ++lk) {
575  Real z_lo = 0.25 * ( z_phys_arr(i_new,j_new ,lk ) + z_phys_arr(i_new+1,j_new ,lk )
576  + z_phys_arr(i_new,j_new+1,lk ) + z_phys_arr(i_new+1,j_new+1,lk ) );
577  Real z_hi = 0.25 * ( z_phys_arr(i_new,j_new ,lk+1) + z_phys_arr(i_new+1,j_new ,lk+1)
578  + z_phys_arr(i_new,j_new+1,lk+1) + z_phys_arr(i_new+1,j_new+1,lk+1) );
579  if (z_target > z_lo && z_target < z_hi){
580  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(lk >= d_radius,
581  "K index must be larger than averaging radius!");
582  amrex::ignore_unused(d_radius);
583  k_arr(i,j,k) = lk;
584  zref_arr(i,j,k) = 0.5 * (z_hi + z_lo) - z_bot_face;
585  break;
586  }
587  }
588 
589  // Destination cell must be contained on the current process!
590  amrex::ignore_unused(gpbx);
591  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(gpbx.contains(i_arr(i,j,k),j_arr(i,j,k),k_arr(i,j,k)),
592  "Query index outside of proc domain!");
593  });
594  }
595 }
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).

655 {
656  ParmParse pp(m_pp_prefix);
657  Real zref_tmp = zref_default;
658  auto read_zref = pp.query("most.zref",zref_tmp);
659  if (!read_zref) {
660  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
661  }
662 
663  // Capture for device
664  Real d_zref = zref_tmp;
665  const auto plo = m_geom[lev].ProbLoArray();
666 
667  RealVect base;
668  const auto dx = m_geom[lev].CellSizeArray();
669  const auto dxInv = m_geom[lev].InvCellSizeArray();
670  IntVect ng = m_x_pos[lev]->nGrowVect(); ng[2]=0;
671  for (MFIter mfi(*m_x_pos[lev], TileNoZ()); mfi.isValid(); ++mfi) {
672  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
673  Box gpbx = mfi.growntilebox(ng);
674  RealBox grb{gpbx,dx.data(),base.dataPtr()};
675 
676  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
677  auto x_pos_arr = m_x_pos[lev]->array(mfi);
678  auto y_pos_arr = m_y_pos[lev]->array(mfi);
679  auto z_pos_arr = m_z_pos[lev]->array(mfi);
680  auto zref_arr = m_zref[lev]->array(mfi);
681  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
682  {
683  // Elements of normal vector
684  Real met_h_xi = Compute_h_xi_AtCellCenter (i,j,k,dxInv,z_phys_arr);
685  Real met_h_eta = Compute_h_eta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
686  Real imag = 1.0 / std::sqrt(met_h_xi*met_h_xi + met_h_eta*met_h_eta + 1.0);
687 
688  // Unit-normal vector scaled by z_ref
689  Real delta_x = -met_h_xi * imag * d_zref;
690  Real delta_y = -met_h_eta * imag * d_zref;
691  Real delta_z = imag * d_zref;
692 
693  // Position of the current node (indx:0,0,1)
694  Real x0 = plo[0] + ((Real) i + 0.5) * dx[0];
695  Real y0 = plo[1] + ((Real) j + 0.5) * dx[1];
696 
697  // Final position at end of vector
698  x_pos_arr(i,j,k) = x0 + delta_x;
699  y_pos_arr(i,j,k) = y0 + delta_y;
700  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
701  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
702  z_pos_arr(i,j,k) = z_bot_face + delta_z;
703 
704  // NOTE: Normal vector end point can be below the surface for concave regions.
705  // Here we protect against that by augmenting the normal if needed.
706  int i_new = (int) ((x_pos_arr(i,j,k) - plo[0]) / dx[0] - 0.5);
707  int j_new = (int) ((y_pos_arr(i,j,k) - plo[1]) / dx[1] - 0.5);
708  Real z_new_bot_face = 0.25 * ( z_phys_arr(i_new,j_new ,k) + z_phys_arr(i_new+1,j_new ,k)
709  + z_phys_arr(i_new,j_new+1,k) + z_phys_arr(i_new+1,j_new+1,k) );
710  if (z_pos_arr(i,j,k) < z_new_bot_face) {
711  z_pos_arr(i,j,k) = z_new_bot_face + delta_z;
712  }
713 
714  zref_arr(i,j,k) = delta_z;
715 
716  // Destination position must be contained on the current process!
717  Real pos[] = {x_pos_arr(i,j,k)-plo[0],y_pos_arr(i,j,k)-plo[1],0.5*dx[2]};
718  amrex::ignore_unused(pos);
719  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(grb.contains(&pos[0]),
720  "Query point outside of proc domain!");
721  });
722  }
723 }

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[4]->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)

604 {
605  ParmParse pp(m_pp_prefix);
606  Real zref_tmp = zref_default;
607  auto read_zref = pp.query("most.zref",zref_tmp);
608  if (!read_zref) {
609  Print() << "most.zref not specified, query distance default is " << zref_tmp << std::endl;
610  } else {
611  m_zref[lev]->setVal(zref_tmp);
612  }
613 
614  // Capture for device
615  Real d_zref = zref_tmp;
616  const auto plo = m_geom[lev].ProbLoArray();
617 
618  RealVect base;
619  const auto dx = m_geom[lev].CellSizeArray();
620  IntVect ng = m_x_pos[lev]->nGrowVect(); ng[2]=0;
621  for (MFIter mfi(*m_x_pos[lev], TileNoZ()); mfi.isValid(); ++mfi) {
622  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
623  Box gpbx = mfi.growntilebox(ng);
624  RealBox grb{gpbx,dx.data(),base.dataPtr()};
625 
626  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
627  auto x_pos_arr = m_x_pos[lev]->array(mfi);
628  auto y_pos_arr = m_y_pos[lev]->array(mfi);
629  auto z_pos_arr = m_z_pos[lev]->array(mfi);
630  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
631  {
632  // Final position at end of vector
633  x_pos_arr(i,j,k) = plo[0] + ((Real) i + 0.5) * dx[0];
634  y_pos_arr(i,j,k) = plo[1] + ((Real) j + 0.5) * dx[1];
635  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
636  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
637  z_pos_arr(i,j,k) = z_bot_face + d_zref;
638 
639  // Destination position must be contained on the current process!
640  Real pos[] = {x_pos_arr(i,j,k)-plo[0],y_pos_arr(i,j,k)-plo[1],0.5*dx[2]};
641  amrex::ignore_unused(pos);
642  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(grb.contains(&pos[0]),
643  "Query point outside of proc domain!");
644  });
645  }
646 }

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
1520 {
1521  // Peel back the level
1522  auto& averages = m_averages[lev];
1523 
1524  int navg = m_navg - 1;
1525 
1526  std::ofstream ofile;
1527  ofile.open ("MOST_averages.txt");
1528  ofile << "Averages computed via MOSTAverages class:\n";
1529 
1530  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1531  Box bx = mfi.tilebox(); bx.setBig(2,0);
1532  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1533  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1534 
1535  for (int j(jl); j <= ju; ++j) {
1536  for (int i(il); i <= iu; ++i) {
1537  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1538  int k = 0;
1539  for (int iavg(0); iavg <= navg; ++iavg) {
1540  auto mf_arr = averages[iavg]->array(mfi);
1541  ofile << "iavg val: "
1542  << iavg << ' '
1543  << mf_arr(i,j,k) << "\n";
1544  }
1545  ofile << "\n";
1546  }
1547  }
1548  }
1549  ofile.close();
1550 }
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
1401 {
1402  // Peel back the level
1403  auto& averages = m_averages[lev];
1404  auto& k_indx = m_k_indx[lev];
1405 
1406  int navg = m_navg - 1;
1407 
1408  std::ofstream ofile;
1409  ofile.open ("MOST_k_indices.txt");
1410  ofile << "K indices used to compute averages via MOSTAverages class:\n";
1411 
1412  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1413  Box bx = mfi.tilebox(); bx.setBig(2,0);
1414  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1415  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1416 
1417  auto k_arr = k_indx->array(mfi);
1418 
1419  for (int j(jl); j <= ju; ++j) {
1420  for (int i(il); i <= iu; ++i) {
1421  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1422  int k = 0;
1423  ofile << "K_ind: "
1424  << k_arr(i,j,k) << "\n";
1425  ofile << "\n";
1426  }
1427  }
1428  }
1429  ofile.close();
1430 }
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
1440 {
1441  // Peel back the level
1442  auto& averages = m_averages[lev];
1443  auto& k_indx = m_k_indx[lev];
1444  auto& j_indx = m_j_indx[lev];
1445  auto& i_indx = m_i_indx[lev];
1446 
1447  int navg = m_navg - 1;
1448 
1449  std::ofstream ofile;
1450  ofile.open ("MOST_ijk_indices.txt");
1451  ofile << "IJK indices used to compute averages via MOSTAverages class:\n";
1452 
1453  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1454  Box bx = mfi.tilebox(); bx.setBig(2,0);
1455  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1456  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1457 
1458  auto k_arr = k_indx->array(mfi);
1459  auto j_arr = j_indx ? j_indx->array(mfi) : Array4<int> {};
1460  auto i_arr = i_indx ? i_indx->array(mfi) : Array4<int> {};
1461 
1462  for (int j(jl); j <= ju; ++j) {
1463  for (int i(il); i <= iu; ++i) {
1464  ofile << "(I1,J1,K1): " << "(" << i << "," << j << "," << 0 << ")" << "\n";
1465 
1466  int k = 0;
1467  int km = k_arr(i,j,k);
1468  int jm = j_arr ? j_arr(i,j,k) : j;
1469  int im = i_arr ? i_arr(i,j,k) : i;
1470 
1471  ofile << "(I2,J2,K2): "
1472  << "(" << im << "," << jm << "," << km << ")" << "\n";
1473  ofile << "\n";
1474  }
1475  }
1476  }
1477  ofile.close();
1478 }
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
1490 {
1491  // Peel back the level
1492  auto& x_pos_mf = m_x_pos[lev];
1493  auto& z_pos_mf = m_z_pos[lev];
1494 
1495  std::ofstream ofile;
1496  ofile.open ("MOST_xz_positions.txt");
1497 
1498  for (MFIter mfi(*x_pos_mf, TileNoZ()); mfi.isValid(); ++mfi) {
1499  Box bx = mfi.tilebox(); bx.setBig(2,0);
1500  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1501 
1502  auto x_pos_arr = x_pos_mf->array(mfi);
1503  auto z_pos_arr = z_pos_mf->array(mfi);
1504 
1505  int k = 0;
1506  for (int i(il); i <= iu; ++i)
1507  ofile << x_pos_arr(i,j,k) << ' ' << z_pos_arr(i,j,k) << "\n";
1508  }
1509  ofile.close();
1510 }
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

amrex::Vector<amrex::Vector<amrex::MultiFab*> > MOSTAverage::m_fields
protected

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