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