/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
80void postprocess_velocity_mfab_liquid(
81 amrex::MultiFab& vel_mfab,
82 amrex::MultiFab& lvs_mfab,
83 const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM> dx)
84{
85 auto vel = vel_mfab.arrays();
86 const auto phi = lvs_mfab.const_arrays();
87 amrex::ParallelFor(
88 vel_mfab, amrex::IntVect(3),
89 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
90 // Set velocity to zero if no liquid present
91 const amrex::Real cell_length_2D =
92 std::sqrt(dx[0] * dx[0] + dx[2] * dx[2]);
93 if (phi[nbx](i, j, k) + cell_length_2D < 0) {
94 vel[nbx](i, j, k, 0) = 0.0;
95 vel[nbx](i, j, k, 1) = 0.0;
96 vel[nbx](i, j, k, 2) = 0.0;
97 }
98 });
99 amrex::Gpu::streamSynchronize();
100}
101
102void postprocess_velocity_field_liquid(
103 amr_wind::Field& vel_field,
104 amr_wind::Field& lvs_field,
105 amrex::Vector<amrex::Geometry>& geom_all)
106{
107 int nlevels = vel_field.repo().num_active_levels();
108 for (int lev = 0; lev < nlevels; ++lev) {
109 const auto& dx_lev = geom_all[lev].CellSizeArray();
110 postprocess_velocity_mfab_liquid(
111 vel_field(lev), lvs_field(lev), dx_lev);
112 }
113}
114
115int update_offset_timestep(const int ntime, const int n0)
116{
117 // Offending timestep (goes too far): ntime + offset
118 // Farthest back timesteps are permitted to go: n0
119 // Subtract offset by offending timestep, add back lower limit
120 // new offset = offset - (ntime + offset) + n0
121 return (-ntime + n0);
122}
123
124void populate_fields_all_levels(
126 amrex::Vector<amrex::Geometry>& geom_all,
127 amr_wind::Field& lvs_field,
128 amr_wind::Field& vel_field,
129 int ntime_off = 0)
130{
131
132 // Get data from modes
133 bool no_EOF = wdata.rmodes.get_data(
134 wdata.ntime + wdata.n_offset + ntime_off, wdata.mX, wdata.mY, wdata.mZ,
135 wdata.mFS);
136 // Navigate when end of file is reached
137 if (!no_EOF) {
138 // End of file detected, reset reading
139 wdata.n_offset = update_offset_timestep(wdata.ntime, wdata.n_winit);
140 // Print warning to screen
141 amrex::Print() << "WARNING (waves2amr_ops): end of mode data file "
142 "detected, resetting to beginning of mode data.\n";
143 // Read data again, now from a valid timestep
144 no_EOF = wdata.rmodes.get_data(
145 wdata.ntime + wdata.n_offset + ntime_off, wdata.mX, wdata.mY,
146 wdata.mZ, wdata.mFS);
147 // If no valid data is detected at this point, abort
148 if (!no_EOF) {
149 amrex::Abort(
150 "waves2amr_ops: end of mode data file detected after "
151 "resetting to beginning; please evaluate HOS_init_time "
152 "or HOS_init_timestep and check the length of the mode "
153 "file.");
154 }
155 }
156
157 // Convert to spatial data in vectors
158 modes_hosgrid::copy_complex(wdata.n0, wdata.n1, wdata.mFS, wdata.eta_mptr);
159 modes_hosgrid::populate_hos_eta(
160 wdata.rmodes, wdata.plan, wdata.eta_mptr, wdata.sp_eta_vec);
161
162 for (int iht = 0; iht < wdata.indvec.size(); ++iht) {
163 // Get sample height
164 amrex::Real ht = wdata.hvec[wdata.indvec[iht]];
165 // Sample velocity
166 modes_hosgrid::populate_hos_vel(
167 wdata.rmodes, ht, wdata.mX, wdata.mY, wdata.mZ, wdata.plan,
168 wdata.u_mptr, wdata.v_mptr, wdata.w_mptr, wdata.sp_u_vec,
169 wdata.sp_v_vec, wdata.sp_w_vec, iht * wdata.n0 * wdata.n1);
170 }
171
172 // Interpolate to fields (vector of MultiFabs)
173 interp_to_mfab::interp_eta_to_levelset_field(
174 wdata.n0, wdata.n1, wdata.dx0, wdata.dx1, wdata.zsl, wdata.sp_eta_vec,
175 lvs_field.vec_ptrs(), geom_all);
176 interp_to_mfab::interp_velocity_to_field(
177 wdata.n0, wdata.n1, wdata.dx0, wdata.dx1, wdata.indvec, wdata.hvec,
178 wdata.sp_u_vec, wdata.sp_v_vec, wdata.sp_w_vec, vel_field.vec_ptrs(),
179 geom_all);
180
181 // Zero velocity in pure gas cells
182 postprocess_velocity_field_liquid(vel_field, lvs_field, geom_all);
183}
184
185} // namespace
186#endif
187
189
190template <>
192{
194 // cppcheck-suppress constParameterReference
195 W2AWaves::DataType& data,
196 const ::amr_wind::utils::MultiParser& pp)
197 {
198// Check for W2A initialization
199#ifndef AMR_WIND_USE_W2A
200 // Assert Waves2AMR must be used for initial condition file
201 amrex::Abort(
202 "ocean_waves/W2AWaves: AMR-Wind was not built with Waves2AMR "
203 "support; associated wave data cannot be processed for relaxation "
204 "zones.");
205
206 amrex::ignore_unused(data, pp);
207#else
208 auto& wdata = data.meta();
209 auto& info = data.info();
210 relaxation_zones::read_inputs(wdata, info, pp);
211
212 if (wdata.current > constants::TIGHT_TOL) {
213 amrex::Abort(
214 "Current is specified as nonzero, but current is not yet "
215 "implemented for W2A Waves.");
216 }
217
218 pp.get("HOS_modes_filename", wdata.modes_file);
219 pp.query("HOS_init_timestep", wdata.ntime);
220 if (!pp.contains("HOS_init_timestep")) {
221 pp.query("HOS_init_time", wdata.t_winit);
222 }
223
224 // Default fftw_plan is deterministic
225 std::string fftw_planner_flag{"estimate"};
226 pp.query("fftw_planner_flag", fftw_planner_flag);
227
228 amrex::Vector<amrex::Real> prob_lo_input(AMREX_SPACEDIM);
229 amrex::ParmParse pp_geom("geometry");
230 pp_geom.getarr("prob_lo", prob_lo_input);
231
232 // Read user inputs to height vector
233 int nheights = 0;
234 int nh_above = 1;
235 amrex::Real dz0 = 0.;
236 pp.get("number_interp_points_in_z", nheights);
237 pp.get("interp_spacing_at_surface", dz0);
238 pp.query("number_interp_above_surface", nh_above);
239
240 // Initialize mode reader
241 bool file_exists = wdata.rmodes.initialize(wdata.modes_file);
242
243 // Abort if file does not exist
244 if (!file_exists) {
245 amrex::Abort(
246 "Waves2AMR ReadInputsOp: modes file requested does not exist");
247 }
248
249 // Get dt of HOS data
250 wdata.dt_modes = wdata.rmodes.get_dtout();
251
252 // Get initial time and timestep synced
253 if (wdata.t_winit > 0.0) {
254 // If initial time was specified
255 // Get time index near requested time
256 wdata.ntime = wdata.rmodes.time2step(wdata.t_winit, wdata.ntime);
257 // Sync time to time index
258 wdata.t_winit = wdata.dt_modes * wdata.ntime;
259 // Save first timestep
260 wdata.n_winit = wdata.ntime;
261 } else {
262 // If initial timestep is being used
263 wdata.t_winit = wdata.dt_modes * wdata.ntime;
264 // Save first timestep
265 wdata.n_winit = wdata.ntime;
266 }
267
268 // Initialize variables to store modes
269 int vsize = wdata.rmodes.get_vector_size();
270 double initval = 0.0;
271 wdata.mX.resize(vsize, initval);
272 wdata.mY.resize(vsize, initval);
273 wdata.mZ.resize(vsize, initval);
274 wdata.mFS.resize(vsize, initval);
275
276 // Get dimensions of data
277 wdata.n0 = wdata.rmodes.get_first_dimension();
278 wdata.n1 = wdata.rmodes.get_second_dimension();
279 // Get resolution
280 wdata.dx0 = wdata.rmodes.get_xlen() / wdata.n0;
281 wdata.dx1 = wdata.rmodes.get_ylen() / wdata.n1;
282 // Get depth
283 const amrex::Real depth = wdata.rmodes.get_depth();
284 // Get dimensional length
285 wdata.dimL = wdata.rmodes.get_L();
286 // Get nominal last timestep of data
287 wdata.n_wstop =
288 (int)((wdata.rmodes.get_Tstop() + 1e-8) / wdata.dt_modes);
289
290 // Check if stop time is exceeded, introduce offset to ntime
291 if (wdata.ntime + wdata.n_offset > wdata.n_wstop) {
292 // If exceeding stop step, calculate new offset
293 wdata.n_offset = update_offset_timestep(wdata.ntime, wdata.n_winit);
294 // Print warning to screen
295 amrex::Print()
296 << "WARNING (waves2amr_ops): available mode data exceeded, "
297 "resetting to beginning of mode data.\n";
298 }
299
300 // Warning if depth does not correspond to simulation
301 if (std::abs(depth - (wdata.zsl - prob_lo_input[2])) > 1e-3 * depth) {
302 amrex::Print()
303 << "WARNING: Mismatch between water depths from AMR-Wind "
304 "domain and HOS data interpreted by Waves2AMR. \n ^This "
305 "warning is not a concern when using waves as terrain.\n";
306 }
307
308 // Allocate pointers for FFTW
309 wdata.eta_mptr = modes_hosgrid::allocate_complex(wdata.n0, wdata.n1);
310 wdata.u_mptr = modes_hosgrid::allocate_complex(wdata.n0, wdata.n1);
311 wdata.v_mptr = modes_hosgrid::allocate_complex(wdata.n0, wdata.n1);
312 wdata.w_mptr = modes_hosgrid::allocate_complex(wdata.n0, wdata.n1);
313
314 // Set up planner flag based on input
315 auto plan_f = modes_hosgrid::planner_flags::estimate;
316 if (fftw_planner_flag == "patient") {
317 plan_f = modes_hosgrid::planner_flags::patient;
318 } else if (fftw_planner_flag == "exhaustive") {
319 plan_f = modes_hosgrid::planner_flags::exhaustive;
320 } else if (fftw_planner_flag == "measure") {
321 plan_f = modes_hosgrid::planner_flags::measure;
322 } else if (!(fftw_planner_flag == "estimate")) {
323 amrex::Print()
324 << "WARNING (waves2amr_ops): invalid fftw_planner_flag "
325 "specified; defaulting to estimate (FFTW_ESTIMATE).\n";
326 }
327 // Set up plan for FFTW
328 wdata.plan = modes_hosgrid::plan_ifftw(
329 wdata.n0, wdata.n1, wdata.eta_mptr, plan_f);
330
331 // Create height vector for velocity mode conversion before
332 // interpolation, with prob_lo as bottom
333 int flag = interp_to_mfab::create_height_vector(
334 wdata.hvec, nheights, dz0, wdata.zsl, prob_lo_input[2], nh_above);
335 // Fail if flag indicates it should
336 if (flag > 0) {
337 amrex::Abort(
338 "Waves2AMR ReadInputsOp: create_height_vector error, failure "
339 "code " +
340 std::to_string(flag));
341 }
342
343 // Always get modes on every processor for initialization
344 bool no_EOF = wdata.rmodes.get_data(
345 wdata.ntime + wdata.n_offset, wdata.mX, wdata.mY, wdata.mZ,
346 wdata.mFS);
347 if (!no_EOF) {
348 // End of file detected, reset reading
349 wdata.n_offset = update_offset_timestep(wdata.ntime, wdata.n_winit);
350 // Print warning to screen
351 amrex::Print()
352 << "WARNING (waves2amr_ops): end of mode data file "
353 "detected, resetting to beginning of mode data.\n";
354 // Read data again, now from a valid timestep
355 no_EOF = wdata.rmodes.get_data(
356 wdata.ntime + wdata.n_offset, wdata.mX, wdata.mY, wdata.mZ,
357 wdata.mFS);
358 // If no valid data is detected at this point, abort
359 if (!no_EOF) {
360 amrex::Abort(
361 "waves2amr_ops: end of mode data file detected after "
362 "resetting to beginning; please evaluate HOS_init_time "
363 "or HOS_init_timestep and check the length of the mode "
364 "file.");
365 }
366 }
367
368 // Convert modes to spatial data
369 modes_hosgrid::copy_complex(
370 wdata.n0, wdata.n1, wdata.mFS, wdata.eta_mptr);
371 wdata.sp_eta_vec.resize(
372 static_cast<size_t>(wdata.n0) * static_cast<size_t>(wdata.n1), 0.0);
373 modes_hosgrid::populate_hos_eta(
374 wdata.rmodes, wdata.plan, wdata.eta_mptr, wdata.sp_eta_vec);
375 // Mesh is not yet created, so get data at every height
376 const auto n_hts = wdata.hvec.size();
377 wdata.sp_u_vec.resize(static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
378 wdata.sp_v_vec.resize(static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
379 wdata.sp_w_vec.resize(static_cast<size_t>(wdata.n0 * wdata.n1) * n_hts);
380 for (int iht = 0; iht < static_cast<int>(n_hts); ++iht) {
381 // Get sample height
382 amrex::Real ht = wdata.hvec[iht];
383 // Sample velocity
384 modes_hosgrid::populate_hos_vel(
385 wdata.rmodes, ht, wdata.mX, wdata.mY, wdata.mZ, wdata.plan,
386 wdata.u_mptr, wdata.v_mptr, wdata.w_mptr, wdata.sp_u_vec,
387 wdata.sp_v_vec, wdata.sp_w_vec, iht * wdata.n0 * wdata.n1);
388 }
389
390 // Declare fields for HOS
391 auto& w2a_levelset =
392 data.sim().repo().declare_field("w2a_levelset", 1, 3, 1);
393 auto& w2a_velocity = data.sim().repo().declare_field(
394 "w2a_velocity", AMREX_SPACEDIM, 3, 1);
395
396 // Extrapolation can work well when finer data is available
397 w2a_levelset.set_default_fillpatch_bc(data.sim().time());
398 w2a_velocity.set_default_fillpatch_bc(data.sim().time());
399#endif
400 }
401}; // namespace ops
402
403template <>
405{
406 void
407 // cppcheck-suppress constParameterReference
410 data, int level, const amrex::Geometry & geom, bool multiphase_mode)
411 {
412
413#ifdef AMR_WIND_USE_W2A
414 auto& wdata = data.meta();
415
416 auto& sim = data.sim();
417
418 // Fill ow fields, then populate flow fields according to setup
419 auto& ow_levelset = sim.repo().get_field("ow_levelset");
420 auto& ow_velocity = sim.repo().get_field("ow_velocity");
421
422 auto& velocity = sim.repo().get_field("velocity");
423 // Set w2a fields to default values to prep for updates
424 auto& w2a_levelset = sim.repo().get_field("w2a_levelset");
425 auto& w2a_velocity = sim.repo().get_field("w2a_velocity");
426
427 Field* levelset{nullptr};
428 if (multiphase_mode) {
429 levelset = &sim.repo().get_field("levelset");
430 }
431
432 const auto& problo = geom.ProbLoArray();
433 const auto& probhi = geom.ProbHiArray();
434 const auto& dx = geom.CellSizeArray();
435
436 // Set t_last to 0.0 to signify information read in
437 wdata.t_last = 0.0;
438
439 // indvec is complete upon initialization (all heights every proc)
440 amrex::Vector<int> indvec(wdata.hvec.size());
441 std::iota(indvec.begin(), indvec.end(), 0);
442 // Interpolate to MultiFabs (one level at a time)
443 interp_to_mfab::interp_eta_to_levelset_multifab(
444 wdata.n0, wdata.n1, wdata.dx0, wdata.dx1, wdata.zsl,
445 wdata.sp_eta_vec, ow_levelset(level), problo, dx);
446 interp_to_mfab::interp_velocity_to_multifab(
447 wdata.n0, wdata.n1, wdata.dx0, wdata.dx1, indvec, wdata.hvec,
448 wdata.sp_u_vec, wdata.sp_v_vec, wdata.sp_w_vec, ow_velocity(level),
449 problo, dx);
450 // Zero velocity in pure air cells
451 postprocess_velocity_mfab_liquid(
452 ow_velocity(level), ow_levelset(level), dx);
453
454 // Populate flow fields according to intended forcing and init setup
455 const auto& ow_phi = ow_levelset(level).const_arrays();
456 const auto& ow_vel = ow_velocity(level).const_arrays();
457 const auto& vel = velocity(level).arrays();
458 const auto& phi_arrs = multiphase_mode
459 ? (*levelset)(level).arrays()
460 : amrex::MultiArray4<amrex::Real>();
461
462 const auto& w2a_phi = w2a_levelset(level).arrays();
463 const auto& w2a_vel = w2a_velocity(level).arrays();
464
465 const amrex::Real gen_length = wdata.gen_length;
466 const amrex::Real beach_length = wdata.beach_length;
467 const amrex::Real zero_sea_level = wdata.zsl;
468
469 const bool has_beach = wdata.has_beach && multiphase_mode;
470 const bool init_wave_field = wdata.init_wave_field || !multiphase_mode;
471
472 amrex::ParallelFor(
473 velocity(level), amrex::IntVect(3),
474 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
475 const amrex::Real x = problo[0] + (i + 0.5) * dx[0];
476 const amrex::Real z = problo[2] + (k + 0.5) * dx[2];
477
478 // Wave profile
479 const utils::WaveVec wave_sol{
480 ow_vel[nbx](i, j, k, 0), ow_vel[nbx](i, j, k, 1),
481 ow_vel[nbx](i, j, k, 2), ow_phi[nbx](i, j, k) + z};
482 // Quiescent profile
483 const utils::WaveVec quiescent{0.0, 0.0, 0.0, zero_sea_level};
484
485 // Specify initial state for each region of domain
486 const auto bulk = init_wave_field ? wave_sol : quiescent;
487 const auto outlet = has_beach ? quiescent : wave_sol;
488
489 const auto local_profile = utils::harmonize_profiles_1d(
490 x, problo[0], gen_length, probhi[0], beach_length, wave_sol,
491 bulk, outlet);
492
493 const amrex::Real phi = local_profile[3] - z;
494 const amrex::Real cell_length_2D =
495 std::sqrt(dx[0] * dx[0] + dx[2] * dx[2]);
496 if (phi + cell_length_2D >= 0) {
497 vel[nbx](i, j, k, 0) = local_profile[0];
498 vel[nbx](i, j, k, 1) = local_profile[1];
499 vel[nbx](i, j, k, 2) = local_profile[2];
500 }
501
502 if (multiphase_mode) {
503 phi_arrs[nbx](i, j, k) = phi;
504 }
505
506 // Default w2a values matter for where no updates happen
507 w2a_phi[nbx](i, j, k) = quiescent[3] - z;
508 w2a_vel[nbx](i, j, k, 0) = quiescent[0];
509 w2a_vel[nbx](i, j, k, 1) = quiescent[0];
510 w2a_vel[nbx](i, j, k, 2) = quiescent[0];
511 });
512 amrex::Gpu::streamSynchronize();
513#else
514 amrex::ignore_unused(data, level, geom, multiphase_mode);
515#endif
516 }
517}; // namespace ocean_waves
518
519template <>
521{
522 // cppcheck-suppress constParameterReference
523 void operator()(W2AWaves::DataType& data, const amrex::Real time)
524 {
525
526#ifdef AMR_WIND_USE_W2A
527 auto& wdata = data.meta();
528 auto& sim = data.sim();
529
530 // Update ow fields every time
531 auto& ow_levelset = sim.repo().get_field("ow_levelset");
532 auto& ow_velocity = sim.repo().get_field("ow_velocity");
533 // Update HOS fields when necessary
534 auto& w2a_levelset = sim.repo().get_field("w2a_levelset");
535 auto& w2a_velocity = sim.repo().get_field("w2a_velocity");
536
537 // Proxy for multiphase mode of ocean waves
538 bool vof_exists = sim.repo().field_exists("vof");
539
540 auto nlevels = sim.repo().num_active_levels();
541 auto geom = sim.mesh().Geom();
542
543 // Get value for time interpolation
544 amrex::Real t_last = wdata.t_last;
545
546 // Check if new HOS data needs to be read
547 bool read_flag = false;
548 // Check if time indicates reading must take place
549 int new_ntime =
550 wdata.rmodes.time2step(time + wdata.t_winit, wdata.ntime);
551 int double_data = evaluate_read_resize(
552 wdata.ntime, read_flag, wdata.resize_flag, wdata.t, wdata.t_last,
553 new_ntime, wdata.t_winit, wdata.dt_modes, time);
554 // Check if stop time is exceeded, introduce offset to ntime
555 if (read_flag) {
556 // Need to only check when reading is occurring
557 if (wdata.ntime + wdata.n_offset > wdata.n_wstop) {
558 // If exceeding stop step, calculate new offset
559 wdata.n_offset =
560 update_offset_timestep(wdata.ntime, wdata.n_winit);
561 // Print warning to screen
562 amrex::Print()
563 << "WARNING (waves2amr_ops): available mode data exceeded, "
564 "resetting to beginning of mode data.\n";
565 }
566 }
567 // Resizing (assuming reading is taking place) must happen after regrid
568 if (wdata.regrid_occurred) {
569 // resize_flag remains true until resizing occurs, but
570 // regrid_occurred resets every timestep
571 wdata.resize_flag = true;
572 }
573
574 // Read HOS data if necessary based on time
575 if (read_flag) {
576
577 if (wdata.resize_flag) {
578 // Reset flag
579 wdata.resize_flag = false;
580 // Flags for indicating overlap, assume none at first
581 bool flag_z = false;
582 bool flag_xlo = false;
583 bool flag_xhi = false;
584 // Get heights for this processor, check overlap in z
585 flag_z =
586 (interp_to_mfab::get_local_height_indices(
587 wdata.indvec, wdata.hvec, ow_velocity.vec_ptrs(),
588 geom) == 1);
589 // No overlap from heights definitely means no interp
590
591 // Check lateral bounds (in x)
592 const int dir = 0;
593 flag_xlo =
594 (interp_to_mfab::check_lateral_overlap_lo(
595 wdata.gen_length, dir, ow_velocity.vec_ptrs(), geom) ==
596 1);
597 // No overlap with gen region means no interp, unless ...
598 if (wdata.has_outprofile) {
599 // ... if overlap exists here, needing interp
600 flag_xhi =
601 (interp_to_mfab::check_lateral_overlap_hi(
602 wdata.beach_length, dir, ow_velocity.vec_ptrs(),
603 geom) == 1);
604 }
605
606 // Wave generation zones are irrelevant for single-phase mode,
607 // indicated by nonexistence of vof
608 if (!vof_exists) {
609 flag_xlo = true;
610 flag_xhi = true;
611 }
612
613 if (flag_z && (flag_xlo || flag_xhi)) {
614 // Interpolation is needed
615 wdata.do_interp = true;
616 // Do resizing
617 wdata.sp_eta_vec.resize(
618 static_cast<size_t>(wdata.n0) *
619 static_cast<size_t>(wdata.n1),
620 0.0);
621 wdata.sp_u_vec.resize(
622 static_cast<size_t>(wdata.n0 * wdata.n1) *
623 wdata.indvec.size());
624 wdata.sp_v_vec.resize(
625 static_cast<size_t>(wdata.n0 * wdata.n1) *
626 wdata.indvec.size());
627 wdata.sp_w_vec.resize(
628 static_cast<size_t>(wdata.n0 * wdata.n1) *
629 wdata.indvec.size());
630 // Sizes will remain constant and need for interpolation
631 // will remain until a regrid occurs
632 } else {
633 // No overlapping with spatial data or no overlapping with
634 // relaxation zones, interpolation can be skipped
635 wdata.do_interp = false;
636 }
637 }
638 // Only perform reading where needed, communicate offset though
639 amrex::ParallelDescriptor::ReduceIntMax(wdata.n_offset);
640
641 // If double read is required, then copy older wave data to ow_
642 // fields and modify interpolation parameters to get things right
643 if (double_data == 1) {
644 if (wdata.do_interp) {
645 populate_fields_all_levels(
646 wdata, geom, ow_levelset, ow_velocity, -1);
647 }
648
649 // Average down to get fine information on coarse grid where
650 // possible (may be unnecessary)
651 for (int lev = nlevels - 1; lev > 0; --lev) {
652 amrex::average_down(
653 ow_velocity(lev), ow_velocity(lev - 1), 0,
654 AMREX_SPACEDIM, sim.mesh().refRatio(lev - 1));
655 amrex::average_down(
656 ow_levelset(lev), ow_levelset(lev - 1), 0, 1,
657 sim.mesh().refRatio(lev - 1));
658 }
659 // Fill patch to get correct ghost cells after average down
660 ow_velocity.fillpatch(sim.time().new_time());
661 ow_levelset.fillpatch(sim.time().new_time());
662
663 // Prior t_last (corresponding to ow fields)
664 t_last = (wdata.ntime - 1) * wdata.dt_modes;
665 } else if (double_data == 2) {
666 // Restarting simulation or taking a big step, new time at ntime
667 // Initialize ow fields to 0 for time interp, will be replaced
668 if (wdata.do_interp) {
669 ow_levelset.setVal(0.0, ow_levelset.num_grow()[0]);
670 ow_velocity.setVal(0.0, ow_levelset.num_grow()[0]);
671 }
672 // No modification needed for t_last, leads to interp factor = 1
673 }
674
675 // After possible prior read, now read data for this ntime
676 if (wdata.do_interp) {
677 populate_fields_all_levels(
678 wdata, geom, w2a_levelset, w2a_velocity);
679 }
680
681 // Average down to get fine information on coarse grid where
682 // possible and update ghost cells (may be unnecessary)
683 for (int lev = nlevels - 1; lev > 0; --lev) {
684 amrex::average_down(
685 w2a_velocity(lev), w2a_velocity(lev - 1), 0, AMREX_SPACEDIM,
686 sim.mesh().refRatio(lev - 1));
687 w2a_velocity(lev - 1).FillBoundary(geom[lev - 1].periodicity());
688 amrex::average_down(
689 w2a_levelset(lev), w2a_levelset(lev - 1), 0, 1,
690 sim.mesh().refRatio(lev - 1));
691 w2a_levelset(lev - 1).FillBoundary(geom[lev - 1].periodicity());
692 }
693 }
694
695 // Temporally interpolate at every timestep to get target solution
696 for (int lev = 0; lev < nlevels; ++lev) {
697 auto phi = ow_levelset(lev).arrays();
698 auto vel = ow_velocity(lev).arrays();
699 const auto W2A_phi = w2a_levelset(lev).const_arrays();
700 const auto W2A_vel = w2a_velocity(lev).const_arrays();
701
702 const amrex::Real W2A_t = wdata.t;
703 amrex::ParallelFor(
704 ow_levelset(lev), amrex::IntVect(3),
705 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
706 // Interpolate temporally every time
707 phi[nbx](i, j, k) +=
708 (W2A_phi[nbx](i, j, k) - phi[nbx](i, j, k)) *
709 (time - t_last) / (W2A_t - t_last + 1e-16);
710 vel[nbx](i, j, k, 0) +=
711 (W2A_vel[nbx](i, j, k, 0) - vel[nbx](i, j, k, 0)) *
712 (time - t_last) / (W2A_t - t_last + 1e-16);
713 vel[nbx](i, j, k, 1) +=
714 (W2A_vel[nbx](i, j, k, 1) - vel[nbx](i, j, k, 1)) *
715 (time - t_last) / (W2A_t - t_last + 1e-16);
716 vel[nbx](i, j, k, 2) +=
717 (W2A_vel[nbx](i, j, k, 2) - vel[nbx](i, j, k, 2)) *
718 (time - t_last) / (W2A_t - t_last + 1e-16);
719 });
720 }
721 amrex::Gpu::streamSynchronize();
722#else
723 amrex::ignore_unused(data, time);
724#endif
725 }
726};
727
728} // namespace amr_wind::ocean_waves::ops
729
730#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
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:377
amrex::Vector< amrex::MultiFab * > vec_ptrs() noexcept
Return a vector of MultiFab pointers for all levels.
Definition Field.cpp:142
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
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
bool field_exists(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:215
Definition OceanWavesTypes.H:59
OceanWavesTrait::MetaType & meta()
Definition OceanWavesTypes.H:89
CFDSim & sim()
Definition OceanWavesTypes.H:83
OceanWavesTrait::InfoType & info()
Definition OceanWavesTypes.H:86
static constexpr amrex::Real TIGHT_TOL
A tight tolerance.
Definition constants.H:19
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:36
bool has_outprofile
Definition RelaxationZones.H:37
amrex::Real beach_length
Definition RelaxationZones.H:27
bool init_wave_field
Definition RelaxationZones.H:32
amrex::Real current
Definition RelaxationZones.H:30
amrex::Real gen_length
Definition RelaxationZones.H:25
bool regrid_occurred
Definition RelaxationZones.H:41
amrex::Real zsl
Definition RelaxationZones.H:17
Definition W2AWaves.H:14
amrex::Real dimL
Definition W2AWaves.H:24
amrex::Real dt_modes
Definition W2AWaves.H:26
amrex::Gpu::DeviceVector< amrex::Real > sp_v_vec
Definition W2AWaves.H:69
std::vector< std::complex< double > > mFS
Definition W2AWaves.H:47
amrex::Vector< int > indvec
Definition W2AWaves.H:63
int n1
Definition W2AWaves.H:19
amrex::Gpu::DeviceVector< amrex::Real > sp_u_vec
Definition W2AWaves.H:69
int n_winit
Definition W2AWaves.H:34
amrex::Vector< amrex::Real > hvec
Definition W2AWaves.H:61
int n_offset
Definition W2AWaves.H:38
amrex::Real dx0
Definition W2AWaves.H:21
int ntime
Definition W2AWaves.H:32
std::vector< std::complex< double > > mX
Definition W2AWaves.H:47
amrex::Real t_last
Definition W2AWaves.H:42
amrex::Gpu::DeviceVector< amrex::Real > sp_w_vec
Definition W2AWaves.H:70
int n0
Definition W2AWaves.H:18
std::vector< std::complex< double > > mY
Definition W2AWaves.H:47
amrex::Real dx1
Definition W2AWaves.H:22
amrex::Real t_winit
Definition W2AWaves.H:30
int n_wstop
Definition W2AWaves.H:36
std::string modes_file
Definition W2AWaves.H:16
std::vector< std::complex< double > > mZ
Definition W2AWaves.H:47
amrex::Gpu::DeviceVector< amrex::Real > sp_eta_vec
Definition W2AWaves.H:69
bool do_interp
Definition W2AWaves.H:65
bool resize_flag
Definition W2AWaves.H:67
amrex::Real t
Definition W2AWaves.H:44
Definition W2AWaves.H:95
void operator()(W2AWaves::DataType &data, int level, const amrex::Geometry &geom, bool multiphase_mode)
Definition waves2amr_ops.H:408
Definition OceanWavesOps.H:14
void operator()(W2AWaves::DataType &data, const ::amr_wind::utils::MultiParser &pp)
Definition waves2amr_ops.H:193
Definition OceanWavesOps.H:11
void operator()(W2AWaves::DataType &data, const amrex::Real time)
Definition waves2amr_ops.H:523