ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
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 #include "ERF_EB.H"
14 
15 class MOSTAverage {
16 public:
17  explicit MOSTAverage (amrex::Vector<amrex::Geometry> geom,
18  const bool& has_zphys,
19  std::string a_pp_prefix,
20  const MeshType& m_mesh_type,
21  const TerrainType& m_terrain_type,
22  const amrex::Vector<const eb_*>& eb_vec = {});
23 
24  // MOSTAverage() = default;
26  {}
27 
28  // Declare a default move constructor so we ensure the destructor is
29  // not called when we return an object of this class by value
30  MOSTAverage (MOSTAverage&&) noexcept = default;
31 
32  // Delete the move assignment operator
33  MOSTAverage& operator=(MOSTAverage&& other) noexcept = delete;
34 
35  // Delete the copy constructor
36  MOSTAverage (const MOSTAverage& other) = delete;
37 
38  // Delete the copy assignment operator
39  MOSTAverage& operator=(const MOSTAverage& other) = delete;
40 
41  // Make data structures at level
42  void make_MOSTAverage_at_level (const int& lev,
43  const amrex::Vector<amrex::MultiFab*>& vars_old,
44  std::unique_ptr<amrex::MultiFab>& Theta_prim,
45  std::unique_ptr<amrex::MultiFab>& Qv_prim,
46  std::unique_ptr<amrex::MultiFab>& Qr_prim,
47  std::unique_ptr<amrex::MultiFab>& z_phys_nd);
48 
49  // Reset the pointers to field MFs
50  void update_field_ptrs (const int& lev,
51  amrex::Vector<amrex::Vector<amrex::MultiFab>>& vars_old,
52  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Theta_prim,
53  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Qv_prim,
54  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& Qr_prim);
55 
56  // Update the rotated fields
57  void set_rotated_fields (const int& lev);
58 
59  // Compute ncells per plane
60  void set_plane_normalization (const int& lev);
61 
62  // Compute total boundary area over EB
63  void set_eb_normalization (const int& lev);
64 
65  // Compute ncells per plane
66  void set_region_normalization (const int& /*lev*/)
67  {m_ncell_region = (2 * m_radius + 1) * (2 * m_radius + 1) * (2 * m_radius + 1);}
68 
69  // Populate a 2D iMF k_indx (w/o terrain)
70  void set_k_indices_N (const int& lev);
71 
72  // Populate a 2D iMF k_indx (w/ terrain)
73  void set_k_indices_T (const int& lev);
74 
75  // Populate all 2D iMFs ijk_indx (w/ terrain)
76  void set_norm_indices_T (const int& lev);
77 
78  // Populate positions (w/ terrain & norm vector & interpolation)
79  void set_z_positions_T (const int& lev);
80 
81  // Populate positions (for EB)
82  void set_z_positions_EB (const int& lev);
83 
84  // Populate positions (w/ terrain & norm vector & interpolation)
85  void set_norm_positions_T (const int& lev);
86 
87  // Driver for the different average policies
88  void compute_averages (const int& lev);
89 
90  // Fill averages for policy::plane
91  void compute_plane_averages (const int& lev);
92 
93  // Fill averages for policy::point
94  void compute_region_averages (const int& lev);
95 
96  // Fill averages for policy::EB
97  void compute_eb_averages (const int& lev);
98 
99  // Write k indices
100  void write_k_indices (const int& lev);
101 
102  // Write ijk indices
103  void write_norm_indices (const int& lev);
104 
105  // Write XZ planar positions
106  void write_xz_positions (const int& lev,
107  const int& j);
108 
109  // Write averages on 2D mf
110  void write_averages (const int& lev);
111 
112  // Get pointer to the 2D mf of averages
113  [[nodiscard]] const amrex::MultiFab* get_average (const int& lev, const int& comp) const { return m_averages[lev][comp].get(); }
114 
115  // Get z_ref (may be computed from specified k_indx)
116  [[nodiscard]] amrex::MultiFab* get_zref (const int& lev) const { return m_zref[lev].get(); }
117 
118  // Get pointer to the k-index mf
119  [[nodiscard]] const amrex::iMultiFab* get_k_indices (const int& lev) const { return m_k_indx[lev].get(); }
120 
121  /**
122  * Function to compute trilinear interpolation with terrain.
123  *
124  * @param[in] xp X-position
125  * @param[in] yp Y-position
126  * @param[out] interp_vals Values interpolated
127  * @param[in] interp_array Array to interpolate on
128  * @param[in] z_arr Physical heights
129  * @param[in] plo Problem lower bounds
130  * @param[in] dxi Inverse cell size array
131  * @param[in] interp_comp Number of components to interpolate
132  */
133  AMREX_GPU_HOST_DEVICE AMREX_INLINE
134  static void trilinear_interp_T (const amrex::Real& xp,
135  const amrex::Real& yp,
136  const amrex::Real& zp,
137  amrex::Real* interp_vals,
138  amrex::Array4<amrex::Real const> const& interp_array,
139  amrex::Array4<amrex::Real const> const& z_arr,
140  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& plo,
141  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& dxi,
142  const int interp_comp)
143  {
144  // Search to get z/k
145  bool found = false;
146  int kmax = ubound(z_arr).z;
147  amrex::Real zval = zero;
148  amrex::Real z_target = zp;
149 
150  // Map position to i,j (must be same mapping in cpp file)
151  amrex::Real ireal = (xp - plo[0]) * dxi[0];
152  amrex::Real jreal = (yp - plo[1]) * dxi[1];
153  int i_new = (int) (ireal - myhalf);
154  int j_new = (int) (jreal - myhalf);
155 
156  for (int lk(0); lk<kmax; ++lk) {
157  amrex::Real z_lo = fourth * ( z_arr(i_new,j_new ,lk ) + z_arr(i_new+1,j_new ,lk )
158  + z_arr(i_new,j_new+1,lk ) + z_arr(i_new+1,j_new+1,lk ) );
159  amrex::Real z_hi = fourth * ( z_arr(i_new,j_new ,lk+1) + z_arr(i_new+1,j_new ,lk+1)
160  + z_arr(i_new,j_new+1,lk+1) + z_arr(i_new+1,j_new+1,lk+1) );
161  if (z_target > z_lo && z_target < z_hi){
162  found = true;
163  zval = (amrex::Real) lk + ((z_target - z_lo) / (z_hi - z_lo)) + myhalf;
164  break;
165  }
166  }
167 
168  amrex::ignore_unused(found);
169  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(found, "MOSTAverage: Height above terrain not found, try increasing z_ref!");
170 
171  // NOTE: This is the point ahead of the current i,j (e.g. i/j_new + 1)
172  const amrex::RealVect lx(ireal + myhalf, jreal + myhalf, zval);
173 
174  const amrex::IntVect ijk = lx.floor();
175 
176  int i = ijk[0]; int j = ijk[1]; int k = ijk[2];
177 
178  // Convert ijk (IntVect) to a RealVect explicitly
179  amrex::RealVect ijk_r(static_cast<amrex::Real>(ijk[0]),
180  static_cast<amrex::Real>(ijk[1]),
181  static_cast<amrex::Real>(ijk[2]));
182 
183  // Weights
184  const amrex::RealVect sx_hi = lx - ijk_r;
185  const amrex::RealVect sx_lo = one - sx_hi;
186 
187  for (int n = 0; n < interp_comp; n++) {
188  interp_vals[n] = sx_lo[0]*sx_lo[1]*sx_lo[2]*interp_array(i-1, j-1, k-1,n) +
189  sx_lo[0]*sx_lo[1]*sx_hi[2]*interp_array(i-1, j-1, k ,n) +
190  sx_lo[0]*sx_hi[1]*sx_lo[2]*interp_array(i-1, j , k-1,n) +
191  sx_lo[0]*sx_hi[1]*sx_hi[2]*interp_array(i-1, j , k ,n) +
192  sx_hi[0]*sx_lo[1]*sx_lo[2]*interp_array(i , j-1, k-1,n) +
193  sx_hi[0]*sx_lo[1]*sx_hi[2]*interp_array(i , j-1, k ,n) +
194  sx_hi[0]*sx_hi[1]*sx_lo[2]*interp_array(i , j , k-1,n) +
195  sx_hi[0]*sx_hi[1]*sx_hi[2]*interp_array(i , j , k ,n);
196  }
197  }
198 
199 protected:
200 
201  // Passed through constructor
202  //--------------------------------------------
203  const amrex::Vector<amrex::Geometry> m_geom; // Geometry at each level
204  amrex::Vector<amrex::Vector<amrex::MultiFab*>> m_fields; // Ptr to fields to be averaged
205  amrex::Vector<amrex::MultiFab*> m_z_phys_nd; // Ptr to terrain height coords
206  std::string m_pp_prefix; // ParmParse prefix
207  MeshType m_mesh_type; // Mesh type
208  TerrainType m_terrain_type; // Terrain type
209 
210  // General vars for multiple or all policies
211  //--------------------------------------------
212  int m_nvar{6}; // 6 fields for U/V/T/Qv/Qr/W
213  int m_navg{6}; // 6 averages for U/V/T/Qv/Tv/Umag
214  int m_maxlev{0}; // Total number of levels
215  int m_policy{0}; // Policy for type of averaging
216  bool m_rotate{false}; // Do vector rotations for terrain?
217  amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_zref; // Height above surface for MOST BC
218  amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_x_pos; // Ptr to 2D mf to hold x position (maxlev)
219  amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_y_pos; // Ptr to 2D mf to hold y position (maxlev)
220  amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_z_pos; // Ptr to 2D mf to hold z position (maxlev)
221  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> m_i_indx; // Ptr to 2D imf to hold i indices (maxlev)
222  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> m_j_indx; // Ptr to 2D imf to hold j indices (maxlev)
223  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> m_k_indx; // Ptr to 2D imf to hold k indices (maxlev)
224  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> m_averages; // Ptr to 2D mf to hold averages (maxlev,navg)
225  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> m_rot_fields; // Rotated field data
226 
227  // Vars for planar average policy
228  //--------------------------------------------
229  amrex::Vector<amrex::Vector<int>> m_ncell_plane; // Number of cells in plane (maxlev,navg)
230  amrex::Vector<amrex::Vector<amrex::Real>> m_plane_average; // Plane avgs (maxlev,navg)
231 
232  // Vars for point/region average policy
233  //--------------------------------------------
234  int m_radius{0}; // Radius around k_index
235  int m_ncell_region{1}; // Number of cells in local region
236  amrex::Vector<int> m_k_in; // Specified k_index for region avg (maxlev)
237 
238  // Vars for normal vector policy
239  //--------------------------------------------
240  bool m_interp{false}; // Do interpolation on destination?
241  bool m_norm_vec{false}; // Use normal vector to find IJK?
242 
243  // Vars for EB averaging policy
244  //--------------------------------------------
245  amrex::Vector<const eb_*> m_eb_vec; // Vector of pointers to EB objects (one per level)
246  amrex::Vector<amrex::Vector<amrex::Real>> m_total_bndry_area; // Total EB surface areas (maxlev, navg)
247 
248  // Time average w/ exponential filter fun
249  //--------------------------------------------
250  bool m_t_avg{false}; // Flag to do moving average in time
251  amrex::Vector<int> m_t_init; // Flag to specify if averages are initialized
252  amrex::Real m_time_window{amrex::Real(1.0e-16)}; // Width of the exp filter function
253  amrex::Real m_fact_new, m_fact_old; // Time average factors for new and old means
254 
255  // Surface velocity correction
256  //--------------------------------------------
257  bool include_subgrid_vel = false;
258  amrex::Vector<amrex::Real> m_Vsg; // Subgrid velocity scale (Mahrt & Sun 1995 MWR)
259 
260  // Default values
261  //--------------------------------------------
263 };
264 #endif
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_MOSTAverage.H:15
int m_navg
Definition: ERF_MOSTAverage.H:213
bool m_t_avg
Definition: ERF_MOSTAverage.H:250
void write_xz_positions(const int &lev, const int &j)
Definition: ERF_MOSTAverage.cpp:1997
void compute_plane_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:929
int m_policy
Definition: ERF_MOSTAverage.H:215
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_averages
Definition: ERF_MOSTAverage.H:224
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_y_pos
Definition: ERF_MOSTAverage.H:219
int m_radius
Definition: ERF_MOSTAverage.H:234
void set_z_positions_EB(const int &lev)
Definition: ERF_MOSTAverage.cpp:558
void write_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:2039
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_i_indx
Definition: ERF_MOSTAverage.H:221
void compute_region_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:1228
amrex::Vector< amrex::MultiFab * > m_z_phys_nd
Definition: ERF_MOSTAverage.H:205
amrex::Vector< int > m_t_init
Definition: ERF_MOSTAverage.H:251
void compute_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:900
void set_rotated_fields(const int &lev)
Definition: ERF_MOSTAverage.cpp:298
amrex::Vector< amrex::Vector< amrex::Real > > m_total_bndry_area
Definition: ERF_MOSTAverage.H:246
void set_region_normalization(const int &)
Definition: ERF_MOSTAverage.H:66
void set_z_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:748
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_x_pos
Definition: ERF_MOSTAverage.H:218
void set_norm_positions_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:811
amrex::Vector< amrex::Real > m_Vsg
Definition: ERF_MOSTAverage.H:258
void set_k_indices_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:580
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > m_rot_fields
Definition: ERF_MOSTAverage.H:225
bool m_rotate
Definition: ERF_MOSTAverage.H:216
amrex::MultiFab * get_zref(const int &lev) const
Definition: ERF_MOSTAverage.H:116
amrex::Vector< amrex::Vector< amrex::Real > > m_plane_average
Definition: ERF_MOSTAverage.H:230
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_z_pos
Definition: ERF_MOSTAverage.H:220
amrex::Vector< amrex::Vector< int > > m_ncell_plane
Definition: ERF_MOSTAverage.H:229
void set_plane_normalization(const int &lev)
Definition: ERF_MOSTAverage.cpp:357
void write_k_indices(const int &lev)
Definition: ERF_MOSTAverage.cpp:1893
TerrainType m_terrain_type
Definition: ERF_MOSTAverage.H:208
std::string m_pp_prefix
Definition: ERF_MOSTAverage.H:206
bool m_norm_vec
Definition: ERF_MOSTAverage.H:241
int m_nvar
Definition: ERF_MOSTAverage.H:212
amrex::Real m_fact_new
Definition: ERF_MOSTAverage.H:253
const amrex::iMultiFab * get_k_indices(const int &lev) const
Definition: ERF_MOSTAverage.H:119
int m_ncell_region
Definition: ERF_MOSTAverage.H:235
void set_norm_indices_T(const int &lev)
Definition: ERF_MOSTAverage.cpp:660
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_zref
Definition: ERF_MOSTAverage.H:217
int m_maxlev
Definition: ERF_MOSTAverage.H:214
const amrex::MultiFab * get_average(const int &lev, const int &comp) const
Definition: ERF_MOSTAverage.H:113
void set_eb_normalization(const int &lev)
Definition: ERF_MOSTAverage.cpp:414
bool include_subgrid_vel
Definition: ERF_MOSTAverage.H:257
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:279
~MOSTAverage()
Definition: ERF_MOSTAverage.H:25
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, const amrex::Vector< const eb_ * > &eb_vec={})
amrex::Real m_time_window
Definition: ERF_MOSTAverage.H:252
MeshType m_mesh_type
Definition: ERF_MOSTAverage.H:207
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_j_indx
Definition: ERF_MOSTAverage.H:222
void set_k_indices_N(const int &lev)
Definition: ERF_MOSTAverage.cpp:502
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:90
MOSTAverage(MOSTAverage &&) noexcept=default
void write_norm_indices(const int &lev)
Definition: ERF_MOSTAverage.cpp:1940
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_fields
Definition: ERF_MOSTAverage.H:204
amrex::Vector< int > m_k_in
Definition: ERF_MOSTAverage.H:236
amrex::Vector< const eb_ * > m_eb_vec
Definition: ERF_MOSTAverage.H:245
const amrex::Real zref_default
Definition: ERF_MOSTAverage.H:262
amrex::Real m_fact_old
Definition: ERF_MOSTAverage.H:253
bool m_interp
Definition: ERF_MOSTAverage.H:240
const amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_MOSTAverage.H:203
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:134
void compute_eb_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:1605
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > m_k_indx
Definition: ERF_MOSTAverage.H:223
Definition: ERF_ConsoleIO.cpp:12