1 #ifndef ERF_HURRICANE_DIAGNOSTICS_H_
2 #define ERF_HURRICANE_DIAGNOSTICS_H_
5 #include <AMReX_MultiFab.H>
6 #include <AMReX_ParallelReduce.H>
16 #define M_PI 3.14159265358979323846
26 const amrex::Vector<amrex::MultiFab>& S_data,
27 const amrex::Vector<amrex::MultiFab>* forecast_state_at_lev,
34 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_xy,
35 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_latlon)
40 amrex::Gpu::copy(amrex::Gpu::deviceToHost, d_val_min_ptr, d_val_min_ptr + 1, &h_val_min);
41 amrex::Gpu::copy(amrex::Gpu::deviceToHost, d_i_min_ptr, d_i_min_ptr + 1, &h_i_min);
42 amrex::Gpu::copy(amrex::Gpu::deviceToHost, d_j_min_ptr, d_j_min_ptr + 1, &h_j_min);
43 amrex::Gpu::synchronize();
46 int local_i_min = h_i_min;
47 int local_j_min = h_j_min;
49 int rank = amrex::ParallelDescriptor::MyProc();
51 in.value = local_val_min;
55 MPI_Allreduce(&
in, &
out, 1, MPI_DOUBLE_INT, MPI_MINLOC, MPI_COMM_WORLD);
60 global_val_min =
out.value;
61 int owner_rank =
out.rank;
64 global_i_min = local_i_min;
65 global_j_min = local_j_min;
67 amrex::ParallelDescriptor::Bcast(&global_i_min, 1, owner_rank);
68 amrex::ParallelDescriptor::Bcast(&global_j_min, 1, owner_rank);
71 amrex::Print() <<
"Global minimum distance to hurricane eye (k=0): "
72 << global_val_min <<
" at (i,j) = ("
73 << global_i_min <<
", " << global_j_min <<
")\n";
76 amrex::Gpu::DeviceScalar<amrex::Real> d_eye_lat(0.0), d_eye_lon(0.0);
82 if (
rank == owner_rank) {
83 for (amrex::MFIter mfi(S_data[
IntVars::cons]); mfi.isValid(); ++mfi) {
84 const amrex::Box& box = mfi.validbox();
85 const auto latlon_arr = (*forecast_state_at_lev)[4].array(mfi);
86 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
87 if (i == global_i_min && j == global_j_min && k == 0) {
88 *d_eye_lat_ptr = latlon_arr(i,j,k,0);
89 *d_eye_lon_ptr = latlon_arr(i,j,k,1);
99 amrex::Gpu::synchronize();
101 amrex::ParallelDescriptor::Bcast(&eye_lat, 1, owner_rank);
102 amrex::ParallelDescriptor::Bcast(&eye_lon, 1, owner_rank);
104 const auto dx = geom.CellSizeArray();
105 const auto prob_lo = geom.ProbLoArray();
107 amrex::Real eye_x = prob_lo[0] + (global_i_min+0.5)*dx[0];
108 amrex::Real eye_y = prob_lo[1] + (global_j_min+0.5)*dx[1];
110 hurricane_eye_track_xy.push_back({eye_x, eye_y});
111 hurricane_eye_track_latlon.push_back({eye_lon, eye_lat});
118 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_tracker_circle)
121 if (hurricane_eye_track_xy.empty())
return;
124 const auto [x_last, y_last] = hurricane_eye_track_xy.back();
127 const int n_points = 100;
131 hurricane_tracker_circle.clear();
132 hurricane_tracker_circle.reserve(n_points);
135 for (
int i = 0; i < n_points; ++i) {
139 hurricane_tracker_circle.push_back({
x,
y});
146 const amrex::Vector<amrex::MultiFab>& S_data,
147 const amrex::Vector<amrex::MultiFab>* forecast_state_at_lev,
148 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_xy,
149 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_latlon,
153 amrex::Gpu::DeviceScalar<amrex::Real> d_val_min(1e10);
154 amrex::Gpu::DeviceScalar<int> d_i_min(-1), d_j_min(-1);
157 int* d_i_min_ptr = d_i_min.dataPtr();
158 int* d_j_min_ptr = d_j_min.dataPtr();
160 for (amrex::MFIter mfi(S_data[
IntVars::cons]); mfi.isValid(); ++mfi) {
161 const amrex::Box& box = mfi.validbox();
162 const auto latlon_arr = (*forecast_state_at_lev)[4].array(mfi);
164 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
167 amrex::Real dlat = latlon_arr(i,j,k,0) - hurricane_eye_latitude;
168 amrex::Real dlon = latlon_arr(i,j,k,1) - hurricane_eye_longitude;
169 amrex::Real dist = std::sqrt(dlat*dlat + dlon*dlon);
171 amrex::Real old = amrex::Gpu::Atomic::Min(&d_val_min_ptr[0], dist);
183 int global_i_min, global_j_min;
186 d_val_min_ptr, d_i_min_ptr, d_j_min_ptr,
187 global_val_min, global_i_min, global_j_min,
188 hurricane_eye_track_xy,
189 hurricane_eye_track_latlon);
195 const amrex::Vector<amrex::MultiFab>& S_data,
196 MoistureType moisture_type,
197 const amrex::Vector<amrex::MultiFab>* forecast_state_at_lev,
198 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_xy,
199 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_latlon)
202 if (hurricane_eye_track_xy.empty()) {
203 amrex::Print() <<
"Error: hurricane_eye_track_xy is empty!\n";
204 amrex::Abort(
"Attempted to access hurricane_eye_track_xy[0]");
207 amrex::Real tmp_x_eye = hurricane_eye_track_xy.back()[0];
208 amrex::Real tmp_y_eye = hurricane_eye_track_xy.back()[1];
210 if(amrex::ParallelDescriptor::IOProcessor()){
211 std::cout <<
"The value of x y are " << tmp_x_eye <<
" " << tmp_y_eye << std::endl;
214 amrex::Gpu::DeviceScalar<amrex::Real> d_val_min(1e10);
215 amrex::Gpu::DeviceScalar<int> d_i_min(-1), d_j_min(-1);
218 int* d_i_min_ptr = d_i_min.dataPtr();
219 int* d_j_min_ptr = d_j_min.dataPtr();
221 bool use_moisture = (moisture_type != MoistureType::None);
224 const auto dx = geom.CellSizeArray();
225 const auto prob_lo = geom.ProbLoArray();
227 for (amrex::MFIter mfi(S_data[
IntVars::cons]); mfi.isValid(); ++mfi) {
228 const amrex::Box& box = mfi.validbox();
229 const amrex::Array4<amrex::Real const>& S_arr = S_data[
IntVars::cons].const_array(mfi);
231 amrex::ParallelFor(box,[=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
235 amrex::Real dist = std::sqrt((
x-tmp_x_eye)*(
x-tmp_x_eye) + (
y-tmp_y_eye)*(
y-tmp_y_eye));
240 amrex::Real old = amrex::Gpu::Atomic::Min(&d_val_min_ptr[0], pressure);
242 if (old > pressure) {
253 int global_i_min, global_j_min;
256 d_val_min_ptr, d_i_min_ptr, d_j_min_ptr,
257 global_val_min, global_i_min, global_j_min,
258 hurricane_eye_track_xy,
259 hurricane_eye_track_latlon);
265 const amrex::Vector<amrex::MultiFab>& S_data,
266 MoistureType moisture_type,
267 const amrex::Vector<amrex::MultiFab>* forecast_state_at_lev,
270 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_xy,
271 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_latlon,
272 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_tracker_circle)
274 static bool is_start =
true;
277 forecast_state_at_lev, hurricane_eye_track_xy,
278 hurricane_eye_track_latlon,
279 hurricane_eye_latitude, hurricane_eye_longitude);
283 forecast_state_at_lev, hurricane_eye_track_xy,
284 hurricane_eye_track_latlon);
292 const amrex::MultiFab& mf_cc_vel,
294 const amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_xy,
295 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_maxvel_vs_time)
297 const int ncomp = AMREX_SPACEDIM;
300 amrex::Gpu::DeviceVector<amrex::Real> d_val_max(1, -1.0e30);
301 d_val_max_ptr = d_val_max.data();
303 const auto [x_last, y_last] = hurricane_eye_track_xy.back();
304 const auto dx = geom.CellSizeArray();
305 const auto prob_lo = geom.ProbLoArray();
310 for (amrex::MFIter mfi(mf_cc_vel); mfi.isValid(); ++mfi) {
311 const amrex::Box& box = mfi.validbox();
312 const auto& vel_arr = mf_cc_vel.const_array(mfi);
314 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
318 (
y-y_eye)*(
y-y_eye));
319 if(k==1 && dist < 200e3) {
321 for (
int comp = 0; comp < ncomp; ++comp) {
325 velmag = std::sqrt(velmag)*3.6;
326 amrex::Gpu::Atomic::Max(&d_val_max_ptr[0], velmag);
331 amrex::Gpu::synchronize();
334 amrex::Gpu::copy(amrex::Gpu::deviceToHost, d_val_max.begin(), d_val_max.end(), &h_val_max_local);
338 MPI_Allreduce(&h_val_max_local, &h_val_max_global, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
340 h_val_max_global = h_val_max_local;
344 hurricane_maxvel_vs_time.push_back({time_in_hrs, h_val_max_global});
350 const amrex::Geometry& geom,
351 const amrex::MultiFab& mf_cons_var,
353 const amrex::Vector<std::array<amrex::Real, 2>>& hurricane_eye_track_xy,
354 amrex::Vector<std::array<amrex::Real, 2>>& hurricane_minpressure_vs_time)
358 amrex::Gpu::DeviceVector<amrex::Real> d_val_min(1, 1.0e30);
359 d_val_min_ptr = d_val_min.data();
361 const amrex::Real x_last = hurricane_eye_track_xy.back()[0];
362 const amrex::Real y_last = hurricane_eye_track_xy.back()[1];
363 const auto dx = geom.CellSizeArray();
364 const auto prob_lo = geom.ProbLoArray();
366 const int ncomp = mf_cons_var.nComp();
367 bool use_moisture = (moisture_type != MoistureType::None);
369 for (amrex::MFIter mfi(mf_cons_var); mfi.isValid(); ++mfi) {
370 const amrex::Box& box = mfi.validbox();
371 const auto& S_arr = mf_cons_var.const_array(mfi);
373 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
377 (
y-y_last)*(
y-y_last);
378 if(k==1 && dist2 < 200e3*200e3) {
382 amrex::Gpu::Atomic::Min(&d_val_min_ptr[0], pressure);
387 amrex::Gpu::synchronize();
390 amrex::Gpu::copy(amrex::Gpu::deviceToHost, d_val_min.begin(), d_val_min.end(), &h_val_min_local);
394 MPI_Allreduce(&h_val_min_local, &h_val_min_global, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
396 h_val_min_global = h_val_min_local;
400 hurricane_minpressure_vs_time.push_back({time_in_hrs, h_val_min_global});
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
AMREX_FORCE_INLINE void ComputeGlobalMinLocation(const amrex::Geometry &geom, const amrex::Vector< amrex::MultiFab > &S_data, const amrex::Vector< amrex::MultiFab > *forecast_state_at_lev, amrex::Real *d_val_min_ptr, int *d_i_min_ptr, int *d_j_min_ptr, amrex::Real &global_val_min, int &global_i_min, int &global_j_min, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_xy, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_latlon)
Definition: ERF_HurricaneDiagnostics.H:25
AMREX_FORCE_INLINE void HurricaneMaxVelTracker(const amrex::Geometry &geom, const amrex::MultiFab &mf_cc_vel, const amrex::Real &time, const amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_xy, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_maxvel_vs_time)
Definition: ERF_HurricaneDiagnostics.H:291
int rank
Definition: ERF_HurricaneDiagnostics.H:21
amrex::Real value
Definition: ERF_HurricaneDiagnostics.H:20
AMREX_FORCE_INLINE void HurricaneEyeTrackerInitial(const amrex::Geometry &geom, const amrex::Vector< amrex::MultiFab > &S_data, const amrex::Vector< amrex::MultiFab > *forecast_state_at_lev, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_xy, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_latlon, const amrex::Real &hurricane_eye_latitude, const amrex::Real &hurricane_eye_longitude)
Definition: ERF_HurricaneDiagnostics.H:145
AMREX_FORCE_INLINE void HurricaneEyeTrackerNotInitial(const amrex::Geometry &geom, const amrex::Vector< amrex::MultiFab > &S_data, MoistureType moisture_type, const amrex::Vector< amrex::MultiFab > *forecast_state_at_lev, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_xy, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_latlon)
Definition: ERF_HurricaneDiagnostics.H:194
AMREX_FORCE_INLINE void HurricaneTrackerCircle(const amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_xy, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_tracker_circle)
Definition: ERF_HurricaneDiagnostics.H:117
AMREX_FORCE_INLINE void HurricaneEyeTracker(const amrex::Geometry &geom, const amrex::Vector< amrex::MultiFab > &S_data, MoistureType moisture_type, const amrex::Vector< amrex::MultiFab > *forecast_state_at_lev, const amrex::Real &hurricane_eye_latitude, const amrex::Real &hurricane_eye_longitude, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_xy, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_latlon, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_tracker_circle)
Definition: ERF_HurricaneDiagnostics.H:264
#define M_PI
Definition: ERF_HurricaneDiagnostics.H:16
AMREX_FORCE_INLINE void HurricaneMinPressureTracker(MoistureType moisture_type, const amrex::Geometry &geom, const amrex::MultiFab &mf_cons_var, const amrex::Real &time, const amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_eye_track_xy, amrex::Vector< std::array< amrex::Real, 2 >> &hurricane_minpressure_vs_time)
Definition: ERF_HurricaneDiagnostics.H:349
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ cons
Definition: ERF_IndexDefines.H:158
@ theta
Definition: ERF_MM5.H:20