ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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 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
 
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::Real > m_Vsg
 

Constructor & Destructor Documentation

◆ MOSTAverage() [1/3]

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

◆ ~MOSTAverage()

MOSTAverage::~MOSTAverage ( )
inline
23  {}

◆ 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
680 {
681  if (m_rotate) set_rotated_fields(lev);
682 
683  switch(m_policy) {
684  case 0: // Standard plane average
686  break;
687  case 1: // Local region/point
689  break;
690  default:
691  AMREX_ASSERT_WITH_MESSAGE(false, "Unknown policy for MOSTAverage!");
692  }
693 
694  // We have initialized the averages
695  if (m_t_avg) m_t_init[lev] = 1;
696 }
bool m_t_avg
Definition: ERF_MOSTAverage.H:229
void compute_plane_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:705
int m_policy
Definition: ERF_MOSTAverage.H:199
void compute_region_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:990
amrex::Vector< int > m_t_init
Definition: ERF_MOSTAverage.H:230
void set_rotated_fields(const int &lev)
Definition: ERF_MOSTAverage.cpp:265
bool m_rotate
Definition: ERF_MOSTAverage.H:200
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
706 {
707  // Peel back the level
708  auto& fields = m_fields[lev];
709  auto& rot_fields = m_rot_fields[lev];
710  auto& averages = m_averages[lev];
711  const auto & geom = m_geom[lev];
712 
713  auto& z_phys = m_z_phys_nd[lev];
714  auto& x_pos = m_x_pos[lev];
715  auto& y_pos = m_y_pos[lev];
716  auto& z_pos = m_z_pos[lev];
717 
718  auto& i_indx = m_i_indx[lev];
719  auto& j_indx = m_j_indx[lev];
720  auto& k_indx = m_k_indx[lev];
721 
722  auto& ncell_plane = m_ncell_plane[lev];
723  auto& plane_average = m_plane_average[lev];
724 
725  // Set factors for time averaging
726  Real d_fact_new, d_fact_old;
727  if (m_t_avg && m_t_init[lev]) {
728  d_fact_new = m_fact_new;
729  d_fact_old = m_fact_old;
730  } else {
731  d_fact_new = 1.0;
732  d_fact_old = 0.0;
733  }
734 
735  // GPU array to accumulate averages into
736  Gpu::DeviceVector<Real> pavg(plane_average.size(), 0.0);
737  Real* plane_avg = pavg.data();
738 
739  // Vectors for normalization and buffer storage
740  Vector<Real> denom(plane_average.size(),0.0);
741  Vector<Real> val_old(plane_average.size(),0.0);
742 
743  //
744  //----------------------------------------------------------
745  // Averages over all the fields
746  //----------------------------------------------------------
747  //
748  Box domain = geom.Domain();
749 
750  Array<int,AMREX_SPACEDIM> is_per = {0,0,0};
751  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
752  if (geom.isPeriodic(idim)) is_per[idim] = 1;
753  }
754 
755  // Averages for U,V,T,Qv (not Qr or W)
756  for (int imf(0); imf < 4; ++imf) {
757 
758  // Continue if no valid Qv pointer
759  if (!fields[imf]) continue;
760 
761  denom[imf] = 1.0 / (Real)ncell_plane[imf];
762  val_old[imf] = plane_average[imf]*d_fact_old;
763 
764 #ifdef _OPENMP
765 #pragma omp parallel if (Gpu::notInLaunchRegion())
766 #endif
767  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
768  Box vbx = mfi.validbox(); // This is the grid (not tile)
769  Box pbx = mfi.tilebox(); // This is the tile (not grid)
770  pbx.setSmall(2,0); pbx.setBig(2,0);
771 
772  // Avoid double counting nodal data by changing the high end when we are
773  // at the high side of the grid (not just of the tile)
774  IndexType ixt = averages[imf]->boxArray().ixType();
775  for (int idim(0); idim < AMREX_SPACEDIM-1; ++idim) {
776  if ( ixt.nodeCentered(idim) && (pbx.bigEnd(idim) == vbx.bigEnd(idim)) ) {
777  int dom_hi = domain.bigEnd(idim)+1;
778  if (pbx.bigEnd(idim) < dom_hi || is_per[idim]) {
779  pbx.growHi(idim,-1);
780  }
781  }
782  }
783 
784  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
785  fields[imf]->const_array(mfi);
786 
787  if (m_interp) {
788  const auto plo = geom.ProbLoArray();
789  const auto dxInv = geom.InvCellSizeArray();
790  const auto z_phys_arr = z_phys->const_array(mfi);
791  auto x_pos_arr = x_pos->array(mfi);
792  auto y_pos_arr = y_pos->array(mfi);
793  auto z_pos_arr = z_pos->array(mfi);
794  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
795  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
796  {
797  Real interp{0};
798  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
799  &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
800  Real val = interp;
801  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
802  });
803  } else {
804  auto k_arr = k_indx->const_array(mfi);
805  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
806  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
807  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
808  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
809  {
810  int mk = k_arr(i,j,k);
811  int mj = j_arr ? j_arr(i,j,k) : j;
812  int mi = i_arr ? i_arr(i,j,k) : i;
813  Real val = mf_arr(mi,mj,mk);
814  Gpu::deviceReduceSum(&plane_avg[imf], val, handler);
815  });
816  }
817  }
818  }
819 
820  //
821  //------------------------------------------------------------------------
822  // Averages for virtual potential temperature
823  // (This is cell-centered so we don't need to worry about double-counting)
824  //------------------------------------------------------------------------
825  //
826  if (fields[3]) // We have water vapor
827  {
828  int iavg = 4;
829  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
830  val_old[iavg] = plane_average[iavg]*d_fact_old;
831 
832 #ifdef _OPENMP
833 #pragma omp parallel if (Gpu::notInLaunchRegion())
834 #endif
835  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi)
836  {
837  Box pbx = mfi.tilebox();
838  pbx.setSmall(2,0); pbx.setBig(2,0);
839 
840  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
841  const Array4<Real const>& qv_mf_arr = (fields[3])? fields[3]->const_array(mfi) : Array4<const Real>{};
842  const Array4<Real const>& qr_mf_arr = (fields[4])? fields[4]->const_array(mfi) : Array4<const Real>{};
843 
844  if (m_interp) {
845  const auto plo = m_geom[lev].ProbLoArray();
846  const auto dxInv = m_geom[lev].InvCellSizeArray();
847  const auto z_phys_arr = z_phys->const_array(mfi);
848  auto x_pos_arr = x_pos->array(mfi);
849  auto y_pos_arr = y_pos->array(mfi);
850  auto z_pos_arr = z_pos->array(mfi);
851  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
852  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
853  {
854  Real T_interp{0};
855  Real qv_interp{0};
856  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
857  &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
858  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
859  &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
860  Real vfac;
861  if (qr_mf_arr) {
862  // We also have liquid water
863  Real qr_interp{0};
864  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
865  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
866  vfac = 1.0 + 0.61*qv_interp - qr_interp;
867  } else {
868  vfac = 1.0 + 0.61*qv_interp;
869  }
870  const Real val = T_interp * vfac;
871  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
872  });
873  } else {
874  auto k_arr = k_indx->const_array(mfi);
875  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
876  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
877  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
878  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
879  {
880  int mk = k_arr(i,j,k);
881  int mj = j_arr ? j_arr(i,j,k) : j;
882  int mi = i_arr ? i_arr(i,j,k) : i;
883  Real vfac;
884  if (qr_mf_arr) {
885  // We also have liquid water
886  vfac = 1.0 + 0.61*qv_mf_arr(mi,mj,mk) - qr_mf_arr(mi,mj,mk);
887  } else {
888  vfac = 1.0 + 0.61*qv_mf_arr(mi,mj,mk);
889  }
890  const Real val = T_mf_arr(mi,mj,mk) * vfac;
891  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
892  });
893  }
894  }
895  }
896  else // copy temperature
897  {
898  int iavg = m_navg - 2;
899  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
900  // plane_avg[iavg] = plane_avg[2]
901  Gpu::copy(Gpu::deviceToDevice, pavg.begin() + 2, pavg.begin() + 3,
902  pavg.begin() + iavg);
903  }
904 
905  //
906  //------------------------------------------------------------------------
907  // Averages for the tangential velocity magnitude
908  // (This is cell-centered so we don't need to worry about double-counting)
909  //------------------------------------------------------------------------
910  //
911  {
912  int imf = 0;
913  int iavg = m_navg - 1;
914  denom[iavg] = 1.0 / (Real)ncell_plane[iavg];
915  val_old[iavg] = plane_average[iavg]*d_fact_old;
916 
917  const Real Vsg = m_Vsg[lev];
918 
919 #ifdef _OPENMP
920 #pragma omp parallel if (Gpu::notInLaunchRegion())
921 #endif
922  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi)
923  {
924  Box pbx = mfi.tilebox();
925  pbx.setSmall(2,0); pbx.setBig(2,0);
926 
927  // Last element is Umag and always cell centered
928  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
929  fields[imf ]->const_array(mfi);
930  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
931  fields[imf+1]->const_array(mfi);
932 
933  if (m_interp) {
934  const auto plo = m_geom[lev].ProbLoArray();
935  const auto dxInv = m_geom[lev].InvCellSizeArray();
936  const auto z_phys_arr = z_phys->const_array(mfi);
937  auto x_pos_arr = x_pos->array(mfi);
938  auto y_pos_arr = y_pos->array(mfi);
939  auto z_pos_arr = z_pos->array(mfi);
940  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
941  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
942  {
943  Real u_interp{0};
944  Real v_interp{0};
945  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
946  &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
947  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
948  &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
949  const Real val = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
950  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
951  });
952  } else {
953  auto k_arr = k_indx->const_array(mfi);
954  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
955  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
956  ParallelFor(Gpu::KernelInfo().setReduction(true), pbx, [=]
957  AMREX_GPU_DEVICE(int i, int j, int k, Gpu::Handler const& handler) noexcept
958  {
959  int mk = k_arr(i,j,k);
960  int mj = j_arr ? j_arr(i,j,k) : j;
961  int mi = i_arr ? i_arr(i,j,k) : i;
962  const Real u_val = 0.5 * (u_mf_arr(mi,mj,mk) + u_mf_arr(mi+1,mj ,mk));
963  const Real v_val = 0.5 * (v_mf_arr(mi,mj,mk) + v_mf_arr(mi ,mj+1,mk));
964  const Real val = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
965  Gpu::deviceReduceSum(&plane_avg[iavg], val, handler);
966  });
967  }
968  }
969  }
970 
971  // Copy to host and sum across procs
972  Gpu::copy(Gpu::deviceToHost, pavg.begin(), pavg.end(), plane_average.begin());
973  ParallelDescriptor::ReduceRealSum(plane_average.data(), plane_average.size());
974 
975  // No spatial variation with plane averages
976  for (int iavg(0); iavg < m_navg; ++iavg){
977  plane_average[iavg] *= denom[iavg]*d_fact_new;
978  plane_average[iavg] += val_old[iavg];
979  averages[iavg]->setVal(plane_average[iavg]);
980  }
981 }
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
int m_navg
Definition: ERF_MOSTAverage.H:197
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_averages
Definition: ERF_MOSTAverage.H:208
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_y_pos
Definition: ERF_MOSTAverage.H:203
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_i_indx
Definition: ERF_MOSTAverage.H:205
amrex::Vector< amrex::MultiFab * > m_z_phys_nd
Definition: ERF_MOSTAverage.H:190
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_x_pos
Definition: ERF_MOSTAverage.H:202
amrex::Vector< amrex::Real > m_Vsg
Definition: ERF_MOSTAverage.H:237
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_rot_fields
Definition: ERF_MOSTAverage.H:209
amrex::Vector< amrex::Vector< amrex::Real > > m_plane_average
Definition: ERF_MOSTAverage.H:214
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_z_pos
Definition: ERF_MOSTAverage.H:204
amrex::Vector< amrex::Vector< int > > m_ncell_plane
Definition: ERF_MOSTAverage.H:213
amrex::Real m_fact_new
Definition: ERF_MOSTAverage.H:232
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_j_indx
Definition: ERF_MOSTAverage.H:206
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_fields
Definition: ERF_MOSTAverage.H:189
amrex::Real m_fact_old
Definition: ERF_MOSTAverage.H:232
bool m_interp
Definition: ERF_MOSTAverage.H:224
const amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_MOSTAverage.H:188
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:119
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_k_indx
Definition: ERF_MOSTAverage.H:207

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
991 {
992  // Peel back the level
993  auto& fields = m_fields[lev];
994  auto& rot_fields = m_rot_fields[lev];
995  auto& averages = m_averages[lev];
996  const auto & geom = m_geom[lev];
997 
998  auto& z_phys = m_z_phys_nd[lev];
999  auto& x_pos = m_x_pos[lev];
1000  auto& y_pos = m_y_pos[lev];
1001  auto& z_pos = m_z_pos[lev];
1002 
1003  auto& i_indx = m_i_indx[lev];
1004  auto& j_indx = m_j_indx[lev];
1005  auto& k_indx = m_k_indx[lev];
1006 
1007  // Set factors for time averaging
1008  Real d_fact_new, d_fact_old;
1009  if (m_t_avg && m_t_init[lev]) {
1010  d_fact_new = m_fact_new;
1011  d_fact_old = m_fact_old;
1012  } else {
1013  d_fact_new = 1.0;
1014  d_fact_old = 0.0;
1015  }
1016 
1017  // Number of cells contained in the local average
1018  const Real denom = 1.0 / (Real) m_ncell_region;
1019 
1020  // Capture radius for device
1021  int d_radius = m_radius;
1022 
1023  //
1024  //----------------------------------------------------------
1025  // Averages for U,V,T,Qv
1026  //----------------------------------------------------------
1027  //
1028  for (int imf(0); imf < 4; ++imf) {
1029 
1030  // Continue if no valid Qv pointer
1031  if (!fields[imf]) continue;
1032 
1033 #ifdef _OPENMP
1034 #pragma omp parallel if (Gpu::notInLaunchRegion())
1035 #endif
1036  for (MFIter mfi(*fields[imf], TileNoZ()); mfi.isValid(); ++mfi) {
1037  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1038 
1039  auto mf_arr = (m_rotate) ? rot_fields[imf]->const_array(mfi) :
1040  fields[imf]->const_array(mfi);
1041  auto ma_arr = averages[imf]->array(mfi);
1042 
1043  if (m_interp) {
1044  const auto plo = geom.ProbLoArray();
1045  const auto dx = geom.CellSizeArray();
1046  const auto dxInv = geom.InvCellSizeArray();
1047  const auto z_phys_arr = z_phys->const_array(mfi);
1048  auto x_pos_arr = x_pos->array(mfi);
1049  auto y_pos_arr = y_pos->array(mfi);
1050  auto z_pos_arr = z_pos->array(mfi);
1051  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1052  {
1053  ma_arr(i,j,k) *= d_fact_old;
1054 
1055  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1056  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1057  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1058  for (int li(-d_radius); li <= (d_radius); ++li) {
1059  Real interp{0};
1060  Real xp = x_pos_arr(i+li,j+lj,k);
1061  Real yp = y_pos_arr(i+li,j+lj,k);
1062  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1063  trilinear_interp_T(xp, yp, zp, &interp, mf_arr, z_phys_arr, plo, dxInv, 1);
1064  Real val = denom * interp * d_fact_new;
1065  ma_arr(i,j,k) += val;
1066  }
1067  }
1068  }
1069  });
1070  } else {
1071  auto k_arr = k_indx->const_array(mfi);
1072  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1073  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1074  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1075  {
1076  ma_arr(i,j,k) *= d_fact_old;
1077 
1078  int mk = k_arr(i,j,k);
1079  int mj = j_arr ? j_arr(i,j,k) : j;
1080  int mi = i_arr ? i_arr(i,j,k) : i;
1081  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1082  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1083  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1084  Real val = denom * mf_arr(li, lj, lk) * d_fact_new;
1085  ma_arr(i,j,k) += val;
1086  }
1087  }
1088  }
1089  });
1090  }
1091  } // MFiter
1092 
1093  // Fill interior ghost cells and any ghost cells outside a periodic domain
1094  //***********************************************************************************
1095  averages[imf]->FillBoundary(geom.periodicity());
1096 
1097  } // imf
1098 
1099  //
1100  //----------------------------------------------------------
1101  // Averages for virtual potential temperature
1102  //----------------------------------------------------------
1103  //
1104  if (fields[3]) // We have water vapor
1105  {
1106  int iavg = 4;
1107 
1108 #ifdef _OPENMP
1109 #pragma omp parallel if (Gpu::notInLaunchRegion())
1110 #endif
1111  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1112  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1113 
1114  const Array4<Real const>& T_mf_arr = fields[2]->const_array(mfi);
1115  const Array4<Real const>& qv_mf_arr = (fields[3])? fields[3]->const_array(mfi) : Array4<const Real>{};
1116  const Array4<Real const>& qr_mf_arr = (fields[4])? fields[4]->const_array(mfi) : Array4<const Real>{};
1117  auto ma_arr = averages[iavg]->array(mfi);
1118 
1119  if (m_interp) {
1120  const auto plo = geom.ProbLoArray();
1121  const auto dx = geom.CellSizeArray();
1122  const auto dxInv = geom.InvCellSizeArray();
1123  const auto z_phys_arr = z_phys->const_array(mfi);
1124  auto x_pos_arr = x_pos->array(mfi);
1125  auto y_pos_arr = y_pos->array(mfi);
1126  auto z_pos_arr = z_pos->array(mfi);
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  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1132  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1133  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1134  for (int li(-d_radius); li <= (d_radius); ++li) {
1135  Real T_interp{0};
1136  Real qv_interp{0};
1137  Real xp = x_pos_arr(i+li,j+lj,k);
1138  Real yp = y_pos_arr(i+li,j+lj,k);
1139  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1140  trilinear_interp_T(xp, yp, zp, &T_interp, T_mf_arr, z_phys_arr, plo, dxInv, 1);
1141  trilinear_interp_T(xp, yp, zp, &qv_interp, qv_mf_arr, z_phys_arr, plo, dxInv, 1);
1142  Real vfac;
1143  if (qr_mf_arr) {
1144  // We also have liquid water
1145  Real qr_interp{0};
1146  trilinear_interp_T(x_pos_arr(i,j,k), y_pos_arr(i,j,k), z_pos_arr(i,j,k),
1147  &qr_interp, qr_mf_arr, z_phys_arr, plo, dxInv, 1);
1148  vfac = 1.0 + 0.61*qv_interp - qr_interp;
1149  } else {
1150  vfac = 1.0 + 0.61*qv_interp;
1151  }
1152  const Real mag = T_interp * vfac;
1153  const Real val = denom * mag * d_fact_new;
1154  ma_arr(i,j,k) += val;
1155  }
1156  }
1157  }
1158  });
1159  } else {
1160  auto k_arr = k_indx->const_array(mfi);
1161  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1162  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1163  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1164  {
1165  ma_arr(i,j,k) *= d_fact_old;
1166 
1167  int mk = k_arr(i,j,k);
1168  int mj = j_arr ? j_arr(i,j,k) : j;
1169  int mi = i_arr ? i_arr(i,j,k) : i;
1170  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1171  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1172  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1173  Real vfac;
1174  if (qr_mf_arr) {
1175  // We also have liquid water
1176  vfac = 1.0 + 0.61*qv_mf_arr(li,lj,lk) - qr_mf_arr(li,lj,lk);
1177  } else {
1178  vfac = 1.0 + 0.61*qv_mf_arr(li,lj,lk);
1179  }
1180  const Real mag = T_mf_arr(li,lj,lk) * vfac;
1181  const Real val = denom * mag * d_fact_new;
1182  ma_arr(i,j,k) += val;
1183  }
1184  }
1185  }
1186  });
1187  }
1188  } // MFiter
1189 
1190  // Fill interior ghost cells and any ghost cells outside a periodic domain
1191  //***********************************************************************************
1192  averages[iavg]->FillBoundary(geom.periodicity());
1193 
1194  }
1195  else // copy temperature
1196  {
1197  int iavg = m_navg - 2;
1198  IntVect ng = averages[iavg]->nGrowVect();
1199  MultiFab::Copy(*(averages[iavg]),*(averages[2]),0,0,1,ng);
1200  }
1201 
1202  //
1203  //----------------------------------------------------------
1204  // Averages for the tangential velocity magnitude
1205  //----------------------------------------------------------
1206  //
1207  {
1208  int imf = 0;
1209  int iavg = m_navg - 1;
1210 
1211  const Real Vsg = m_Vsg[lev];
1212 
1213 #ifdef _OPENMP
1214 #pragma omp parallel if (Gpu::notInLaunchRegion())
1215 #endif
1216  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1217  Box pbx = mfi.tilebox(); pbx.setSmall(2,0); pbx.setBig(2,0);
1218 
1219  auto u_mf_arr = (m_rotate) ? rot_fields[imf ]->const_array(mfi) :
1220  fields[imf ]->const_array(mfi);
1221  auto v_mf_arr = (m_rotate) ? rot_fields[imf+1]->const_array(mfi) :
1222  fields[imf+1]->const_array(mfi);
1223  auto ma_arr = averages[iavg]->array(mfi);
1224 
1225  if (m_interp) {
1226  const auto plo = geom.ProbLoArray();
1227  const auto dx = geom.CellSizeArray();
1228  const auto dxInv = geom.InvCellSizeArray();
1229  const auto z_phys_arr = z_phys->const_array(mfi);
1230  auto x_pos_arr = x_pos->array(mfi);
1231  auto y_pos_arr = y_pos->array(mfi);
1232  auto z_pos_arr = z_pos->array(mfi);
1233  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1234  {
1235  ma_arr(i,j,k) *= d_fact_old;
1236 
1237  Real met_h_zeta = Compute_h_zeta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
1238  for (int lk(-d_radius); lk <= (d_radius); ++lk) {
1239  for (int lj(-d_radius); lj <= (d_radius); ++lj) {
1240  for (int li(-d_radius); li <= (d_radius); ++li) {
1241  Real u_interp{0};
1242  Real v_interp{0};
1243  Real xp = x_pos_arr(i+li,j+lj,k);
1244  Real yp = y_pos_arr(i+li,j+lj,k);
1245  Real zp = z_pos_arr(i+li,j+lj,k) + met_h_zeta*lk*dx[2];
1246  trilinear_interp_T(xp, yp, zp, &u_interp, u_mf_arr, z_phys_arr, plo, dxInv, 1);
1247  trilinear_interp_T(xp, yp, zp, &v_interp, v_mf_arr, z_phys_arr, plo, dxInv, 1);
1248  const Real mag = std::sqrt(u_interp*u_interp + v_interp*v_interp + Vsg*Vsg);
1249  Real val = denom * mag * d_fact_new;
1250  ma_arr(i,j,k) += val;
1251  }
1252  }
1253  }
1254  });
1255  } else {
1256  auto k_arr = k_indx->const_array(mfi);
1257  auto j_arr = j_indx ? j_indx->const_array(mfi) : Array4<const int> {};
1258  auto i_arr = i_indx ? i_indx->const_array(mfi) : Array4<const int> {};
1259  ParallelFor(pbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1260  {
1261  ma_arr(i,j,k) *= d_fact_old;
1262 
1263  int mk = k_arr(i,j,k);
1264  int mj = j_arr ? j_arr(i,j,k) : j;
1265  int mi = i_arr ? i_arr(i,j,k) : i;
1266  for (int lk(mk-d_radius); lk <= (mk+d_radius); ++lk) {
1267  for (int lj(mj-d_radius); lj <= (mj+d_radius); ++lj) {
1268  for (int li(mi-d_radius); li <= (mi+d_radius); ++li) {
1269  const Real u_val = 0.5 * (u_mf_arr(li,lj,lk) + u_mf_arr(li+1,lj ,lk));
1270  const Real v_val = 0.5 * (v_mf_arr(li,lj,lk) + v_mf_arr(li ,lj+1,lk));
1271  const Real mag = std::sqrt(u_val*u_val + v_val*v_val + Vsg*Vsg);
1272  Real val = denom * mag * d_fact_new;
1273  ma_arr(i,j,k) += val;
1274  }
1275  }
1276  }
1277  });
1278  }
1279  } // MFiter
1280 
1281  // Fill interior ghost cells and any ghost cells outside a periodic domain
1282  //***********************************************************************************
1283  averages[iavg]->FillBoundary(geom.periodicity());
1284 
1285  }
1286 
1287  // NOTE: Checking periodicity with the geom structure is not
1288  // sufficient at higher levels. The BA may be contained
1289  // within the domain and it's exterior ghost cells filled
1290  // from interpolation; yet the domain BCs are periodic.
1291 
1292  // Need to fill ghost cells outside the domain if not periodic
1293  bool not_per_x = !(geom.periodicity().isPeriodic(0));
1294  bool not_per_y = !(geom.periodicity().isPeriodic(1));
1295  Box cc_bnd_bx = (m_fields[lev][2]->boxArray()).minimalBox();
1296  Box domain = geom.Domain();
1297  if (domain.contains(cc_bnd_bx) || (not_per_x || not_per_y)) {
1298  for (int iavg(0); iavg < m_navg; ++iavg) {
1299  IntVect ng = averages[iavg]->nGrowVect(); ng[2]=0;
1300 
1301  // NOTE: Level 0 spans the whole domain, but finer
1302  // levels do not have such a restriction.
1303  // For now, use the bounding box of the boxArray.
1304 
1305  // NOTE2: The fields and averages have different indexing.
1306  // The averages are: U/V/T/Qv/Tv/Umag
1307  // The fields are: U/V/T/Qv/Qr/W
1308  // We clip iavg at 2 since all the remaining data is CC
1309 
1310  // Bounded box of CC data used for normalization
1311  int imf = min(iavg,2);
1312  Box bnd_bx = (m_fields[lev][imf]->boxArray()).minimalBox();
1313 #ifdef _OPENMP
1314 #pragma omp parallel if (Gpu::notInLaunchRegion())
1315 #endif
1316  for (MFIter mfi(*averages[iavg], TileNoZ()); mfi.isValid(); ++mfi) {
1317  Box gpbx = mfi.growntilebox(ng); gpbx.setSmall(2,0); gpbx.setBig(2,0);
1318 
1319  if (bnd_bx.contains(gpbx)) continue;
1320 
1321  auto ma_arr = averages[iavg]->array(mfi);
1322 
1323  int i_lo = bnd_bx.smallEnd(0); int i_hi = bnd_bx.bigEnd(0);
1324  int j_lo = bnd_bx.smallEnd(1); int j_hi = bnd_bx.bigEnd(1);
1325  ParallelFor(gpbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1326  {
1327  int li, lj;
1328  li = i < i_lo ? i_lo : i;
1329  li = li > i_hi ? i_hi : li;
1330  lj = j < j_lo ? j_lo : j;
1331  lj = lj > j_hi ? j_hi : lj;
1332 
1333  ma_arr(i,j,k) = ma_arr(li,lj,k);
1334  });
1335  } // MFiter
1336  } // iavg
1337  } // Not periodic
1338 }
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:218
int m_ncell_region
Definition: ERF_MOSTAverage.H:219
@ 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
101 { 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
104 { return m_zref; }
amrex::Real m_zref
Definition: ERF_MOSTAverage.H:201

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

383 {
384  ParmParse pp(m_pp_prefix);
385  auto read_z = pp.query("most.zref",m_zref);
386  auto read_k = pp.queryarr("most.k_arr_in",m_k_in);
387 
388  // Default behavior is to use the first cell center
389  if (!read_z && !read_k) {
390  Real m_zlo = m_geom[0].ProbLo(2);
391  Real m_dz = m_geom[0].CellSize(2);
392  m_zref = m_zlo + 0.5 * m_dz;
393  Print() << "Reference height for MOST set to " << m_zref << std::endl;
394  read_z = true;
395  }
396 
397  // Specify z_ref & compute k_indx (z_ref takes precedence)
398  if (read_z) {
399  Real m_zlo = m_geom[lev].ProbLo(2);
400  Real m_zhi = m_geom[lev].ProbHi(2);
401  Real m_dz = m_geom[lev].CellSize(2);
402 
403  amrex::ignore_unused(m_zhi);
404 
405  AMREX_ASSERT_WITH_MESSAGE(m_zref >= m_zlo + 0.5 * m_dz,
406  "Query point must be past first z-cell!");
407 
408  AMREX_ASSERT_WITH_MESSAGE(m_zref <= m_zhi - 0.5 * m_dz,
409  "Query point must be below the last z-cell!");
410 
411  int lk = static_cast<int>(floor((m_zref - m_zlo) / m_dz - 0.5));
412 
413  m_zref = (lk + 0.5) * m_dz + m_zlo;
414 
415  AMREX_ALWAYS_ASSERT(lk >= m_radius);
416 
417  m_k_indx[lev]->setVal(lk);
418  // Specified k_indx & compute z_ref
419  } else if (read_k) {
420  AMREX_ASSERT_WITH_MESSAGE(m_k_in[lev] >= m_radius,
421  "K index must be larger than averaging radius!");
422  m_k_indx[lev]->setVal(m_k_in[lev]);
423 
424  // TODO: check that z_ref is constant across levels
425  Real m_zlo = m_geom[0].ProbLo(2);
426  Real m_dz = m_geom[0].CellSize(2);
427  m_zref = ((Real)m_k_in[0] + 0.5) * m_dz + m_zlo;
428  }
429 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:219
std::string m_pp_prefix
Definition: ERF_MOSTAverage.H:191
amrex::Vector< int > m_k_in
Definition: ERF_MOSTAverage.H:220

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.

438 {
439  ParmParse pp(m_pp_prefix);
440  auto read_z = pp.query("most.zref",m_zref);
441  auto read_k = pp.queryarr("most.k_arr_in",m_k_in);
442 
443  // No default behavior with terrain (we can't tell the difference between
444  // vertical grid stretching and true terrain)
445  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(read_z != read_k,
446  "Need to specify zref or k_arr_in for MOST");
447 
448  // Capture for device
449  Real d_zref = m_zref;
450  Real d_radius = m_radius;
451  amrex::ignore_unused(d_radius);
452 
453  // Specify z_ref & compute k_indx (z_ref takes precedence)
454  if (read_z) {
455  int kmax = m_geom[lev].Domain().bigEnd(2);
456  for (MFIter mfi(*m_k_indx[lev], TileNoZ()); mfi.isValid(); ++mfi) {
457  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
458  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
459  auto k_arr = m_k_indx[lev]->array(mfi);
460  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
461  {
462  k_arr(i,j,k) = 0;
463  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
464  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
465  Real z_target = z_bot_face + d_zref;
466  for (int lk(0); lk<=kmax; ++lk) {
467  Real z_lo = 0.25 * ( z_phys_arr(i,j ,lk ) + z_phys_arr(i+1,j ,lk )
468  + z_phys_arr(i,j+1,lk ) + z_phys_arr(i+1,j+1,lk ) );
469  Real z_hi = 0.25 * ( z_phys_arr(i,j ,lk+1) + z_phys_arr(i+1,j ,lk+1)
470  + z_phys_arr(i,j+1,lk+1) + z_phys_arr(i+1,j+1,lk+1) );
471  if (z_target > z_lo && z_target < z_hi){
472  AMREX_ASSERT_WITH_MESSAGE(lk >= d_radius,
473  "K index must be larger than averaging radius!");
474  k_arr(i,j,k) = lk;
475  break;
476  }
477  }
478  });
479  }
480  // Specified k_indx & compute z_ref
481  } else if (read_k) {
482  AMREX_ASSERT_WITH_MESSAGE(false, "Specified k-indx with terrain not implemented!");
483  }
484 }

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.

493 {
494  ParmParse pp(m_pp_prefix);
495  pp.get("most.zref",m_zref);
496 
497  // Capture for device
498  Real d_zref = m_zref;
499  Real d_radius = m_radius;
500 
501  int kmax = m_geom[lev].Domain().bigEnd(2);
502  const auto dxInv = m_geom[lev].InvCellSizeArray();
503  IntVect ng = m_k_indx[lev]->nGrowVect(); ng[2]=0;
504  for (MFIter mfi(*m_k_indx[lev], TileNoZ()); mfi.isValid(); ++mfi) {
505  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
506  Box gpbx = mfi.growntilebox(ng);
507  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
508  auto i_arr = m_i_indx[lev]->array(mfi);
509  auto j_arr = m_j_indx[lev]->array(mfi);
510  auto k_arr = m_k_indx[lev]->array(mfi);
511  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
512  {
513  // Elements of normal vector
514  Real met_h_xi = Compute_h_xi_AtCellCenter (i,j,k,dxInv,z_phys_arr);
515  Real met_h_eta = Compute_h_eta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
516  Real mag = std::sqrt(met_h_xi*met_h_xi + met_h_eta*met_h_eta + 1.0);
517 
518  // Unit-normal vector scaled by z_ref
519  Real delta_x = -met_h_xi/mag * d_zref;
520  Real delta_y = -met_h_eta/mag * d_zref;
521  Real delta_z = 1.0/mag * d_zref;
522 
523  // Compute i & j as displacements (no grid stretching)
524  int delta_i = static_cast<int>(std::round(delta_x*dxInv[0]));
525  int delta_j = static_cast<int>(std::round(delta_y*dxInv[1]));
526  int i_new = i + delta_i;
527  int j_new = j + delta_j;
528  i_arr(i,j,k) = i_new;
529  j_arr(i,j,k) = j_new;
530 
531  // Search for k (grid is stretched in z)
532  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
533  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
534  Real z_target = z_bot_face + delta_z;
535  for (int lk(0); lk<=kmax; ++lk) {
536  Real z_lo = 0.25 * ( z_phys_arr(i_new,j_new ,lk ) + z_phys_arr(i_new+1,j_new ,lk )
537  + z_phys_arr(i_new,j_new+1,lk ) + z_phys_arr(i_new+1,j_new+1,lk ) );
538  Real z_hi = 0.25 * ( z_phys_arr(i_new,j_new ,lk+1) + z_phys_arr(i_new+1,j_new ,lk+1)
539  + z_phys_arr(i_new,j_new+1,lk+1) + z_phys_arr(i_new+1,j_new+1,lk+1) );
540  if (z_target > z_lo && z_target < z_hi){
541  AMREX_ASSERT_WITH_MESSAGE(lk >= d_radius,
542  "K index must be larger than averaging radius!");
543  amrex::ignore_unused(d_radius);
544  k_arr(i,j,k) = lk;
545  break;
546  }
547  }
548 
549  // Destination cell must be contained on the current process!
550  amrex::ignore_unused(gpbx);
551  AMREX_ASSERT_WITH_MESSAGE(gpbx.contains(i_arr(i,j,k),j_arr(i,j,k),k_arr(i,j,k)),
552  "Query index outside of proc domain!");
553  });
554  }
555 }
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.

