1 #ifndef ERF_MOSTAverage_H
2 #define ERF_MOSTAverage_H
5 #include "AMReX_FArrayBox.H"
6 #include "AMReX_MultiFab.H"
7 #include "AMReX_iMultiFab.H"
8 #include "AMReX_ParmParse.H"
17 const bool& has_zphys,
18 std::string a_pp_prefix,
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);
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);
101 [[nodiscard]]
const amrex::MultiFab*
get_average (
const int& lev,
const int& comp)
const {
return m_averages[lev][comp].get(); }
118 AMREX_GPU_HOST_DEVICE AMREX_INLINE
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)
131 int kmax = ubound(z_arr).z;
132 amrex::Real zval = 0.0;
133 amrex::Real z_target = zp;
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);
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){
148 zval = (amrex::Real) lk + ((z_target - z_lo) / (z_hi - z_lo)) + 0.5;
153 amrex::ignore_unused(found);
154 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(found,
"MOSTAverage: Height above terrain not found, try increasing z_ref!");
157 const amrex::RealVect lx(ireal + 0.5, jreal + 0.5, zval);
159 const amrex::IntVect ijk = lx.floor();
161 int i = ijk[0];
int j = ijk[1];
int k = ijk[2];
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]));
169 const amrex::RealVect sx_hi = lx - ijk_r;
170 const amrex::RealVect sx_lo = 1.0 - sx_hi;
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);
188 const amrex::Vector<amrex::Geometry>
m_geom;
189 amrex::Vector<amrex::Vector<amrex::MultiFab*>>
m_fields;
202 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
m_x_pos;
203 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
m_y_pos;
204 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
m_z_pos;
205 amrex::Vector<std::unique_ptr<amrex::iMultiFab>>
m_i_indx;
206 amrex::Vector<std::unique_ptr<amrex::iMultiFab>>
m_j_indx;
207 amrex::Vector<std::unique_ptr<amrex::iMultiFab>>
m_k_indx;
208 amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>>
m_averages;
209 amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>>
m_rot_fields;
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