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::Real get_zref () 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::Real m_zref {10.0}
 
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
 

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
709 {
710  if (m_rotate) set_rotated_fields(lev);
711 
712  switch(m_policy) {
713  case 0: // Standard plane average
715  break;
716  case 1: // Local region/point
718  break;
719  default:
720  AMREX_ASSERT_WITH_MESSAGE(false, "Unknown policy for MOSTAverage!");
721  }
722 
723  // We have initialized the averages
724  if (m_t_avg) m_t_init[lev] = 1;
725 }
bool m_t_avg
Definition: ERF_MOSTAverage.H:231
void compute_plane_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:734
int m_policy
Definition: ERF_MOSTAverage.H:201
void compute_region_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:1019
amrex::Vector< int > m_t_init
Definition: ERF_MOSTAverage.H:232
void set_rotated_fields(const int &lev)
Definition: ERF_MOSTAverage.cpp:279
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
735 {
736  // Peel back the level
737  auto& fields = m_fields[lev];
738  auto& rot_fields = m_rot_fields[lev];
739  auto& averages = m_averages[lev];
740  const auto & geom = m_geom[lev];
741 
742  auto& z_phys = m_z_phys_nd[lev];
743  auto& x_pos = m_x_pos[lev];
744  auto& y_pos = m_y_pos[lev];
745  auto& z_pos = m_z_pos[lev];
746 
747  auto& i_indx = m_i_indx[lev];
748  auto& j_indx = m_j_indx[lev];
749  auto& k_indx = m_k_indx[lev];
750 
751  auto& ncell_plane = m_ncell_plane[lev];
752  auto& plane_average = m_plane_average[lev];
753 
754  // Set factors for time averaging
755  Real d_fact_new, d_fact_old;
756  if (m_t_avg && m_t_init[lev]) {
757  d_fact_new = m_fact_new;
758  d_fact_old = m_fact_old;
759  } else {
760  d_fact_new = 1.0;
761  d_fact_old = 0.0;
762  }
763 
764  // GPU array to accumulate averages into
765  Gpu::DeviceVector<Real> pavg(plane_average.size(), 0.0);
766  Real* plane_avg = pavg.data();
767 
768  // Vectors for normalization and buffer storage
769  Vector<Real> denom(plane_average.size(),0.0);
770  Vector<Real> val_old(plane_average.size(),0.0);
771 
772  //
773  //----------------------------------------------------------
774  // Averages over all the fields
775  //----------------------------------------------------------
776  //
777  Box domain = geom.Domain();
778 
779  Array<int,AMREX_SPACEDIM> is_per = {0,0,0};
780  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
781  if (geom.isPeriodic(idim)) is_per[idim] = 1;
782  }
783 
784  // Averages for U,V,T,Qv (not Qr or W)
785  for (int imf(0); imf < 4; ++imf) {
786 
787  // Continue if no valid Qv pointer
788  if (!fields[imf]) continue;
789 
790  denom[imf] = 1.0 / (Real)ncell_plane[imf];
791  val_old[imf] = plane_average[imf]*d_fact_old;
792 
793 #ifdef _OPENMP
794 #pragma omp parallel if (Gpu::notInLaunchRegion())
795 #endif
796  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
797  Box vbx = mfi.validbox(); // This is the grid (not tile)
798  Box pbx = mfi.tilebox(); // This is the tile (not grid)
799  pbx.setSmall(2,0); pbx.setBig(2,0);
800 
801  // Avoid double counting nodal data by changing the high end when we are
802  // at the high side of the grid (not just of the tile)
803  IndexType ixt = averages[imf]->boxArray().ixType();
804  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
805  if ( ixt.nodeCentered(idim) && (pbx.bigEnd(idim) == vbx.bigEnd(idim)) ) {
806  int dom_hi = domain.bigEnd(idim)+1;
807  if (pbx.bigEnd(idim) < dom_hi || is_per[idim]) {
808  pbx.growHi(idim,-1);
809  }
810  }
811  }
812 
813  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
814  fields[imf]->const_array(mfi);
815 
816  if (m_interp) {
817  const auto plo = geom.ProbLoArray();
818  const auto dxInv = geom.InvCellSizeArray();
819  const auto z_phys_arr = z_phys->const_array(mfi);
820  auto x_pos_arr = x_pos->array(mfi);
821  auto y_pos_arr = y_pos->array(mfi);
822  auto z_pos_arr = z_pos->array(mfi);
823  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
824  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
825  {
826  Real interp{0};
827  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
828  &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
829  Real val = interp;
830  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
831  });
832  } else {
833  auto k_arr = k_indx->const_array(mfi);
834  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
835  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
836  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
837  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
838  {
839  int mk = k_arr(i,j,k);
840  int mj = j_arr ? j_arr(i,j,k) : j;
841  int mi = i_arr ? i_arr(i,j,k) : i;
842  Real val = mf_arr(mi,mj,mk);
843  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
844  });
845  }
846  }
847  }
848 
849  //
850  //------------------------------------------------------------------------
851  // Averages for virtual potential temperature
852  // (This is cell-centered so we don't need to worry about double-counting)
853  //------------------------------------------------------------------------
854  //
855  if (fields[3]) // We have water vapor
856  {
857  int iavg = 4;
858  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
859  val_old[iavg] = plane_average[iavg]*d_fact_old;
860 
861 #ifdef _OPENMP
862 #pragma omp parallel if (Gpu::notInLaunchRegion())
863 #endif
864  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi)
865  {
866  Box pbx = mfi.tilebox();
867  pbx.setSmall(2,0); pbx.setBig(2,0);
868 
869  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
870  const Array4<Real const>& qv_mf_arr = (fields[3])? fields[3]->const_array(mfi) : Array4<const Real>{};
871  const Array4<Real const>& qr_mf_arr = (fields[4])? fields[4]->const_array(mfi) : Array4<const Real>{};
872 
873  if (m_interp) {
874  const auto plo = m_geom[lev].ProbLoArray();
875  const auto dxInv = m_geom[lev].InvCellSizeArray();
876  const auto z_phys_arr = z_phys->const_array(mfi);
877  auto x_pos_arr = x_pos->array(mfi);
878  auto y_pos_arr = y_pos->array(mfi);
879  auto z_pos_arr = z_pos->array(mfi);
880  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
881  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
882  {
883  Real T_interp{0};
884  Real qv_interp{0};
885  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
886  &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
887  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
888  &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
889  Real vfac;
890  if (qr_mf_arr) {
891  // We also have liquid water
892  Real qr_interp{0};
893  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
894  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
895  vfac = 1.0 + 0.61*qv_interp - qr_interp;
896  } else {
897  vfac = 1.0 + 0.61*qv_interp;
898  }
899  const Real val = T_interp * vfac;
900  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
901  });
902  } else {
903  auto k_arr = k_indx->const_array(mfi);
904  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
905  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
906  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
907  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
908  {
909  int mk = k_arr(i,j,k);
910  int mj = j_arr ? j_arr(i,j,k) : j;
911  int mi = i_arr ? i_arr(i,j,k) : i;
912  Real vfac;
913  if (qr_mf_arr) {
914  // We also have liquid water
915  vfac = 1.0 + 0.61*qv_mf_arr(mi,mj,mk) - qr_mf_arr(mi,mj,mk);
916  } else {
917  vfac = 1.0 + 0.61*qv_mf_arr(mi,mj,mk);
918  }
919  const Real val = T_mf_arr(mi,mj,mk) * vfac;
920  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
921  });
922  }
923  }
924  }
925  else // copy temperature
926  {
927  int iavg = m_navg - 2;
928  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
929  // plane_avg[iavg] = plane_avg[2]
930  Gpu::copy(Gpu::deviceToDevice, pavg.begin() + 2, pavg.begin() + 3,
931  pavg.begin() + iavg);
932  }
933 
934  //
935  //------------------------------------------------------------------------
936  // Averages for the tangential velocity magnitude
937  // (This is cell-centered so we don't need to worry about double-counting)
938  //------------------------------------------------------------------------
939  //
940  {
941  int imf = 0;
942  int iavg = m_navg - 1;
943  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
944  val_old[iavg] = plane_average[iavg]*d_fact_old;
945 
946  const Real Vsg = m_Vsg[lev];
947 
948 #ifdef _OPENMP
949 #pragma omp parallel if (Gpu::notInLaunchRegion())
950 #endif
951  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi)
952  {
953  Box pbx = mfi.tilebox();
954  pbx.setSmall(2,0); pbx.setBig(2,0);
955 
956  // Last element is Umag and always cell centered
957  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
958  fields[imf ]->const_array(mfi);
959  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
960  fields[imf+1]->const_array(mfi);
961 
962  if (m_interp) {
963  const auto plo = m_geom[lev].ProbLoArray();
964  const auto dxInv = m_geom[lev].InvCellSizeArray();
965  const auto z_phys_arr = z_phys->const_array(mfi);
966  auto x_pos_arr = x_pos->array(mfi);
967  auto y_pos_arr = y_pos->array(mfi);
968  auto z_pos_arr = z_pos->array(mfi);
969  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
970  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
971  {
972  Real u_interp{0};
973  Real v_interp{0};
974  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
975  &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
976  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
977  &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
978  const Real val = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
979  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
980  });
981  } else {
982  auto k_arr = k_indx->const_array(mfi);
983  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
984  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
985  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
986  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
987  {
988  int mk = k_arr(i,j,k);
989  int mj = j_arr ? j_arr(i,j,k) : j;
990  int mi = i_arr ? i_arr(i,j,k) : i;
991  const Real u_val = 0.5 * (u_mf_arr(mi,mj,mk) + u_mf_arr(mi+1,mj ,mk));
992  const Real v_val = 0.5 * (v_mf_arr(mi,mj,mk) + v_mf_arr(mi ,mj+1,mk));
993  const Real val = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
994  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
995  });
996  }
997  }
998  }
999 
1000  // Copy to host and sum across procs
1001  Gpu::copy(Gpu::deviceToHost, pavg.begin(), pavg.end(), plane_average.begin());
1002  ParallelDescriptor::ReduceRealSum(plane_average.data(), plane_average.size());
1003 
1004  // No spatial variation with plane averages
1005  for (int iavg(0); iavg < m_navg; ++iavg){
1006  plane_average[iavg] *= denom[iavg]*d_fact_new;
1007  plane_average[iavg] += val_old[iavg];
1008  averages[iavg]->setVal(plane_average[iavg]);
1009  }
1010 }
const auto & dom_hi
Definition: ERF_DiffSetup.H:10
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
1020 {
1021  // Peel back the level
1022  auto& fields = m_fields[lev];
1023  auto& rot_fields = m_rot_fields[lev];
1024  auto& averages = m_averages[lev];
1025  const auto & geom = m_geom[lev];
1026 
1027  auto& z_phys = m_z_phys_nd[lev];
1028  auto& x_pos = m_x_pos[lev];
1029  auto& y_pos = m_y_pos[lev];
1030  auto& z_pos = m_z_pos[lev];
1031 
1032  auto& i_indx = m_i_indx[lev];
1033  auto& j_indx = m_j_indx[lev];
1034  auto& k_indx = m_k_indx[lev];
1035 
1036  // Set factors for time averaging
1037  Real d_fact_new, d_fact_old;
1038  if (m_t_avg && m_t_init[lev]) {
1039  d_fact_new = m_fact_new;
1040  d_fact_old = m_fact_old;
1041  } else {
1042  d_fact_new = 1.0;
1043  d_fact_old = 0.0;
1044  }
1045 
1046  // Number of cells contained in the local average
1047  const Real denom = 1.0 / (Real) m_ncell_region;
1048 
1049  // Capture radius for device
1050  int d_radius = m_radius;
1051 
1052  //
1053  //----------------------------------------------------------
1054  // Averages for U,V,T,Qv
1055  //----------------------------------------------------------
1056  //
1057  for (int imf(0); imf < 4; ++imf) {
1058 
1059  // Continue if no valid Qv pointer
1060  if (!fields[imf]) continue;
1061 
1062 #ifdef _OPENMP
1063 #pragma omp parallel if (Gpu::notInLaunchRegion())
1064 #endif
1065  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
1066  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1067 
1068  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
1069  fields[imf]->const_array(mfi);
1070  auto ma_arr = averages[imf]->array(mfi);
1071 
1072  if (m_interp) {
1073  const auto plo = geom.ProbLoArray();
1074  const auto dx = geom.CellSizeArray();
1075  const auto dxInv = geom.InvCellSizeArray();
1076  const auto z_phys_arr = z_phys->const_array(mfi);
1077  auto x_pos_arr = x_pos->array(mfi);
1078  auto y_pos_arr = y_pos->array(mfi);
1079  auto z_pos_arr = z_pos->array(mfi);
1080  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1081  {
1082  ma_arr(i,j,k) *= d_fact_old;
1083 
1084  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1085  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1086  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1087  for (int li(-d_radius); li <= (d_radius); ++li) {
1088  Real interp{0};
1089  Real xp = x_pos_arr(i+li,j+lj,k);
1090  Real yp = y_pos_arr(i+li,j+lj,k);
1091  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1092  trilinear_interp_T(xp, yp, zp, &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
1093  Real val = denom * interp * d_fact_new;
1094  ma_arr(i,j,k) += val;
1095  }
1096  }
1097  }
1098  });
1099  } else {
1100  auto k_arr = k_indx->const_array(mfi);
1101  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1102  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1103  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1104  {
1105  ma_arr(i,j,k) *= d_fact_old;
1106 
1107  int mk = k_arr(i,j,k);
1108  int mj = j_arr ? j_arr(i,j,k) : j;
1109  int mi = i_arr ? i_arr(i,j,k) : i;
1110  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1111  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1112  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1113  Real val = denom * mf_arr(li, lj, lk) * d_fact_new;
1114  ma_arr(i,j,k) += val;
1115  }
1116  }
1117  }
1118  });
1119  }
1120  } // MFiter
1121 
1122  // Fill interior ghost cells and any ghost cells outside a periodic domain
1123  //***********************************************************************************
1124  averages[imf]->FillBoundary(geom.periodicity());
1125 
1126  } // imf
1127 
1128  //
1129  //----------------------------------------------------------
1130  // Averages for virtual potential temperature
1131  //----------------------------------------------------------
1132  //
1133  if (fields[3]) // We have water vapor
1134  {
1135  int iavg = 4;
1136 
1137 #ifdef _OPENMP
1138 #pragma omp parallel if (Gpu::notInLaunchRegion())
1139 #endif
1140  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1141  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1142 
1143  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
1144  const Array4<Real const>& qv_mf_arr = (fields[3])? fields[3]->const_array(mfi) : Array4<const Real>{};
1145  const Array4<Real const>& qr_mf_arr = (fields[4])? fields[4]->const_array(mfi) : Array4<const Real>{};
1146  auto ma_arr = averages[iavg]->array(mfi);
1147 
1148  if (m_interp) {
1149  const auto plo = geom.ProbLoArray();
1150  const auto dx = geom.CellSizeArray();
1151  const auto dxInv = geom.InvCellSizeArray();
1152  const auto z_phys_arr = z_phys->const_array(mfi);
1153  auto x_pos_arr = x_pos->array(mfi);
1154  auto y_pos_arr = y_pos->array(mfi);
1155  auto z_pos_arr = z_pos->array(mfi);
1156  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1157  {
1158  ma_arr(i,j,k) *= d_fact_old;
1159 
1160  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1161  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1162  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1163  for (int li(-d_radius); li <= (d_radius); ++li) {
1164  Real T_interp{0};
1165  Real qv_interp{0};
1166  Real xp = x_pos_arr(i+li,j+lj,k);
1167  Real yp = y_pos_arr(i+li,j+lj,k);
1168  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1169  trilinear_interp_T(xp, yp, zp, &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
1170  trilinear_interp_T(xp, yp, zp, &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
1171  Real vfac;
1172  if (qr_mf_arr) {
1173  // We also have liquid water
1174  Real qr_interp{0};
1175  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
1176  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
1177  vfac = 1.0 + 0.61*qv_interp - qr_interp;
1178  } else {
1179  vfac = 1.0 + 0.61*qv_interp;
1180  }
1181  const Real mag = T_interp * vfac;
1182  const Real val = denom * mag * d_fact_new;
1183  ma_arr(i,j,k) += val;
1184  }
1185  }
1186  }
1187  });
1188  } else {
1189  auto k_arr = k_indx->const_array(mfi);
1190  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1191  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1192  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1193  {
1194  ma_arr(i,j,k) *= d_fact_old;
1195 
1196  int mk = k_arr(i,j,k);
1197  int mj = j_arr ? j_arr(i,j,k) : j;
1198  int mi = i_arr ? i_arr(i,j,k) : i;
1199  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1200  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1201  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1202  Real vfac;
1203  if (qr_mf_arr) {
1204  // We also have liquid water
1205  vfac = 1.0 + 0.61*qv_mf_arr(li,lj,lk) - qr_mf_arr(li,lj,lk);
1206  } else {
1207  vfac = 1.0 + 0.61*qv_mf_arr(li,lj,lk);
1208  }
1209  const Real mag = T_mf_arr(li,lj,lk) * vfac;
1210  const Real val = denom * mag * d_fact_new;
1211  ma_arr(i,j,k) += val;
1212  }
1213  }
1214  }
1215  });
1216  }
1217  } // MFiter
1218 
1219  // Fill interior ghost cells and any ghost cells outside a periodic domain
1220  //***********************************************************************************
1221  averages[iavg]->FillBoundary(geom.periodicity());
1222 
1223  }
1224  else // copy temperature
1225  {
1226  int iavg = m_navg - 2;
1227  IntVect ng = averages[iavg]->nGrowVect();
1228  MultiFab::Copy(*(averages[iavg]),*(averages[2]),0,0,1,ng);
1229  }
1230 
1231  //
1232  //----------------------------------------------------------
1233  // Averages for the tangential velocity magnitude
1234  //----------------------------------------------------------
1235  //
1236  {
1237  int imf = 0;
1238  int iavg = m_navg - 1;
1239 
1240  const Real Vsg = m_Vsg[lev];
1241 
1242 #ifdef _OPENMP
1243 #pragma omp parallel if (Gpu::notInLaunchRegion())
1244 #endif
1245  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1246  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1247 
1248  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
1249  fields[imf ]->const_array(mfi);
1250  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
1251  fields[imf+1]->const_array(mfi);
1252  auto ma_arr = averages[iavg]->array(mfi);
1253 
1254  if (m_interp) {
1255  const auto plo = geom.ProbLoArray();
1256  const auto dx = geom.CellSizeArray();
1257  const auto dxInv = geom.InvCellSizeArray();
1258  const auto z_phys_arr = z_phys->const_array(mfi);
1259  auto x_pos_arr = x_pos->array(mfi);
1260  auto y_pos_arr = y_pos->array(mfi);
1261  auto z_pos_arr = z_pos->array(mfi);
1262  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1263  {
1264  ma_arr(i,j,k) *= d_fact_old;
1265 
1266  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1267  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1268  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1269  for (int li(-d_radius); li <= (d_radius); ++li) {
1270  Real u_interp{0};
1271  Real v_interp{0};
1272  Real xp = x_pos_arr(i+li,j+lj,k);
1273  Real yp = y_pos_arr(i+li,j+lj,k);
1274  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1275  trilinear_interp_T(xp, yp, zp, &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
1276  trilinear_interp_T(xp, yp, zp, &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
1277  const Real mag = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
1278  Real val = denom * mag * d_fact_new;
1279  ma_arr(i,j,k) += 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 k) noexcept
1289  {
1290  ma_arr(i,j,k) *= d_fact_old;
1291 
1292  int mk = k_arr(i,j,k);
1293  int mj = j_arr ? j_arr(i,j,k) : j;
1294  int mi = i_arr ? i_arr(i,j,k) : 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  const Real u_val = 0.5 * (u_mf_arr(li,lj,lk) + u_mf_arr(li+1,lj ,lk));
1299  const Real v_val = 0.5 * (v_mf_arr(li,lj,lk) + v_mf_arr(li ,lj+1,lk));
1300  const Real mag = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
1301  Real val = denom * mag * d_fact_new;
1302  ma_arr(i,j,k) += val;
1303  }
1304  }
1305  }
1306  });
1307  }
1308  } // MFiter
1309 
1310  // Fill interior ghost cells and any ghost cells outside a periodic domain
1311  //***********************************************************************************
1312  averages[iavg]->FillBoundary(geom.periodicity());
1313 
1314  }
1315 
1316  // NOTE: Checking periodicity with the geom structure is not
1317  // sufficient at higher levels. The BA may be contained
1318  // within the domain and it's exterior ghost cells filled
1319  // from interpolation; yet the domain BCs are periodic.
1320 
1321  // Need to fill ghost cells outside the domain if not periodic
1322  bool not_per_x = !(geom.periodicity().isPeriodic(0));
1323  bool not_per_y = !(geom.periodicity().isPeriodic(1));
1324  Box cc_bnd_bx = (m_fields[lev][2]->boxArray()).minimalBox();
1325  Box domain = geom.Domain();
1326  if (domain.contains(cc_bnd_bx) || (not_per_x || not_per_y)) {
1327  for (int iavg(0); iavg < m_navg; ++iavg) {
1328  IntVect ng = averages[iavg]->nGrowVect(); ng[2]=0;
1329 
1330  // NOTE: Level 0 spans the whole domain, but finer
1331  // levels do not have such a restriction.
1332  // For now, use the bounding box of the boxArray.
1333 
1334  // NOTE2: The fields and averages have different indexing.
1335  // The averages are: U/V/T/Qv/Tv/Umag
1336  // The fields are: U/V/T/Qv/Qr/W
1337  // We clip iavg at 2 since all the remaining data is CC
1338 
1339  // Bounded box of CC data used for normalization
1340  int imf = min(iavg,2);
1341  Box bnd_bx = (m_fields[lev][imf]->boxArray()).minimalBox();
1342 #ifdef _OPENMP
1343 #pragma omp parallel if (Gpu::notInLaunchRegion())
1344 #endif
1345  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1346  Box gpbx = mfi.growntilebox(ng); gpbx.setSmall(2,0); gpbx.setBig(2,0);
1347 
1348  if (bnd_bx.contains(gpbx)) continue;
1349 
1350  auto ma_arr = averages[iavg]->array(mfi);
1351 
1352  int i_lo = bnd_bx.smallEnd(0); int i_hi = bnd_bx.bigEnd(0);
1353  int j_lo = bnd_bx.smallEnd(1); int j_hi = bnd_bx.bigEnd(1);
1354  ParallelFor(gpbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1355  {
1356  int li, lj;
1357  li = i < i_lo ? i_lo : i;
1358  li = li > i_hi ? i_hi : li;
1359  lj = j < j_lo ? j_lo : j;
1360  lj = lj > j_hi ? j_hi : lj;
1361 
1362  ma_arr(i,j,k) = ma_arr(li,lj,k);
1363  });
1364  } // MFiter
1365  } // iavg
1366  } // Not periodic
1367 }
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:47
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::Real MOSTAverage::get_zref ( ) const
inline
105 { return m_zref; }
amrex::Real 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 
)
90 {
91  m_fields[lev].resize(m_nvar);
92  m_rot_fields[lev].resize(m_nvar-1);
93  m_averages[lev].resize(m_navg);
94  m_z_phys_nd[lev] = z_phys_nd.get();
95 
96  bool use_terrain_fitted_coords = ( (m_terrain_type == TerrainType::StaticFittedMesh) ||
97  (m_terrain_type == TerrainType::MovingFittedMesh) );
98 
99  { // Nodal in x
100  auto& mf = *vars_old[Vars::xvel];
101  // Create a 2D ba, dm, & ghost cells
102  const BoxArray& ba = mf.boxArray();
103  BoxList bl2d = ba.boxList();
104  for (auto& b : bl2d) b.setRange(2,0);
105  BoxArray ba2d(std::move(bl2d));
106  const DistributionMapping& dm = mf.DistributionMap();
107  const int ncomp = 1;
108  IntVect ng = mf.nGrowVect(); ng[2]=0;
109 
110  m_fields[lev][0] = vars_old[Vars::xvel];
111  m_averages[lev][0] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
112  m_averages[lev][0]->setVal(1.E34);
113  if (m_rotate) {
114  m_rot_fields[lev][0] = std::make_unique<MultiFab>(ba,dm,ncomp,ng);
115  MultiFab::Copy(*m_rot_fields[lev][0],mf,0,0,1,ng);
116  } else {
117  m_rot_fields[lev][0] = nullptr;
118  }
119  }
120  { // Nodal in y
121  auto& mf = *vars_old[Vars::yvel];
122  // Create a 2D ba, dm, & ghost cells
123  const BoxArray& ba = mf.boxArray();
124  BoxList bl2d = ba.boxList();
125  for (auto& b : bl2d) b.setRange(2,0);
126  BoxArray ba2d(std::move(bl2d));
127  const DistributionMapping& dm = mf.DistributionMap();
128  const int ncomp = 1;
129  IntVect ng = mf.nGrowVect(); ng[2]=0;
130 
131  m_fields[lev][1] = vars_old[Vars::yvel];
132  m_averages[lev][1] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
133  m_averages[lev][1]->setVal(1.E34);
134  if (m_rotate) {
135  m_rot_fields[lev][1] = std::make_unique<MultiFab>(ba,dm,ncomp,ng);
136  MultiFab::Copy(*m_rot_fields[lev][1],mf,0,0,1,ng);
137  } else {
138  m_rot_fields[lev][1] = nullptr;
139  }
140  }
141  { // CC vars
142  auto& mf = *Theta_prim;
143  // Create a 2D ba, dm, & ghost cells
144  const BoxArray& ba = mf.boxArray();
145  BoxList bl2d = ba.boxList();
146  for (auto& b : bl2d) b.setRange(2,0);
147  BoxArray ba2d(std::move(bl2d));
148  const DistributionMapping& dm = mf.DistributionMap();
149  const int ncomp = 1;
150  const int incomp = 1;
151  IntVect ng = mf.nGrowVect(); ng[2]=0;
152 
153  // Get field pointers
154  m_fields[lev][2] = Theta_prim.get();
155  m_fields[lev][3] = Qv_prim.get();
156  m_fields[lev][4] = Qr_prim.get();
157 
158  // Initialize remaining multifabs
159  for (int iavg(2); iavg < m_navg; ++iavg) {
160  m_averages[lev][iavg] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
161  m_averages[lev][iavg]->setVal(1.E34);
162  }
163 
164  // Default to dry
165  m_averages[lev][3]->setVal(0.0);
166 
167  if (m_rotate) {
168  m_rot_fields[lev][2] = std::make_unique<MultiFab>(ba,dm,ncomp,ng);
169  m_rot_fields[lev][3] = std::make_unique<MultiFab>(ba,dm,ncomp,ng);
170  MultiFab::Copy(*m_rot_fields[lev][2],*Theta_prim,0,0,1,ng);
171  if (Qv_prim) MultiFab::Copy(*m_rot_fields[lev][3],*Qv_prim,0,0,1,ng);
172  } else {
173  m_rot_fields[lev][2] = nullptr;
174  m_rot_fields[lev][3] = nullptr;
175  }
176 
177  if (use_terrain_fitted_coords && m_norm_vec && m_interp) {
178  m_x_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
179  m_y_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
180  m_z_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
181  } else if (use_terrain_fitted_coords && m_interp) {
182  m_x_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
183  m_y_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
184  m_z_pos[lev] = std::make_unique<MultiFab>(ba2d,dm,ncomp,ng);
185  } else if (use_terrain_fitted_coords && m_norm_vec) {
186  m_i_indx[lev] = std::make_unique<iMultiFab>(ba2d,dm,incomp,ng);
187  m_j_indx[lev] = std::make_unique<iMultiFab>(ba2d,dm,incomp,ng);
188  m_k_indx[lev] = std::make_unique<iMultiFab>(ba2d,dm,incomp,ng);
189  } else {
190  m_k_indx[lev] = std::make_unique<iMultiFab>(ba2d,dm,incomp,ng);
191  }
192  }
193  // Nodal in z (only used with terrain stress rotations)
194  m_fields[lev][4] = vars_old[Vars::zvel];
195 
196  // Setup auxiliary data for spatial configuration & policy
197  //--------------------------------------------------------
198  if (use_terrain_fitted_coords && m_norm_vec && m_interp) { // Terrain w/ norm & w/ interpolation
200  } else if (use_terrain_fitted_coords && m_interp) { // Terrain w/ interpolation
201  set_z_positions_T(lev);
202  } else if (use_terrain_fitted_coords && m_norm_vec) { // Terrain w/ norm & w/o interpolation
203  set_norm_indices_T(lev);
204  } else if (use_terrain_fitted_coords) { // Terrain
205  set_k_indices_T(lev);
206  } else { // No Terrain
207  set_k_indices_N(lev);
208  }
209 
210  // Setup normalization data for the chosen policy
211  //--------------------------------------------------------
212  switch(m_policy) {
213  case 0: // Plane average
215  break;
216  case 1: // Local region/point
218  break;
219  default:
220  AMREX_ASSERT_WITH_MESSAGE(false, "Unknown policy for MOSTAverage!");
221  }
222 
223  // Set up the exponential time filtering
224  //--------------------------------------------------------
225  if (m_t_avg) {
226  // Exponential filter function
227  m_fact_old = std::exp(-1.0 / m_time_window);
228 
229  // Enforce discrete normalization: (mfn*val_new + mfo*val_old)
230  m_fact_new = 1.0 - m_fact_old;
231 
232  // None of the averages are initialized
233  m_t_init.resize(m_maxlev,0);
234  }
235 
236  // Corrections to the mean surface velocity
237  m_Vsg = Vector<Real>(m_maxlev, 0.0);
238  if (include_subgrid_vel) {
239  if (include_subgrid_vel) {
240  Print() << "Subgrid velocity scale correction at level : " << lev << ' ';
241  const auto dxArr = m_geom[lev].CellSizeArray();
242  Real dx = std::sqrt(dxArr[0]*dxArr[1]);
243  if (dx > 5000.) {
244  m_Vsg[lev] = 0.32 * std::pow(dx/5000.-1, 0.33);
245  }
246  Print() << m_Vsg[lev] << std::endl;
247  }
248  }
249 }
void set_region_normalization(const int &)
Definition: ERF_MOSTAverage.H:61
void set_z_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:586
void set_norm_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:634
void set_k_indices_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:451
void set_plane_normalization(const int &lev)
Definition: ERF_MOSTAverage.cpp:338
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:512
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:396
@ 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.

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