609 {
610  ParmParse pp(m_pp_prefix);
611  pp.get("most.zref",m_zref);
612 
613  // Capture for device
614  Real d_zref = m_zref;
615  const auto plo = m_geom[lev].ProbLoArray();
616 
617  RealVect base;
618  const auto dx = m_geom[lev].CellSizeArray();
619  const auto dxInv = m_geom[lev].InvCellSizeArray();
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  // Elements of normal vector
633  Real met_h_xi = Compute_h_xi_AtCellCenter (i,j,k,dxInv,z_phys_arr);
634  Real met_h_eta = Compute_h_eta_AtCellCenter(i,j,k,dxInv,z_phys_arr);
635  Real imag = 1.0 / std::sqrt(met_h_xi*met_h_xi + met_h_eta*met_h_eta + 1.0);
636 
637  // Unit-normal vector scaled by z_ref
638  Real delta_x = -met_h_xi * imag * d_zref;
639  Real delta_y = -met_h_eta * imag * d_zref;
640  Real delta_z = imag * d_zref;
641 
642  // Position of the current node (indx:0,0,1)
643  Real x0 = plo[0] + ((Real) i + 0.5) * dx[0];
644  Real y0 = plo[1] + ((Real) j + 0.5) * dx[1];
645 
646  // Final position at end of vector
647  x_pos_arr(i,j,k) = x0 + delta_x;
648  y_pos_arr(i,j,k) = y0 + delta_y;
649  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
650  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
651  z_pos_arr(i,j,k) = z_bot_face + delta_z;
652 
653  // NOTE: Normal vector end point can be below the surface for concave regions.
654  // Here we protect against that by augmenting the normal if needed.
655  int i_new = (int) ((x_pos_arr(i,j,k) - plo[0]) / dx[0] - 0.5);
656  int j_new = (int) ((y_pos_arr(i,j,k) - plo[1]) / dx[1] - 0.5);
657  Real z_new_bot_face = 0.25 * ( z_phys_arr(i_new,j_new ,k) + z_phys_arr(i_new+1,j_new ,k)
658  + z_phys_arr(i_new,j_new+1,k) + z_phys_arr(i_new+1,j_new+1,k) );
659  if (z_pos_arr(i,j,k) < z_new_bot_face) {
660  z_pos_arr(i,j,k) = z_new_bot_face + delta_z;
661  }
662 
663  // Destination position must be contained on the current process!
664  Real pos[] = {x_pos_arr(i,j,k)-plo[0],y_pos_arr(i,j,k)-plo[1],0.5*dx[2]};
665  amrex::ignore_unused(pos);
666  AMREX_ASSERT_WITH_MESSAGE( grb.contains(&pos[0]),
667  "Query point outside of proc domain!");
668  });
669  }
670 }

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.

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

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
61  {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.

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

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.

564 {
565  ParmParse pp(m_pp_prefix);
566  pp.get("most.zref",m_zref);
567 
568  // Capture for device
569  Real d_zref = m_zref;
570  const auto plo = m_geom[lev].ProbLoArray();
571 
572  RealVect base;
573  const auto dx = m_geom[lev].CellSizeArray();
574  IntVect ng = m_x_pos[lev]->nGrowVect(); ng[2]=0;
575  for (MFIter mfi(*m_x_pos[lev], TileNoZ()); mfi.isValid(); ++mfi) {
576  Box npbx = mfi.tilebox(IntVect(1,1,0),IntVect(1,1,0));
577  Box gpbx = mfi.growntilebox(ng);
578  RealBox grb{gpbx,dx.data(),base.dataPtr()};
579 
580  const auto z_phys_arr = m_z_phys_nd[lev]->const_array(mfi);
581  auto x_pos_arr = m_x_pos[lev]->array(mfi);
582  auto y_pos_arr = m_y_pos[lev]->array(mfi);
583  auto z_pos_arr = m_z_pos[lev]->array(mfi);
584  ParallelFor(npbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
585  {
586  // Final position at end of vector
587  x_pos_arr(i,j,k) = plo[0] + ((Real) i + 0.5) * dx[0];
588  y_pos_arr(i,j,k) = plo[1] + ((Real) j + 0.5) * dx[1];
589  Real z_bot_face = 0.25 * ( z_phys_arr(i ,j ,k) + z_phys_arr(i+1,j ,k)
590  + z_phys_arr(i ,j+1,k) + z_phys_arr(i+1,j+1,k) );
591  z_pos_arr(i,j,k) = z_bot_face + d_zref;
592 
593  // Destination position must be contained on the current process!
594  Real pos[] = {x_pos_arr(i,j,k)-plo[0],y_pos_arr(i,j,k)-plo[1],0.5*dx[2]};
595  amrex::ignore_unused(pos);
596  AMREX_ASSERT_WITH_MESSAGE( grb.contains(&pos[0]),
597  "Query point outside of proc domain!");
598  });
599  }
600 }

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

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
251 {
252  m_fields[lev][0] = &vars_old[lev][Vars::xvel];
253  m_fields[lev][1] = &vars_old[lev][Vars::yvel];
254  m_fields[lev][2] = Theta_prim[lev].get();
255  m_fields[lev][3] = Qv_prim[lev].get();
256  m_fields[lev][4] = Qr_prim[lev].get();
257  m_fields[lev][5] = &vars_old[lev][Vars::zvel];
258 }

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
1467 {
1468  // Peel back the level
1469  auto& averages = m_averages[lev];
1470 
1471  int navg = m_navg - 1;
1472 
1473  std::ofstream ofile;
1474  ofile.open ("MOST_averages.txt");
1475  ofile << "Averages computed via MOSTAverages class:\n";
1476 
1477  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1478  Box bx = mfi.tilebox(); bx.setBig(2,0);
1479  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1480  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1481 
1482  for (int j(jl); j <= ju; ++j) {
1483  for (int i(il); i <= iu; ++i) {
1484  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1485  int k = 0;
1486  for (int iavg(0); iavg <= navg; ++iavg) {
1487  auto mf_arr = averages[iavg]->array(mfi);
1488  ofile << "iavg val: "
1489  << iavg << ' '
1490  << mf_arr(i,j,k) << "\n";
1491  }
1492  ofile << "\n";
1493  }
1494  }
1495  }
1496  ofile.close();
1497 }
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
1348 {
1349  // Peel back the level
1350  auto& averages = m_averages[lev];
1351  auto& k_indx = m_k_indx[lev];
1352 
1353  int navg = m_navg - 1;
1354 
1355  std::ofstream ofile;
1356  ofile.open ("MOST_k_indices.txt");
1357  ofile << "K indices used to compute averages via MOSTAverages class:\n";
1358 
1359  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1360  Box bx = mfi.tilebox(); bx.setBig(2,0);
1361  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1362  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1363 
1364  auto k_arr = k_indx->array(mfi);
1365 
1366  for (int j(jl); j <= ju; ++j) {
1367  for (int i(il); i <= iu; ++i) {
1368  ofile << "(I,J): " << "(" << i << "," << j << ")" << "\n";
1369  int k = 0;
1370  ofile << "K_ind: "
1371  << k_arr(i,j,k) << "\n";
1372  ofile << "\n";
1373  }
1374  }
1375  }
1376  ofile.close();
1377 }
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
1387 {
1388  // Peel back the level
1389  auto& averages = m_averages[lev];
1390  auto& k_indx = m_k_indx[lev];
1391  auto& j_indx = m_j_indx[lev];
1392  auto& i_indx = m_i_indx[lev];
1393 
1394  int navg = m_navg - 1;
1395 
1396  std::ofstream ofile;
1397  ofile.open ("MOST_ijk_indices.txt");
1398  ofile << "IJK indices used to compute averages via MOSTAverages class:\n";
1399 
1400  for (MFIter mfi(*averages[navg], TileNoZ()); mfi.isValid(); ++mfi) {
1401  Box bx = mfi.tilebox(); bx.setBig(2,0);
1402  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1403  int jl = bx.smallEnd(1); int ju = bx.bigEnd(1);
1404 
1405  auto k_arr = k_indx->array(mfi);
1406  auto j_arr = j_indx ? j_indx->array(mfi) : Array4<int> {};
1407  auto i_arr = i_indx ? i_indx->array(mfi) : Array4<int> {};
1408 
1409  for (int j(jl); j <= ju; ++j) {
1410  for (int i(il); i <= iu; ++i) {
1411  ofile << "(I1,J1,K1): " << "(" << i << "," << j << "," << 0 << ")" << "\n";
1412 
1413  int k = 0;
1414  int km = k_arr(i,j,k);
1415  int jm = j_arr ? j_arr(i,j,k) : j;
1416  int im = i_arr ? i_arr(i,j,k) : i;
1417 
1418  ofile << "(I2,J2,K2): "
1419  << "(" << im << "," << jm << "," << km << ")" << "\n";
1420  ofile << "\n";
1421  }
1422  }
1423  }
1424  ofile.close();
1425 }
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
1437 {
1438  // Peel back the level
1439  auto& x_pos_mf = m_x_pos[lev];
1440  auto& z_pos_mf = m_z_pos[lev];
1441 
1442  std::ofstream ofile;
1443  ofile.open ("MOST_xz_positions.txt");
1444 
1445  for (MFIter mfi(*x_pos_mf, TileNoZ()); mfi.isValid(); ++mfi) {
1446  Box bx = mfi.tilebox(); bx.setBig(2,0);
1447  int il = bx.smallEnd(0); int iu = bx.bigEnd(0);
1448 
1449  auto x_pos_arr = x_pos_mf->array(mfi);
1450  auto z_pos_arr = z_pos_mf->array(mfi);
1451 
1452  int k = 0;
1453  for (int i(il); i <= iu; ++i)
1454  ofile << x_pos_arr(i,j,k) << ' ' << z_pos_arr(i,j,k) << "\n";
1455  }
1456  ofile.close();
1457 }
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_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::Real MOSTAverage::m_zref {10.0}
protected

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