ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_MOSTAverage.H
Go to the documentation of this file.
1 #ifndef ERF_MOSTAverage_H
2 #define ERF_MOSTAverage_H
3 
4 #include "AMReX_Gpu.H"
5 #include "AMReX_FArrayBox.H"
6 #include "AMReX_MultiFab.H"
7 #include "AMReX_iMultiFab.H"
8 #include "AMReX_ParmParse.H"
9 
10 #include "ERF_IndexDefines.H"
11 #include "ERF_TerrainMetrics.H"
12 #include "ERF_DataStruct.H"
13 
14 class MOSTAverage {
15 public:
16  explicit MOSTAverage (amrex::Vector<amrex::Geometry> geom,
17  const bool& has_zphys,
18  std::string a_pp_prefix,
19  const TerrainType& m_terrain_type);
20 
21  // MOSTAverage() = default;
23  {}
24 
25  // Declare a default move constructor so we ensure the destructor is
26  // not called when we return an object of this class by value
27  MOSTAverage (MOSTAverage&&) noexcept = default;
28 
29  // Delete the move assignment operator
30  MOSTAverage& operator=(MOSTAverage&& other) noexcept = delete;
31 
32  // Delete the copy constructor
33  MOSTAverage (const MOSTAverage& other) = delete;
34 
35  // Delete the copy assignment operator
36  MOSTAverage& operator=(const MOSTAverage& other) = delete;
37 
38  // Make data structures at level
39  void make_MOSTAverage_at_level (const int& lev,
40  const amrex::Vector<amrex::MultiFab*>& vars_old,
41  std::unique_ptr<amrex::MultiFab>& Theta_prim,
42  std::unique_ptr<amrex::MultiFab>& Qv_prim,
43  std::unique_ptr<amrex::MultiFab>& Qr_prim,
44  std::unique_ptr<amrex::MultiFab>& z_phys_nd);
45 
46  // Reset the pointers to field MFs
47  void update_field_ptrs (const int& lev,
48  amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_old,
49  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Theta_prim,
50  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Qv_prim,
51  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Qr_prim);
52 
53  // Update the rotated fields
54  void set_rotated_fields (const int& lev);
55 
56  // Compute ncells per plane
57  void set_plane_normalization (const int& lev);
58 
59  // Compute ncells per plane
60  void set_region_normalization (const int& /*lev*/)
61  {m_ncell_region = (2 * m_radius + 1) * (2 * m_radius + 1) * (2 * m_radius + 1);}
62 
63  // Populate a 2D iMF k_indx (w/o terrain)
64  void set_k_indices_N (const int& lev);
65 
66  // Populate a 2D iMF k_indx (w/ terrain)
67  void set_k_indices_T (const int& lev);
68 
69  // Populate all 2D iMFs ijk_indx (w/ terrain)
70  void set_norm_indices_T (const int& lev);
71 
72  // Populate positions (w/ terrain & norm vector & interpolation)
73  void set_z_positions_T (const int& lev);
74 
75  // Populate positions (w/ terrain & norm vector & interpolation)
76  void set_norm_positions_T (const int& lev);
77 
78  // Driver for the different average policies
79  void compute_averages (const int& lev);
80 
81  // Fill averages for policy::plane
82  void compute_plane_averages (const int& lev);
83 
84  // Fill averages for policy::point
85  void compute_region_averages (const int& lev);
86 
87  // Write k indices
88  void write_k_indices (const int& lev);
89 
90  // Write ijk indices
91  void write_norm_indices (const int& lev);
92 
93  // Write XZ planar positions
94  void write_xz_positions (const int& lev,
95  const int& j);
96 
97  // Write averages on 2D mf
98  void write_averages (const int& lev);
99 
100  // Get pointer to the 2D mf of averages
101  [[nodiscard]] const amrex::MultiFab* get_average (const int& lev, const int& comp) const { return m_averages[lev][comp].get(); }
102 
103  // Get z_ref (may be computed from specified k_indx)
104  [[nodiscard]] amrex::Real get_zref () const { return m_zref; }
105 
106  /**
107  * Function to compute trilinear interpolation with terrain.
108  *
109  * @param[in] xp X-position
110  * @param[in] yp Y-position
111  * @param[out] interp_vals Values interpolated
112  * @param[in] interp_array Array to interpolate on
113  * @param[in] z_arr Physical heights
114  * @param[in] plo Problem lower bounds
115  * @param[in] dxi Inverse cell size array
116  * @param[in] interp_comp Number of components to interpolate
117  */
118  AMREX_GPU_HOST_DEVICE AMREX_INLINE
119  static void trilinear_interp_T (const amrex::Real& xp,
120  const amrex::Real& yp,
121  const amrex::Real& zp,
122  amrex::Real* interp_vals,
123  amrex::Array4<amrex::Real const> const& interp_array,
124  amrex::Array4<amrex::Real const> const& z_arr,
125  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& plo,
126  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& dxi,
127  const int interp_comp)
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  }
183 
184 protected:
185 
186  // Passed through constructor
187  //--------------------------------------------
188  const amrex::Vector<amrex::Geometry> m_geom; // Geometry at each level
189  amrex::Vector<amrex::Vector<amrex::MultiFab*>> m_fields; // Ptr to fields to be averaged
190  amrex::Vector<amrex::MultiFab*> m_z_phys_nd; // Ptr to terrain height coords
191  std::string m_pp_prefix; // ParmParse prefix
192  TerrainType m_terrain_type; // Terrain type
193 
194  // General vars for multiple or all policies
195  //--------------------------------------------
196  int m_nvar{6}; // 6 fields for U/V/T/Qv/Qr/W
197  int m_navg{6}; // 6 averages for U/V/T/Qv/Tv/Umag
198  int m_maxlev{0}; // Total number of levels
199  int m_policy{0}; // Policy for type of averaging
200  bool m_rotate{false}; // Do vector rotations for terrain?
201  amrex::Real m_zref{10.0}; // Height above surface for MOST BC
202  amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_x_pos; // Ptr to 2D mf to hold x position (maxlev)
203  amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_y_pos; // Ptr to 2D mf to hold y position (maxlev)
204  amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_z_pos; // Ptr to 2D mf to hold z position (maxlev)
205  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> m_i_indx; // Ptr to 2D imf to hold i indices (maxlev)
206  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> m_j_indx; // Ptr to 2D imf to hold j indices (maxlev)
207  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> m_k_indx; // Ptr to 2D imf to hold k indices (maxlev)
208  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> m_averages; // Ptr to 2D mf to hold averages (maxlev,navg)
209  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> m_rot_fields; // Rotated field data
210 
211  // Vars for planar average policy
212  //--------------------------------------------
213  amrex::Vector<amrex::Vector<int>> m_ncell_plane; // Number of cells in plane (maxlev,navg)
214  amrex::Vector<amrex::Vector<amrex::Real>> m_plane_average; // Plane avgs (maxlev,navg)
215 
216  // Vars for point/region average policy
217  //--------------------------------------------
218  int m_radius{0}; // Radius around k_index
219  int m_ncell_region{1}; // Number of cells in local region
220  amrex::Vector<int> m_k_in; // Specified k_index for region avg (maxlev)
221 
222  // Vars for normal vector policy
223  //--------------------------------------------
224  bool m_interp{false}; // Do interpolation on destination?
225  bool m_norm_vec{false}; // Use normal vector to find IJK?
226 
227  // Time average w/ exponential filter fun
228  //--------------------------------------------
229  bool m_t_avg{false}; // Flag to do moving average in time
230  amrex::Vector<int> m_t_init; // Flag to specify if averages are initialized
231  amrex::Real m_time_window{1.0e-16}; // Width of the exp filter function
232  amrex::Real m_fact_new, m_fact_old; // Time average factors for new and old means
233 
234  // Surface velocity correction
235  //--------------------------------------------
236  bool include_subgrid_vel = false;
237  amrex::Vector<amrex::Real> m_Vsg; // Subgrid velocity scale (Mahrt & Sun 1995 MWR)
238 };
239 #endif
Definition: ERF_MOSTAverage.H:14
int m_navg
Definition: ERF_MOSTAverage.H:197
bool m_t_avg
Definition: ERF_MOSTAverage.H:229
void write_xz_positions(const int &lev, const int &j)
Definition: ERF_MOSTAverage.cpp:1435
void compute_plane_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:705
int m_policy
Definition: ERF_MOSTAverage.H:199
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
int m_radius
Definition: ERF_MOSTAverage.H:218
void write_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:1466
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_i_indx
Definition: ERF_MOSTAverage.H:205
void compute_region_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:990
amrex::Vector< amrex::MultiFab * > m_z_phys_nd
Definition: ERF_MOSTAverage.H:190
amrex::Vector< int > m_t_init
Definition: ERF_MOSTAverage.H:230
void compute_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:679
void set_rotated_fields(const int &lev)
Definition: ERF_MOSTAverage.cpp:265
void set_region_normalization(const int &)
Definition: ERF_MOSTAverage.H:60
void set_z_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:563
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_x_pos
Definition: ERF_MOSTAverage.H:202
amrex::Real get_zref() const
Definition: ERF_MOSTAverage.H:104
void set_norm_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:608
amrex::Vector< amrex::Real > m_Vsg
Definition: ERF_MOSTAverage.H:237
void set_k_indices_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:437
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_rot_fields
Definition: ERF_MOSTAverage.H:209
bool m_rotate
Definition: ERF_MOSTAverage.H:200
MOSTAverage(amrex::Vector< amrex::Geometry > geom, const bool &has_zphys, std::string a_pp_prefix, const TerrainType &m_terrain_type)
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
void set_plane_normalization(const int &lev)
Definition: ERF_MOSTAverage.cpp:324
void write_k_indices(const int &lev)
Definition: ERF_MOSTAverage.cpp:1347
TerrainType m_terrain_type
Definition: ERF_MOSTAverage.H:192
std::string m_pp_prefix
Definition: ERF_MOSTAverage.H:191
bool m_norm_vec
Definition: ERF_MOSTAverage.H:225
int m_nvar
Definition: ERF_MOSTAverage.H:196
amrex::Real m_fact_new
Definition: ERF_MOSTAverage.H:232
int m_ncell_region
Definition: ERF_MOSTAverage.H:219
void set_norm_indices_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:492
int m_maxlev
Definition: ERF_MOSTAverage.H:198
const amrex::MultiFab * get_average(const int &lev, const int &comp) const
Definition: ERF_MOSTAverage.H:101
bool include_subgrid_vel
Definition: ERF_MOSTAverage.H:236
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)
Definition: ERF_MOSTAverage.cpp:246
~MOSTAverage()
Definition: ERF_MOSTAverage.H:22
amrex::Real m_time_window
Definition: ERF_MOSTAverage.H:231
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_j_indx
Definition: ERF_MOSTAverage.H:206
void set_k_indices_N(const int &lev)
Definition: ERF_MOSTAverage.cpp:382
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)
Definition: ERF_MOSTAverage.cpp:70
MOSTAverage(MOSTAverage &&) noexcept=default
void write_norm_indices(const int &lev)
Definition: ERF_MOSTAverage.cpp:1386
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_fields
Definition: ERF_MOSTAverage.H:189
amrex::Vector< int > m_k_in
Definition: ERF_MOSTAverage.H:220
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::Real m_zref
Definition: ERF_MOSTAverage.H:201
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_k_indx
Definition: ERF_MOSTAverage.H:207
Definition: ERF_ConsoleIO.cpp:12