452 {
453  ParmParse pp(m_pp_prefix);
454  auto read_z = pp.query("most.zref",m_zref);
455  auto read_k = pp.queryarr("most.k_arr_in",m_k_in);
456 
457  // Allow default zref
458  if (!read_z) {
459  Print() << "most.zref not specified, query distance default is " << m_zref << std::endl;
460  read_z = true;
461  }
462 
463  // No default behavior with terrain (we can't tell the difference between
464  // vertical grid stretching and true terrain)
465  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(read_z != read_k,
466  "Need to specify zref or k_arr_in for MOST");
467 
468  // Capture for device
469  Real d_zref = m_zref;
470  Real d_radius = m_radius;
471  amrex::ignore_unused(d_radius);
472 
473  // Specify z_ref & compute k_indx (z_ref takes precedence)
474  if (read_z) {
475  int kmax = m_geom[lev].Domain().bigEnd(2);
476  for (MFIter mfi(*m_k_indx[lev], TileNoZ()); mfi.isValid(); ++mfi) {
477  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
478  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
479  auto k_arr = m_k_indx[lev]->array(mfi);
480  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
481  {
482  k_arr(i,j,k) = 0;
483  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
484  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
485  Real z_target = z_bot_face + d_zref;
486  for (int lk(0); lk<=kmax; ++lk) {
487  Real z_lo = 0.25 * ( z_phys_arr(i,j ,lk ) + z_phys_arr(i+1,j ,lk )
488  + z_phys_arr(i,j+1,lk ) + z_phys_arr(i+1,j+1,lk ) );
489  Real z_hi = 0.25 * ( z_phys_arr(i,j ,lk+1) + z_phys_arr(i+1,j ,lk+1)
490  + z_phys_arr(i,j+1,lk+1) + z_phys_arr(i+1,j+1,lk+1) );
491  if (z_target > z_lo && z_target < z_hi){
492  AMREX_ASSERT_WITH_MESSAGE(lk >= d_radius,
493  "K index must be larger than averaging radius!");
494  k_arr(i,j,k) = lk;
495  break;
496  }
497  }
498  });
499  }
500  // Specified k_indx & compute z_ref
501  } else if (read_k) {
502  AMREX_ASSERT_WITH_MESSAGE(false, "Specified k-indx with terrain not implemented!");
503  }
504 }

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

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

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

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

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.

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

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.

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

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)

