/home/runner/work/amr-wind/amr-wind/amr-wind/ocean_waves/relaxation_zones/waves2amr_ops.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/ocean_waves/relaxation_zones/waves2amr_ops.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
waves2amr_ops.H
Go to the documentation of this file.
1#ifndef W2A_WAVES_OPS_H
2#define W2A_WAVES_OPS_H
3
13#include "AMReX_MultiFabUtil.H"
14#include "AMReX_REAL.H"
15
16using namespace amrex::literals;
17
18#ifdef AMR_WIND_USE_W2A
19namespace {
20int evaluate_read_resize(
21 int& ntime,
22 bool& read_flag,
23 bool& resize_flag,
24 amrex::Real& wtime,
25 amrex::Real& t_last,
26 const int new_ntime,
27 const amrex::Real wtinit,
28 const amrex::Real wdt,
29 const amrex::Real time)
30{
31 // Flag to indicate that data must be read twice for interpolation
32 int double_data = 0;
33 // Check if time indicates reading must take place
34 if (new_ntime != ntime) {
35 // New data is needed - reading should happen
36 read_flag = true;
37 // If time index has proceeded more than 1 step
38 if (new_ntime > ntime + 1) {
39 // Double reading is necessary
40 double_data = 1;
41 }
42 // Time index for reading
43 ntime = new_ntime;
44 // Sim time to go with recorded data
45 wtime = new_ntime * wdt - wtinit;
46 // If double reading is deemed necessary, check for convenience
47 if (double_data == 1 && std::abs(wtime - time) <= 1.0e-10_rt) {
48 // Reading can be done just once, w2a fields replace ow fields
49 double_data = 2;
50 }
51 }
52 // Check if reading must take place for other reasons
53 if (t_last < -1.0e-10_rt) {
54 // Signifies initialization from scratch without waves or a restart
55 read_flag = true;
56 // Resizing needs to happen for the first time
57 resize_flag = true;
58
59 // Confirm that new time is not coincident with modes time step
60 if (std::abs(wtime - time) > 1.0e-10_rt) {
61 // Data must be read before and after wtime
62 double_data = 1;
63 } else {
64 // Data need only be read once, w2a fields will replace ow fields
65 double_data = 2;
66 }
67
68 } else if (std::abs(t_last) < 1.0e-10_rt) {
69 // Signifies initialization with waves
70 read_flag = true;
71 // Resizing needs to happen for the first time
72 resize_flag = true;
73
74 // levelset and velocity fields are up-to-date at t=0
75 // interpolation is ready to go
76 }
77 // Record latest time as 'last' for next timestep
78 t_last = time;
79 // Return flag regarding double reading
80 return double_data;
81}
82
83int update_offset_timestep(const int ntime, const int n0)
84{
85 // Offending timestep (goes too far): ntime + offset
86 // Farthest back timesteps are permitted to go: n0
87 // Subtract offset by offending timestep, add back lower limit
88 // new offset = offset - (ntime + offset) + n0
89 return (-ntime + n0);
90}
91
92void populate_fields_all_levels(
94 amrex::Vector<amrex::Geometry>& geom_all,
95 amr_wind::Field& lvs_field,
96 amr_wind::Field& vel_field,
97 int ntime_off = 0)
98{
99
100 // Get data from modes
101 bool no_EOF =
102 wdata.is_ocean
103 ? wdata.c_rmodes.get_data(
104 wdata.ntime + wdata.n_offset + ntime_off, wdata.c_mX,
105 wdata.c_mY, wdata.c_mZ, wdata.c_mFS)
106 : wdata.r_rmodes.get_data(
107 wdata.ntime + wdata.n_offset + ntime_off, wdata.r_mX,
108 wdata.r_mY, wdata.r_mZ, wdata.r_mFS, wdata.r_mAdd);
109 // Navigate when end of file is reached
110 if (!no_EOF) {
111 // End of file detected, reset reading
112 wdata.n_offset = update_offset_timestep(wdata.ntime, wdata.n_winit);
113 // Print warning to screen
114 amrex::Print() << "WARNING (waves2amr_ops): end of mode data file "
115 "detected, resetting to beginning of mode data.\n";
116 // Read data again, now from a valid timestep
117 no_EOF = wdata.is_ocean
118 ? wdata.c_rmodes.get_data(
119 wdata.ntime + wdata.n_offset + ntime_off, wdata.c_mX,
120 wdata.c_mY, wdata.c_mZ, wdata.c_mFS)
121 : wdata.r_rmodes.get_data(
122 wdata.ntime + wdata.n_offset + ntime_off, wdata.r_mX,
123 wdata.r_mY, wdata.r_mZ, wdata.r_mFS, wdata.r_mAdd);
124 // If no valid data is detected at this point, abort
125 if (!no_EOF) {
126 amrex::Abort(
127 "waves2amr_ops: end of mode data file detected after "
128 "resetting to beginning; please evaluate HOS_init_time "
129 "or HOS_init_timestep and check the length of the mode "
130 "file.");
131 }
132 }
133
134 // Convert to spatial data in vectors
135 if (wdata.is_ocean) {
136 modes_hosgrid::copy_complex(
137 wdata.n0, wdata.n1, wdata.c_mFS, wdata.c_eta_mptr);
138 modes_hosgrid::populate_hos_eta(
139 wdata.c_rmodes, wdata.plan_vector, wdata.c_eta_mptr,
140 wdata.sp_eta_vec);
141 } else {
142 modes_hosgrid::copy_real(
143 wdata.n0, wdata.n1, wdata.r_mFS, wdata.r_eta_mptr);
144 modes_hosgrid::populate_hos_eta(
145 wdata.r_rmodes, wdata.plan_vector, wdata.r_eta_mptr,
146 wdata.sp_eta_vec);
147 }
148
149 for (int iht = 0; iht < wdata.indvec.size(); ++iht) {
150 // Get sample height
151 amrex::Real ht = wdata.hvec[wdata.indvec[iht]];
152 // Sample velocity
153 if (wdata.is_ocean) {
154 modes_hosgrid::populate_hos_vel(
155 wdata.c_rmodes, ht, wdata.zsl, wdata.c_mX, wdata.c_mY,
156 wdata.c_mZ, wdata.plan_vector, wdata.c_u_mptr, wdata.c_v_mptr,
157 wdata.c_w_mptr, wdata.sp_u_vec, wdata.sp_v_vec, wdata.sp_w_vec,
158 iht * wdata.n0 * wdata.n1);
159 } else {
160 modes_hosgrid::populate_hos_vel(
161 wdata.r_rmodes, ht, wdata.zsl, wdata.r_mX, wdata.r_mY,
162 wdata.r_mZ, wdata.r_mAdd, wdata.plan_vector, wdata.r_u_mptr,
163 wdata.r_v_mptr, wdata.r_w_mptr, wdata.au_mptr, wdata.av_mptr,
164 wdata.aw_mptr, wdata.sp_u_vec, wdata.sp_v_vec, wdata.sp_w_vec,
165 iht * wdata.n0 * wdata.n1);
166 }
167 }
168
169 // Interpolate to fields (vector of MultiFabs)
170 interp_to_mfab::interp_eta_to_levelset_field(
171 wdata.n0_sp, wdata.n1_sp, wdata.dx0, wdata.dx1, wdata.xlo0, wdata.xlo1,
172 wdata.zsl, wdata.is_ocean, wdata.sp_eta_vec, lvs_field.vec_ptrs(),
173 geom_all);
174 interp_to_mfab::interp_velocity_to_field(
175 wdata.n0_sp, wdata.n1_sp, wdata.dx0, wdata.dx1, wdata.xlo0, wdata.xlo1,
176 wdata.is_ocean, wdata.indvec, wdata.hvec, wdata.sp_u_vec,
177 wdata.sp_v_vec, wdata.sp_w_vec, vel_field.vec_ptrs(), geom_all);
178}
179
180void time_interpolate_wave_fields(
181 amr_wind::Field& lvs_field,
182 amr_wind::Field& vel_field,
183 const amr_wind::Field& future_lvs_field,
184 const amr_wind::Field& future_vel_field,
185 const amrex::Real t_old,
186 const amrex::Real t_new,
187 const amrex::Real t_future)
188{
189 const amrex::Real b =
190 (t_new - t_old) /
191 (t_future - t_old + std::numeric_limits<amrex::Real>::epsilon());
192 const amrex::Real a = 1.0_rt - b;
194 lvs_field, a, lvs_field, 0, b, future_lvs_field, 0, 0, 1,
195 lvs_field.num_grow());
197 vel_field, a, vel_field, 0, b, future_vel_field, 0, 0, AMREX_SPACEDIM,
198 vel_field.num_grow());
199}
200
201} // namespace
202#endif
203
205
206template <>
208{
210 // cppcheck-suppress constParameterReference
211 W2AWaves::DataType& data,
212 const ::amr_wind::utils::MultiParser& pp)
213 {
214// Check for W2A initialization
215#ifndef AMR_WIND_USE_W2A
216 // Assert Waves2AMR must be used for initial condition file
217 amrex::Abort(
218 "ocean_waves/W2AWaves: AMR-Wind was not built with Waves2AMR "
219 "support; associated wave data cannot be processed for relaxation "
220 "zones.");
221
222 amrex::ignore_unused(data, pp);
223#else
224 auto& wdata = data.meta();
225 auto& info = data.info();
226 relaxation_zones::read_inputs(wdata, info, pp);
227
228 if (wdata.current > constants::TIGHT_TOL) {
229 amrex::Abort(
230 "Current is specified as nonzero, but current is not yet "
231 "implemented for W2A Waves.");
232 }
233
234 pp.get("HOS_modes_filename", wdata.modes_file);
235 pp.query("HOS_simulation_is_ocean", wdata.is_ocean);
236 pp.query("HOS_init_timestep", wdata.ntime);
237 if (!pp.contains("HOS_init_timestep")) {
238 pp.query("HOS_init_time", wdata.t_winit);
239 }
240 pp.query("HOS_domain_offset_x", wdata.xlo0);
241 pp.query("HOS_domain_offset_y", wdata.xlo1);
242
243 // Default fftw_plan is deterministic
244 std::string fftw_planner_flag{"estimate"};
245 pp.query("fftw_planner_flag", fftw_planner_flag);
246
247 amrex::Vector<amrex::Real> prob_lo_input(AMREX_SPACEDIM);
248 amrex::ParmParse pp_geom("geometry");
249 pp_geom.getarr("prob_lo", prob_lo_input);
250
251 // Read user inputs to height vector
252 int nheights = 0;
253 int nh_above = 1;
254 amrex::Real dz0 = 0.0_rt;
255 pp.get("number_interp_points_in_z", nheights);
256 pp.get("interp_spacing_at_surface", dz0);
257 pp.query("number_interp_above_surface", nh_above);
258
259 // Initialize mode reader
260 bool file_exists{false};
261 if (wdata.is_ocean) {
262 file_exists = wdata.c_rmodes.initialize(wdata.modes_file, true);
263 } else {
264 file_exists = wdata.r_rmodes.initialize(wdata.modes_file, false);
265 }
266
267 // Abort if file does not exist
268 if (!file_exists) {
269 amrex::Abort(
270 "Waves2AMR ReadInputsOp: modes file requested does not exist");
271 }
272
273 // Get dt of HOS data
274 wdata.dt_modes = wdata.is_ocean ? wdata.c_rmodes.get_dtout()
275 : wdata.r_rmodes.get_dtout();
276
277 // Get initial time and timestep synced
278 if (wdata.t_winit > 0.0_rt) {
279 // If initial time was specified
280 // Get time index near requested time
281 wdata.ntime =
282 wdata.is_ocean
283 ? wdata.c_rmodes.time2step(wdata.t_winit, wdata.ntime)
284 : wdata.r_rmodes.time2step(wdata.t_winit, wdata.ntime);
285 ;
286 // Sync time to time index
287 wdata.t_winit = wdata.dt_modes * wdata.ntime;
288 // Save first timestep
289 wdata.n_winit = wdata.ntime;
290 } else {
291 // If initial timestep is being used
292 wdata.t_winit = wdata.dt_modes * wdata.ntime;
293 // Save first timestep
294 wdata.n_winit = wdata.ntime;
295 }
296
297 amrex::Print() << "OceanWaves: Initializing Waves2AMR profile at time "
298 << wdata.t_winit << ", step " << wdata.n_winit
299 << std::endl;
300
301 // Initialize variables to store modes
302 amrex::Real depth{0.0_rt};
303 if (wdata.is_ocean) {
304 int vsize = wdata.c_rmodes.get_vector_size();
305 amrex::Real initval = 0.0_rt;
306 wdata.c_mX.resize(vsize, initval);
307 wdata.c_mY.resize(vsize, initval);
308 wdata.c_mZ.resize(vsize, initval);
309 wdata.c_mFS.resize(vsize, initval);
310
311 // Get dimensions of data
312 wdata.n0 = wdata.c_rmodes.get_first_fft_dimension();
313 wdata.n1 = wdata.c_rmodes.get_second_fft_dimension();
314 wdata.n0_sp = wdata.c_rmodes.get_first_spatial_dimension();
315 wdata.n1_sp = wdata.c_rmodes.get_second_spatial_dimension();
316 // Get resolution
317 wdata.dx0 = wdata.c_rmodes.get_xlen() / wdata.n0_sp;
318 wdata.dx1 = wdata.c_rmodes.get_ylen() / wdata.n1_sp;
319 // Get depth
320 depth = wdata.c_rmodes.get_depth();
321 // Get dimensional length
322 wdata.dimL = wdata.c_rmodes.get_L();
323 // Get nominal last timestep of data
324 wdata.n_wstop = static_cast<int>(
325 (wdata.c_rmodes.get_Tstop() + 1.0e-8_rt) / wdata.dt_modes);
326
327 amrex::Print()
328 << "OceanWaves: Waves2AMR parameters from HOS-Ocean file: \n"
329 << " nx " << wdata.n0_sp << " dx " << wdata.dx0 << " Lx "
330 << wdata.c_rmodes.get_xlen() << std::endl
331 << " ny " << wdata.n1_sp << " dy " << wdata.dx1 << " Ly "
332 << wdata.c_rmodes.get_ylen() << std::endl
333 << " depth " << depth << std::endl
334 << " output time interval " << wdata.dt_modes
335 << " nominal last time " << wdata.n_wstop * wdata.dt_modes
336 << std::endl;
337 } else {
338 int vsize = wdata.r_rmodes.get_vector_size();
339 int vasize = wdata.r_rmodes.get_addl_vector_size();
340 amrex::Real initval = 0.0_rt;
341 wdata.r_mX.resize(vsize, initval);
342 wdata.r_mY.resize(vsize, initval);
343 wdata.r_mZ.resize(vsize, initval);
344 wdata.r_mFS.resize(vsize, initval);
345 wdata.r_mAdd.resize(vasize, initval);
346
347 // Get dimensions of data
348 wdata.n0 = wdata.r_rmodes.get_first_fft_dimension();
349 wdata.n1 = wdata.r_rmodes.get_second_fft_dimension();
350 wdata.n2 = wdata.r_rmodes.get_third_dimension();
351 wdata.n0_sp = wdata.r_rmodes.get_first_spatial_dimension();
352 wdata.n1_sp = wdata.r_rmodes.get_second_spatial_dimension();
353 // Get resolution (slightly different for NWT)
354 wdata.dx0 = wdata.r_rmodes.get_xlen() / (wdata.n0_sp - 1);
355 wdata.dx1 =
356 wdata.r_rmodes.get_ylen() / amrex::max(wdata.n1_sp - 1, 1);
357 // Get depth
358 depth = wdata.r_rmodes.get_depth();
359 // Get dimensional length
360 wdata.dimL = wdata.r_rmodes.get_L();
361 // Get nominal last timestep of data
362 wdata.n_wstop = static_cast<int>(
363 (wdata.r_rmodes.get_Tstop() + 1.0e-8_rt) / wdata.dt_modes);
364
365 amrex::Print()
366 << "OceanWaves: Waves2AMR parameters from HOS-NWT file: \n"
367 << " nx " << wdata.n0_sp << " dx " << wdata.dx0 << " Lx "
368 << wdata.r_rmodes.get_xlen() << std::endl
369 << " ny " << wdata.n1_sp << " dy " << wdata.dx1 << " Ly "
370 << wdata.r_rmodes.get_ylen() << std::endl
371 << " nz " << wdata.n2 << " depth " << depth << std::endl
372 << " output time interval " << wdata.dt_modes
373 << " nominal last time " << wdata.n_wstop * wdata.dt_modes
374 << std::endl;
375 }
376
377 // Check if stop time is exceeded, introduce offset to ntime
378 if (wdata.ntime + wdata.n_offset > wdata.n_wstop) {
379 // If exceeding stop step, calculate new offset
380 wdata.n_offset = update_offset_timestep(wdata.ntime, wdata.n_winit);
381 // Print warning to screen
382 amrex::Print()
383 << "WARNING (waves2amr_ops): available mode data exceeded, "
384 "resetting to beginning of mode data.\n";
385 }
386
387 // Warning if depth does not correspond to simulation
388 if (std::abs(depth - (wdata.zsl - prob_lo_input[2])) >
389 1.0e-3_rt * depth) {
390 amrex::Print()
391 << "WARNING: Mismatch between water depths from AMR-Wind "
392 "domain and HOS data interpreted by Waves2AMR. \n ^This "
393 "warning is not a concern when using waves as terrain.\n";
394 }
395
396 // Allocate pointers for FFTW
397 if (wdata.is_ocean) {
398 wdata.c_eta_mptr =
399 modes_hosgrid::allocate_complex(wdata.n0, wdata.n1);
400 wdata.c_u_mptr =
401 modes_hosgrid::allocate_complex(wdata.n0, wdata.n1);
402 wdata.c_v_mptr =
403 modes_hosgrid::allocate_complex(wdata.n0, wdata.n1);
404 wdata.c_w_mptr =
405 modes_hosgrid::allocate_complex(wdata.n0, wdata.n1);
406 } else {
407 wdata.r_eta_mptr = modes_hosgrid::allocate_real(wdata.n0, wdata.n1);
408 wdata.r_u_mptr = modes_hosgrid::allocate_real(wdata.n0, wdata.n1);
409 wdata.r_v_mptr = modes_hosgrid::allocate_real(wdata.n0, wdata.n1);
410 wdata.r_w_mptr = modes_hosgrid::allocate_real(wdata.n0, wdata.n1);
411 wdata.au_mptr = modes_hosgrid::allocate_real(wdata.n0, 1);
412 wdata.av_mptr = modes_hosgrid::allocate_real(wdata.n0, 1);
413 wdata.aw_mptr = modes_hosgrid::allocate_real(wdata.n0, 1);
414 }
415 wdata.plan_out = modes_hosgrid::allocate_real(wdata.n0, wdata.n1);
416
417 // Set up planner flag based on input
418 auto plan_f = modes_hosgrid::planner_flags::estimate;
419 if (fftw_planner_flag == "patient") {
420 plan_f = modes_hosgrid::planner_flags::patient;
421 } else if (fftw_planner_flag == "exhaustive") {
422 plan_f = modes_hosgrid::planner_flags::exhaustive;
423 } else if (fftw_planner_flag == "measure") {
424 plan_f = modes_hosgrid::planner_flags::measure;
425 } else if (!(fftw_planner_flag == "estimate")) {
426 amrex::Print()
427 << "WARNING (waves2amr_ops): invalid fftw_planner_flag "
428 "specified; defaulting to estimate (FFTW_ESTIMATE).\n";
429 }
430 // Set up plan for FFTW
431 if (wdata.is_ocean) {
432 wdata.plan_vector.emplace_back(
433 modes_hosgrid::plan_ifftw(
434 wdata.n0, wdata.n1, wdata.c_eta_mptr, wdata.plan_out,
435 plan_f));
436 } else {
437 modes_hosgrid::plan_ifftw_nwt(
438 wdata.n0, wdata.n1, wdata.plan_vector, wdata.r_eta_mptr,
439 wdata.plan_out, plan_f);
440 }
441
442 // Create height vector for velocity mode conversion before
443 // interpolation, with prob_lo as bottom
444 int flag = interp_to_mfab::create_height_vector(
445 wdata.hvec, nheights, dz0, wdata.zsl, prob_lo_input[2], nh_above);
446 // Fail if flag indicates it should
447 if (flag > 0) {
448 amrex::Abort(
449 "Waves2AMR ReadInputsOp: create_height_vector error, failure "
450 "code " +
451 std::to_string(flag));
452 }
453
454 amrex::Real z_double{0.0_rt};
455 bool found_z_double{false};
456 for (int n = 0; n < wdata.hvec.size() - 1; ++n) {
457 const amrex::Real dz = wdata.hvec[n] - wdata.hvec[n + 1];
458 if (!found_z_double && dz >= 2.0_rt * dz0) {
459 z_double = 0.5_rt * (wdata.hvec[n] + wdata.hvec[n + 1]);
460 found_z_double = true;
461 }
462 }
463
464 amrex::Print()
465 << "OceanWaves: Waves2AMR height vector for interpolation: \n"
466 << " start " << wdata.hvec[0] << " end "
467 << wdata.hvec[wdata.hvec.size() - 1] << " num points "
468 << wdata.hvec.size() << std::endl
469 << " initial dz " << wdata.hvec[0] - wdata.hvec[1]
470 << " first two grown dz "
471 << wdata.hvec[nh_above] - wdata.hvec[nh_above + 1] << " "
472 << wdata.hvec[nh_above + 1] - wdata.hvec[nh_above + 2]
473 << " final dz "
474 << wdata.hvec[wdata.hvec.size() - 2] -
475 wdata.hvec[wdata.hvec.size() - 1]
476 << std::endl;
477
478 if (found_z_double) {
479 amrex::Print()
480 << " interp spacing reaches double that of surface at "
481 << z_double << std::endl;
482 } else {
483 amrex::Print() << " interp spacing remains less than double "
484 "that of surface "
485 << std::endl;
486 }
487
488 // Always get modes on every processor for initialization
489 bool no_EOF =
490 wdata.is_ocean
491 ? wdata.c_rmodes.get_data(
492 wdata.ntime + wdata.n_offset, wdata.c_mX, wdata.c_mY,
493 wdata.c_mZ, wdata.c_mFS)
494 : wdata.r_rmodes.get_data(
495 wdata.ntime + wdata.n_offset, wdata.r_mX, wdata.r_mY,
496 wdata.r_mZ, wdata.r_mFS, wdata.r_mAdd);
497 if (!no_EOF) {
498 // End of file detected, reset reading
499 wdata.n_offset = update_offset_timestep(wdata.ntime, wdata.n_winit);
500 // Print warning to screen
501 amrex::Print()
502 << "WARNING (waves2amr_ops): end of mode data file "
503 "detected, resetting to beginning of mode data.\n";
504 // Read data again, now from a valid timestep
505 no_EOF =
506 wdata.is_ocean
507 ? wdata.c_rmodes.get_data(
508 wdata.ntime + wdata.n_offset, wdata.c_mX, wdata.c_mY,
509 wdata.c_mZ, wdata.c_mFS)
510 : wdata.r_rmodes.get_data(
511 wdata.ntime + wdata.n_offset, wdata.r_mX, wdata.r_mY,
512 wdata.r_mZ, wdata.r_mFS, wdata.r_mAdd);
513 // If no valid data is detected at this point, abort
514 if (!no_EOF) {
515 amrex::Abort(
516 "waves2amr_ops: end of mode data file detected after "
517 "resetting to beginning; please evaluate HOS_init_time "
518 "or HOS_init_timestep and check the length of the mode "
519 "file.");
520 }
521 }
522
523 // Convert modes to spatial data
524 if (wdata.is_ocean) {
525 modes_hosgrid::copy_complex(
526 wdata.n0, wdata.n1, wdata.c_mFS, wdata.c_eta_mptr);
527 wdata.sp_eta_vec.resize(
528 static_cast<size_t>(wdata.n0) * static_cast<size_t>(wdata.n1),
529 0.0_rt);
530 modes_hosgrid::populate_hos_eta(
531 wdata.c_rmodes, wdata.plan_vector, wdata.c_eta_mptr,
532 wdata.sp_eta_vec);
533 // Mesh is not yet created, so get data at every height
534 const auto n_hts = wdata.hvec.size();
535 wdata.sp_u_vec.resize(
536 static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
537 wdata.sp_v_vec.resize(
538 static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
539 wdata.sp_w_vec.resize(
540 static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
541 for (int iht = 0; iht < static_cast<int>(n_hts); ++iht) {
542 // Get sample height
543 amrex::Real ht = wdata.hvec[iht];
544 // Sample velocity
545 modes_hosgrid::populate_hos_vel(
546 wdata.c_rmodes, ht, wdata.zsl, wdata.c_mX, wdata.c_mY,
547 wdata.c_mZ, wdata.plan_vector, wdata.c_u_mptr,
548 wdata.c_v_mptr, wdata.c_w_mptr, wdata.sp_u_vec,
549 wdata.sp_v_vec, wdata.sp_w_vec, iht * wdata.n0 * wdata.n1);
550 }
551 } else {
552 modes_hosgrid::copy_real(
553 wdata.n0, wdata.n1, wdata.r_mFS, wdata.r_eta_mptr);
554 wdata.sp_eta_vec.resize(
555 static_cast<size_t>(wdata.n0) * static_cast<size_t>(wdata.n1),
556 0.0_rt);
557 modes_hosgrid::populate_hos_eta(
558 wdata.r_rmodes, wdata.plan_vector, wdata.r_eta_mptr,
559 wdata.sp_eta_vec);
560 // Mesh is not yet created, so get data at every height
561 const auto n_hts = wdata.hvec.size();
562 wdata.sp_u_vec.resize(
563 static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
564 wdata.sp_v_vec.resize(
565 static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
566 wdata.sp_w_vec.resize(
567 static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
568 for (int iht = 0; iht < static_cast<int>(n_hts); ++iht) {
569 // Get sample height
570 amrex::Real ht = wdata.hvec[iht];
571 // Sample velocity
572 modes_hosgrid::populate_hos_vel(
573 wdata.r_rmodes, ht, wdata.zsl, wdata.r_mX, wdata.r_mY,
574 wdata.r_mZ, wdata.r_mAdd, wdata.plan_vector, wdata.r_u_mptr,
575 wdata.r_v_mptr, wdata.r_w_mptr, wdata.au_mptr,
576 wdata.av_mptr, wdata.aw_mptr, wdata.sp_u_vec,
577 wdata.sp_v_vec, wdata.sp_w_vec, iht * wdata.n0 * wdata.n1);
578 }
579 }
580
581 // Declare fields for HOS
582 auto& w2a_levelset =
583 data.sim().repo().declare_field("w2a_levelset", 1, 3, 2);
584 auto& w2a_velocity = data.sim().repo().declare_field(
585 "w2a_velocity", AMREX_SPACEDIM, 3, 2);
586
587 // Extrapolation can work well when finer data is available
588 w2a_levelset.set_default_fillpatch_bc(data.sim().time());
589 w2a_velocity.set_default_fillpatch_bc(data.sim().time());
590#endif
591 }
592}; // namespace ops
593
594template <>
596{
597 void
598 // cppcheck-suppress constParameterReference
600 W2AWaves::DataType& data,
601 int level,
602 const amrex::Geometry& geom,
603 bool multiphase_mode)
604 {
605
606#ifdef AMR_WIND_USE_W2A
607 auto& wdata = data.meta();
608
609 auto& sim = data.sim();
610
611 // Fill ow fields, then populate flow fields according to setup
612 auto& ow_levelset = sim.repo().get_field("ow_levelset");
613 auto& ow_velocity = sim.repo().get_field("ow_velocity");
614
615 auto& velocity = sim.repo().get_field("velocity");
616 // Set w2a fields to default values to prep for updates
617 auto& w2a_levelset = sim.repo().get_field("w2a_levelset");
618 auto& w2a_velocity = sim.repo().get_field("w2a_velocity");
619 // Copy current fields to old w2a fields for time interp
620 auto& w2a_lvs_old = w2a_levelset.state(amr_wind::FieldState::Old);
621 auto& w2a_vel_old = w2a_velocity.state(amr_wind::FieldState::Old);
622
623 Field* levelset{nullptr};
624 if (multiphase_mode) {
625 levelset = &sim.repo().get_field("levelset");
626 }
627
628 const auto& problo = geom.ProbLoArray();
629 const auto& probhi = geom.ProbHiArray();
630 const auto& dx = geom.CellSizeArray();
631
632 // Set t_last to 0.0_rt to signify information read in
633 wdata.t_last = 0.0_rt;
634
635 // indvec is complete upon initialization (all heights every proc)
636 amrex::Vector<int> indvec(wdata.hvec.size());
637 std::iota(indvec.begin(), indvec.end(), 0);
638 // Interpolate to MultiFabs (one level at a time)
639 interp_to_mfab::interp_eta_to_levelset_multifab(
640 wdata.n0_sp, wdata.n1_sp, wdata.dx0, wdata.dx1, wdata.xlo0,
641 wdata.xlo1, wdata.zsl, wdata.is_ocean, wdata.sp_eta_vec,
642 ow_levelset(level), problo, dx);
643 interp_to_mfab::interp_velocity_to_multifab(
644 wdata.n0_sp, wdata.n1_sp, wdata.dx0, wdata.dx1, wdata.xlo0,
645 wdata.xlo1, wdata.is_ocean, indvec, wdata.hvec, wdata.sp_u_vec,
646 wdata.sp_v_vec, wdata.sp_w_vec, ow_velocity(level), problo, dx);
647
648 // Populate flow fields according to intended forcing and init setup
649 const auto& ow_phi = ow_levelset(level).const_arrays();
650 const auto& ow_vel = ow_velocity(level).const_arrays();
651 const auto& vel = velocity(level).arrays();
652 const auto& phi_arrs = multiphase_mode
653 ? (*levelset)(level).arrays()
654 : amrex::MultiArray4<amrex::Real>();
655
656 const auto& w2a_phi = w2a_levelset(level).arrays();
657 const auto& w2a_vel = w2a_velocity(level).arrays();
658 const auto& w2a_phi_o = w2a_lvs_old(level).arrays();
659 const auto& w2a_vel_o = w2a_vel_old(level).arrays();
660
661 const amrex::Real gen_length = wdata.gen_length;
662 const amrex::Real beach_length = wdata.beach_length;
663 const amrex::Real y_length = wdata.zone_length_y;
664 const amrex::Real zero_sea_level = wdata.zsl;
665
666 const bool has_beach = wdata.has_beach && multiphase_mode;
667 const bool init_wave_field = wdata.init_wave_field || !multiphase_mode;
668
669 amrex::MultiFab phi_tmp_fab(
670 velocity(level).boxArray(), velocity(level).DistributionMap(), 1,
671 3);
672 amrex::MultiFab vel_tmp_fab(
673 velocity(level).boxArray(), velocity(level).DistributionMap(),
674 AMREX_SPACEDIM, 3);
675
676 const auto& phi_tmp = phi_tmp_fab.arrays();
677 const auto& vel_tmp = vel_tmp_fab.arrays();
678
679 amrex::ParallelFor(
680 velocity(level), amrex::IntVect(3),
681 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
682 const amrex::Real x = problo[0] + (i + 0.5_rt) * dx[0];
683 const amrex::Real y = problo[1] + (j + 0.5_rt) * dx[1];
684 const amrex::Real z = problo[2] + (k + 0.5_rt) * dx[2];
685
686 // Wave profile
687 const utils::WaveVec wave_sol{
688 ow_vel[nbx](i, j, k, 0), ow_vel[nbx](i, j, k, 1),
689 ow_vel[nbx](i, j, k, 2), ow_phi[nbx](i, j, k) + z};
690 // Quiescent profile
691 const utils::WaveVec quiescent{
692 0.0_rt, 0.0_rt, 0.0_rt, zero_sea_level};
693
694 // Specify initial state for each region of domain
695 const auto bulk = init_wave_field ? wave_sol : quiescent;
696 const auto outlet = has_beach ? quiescent : wave_sol;
697
698 const auto local_profile = utils::harmonize_profiles_2d(
699 x, y, problo[0], gen_length, probhi[0], beach_length,
700 problo[1], probhi[1], y_length, wave_sol, bulk, outlet);
701
702 phi_tmp[nbx](i, j, k) = local_profile[3] - z;
703 vel_tmp[nbx](i, j, k, 0) = local_profile[0];
704 vel_tmp[nbx](i, j, k, 1) = local_profile[1];
705 vel_tmp[nbx](i, j, k, 2) = local_profile[2];
706
707 if (multiphase_mode) {
708 phi_arrs[nbx](i, j, k) = phi_tmp[nbx](i, j, k);
709 }
710
711 // Start old values at w2a solution (for the sake of time
712 // interpolation)
713 w2a_phi_o[nbx](i, j, k) = wave_sol[3] - z;
714 w2a_vel_o[nbx](i, j, k, 0) = wave_sol[0];
715 w2a_vel_o[nbx](i, j, k, 1) = wave_sol[1];
716 w2a_vel_o[nbx](i, j, k, 2) = wave_sol[2];
717
718 // Default w2a values matter for where no updates happen
719 w2a_phi[nbx](i, j, k) = quiescent[3] - z;
720 w2a_vel[nbx](i, j, k, 0) = quiescent[0];
721 w2a_vel[nbx](i, j, k, 1) = quiescent[1];
722 w2a_vel[nbx](i, j, k, 2) = quiescent[2];
723 });
724 amrex::Gpu::streamSynchronize();
725
726 amrex::ParallelFor(
727 velocity(level), amrex::IntVect(2),
728 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
729 const amrex::Real eps =
730 2.0_rt * std::cbrt(dx[0] * dx[1] * dx[2]);
731 const amrex::Real vof_local =
732 multiphase::levelset_to_vof(i, j, k, eps, phi_tmp[nbx]);
733
734 if (vof_local > 1.0_rt - constants::TIGHT_TOL) {
735 // Entire cell is liquid
736 vel[nbx](i, j, k, 0) = vel_tmp[nbx](i, j, k, 0);
737 vel[nbx](i, j, k, 1) = vel_tmp[nbx](i, j, k, 1);
738 vel[nbx](i, j, k, 2) = vel_tmp[nbx](i, j, k, 2);
739 } else if (vof_local > 1.0e-3_rt) {
740 // Velocity of liquid becomes velocity of entire cell
741 // Interpolate using cell-centered values to liquid centroid
742 // Consider only z-direction to find centroid
743 const amrex::Real z_cell = problo[2] + (k + 0.5_rt) * dx[2];
744 const amrex::Real z_cell_below = z_cell - dx[2];
745 const amrex::Real z_cell_bottom = z_cell - 0.5_rt * dx[2];
746 const amrex::Real z_liq_center =
747 z_cell_bottom + 0.5_rt * vof_local * dx[2];
748 for (int n = 0; n < AMREX_SPACEDIM; ++n) {
749 vel[nbx](i, j, k, n) =
750 vel_tmp[nbx](i, j, k - 1, n) +
751 (vel_tmp[nbx](i, j, k, n) -
752 vel_tmp[nbx](i, j, k - 1, n)) *
753 ((z_liq_center - z_cell_below) /
754 (z_cell - z_cell_below + constants::EPS));
755 }
756 }
757 });
758 amrex::Gpu::streamSynchronize();
759#else
760 amrex::ignore_unused(data, level, geom, multiphase_mode);
761#endif
762 }
763}; // namespace ocean_waves
764
765template <>
767{
768 // cppcheck-suppress constParameterReference
769 void operator()(W2AWaves::DataType& data, const amrex::Real time)
770 {
771
772#ifdef AMR_WIND_USE_W2A
773 auto& wdata = data.meta();
774 auto& sim = data.sim();
775
776 // Update ow fields every time
777 auto& ow_levelset = sim.repo().get_field("ow_levelset");
778 auto& ow_velocity = sim.repo().get_field("ow_velocity");
779 // Update HOS fields when necessary
780 auto& w2a_levelset = sim.repo().get_field("w2a_levelset");
781 auto& w2a_velocity = sim.repo().get_field("w2a_velocity");
782 // Old state of HOS fields used for time interpolation
783 auto& w2a_lvs_old = w2a_levelset.state(amr_wind::FieldState::Old);
784 auto& w2a_vel_old = w2a_velocity.state(amr_wind::FieldState::Old);
785
786 // Proxy for multiphase mode of ocean waves
787 bool vof_exists = sim.repo().field_exists("vof");
788
789 auto nlevels = sim.repo().num_active_levels();
790 auto geom = sim.mesh().Geom();
791
792 // Get value for time interpolation
793 amrex::Real t_last = wdata.t_last;
794
795 // Check if new HOS data needs to be read
796 bool read_flag = false;
797 // Check if time indicates reading must take place
798 int new_ntime =
799 wdata.is_ocean
800 ? wdata.c_rmodes.time2step(time + wdata.t_winit, wdata.ntime)
801 : wdata.r_rmodes.time2step(time + wdata.t_winit, wdata.ntime);
802 int double_data = evaluate_read_resize(
803 wdata.ntime, read_flag, wdata.resize_flag, wdata.t, wdata.t_last,
804 new_ntime, wdata.t_winit, wdata.dt_modes, time);
805 // Check if stop time is exceeded, introduce offset to ntime
806 if (read_flag) {
807 // Need to only check when reading is occurring
808 if (wdata.ntime + wdata.n_offset > wdata.n_wstop) {
809 // If exceeding stop step, calculate new offset
810 wdata.n_offset =
811 update_offset_timestep(wdata.ntime, wdata.n_winit);
812 // Print warning to screen
813 amrex::Print()
814 << "WARNING (waves2amr_ops): available mode data exceeded, "
815 "resetting to beginning of mode data.\n";
816 }
817 }
818 // Resizing (assuming reading is taking place) must happen after regrid
819 if (wdata.regrid_occurred) {
820 // resize_flag remains true until resizing occurs, but
821 // regrid_occurred resets every timestep
822 wdata.resize_flag = true;
823 }
824
825 // Read HOS data if necessary based on time
826 if (read_flag) {
827
828 amrex::Print() << "OceanWaves: Waves2AMR reading from file: step "
829 << wdata.ntime << std::endl
830 << " last interp time " << t_last << std::endl
831 << " simulation time " << time << std::endl
832 << " new data time " << wdata.t << std::endl
833 << " double data read " << double_data
834 << std::endl;
835
836 if (wdata.resize_flag) {
837 // Reset flag
838 wdata.resize_flag = false;
839 // Flags for indicating overlap, assume none at first
840 bool flag_z = false;
841 bool flag_xlo = false;
842 bool flag_xhi = false;
843 bool flag_ylo = false;
844 bool flag_yhi = false;
845 // Get heights for this processor, check overlap in z
846 flag_z =
847 (interp_to_mfab::get_local_height_indices(
848 wdata.indvec, wdata.hvec, ow_velocity.vec_ptrs(),
849 geom) == 1);
850 // No overlap from heights definitely means no interp
851
852 // Check lateral bounds (in x)
853 flag_xlo =
854 (interp_to_mfab::check_lateral_overlap_lo(
855 wdata.gen_length, 0, ow_velocity.vec_ptrs(), geom) ==
856 1);
857 // No overlap with gen region means no interp, unless ...
858 if (!wdata.has_beach) {
859 // ... if overlap exists here, needing interp
860 flag_xhi =
861 (interp_to_mfab::check_lateral_overlap_hi(
862 wdata.beach_length, 0, ow_velocity.vec_ptrs(),
863 geom) == 1);
864 }
865 // Then in y, if y zones are present
866 if (wdata.zone_length_y > constants::EPS) {
867 flag_ylo =
868 (interp_to_mfab::check_lateral_overlap_lo(
869 wdata.zone_length_y, 1, ow_velocity.vec_ptrs(),
870 geom) == 1);
871 flag_yhi =
872 (interp_to_mfab::check_lateral_overlap_hi(
873 wdata.zone_length_y, 1, ow_velocity.vec_ptrs(),
874 geom) == 1);
875 }
876
877 // Wave generation zones are irrelevant for single-phase mode,
878 // indicated by nonexistence of vof
879 if (!vof_exists) {
880 flag_xlo = true;
881 flag_xhi = true;
882 }
883
884 if (flag_z && (flag_xlo || flag_xhi || flag_ylo || flag_yhi)) {
885 // Interpolation is needed
886 wdata.do_interp = true;
887 // Do resizing
888 wdata.sp_eta_vec.resize(
889 static_cast<size_t>(wdata.n0) *
890 static_cast<size_t>(wdata.n1),
891 0.0_rt);
892 wdata.sp_u_vec.resize(
893 static_cast<size_t>(wdata.n0 * wdata.n1) *
894 wdata.indvec.size());
895 wdata.sp_v_vec.resize(
896 static_cast<size_t>(wdata.n0 * wdata.n1) *
897 wdata.indvec.size());
898 wdata.sp_w_vec.resize(
899 static_cast<size_t>(wdata.n0 * wdata.n1) *
900 wdata.indvec.size());
901 // Sizes will remain constant and need for interpolation
902 // will remain until a regrid occurs
903 } else {
904 // No overlapping with spatial data or no overlapping with
905 // relaxation zones, interpolation can be skipped
906 wdata.do_interp = false;
907 }
908 }
909 // Only perform reading where needed, communicate offset though
910 amrex::ParallelDescriptor::ReduceIntMax(wdata.n_offset);
911
912 // Previous W2A data is more recent than last interpolated data
913 t_last = (wdata.ntime - 1) * wdata.dt_modes;
914
915 if (double_data == 0) {
916 // Previous W2A data needs to be copied into old fields
918 w2a_lvs_old, w2a_levelset, 0, 0, 1,
919 w2a_levelset.num_grow());
921 w2a_vel_old, w2a_velocity, 0, 0, AMREX_SPACEDIM,
922 w2a_velocity.num_grow());
923 } else if (double_data == 1) {
924 // Previous W2A data is unknown, read into old fields
925 if (wdata.do_interp) {
926 populate_fields_all_levels(
927 wdata, geom, w2a_lvs_old, w2a_vel_old, -1);
928 }
929
930 // Average down to get fine information on coarse grid where
931 // possible (may be unnecessary)
932 for (int lev = nlevels - 1; lev > 0; --lev) {
933 amrex::average_down(
934 w2a_vel_old(lev), w2a_vel_old(lev - 1), 0,
935 AMREX_SPACEDIM, sim.mesh().refRatio(lev - 1));
936 amrex::average_down(
937 w2a_lvs_old(lev), w2a_lvs_old(lev - 1), 0, 1,
938 sim.mesh().refRatio(lev - 1));
939 }
940 // Fill patch to get correct ghost cells after average down
941 w2a_vel_old.fillpatch(sim.time().new_time());
942 w2a_lvs_old.fillpatch(sim.time().new_time());
943 }
944 // ow fields cancel when double_data == 2, no modification
945
946 // After possible prior read, now read data for this ntime
947 if (wdata.do_interp) {
948 populate_fields_all_levels(
949 wdata, geom, w2a_levelset, w2a_velocity);
950 }
951
952 // Average down to get fine information on coarse grid where
953 // possible and update ghost cells (may be unnecessary)
954 for (int lev = nlevels - 1; lev > 0; --lev) {
955 amrex::average_down(
956 w2a_velocity(lev), w2a_velocity(lev - 1), 0, AMREX_SPACEDIM,
957 sim.mesh().refRatio(lev - 1));
958 w2a_velocity(lev - 1).FillBoundary(geom[lev - 1].periodicity());
959 amrex::average_down(
960 w2a_levelset(lev), w2a_levelset(lev - 1), 0, 1,
961 sim.mesh().refRatio(lev - 1));
962 w2a_levelset(lev - 1).FillBoundary(geom[lev - 1].periodicity());
963 }
964 }
965
966 // Temporally interpolate at every timestep to get target solution
967 time_interpolate_wave_fields(
968 w2a_lvs_old, w2a_vel_old, w2a_levelset, w2a_velocity, t_last, time,
969 wdata.t);
970 // Copy to ow fields, which may be modified for beach
972 ow_levelset, w2a_lvs_old, 0, 0, 1, ow_levelset.num_grow());
974 ow_velocity, w2a_vel_old, 0, 0, AMREX_SPACEDIM,
975 ow_velocity.num_grow());
976#else
977 amrex::ignore_unused(data, time);
978#endif
979 }
980};
981
982} // namespace amr_wind::ocean_waves::ops
983
984#endif /* WAVES2AMR_OPS_H */
FieldRepo & repo()
Return the field repository.
Definition CFDSim.H:75
SimTime & time()
Return simulation time control.
Definition CFDSim.H:65
Definition Field.H:116
const amrex::IntVect & num_grow() const
Ghost cells.
Definition Field.H:137
FieldRepo & repo() const
FieldRepo instance that manages this field.
Definition Field.H:159
void set_default_fillpatch_bc(const SimTime &time, const amrex::BCType::mathematicalBndryTypes bctype=amrex::BCType::hoextrap) noexcept
Definition Field.cpp:381
amrex::Vector< amrex::MultiFab * > vec_ptrs() noexcept
Return a vector of MultiFab pointers for all levels.
Definition Field.cpp:145
Field & state(const FieldState fstate)
Return field at a different time state.
Definition Field.cpp:117
Field & declare_field(const std::string &name, const int ncomp=1, const int ngrow=0, const int nstates=1, const FieldLoc floc=FieldLoc::CELL)
Definition FieldRepo.cpp:86
Field & get_field(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:152
bool field_exists(const std::string &name, const FieldState fstate=FieldState::New) const
Query if field uniquely identified by name and time state exists in repository.
Definition FieldRepo.cpp:218
OceanWavesTrait::MetaType & meta()
Definition OceanWavesTypes.H:89
CFDSim & sim()
Definition OceanWavesTypes.H:83
OceanWavesTrait::InfoType & info()
Definition OceanWavesTypes.H:86
void copy(T1 &dst, const T2 &src, const int srccomp, const int dstcomp, const int numcomp, const amrex::IntVect &nghost)
Definition field_ops.H:157
void lincomb(T1 &dst, const amrex::Real a, const T2 &x, const int xcomp, const amrex::Real b, const T3 &y, const int ycomp, const int dstcomp, const int numcomp, const amrex::IntVect &nghost)
Definition field_ops.H:320
@ Old
Same as FieldState::N.
Definition FieldDescTypes.H:21
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:17
static constexpr amrex::Real EPS
A number very close to zero.
Definition constants.H:13
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real levelset_to_vof(const int i, const int j, const int k, const amrex::Real eps, amrex::Array4< amrex::Real const > const &phi) noexcept
Definition volume_fractions.H:532
Definition OceanWavesOps.H:8
void read_inputs(RelaxZonesBaseData &wdata, OceanWavesInfo &, const ::amr_wind::utils::MultiParser &pp)
Definition relaxation_zones_ops.cpp:16
amrex::GpuArray< amrex::Real, 4 > WaveVec
Definition wave_utils_K.H:13
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE WaveVec harmonize_profiles_2d(const amrex::Real x, const amrex::Real y, const amrex::Real left_bdy, const amrex::Real left_length, const amrex::Real right_bdy, const amrex::Real right_length, const amrex::Real lo_y_bdy, const amrex::Real hi_y_bdy, const amrex::Real y_length, const WaveVec left, const WaveVec bulk, const WaveVec right)
Definition wave_utils_K.H:106
bool has_beach
Definition RelaxationZones.H:41
amrex::Real zone_length_y
Definition RelaxationZones.H:33
amrex::Real beach_length
Definition RelaxationZones.H:30
bool init_wave_field
Definition RelaxationZones.H:37
amrex::Real current
Definition RelaxationZones.H:35
amrex::Real gen_length
Definition RelaxationZones.H:28
bool regrid_occurred
Definition RelaxationZones.H:45
amrex::Real zsl
Definition RelaxationZones.H:20
amrex::Real xlo1
Definition W2AWaves.H:32
amrex::Real dimL
Definition W2AWaves.H:34
amrex::Real dt_modes
Definition W2AWaves.H:36
amrex::Gpu::DeviceVector< amrex::Real > sp_v_vec
Definition W2AWaves.H:87
std::vector< amrex::Real > r_mFS
Definition W2AWaves.H:58
std::vector< std::complex< amrex::Real > > c_mFS
Definition W2AWaves.H:57
amrex::Vector< int > indvec
Definition W2AWaves.H:79
std::vector< amrex::Real > r_mAdd
Definition W2AWaves.H:58
int n1
Definition W2AWaves.H:22
amrex::Gpu::DeviceVector< amrex::Real > sp_u_vec
Definition W2AWaves.H:87
int n_winit
Definition W2AWaves.H:44
amrex::Vector< amrex::Real > hvec
Definition W2AWaves.H:77
int n_offset
Definition W2AWaves.H:48
amrex::Real dx0
Definition W2AWaves.H:28
int ntime
Definition W2AWaves.H:42
int n1_sp
Definition W2AWaves.H:26
std::vector< std::complex< amrex::Real > > c_mZ
Definition W2AWaves.H:57
amrex::Real t_last
Definition W2AWaves.H:52
amrex::Real xlo0
Definition W2AWaves.H:31
amrex::Gpu::DeviceVector< amrex::Real > sp_w_vec
Definition W2AWaves.H:88
std::vector< std::complex< amrex::Real > > c_mX
Definition W2AWaves.H:57
int n0
Definition W2AWaves.H:21
int n0_sp
Definition W2AWaves.H:25
std::vector< amrex::Real > r_mY
Definition W2AWaves.H:58
amrex::Real dx1
Definition W2AWaves.H:29
amrex::Real t_winit
Definition W2AWaves.H:40
int n2
Definition W2AWaves.H:23
int n_wstop
Definition W2AWaves.H:46
std::vector< amrex::Real > r_mZ
Definition W2AWaves.H:58
std::string modes_file
Definition W2AWaves.H:19
std::vector< std::complex< amrex::Real > > c_mY
Definition W2AWaves.H:57
bool is_ocean
Definition W2AWaves.H:81
amrex::Gpu::DeviceVector< amrex::Real > sp_eta_vec
Definition W2AWaves.H:87
bool do_interp
Definition W2AWaves.H:83
bool resize_flag
Definition W2AWaves.H:85
std::vector< amrex::Real > r_mX
Definition W2AWaves.H:58
amrex::Real t
Definition W2AWaves.H:54
Definition W2AWaves.H:138
OceanWavesDataHolder< W2AWaves > DataType
Definition W2AWaves.H:141
W2AWavesData MetaType
Definition W2AWaves.H:140
void operator()(W2AWaves::DataType &data, int level, const amrex::Geometry &geom, bool multiphase_mode)
Definition waves2amr_ops.H:599
Definition OceanWavesOps.H:14
void operator()(W2AWaves::DataType &data, const ::amr_wind::utils::MultiParser &pp)
Definition waves2amr_ops.H:209
Definition OceanWavesOps.H:11
void operator()(W2AWaves::DataType &data, const amrex::Real time)
Definition waves2amr_ops.H:769