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