587 {
588  ParmParse pp(m_pp_prefix);
589  auto read_zref = pp.query("most.zref",m_zref);
590  if (!read_zref) {
591  Print() << "most.zref not specified, query distance default is " << m_zref << std::endl;
592  }
593 
594  // Capture for device
595  Real d_zref = m_zref;
596  const auto plo = m_geom[lev].ProbLoArray();
597 
598  RealVect base;
599  const auto dx = m_geom[lev].CellSizeArray();
600  IntVect ng = m_x_pos[lev]->nGrowVect(); ng[2]=0;
601  for (MFIter mfi(*m_x_pos[lev], TileNoZ()); mfi.isValid(); ++mfi) {
602  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
603  Box gpbx = mfi.growntilebox(ng);
604  RealBox grb{gpbx,dx.data(),base.dataPtr()};
605 
606  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
607  auto x_pos_arr = m_x_pos[lev]->array(mfi);
608  auto y_pos_arr = m_y_pos[lev]->array(mfi);
609  auto z_pos_arr = m_z_pos[lev]->array(mfi);
610  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
611  {
612  // Final position at end of vector
613  x_pos_arr(i,j,k) = plo[0] + ((Real) i + 0.5) * dx[0];
614  y_pos_arr(i,j,k) = plo[1] + ((Real) j + 0.5) * dx[1];
615  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
616  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
617  z_pos_arr(i,j,k) = z_bot_face + d_zref;
618 
619  // Destination position must be contained on the current process!
620  Real pos[] = {x_pos_arr(i,j,k)-plo[0],y_pos_arr(i,j,k)-plo[1],0.5*dx[2]};
621  amrex::ignore_unused(pos);
622  AMREX_ASSERT_WITH_MESSAGE( grb.contains(&pos[0]),
623  "Query point outside of proc domain!");
624  });
625  }
626 }

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
265 {
266  m_fields[lev][0] = &vars_old[lev][Vars::xvel];
267  m_fields[lev][1] = &vars_old[lev][Vars::yvel];
268  m_fields[lev][2] = Theta_prim[lev].get();
269  m_fields[lev][3] = Qv_prim[lev].get();
270  m_fields[lev][4] = Qr_prim[lev].get();
271  m_fields[lev][5] = &vars_old[lev][Vars::zvel];
272 }

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
1496 {
1497  // Peel back the level
1498  auto& averages = m_averages[lev];
1499 
1500  int navg = m_navg - 1;
1501 
1502  std::ofstream ofile;
1503  ofile.open ("MOST_averages.txt");
1504  ofile << "Averages computed via MOSTAverages class:\n";
1505 
1506  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1507  Box bx = mfi.tilebox(); bx.setBig(2,0);
1508  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1509  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1510 
1511  for (int j(jl); j <= ju; ++j) {
1512  for (int i(il); i <= iu; ++i) {
1513  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1514  int k = 0;
1515  for (int iavg(0); iavg <= navg; ++iavg) {
1516  auto mf_arr = averages[iavg]->array(mfi);
1517  ofile << "iavg val: "
1518  << iavg << ' '
1519  << mf_arr(i,j,k) << "\n";
1520  }
1521  ofile << "\n";
1522  }
1523  }
1524  }
1525  ofile.close();
1526 }
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
1377 {
1378  // Peel back the level
1379  auto& averages = m_averages[lev];
1380  auto& k_indx = m_k_indx[lev];
1381 
1382  int navg = m_navg - 1;
1383 
1384  std::ofstream ofile;
1385  ofile.open ("MOST_k_indices.txt");
1386  ofile << "K indices used to compute averages via MOSTAverages class:\n";
1387 
1388  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1389  Box bx = mfi.tilebox(); bx.setBig(2,0);
1390  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1391  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1392 
1393  auto k_arr = k_indx->array(mfi);
1394 
1395  for (int j(jl); j <= ju; ++j) {
1396  for (int i(il); i <= iu; ++i) {
1397  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1398  int k = 0;
1399  ofile << "K_ind: "
1400  << k_arr(i,j,k) << "\n";
1401  ofile << "\n";
1402  }
1403  }
1404  }
1405  ofile.close();
1406 }
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
1416 {
1417  // Peel back the level
1418  auto& averages = m_averages[lev];
1419  auto& k_indx = m_k_indx[lev];
1420  auto& j_indx = m_j_indx[lev];
1421  auto& i_indx = m_i_indx[lev];
1422 
1423  int navg = m_navg - 1;
1424 
1425  std::ofstream ofile;
1426  ofile.open ("MOST_ijk_indices.txt");
1427  ofile << "IJK indices used to compute averages via MOSTAverages class:\n";
1428 
1429  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1430  Box bx = mfi.tilebox(); bx.setBig(2,0);
1431  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1432  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1433 
1434  auto k_arr = k_indx->array(mfi);
1435  auto j_arr = j_indx ? j_indx->array(mfi) : Array4<int> {};
1436  auto i_arr = i_indx ? i_indx->array(mfi) : Array4<int> {};
1437 
1438  for (int j(jl); j <= ju; ++j) {
1439  for (int i(il); i <= iu; ++i) {
1440  ofile << "(I1,J1,K1): " << "(" << i << "," << j << "," << 0 << ")" << "\n";
1441 
1442  int k = 0;
1443  int km = k_arr(i,j,k);
1444  int jm = j_arr ? j_arr(i,j,k) : j;
1445  int im = i_arr ? i_arr(i,j,k) : i;
1446 
1447  ofile << "(I2,J2,K2): "
1448  << "(" << im << "," << jm << "," << km << ")" << "\n";
1449  ofile << "\n";
1450  }
1451  }
1452  }
1453  ofile.close();
1454 }
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
1466 {
1467  // Peel back the level
1468  auto& x_pos_mf = m_x_pos[lev];
1469  auto& z_pos_mf = m_z_pos[lev];
1470 
1471  std::ofstream ofile;
1472  ofile.open ("MOST_xz_positions.txt");
1473 
1474  for (MFIter mfi(*x_pos_mf, TileNoZ()); mfi.isValid(); ++mfi) {
1475  Box bx = mfi.tilebox(); bx.setBig(2,0);
1476  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1477 
1478  auto x_pos_arr = x_pos_mf->array(mfi);
1479  auto z_pos_arr = z_pos_mf->array(mfi);
1480 
1481  int k = 0;
1482  for (int i(il); i <= iu; ++i)
1483  ofile << x_pos_arr(i,j,k) << ' ' << z_pos_arr(i,j,k) << "\n";
1484  }
1485  ofile.close();
1486 }
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


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