13#include "AMReX_MultiFabUtil.H"
14#include "AMReX_REAL.H"
16using namespace amrex::literals;
18#ifdef AMR_WIND_USE_W2A
20int evaluate_read_resize(
27 const amrex::Real wtinit,
28 const amrex::Real wdt,
29 const amrex::Real time)
34 if (new_ntime != ntime) {
38 if (new_ntime > ntime + 1) {
45 wtime = new_ntime * wdt - wtinit;
47 if (double_data == 1 && std::abs(wtime - time) <= 1.0e-10_rt) {
53 if (t_last < -1.0e-10_rt) {
60 if (std::abs(wtime - time) > 1.0e-10_rt) {
68 }
else if (std::abs(t_last) < 1.0e-10_rt) {
83int update_offset_timestep(
const int ntime,
const int n0)
92void populate_fields_all_levels(
94 amrex::Vector<amrex::Geometry>& geom_all,
103 ? wdata.c_rmodes.get_data(
106 : wdata.r_rmodes.get_data(
114 amrex::Print() <<
"WARNING (waves2amr_ops): end of mode data file "
115 "detected, resetting to beginning of mode data.\n";
118 ? wdata.c_rmodes.get_data(
121 : wdata.r_rmodes.get_data(
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 "
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,
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,
149 for (
int iht = 0; iht < wdata.
indvec.size(); ++iht) {
151 amrex::Real ht = wdata.
hvec[wdata.
indvec[iht]];
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,
158 iht * wdata.
n0 * wdata.
n1);
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,
165 iht * wdata.
n0 * wdata.
n1);
170 interp_to_mfab::interp_eta_to_levelset_field(
174 interp_to_mfab::interp_velocity_to_field(
180void time_interpolate_wave_fields(
185 const amrex::Real t_old,
186 const amrex::Real t_new,
187 const amrex::Real t_future)
189 const amrex::Real b =
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,
197 vel_field, a, vel_field, 0, b, future_vel_field, 0, 0, AMREX_SPACEDIM,
212 const ::amr_wind::utils::MultiParser& pp)
215#ifndef AMR_WIND_USE_W2A
218 "ocean_waves/W2AWaves: AMR-Wind was not built with Waves2AMR "
219 "support; associated wave data cannot be processed for relaxation "
222 amrex::ignore_unused(data, pp);
224 auto& wdata = data.
meta();
225 auto& info = data.
info();
230 "Current is specified as nonzero, but current is not yet "
231 "implemented for W2A Waves.");
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);
240 pp.query(
"HOS_domain_offset_x", wdata.
xlo0);
241 pp.query(
"HOS_domain_offset_y", wdata.
xlo1);
244 std::string fftw_planner_flag{
"estimate"};
245 pp.query(
"fftw_planner_flag", fftw_planner_flag);
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);
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);
260 bool file_exists{
false};
262 file_exists = wdata.c_rmodes.initialize(wdata.
modes_file,
true);
264 file_exists = wdata.r_rmodes.initialize(wdata.
modes_file,
false);
270 "Waves2AMR ReadInputsOp: modes file requested does not exist");
275 : wdata.r_rmodes.get_dtout();
297 amrex::Print() <<
"OceanWaves: Initializing Waves2AMR profile at time "
302 amrex::Real depth{0.0_rt};
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);
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();
317 wdata.
dx0 = wdata.c_rmodes.get_xlen() / wdata.
n0_sp;
318 wdata.
dx1 = wdata.c_rmodes.get_ylen() / wdata.
n1_sp;
320 depth = wdata.c_rmodes.get_depth();
322 wdata.
dimL = wdata.c_rmodes.get_L();
324 wdata.
n_wstop =
static_cast<int>(
325 (wdata.c_rmodes.get_Tstop() + 1.0e-8_rt) / wdata.
dt_modes);
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
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);
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();
354 wdata.
dx0 = wdata.r_rmodes.get_xlen() / (wdata.
n0_sp - 1);
356 wdata.r_rmodes.get_ylen() / amrex::max(wdata.
n1_sp - 1, 1);
358 depth = wdata.r_rmodes.get_depth();
360 wdata.
dimL = wdata.r_rmodes.get_L();
362 wdata.
n_wstop =
static_cast<int>(
363 (wdata.r_rmodes.get_Tstop() + 1.0e-8_rt) / wdata.
dt_modes);
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
383 <<
"WARNING (waves2amr_ops): available mode data exceeded, "
384 "resetting to beginning of mode data.\n";
388 if (std::abs(depth - (wdata.
zsl - prob_lo_input[2])) >
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";
399 modes_hosgrid::allocate_complex(wdata.
n0, wdata.
n1);
401 modes_hosgrid::allocate_complex(wdata.
n0, wdata.
n1);
403 modes_hosgrid::allocate_complex(wdata.
n0, wdata.
n1);
405 modes_hosgrid::allocate_complex(wdata.
n0, wdata.
n1);
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);
415 wdata.plan_out = modes_hosgrid::allocate_real(wdata.
n0, wdata.
n1);
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")) {
427 <<
"WARNING (waves2amr_ops): invalid fftw_planner_flag "
428 "specified; defaulting to estimate (FFTW_ESTIMATE).\n";
432 wdata.plan_vector.emplace_back(
433 modes_hosgrid::plan_ifftw(
434 wdata.
n0, wdata.
n1, wdata.c_eta_mptr, wdata.plan_out,
437 modes_hosgrid::plan_ifftw_nwt(
438 wdata.
n0, wdata.
n1, wdata.plan_vector, wdata.r_eta_mptr,
439 wdata.plan_out, plan_f);
444 int flag = interp_to_mfab::create_height_vector(
445 wdata.
hvec, nheights, dz0, wdata.
zsl, prob_lo_input[2], nh_above);
449 "Waves2AMR ReadInputsOp: create_height_vector error, failure "
451 std::to_string(flag));
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;
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]
474 << wdata.
hvec[wdata.
hvec.size() - 2] -
478 if (found_z_double) {
480 <<
" interp spacing reaches double that of surface at "
481 << z_double << std::endl;
483 amrex::Print() <<
" interp spacing remains less than double "
491 ? wdata.c_rmodes.get_data(
494 : wdata.r_rmodes.get_data(
502 <<
"WARNING (waves2amr_ops): end of mode data file "
503 "detected, resetting to beginning of mode data.\n";
507 ? wdata.c_rmodes.get_data(
510 : wdata.r_rmodes.get_data(
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 "
525 modes_hosgrid::copy_complex(
526 wdata.
n0, wdata.
n1, wdata.
c_mFS, wdata.c_eta_mptr);
528 static_cast<size_t>(wdata.
n0) *
static_cast<size_t>(wdata.
n1),
530 modes_hosgrid::populate_hos_eta(
531 wdata.c_rmodes, wdata.plan_vector, wdata.c_eta_mptr,
534 const auto n_hts = wdata.
hvec.size();
536 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
538 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
540 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
541 for (
int iht = 0; iht < static_cast<int>(n_hts); ++iht) {
543 amrex::Real ht = wdata.
hvec[iht];
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,
552 modes_hosgrid::copy_real(
553 wdata.
n0, wdata.
n1, wdata.
r_mFS, wdata.r_eta_mptr);
555 static_cast<size_t>(wdata.
n0) *
static_cast<size_t>(wdata.
n1),
557 modes_hosgrid::populate_hos_eta(
558 wdata.r_rmodes, wdata.plan_vector, wdata.r_eta_mptr,
561 const auto n_hts = wdata.
hvec.size();
563 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
565 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
567 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
568 for (
int iht = 0; iht < static_cast<int>(n_hts); ++iht) {
570 amrex::Real ht = wdata.
hvec[iht];
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,
585 "w2a_velocity", AMREX_SPACEDIM, 3, 2);
589 w2a_velocity.set_default_fillpatch_bc(data.
sim().
time());
602 const amrex::Geometry& geom,
603 bool multiphase_mode)
606#ifdef AMR_WIND_USE_W2A
607 auto& wdata = data.
meta();
609 auto& sim = data.
sim();
623 Field* levelset{
nullptr};
624 if (multiphase_mode) {
628 const auto& problo = geom.ProbLoArray();
629 const auto& probhi = geom.ProbHiArray();
630 const auto& dx = geom.CellSizeArray();
636 amrex::Vector<int> indvec(wdata.
hvec.size());
637 std::iota(indvec.begin(), indvec.end(), 0);
639 interp_to_mfab::interp_eta_to_levelset_multifab(
642 ow_levelset(level), problo, dx);
643 interp_to_mfab::interp_velocity_to_multifab(
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>();
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();
661 const amrex::Real gen_length = wdata.
gen_length;
664 const amrex::Real zero_sea_level = wdata.
zsl;
666 const bool has_beach = wdata.
has_beach && multiphase_mode;
667 const bool init_wave_field = wdata.
init_wave_field || !multiphase_mode;
669 amrex::MultiFab phi_tmp_fab(
670 velocity(level).boxArray(), velocity(level).DistributionMap(), 1,
672 amrex::MultiFab vel_tmp_fab(
673 velocity(level).boxArray(), velocity(level).DistributionMap(),
676 const auto& phi_tmp = phi_tmp_fab.arrays();
677 const auto& vel_tmp = vel_tmp_fab.arrays();
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];
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};
692 0.0_rt, 0.0_rt, 0.0_rt, zero_sea_level};
695 const auto bulk = init_wave_field ? wave_sol : quiescent;
696 const auto outlet = has_beach ? quiescent : wave_sol;
699 x, y, problo[0], gen_length, probhi[0], beach_length,
700 problo[1], probhi[1], y_length, wave_sol, bulk, outlet);
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];
707 if (multiphase_mode) {
708 phi_arrs[nbx](i, j, k) = phi_tmp[nbx](i, j, k);
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];
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];
724 amrex::Gpu::streamSynchronize();
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 =
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) {
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) /
758 amrex::Gpu::streamSynchronize();
760 amrex::ignore_unused(data, level, geom, multiphase_mode);
772#ifdef AMR_WIND_USE_W2A
773 auto& wdata = data.
meta();
774 auto& sim = data.
sim();
789 auto nlevels = sim.repo().num_active_levels();
790 auto geom = sim.mesh().Geom();
793 amrex::Real t_last = wdata.
t_last;
796 bool read_flag =
false;
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(
814 <<
"WARNING (waves2amr_ops): available mode data exceeded, "
815 "resetting to beginning of mode data.\n";
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
841 bool flag_xlo =
false;
842 bool flag_xhi =
false;
843 bool flag_ylo =
false;
844 bool flag_yhi =
false;
847 (interp_to_mfab::get_local_height_indices(
848 wdata.
indvec, wdata.
hvec, ow_velocity.vec_ptrs(),
854 (interp_to_mfab::check_lateral_overlap_lo(
855 wdata.
gen_length, 0, ow_velocity.vec_ptrs(), geom) ==
861 (interp_to_mfab::check_lateral_overlap_hi(
868 (interp_to_mfab::check_lateral_overlap_lo(
872 (interp_to_mfab::check_lateral_overlap_hi(
884 if (flag_z && (flag_xlo || flag_xhi || flag_ylo || flag_yhi)) {
889 static_cast<size_t>(wdata.
n0) *
890 static_cast<size_t>(wdata.
n1),
893 static_cast<size_t>(wdata.
n0 * wdata.
n1) *
896 static_cast<size_t>(wdata.
n0 * wdata.
n1) *
899 static_cast<size_t>(wdata.
n0 * wdata.
n1) *
910 amrex::ParallelDescriptor::ReduceIntMax(wdata.
n_offset);
915 if (double_data == 0) {
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) {
926 populate_fields_all_levels(
927 wdata, geom, w2a_lvs_old, w2a_vel_old, -1);
932 for (
int lev = nlevels - 1; lev > 0; --lev) {
934 w2a_vel_old(lev), w2a_vel_old(lev - 1), 0,
935 AMREX_SPACEDIM, sim.mesh().refRatio(lev - 1));
937 w2a_lvs_old(lev), w2a_lvs_old(lev - 1), 0, 1,
938 sim.mesh().refRatio(lev - 1));
941 w2a_vel_old.fillpatch(sim.time().new_time());
942 w2a_lvs_old.fillpatch(sim.time().new_time());
948 populate_fields_all_levels(
949 wdata, geom, w2a_levelset, w2a_velocity);
954 for (
int lev = nlevels - 1; lev > 0; --lev) {
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());
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());
967 time_interpolate_wave_fields(
968 w2a_lvs_old, w2a_vel_old, w2a_levelset, w2a_velocity, t_last, time,
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());
977 amrex::ignore_unused(data, time);
FieldRepo & repo()
Return the field repository.
Definition CFDSim.H:75
SimTime & time()
Return simulation time control.
Definition CFDSim.H:65
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 amrex::Real time)
Definition waves2amr_ops.H:769
Definition OceanWavesOps.H:17