13#include "AMReX_MultiFabUtil.H"
15#ifdef AMR_WIND_USE_W2A
17int evaluate_read_resize(
24 const amrex::Real wtinit,
25 const amrex::Real wdt,
26 const amrex::Real time)
31 if (new_ntime != ntime) {
35 if (new_ntime > ntime + 1) {
42 wtime = new_ntime * wdt - wtinit;
44 if (double_data == 1 && std::abs(wtime - time) <= 1e-10) {
50 if (t_last < -1e-10) {
57 if (std::abs(wtime - time) > 1e-10) {
65 }
else if (std::abs(t_last) < 1e-10) {
80int update_offset_timestep(
const int ntime,
const int n0)
89void populate_fields_all_levels(
91 amrex::Vector<amrex::Geometry>& geom_all,
100 ? wdata.c_rmodes.get_data(
103 : wdata.r_rmodes.get_data(
111 amrex::Print() <<
"WARNING (waves2amr_ops): end of mode data file "
112 "detected, resetting to beginning of mode data.\n";
115 ? wdata.c_rmodes.get_data(
118 : wdata.r_rmodes.get_data(
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 "
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,
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,
146 for (
int iht = 0; iht < wdata.
indvec.size(); ++iht) {
148 amrex::Real ht = wdata.
hvec[wdata.
indvec[iht]];
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,
155 iht * wdata.
n0 * wdata.
n1);
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,
162 iht * wdata.
n0 * wdata.
n1);
167 interp_to_mfab::interp_eta_to_levelset_field(
171 interp_to_mfab::interp_velocity_to_field(
177void time_interpolate_wave_fields(
182 const amrex::Real t_old,
183 const amrex::Real t_new,
184 const amrex::Real t_future)
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,
192 vel_field, a, vel_field, 0, b, future_vel_field, 0, 0, AMREX_SPACEDIM,
207 const ::amr_wind::utils::MultiParser& pp)
210#ifndef AMR_WIND_USE_W2A
213 "ocean_waves/W2AWaves: AMR-Wind was not built with Waves2AMR "
214 "support; associated wave data cannot be processed for relaxation "
217 amrex::ignore_unused(data, pp);
219 auto& wdata = data.
meta();
220 auto& info = data.
info();
225 "Current is specified as nonzero, but current is not yet "
226 "implemented for W2A Waves.");
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);
235 pp.query(
"HOS_domain_offset_x", wdata.
xlo0);
236 pp.query(
"HOS_domain_offset_y", wdata.
xlo1);
239 std::string fftw_planner_flag{
"estimate"};
240 pp.query(
"fftw_planner_flag", fftw_planner_flag);
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);
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);
255 bool file_exists{
false};
257 file_exists = wdata.c_rmodes.initialize(wdata.
modes_file,
true);
259 file_exists = wdata.r_rmodes.initialize(wdata.
modes_file,
false);
265 "Waves2AMR ReadInputsOp: modes file requested does not exist");
270 : wdata.r_rmodes.get_dtout();
292 amrex::Print() <<
"OceanWaves: Initializing Waves2AMR profile at time "
297 amrex::Real depth{0.0};
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);
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();
312 wdata.
dx0 = wdata.c_rmodes.get_xlen() / wdata.
n0_sp;
313 wdata.
dx1 = wdata.c_rmodes.get_ylen() / wdata.
n1_sp;
315 depth = wdata.c_rmodes.get_depth();
317 wdata.
dimL = wdata.c_rmodes.get_L();
320 (int)((wdata.c_rmodes.get_Tstop() + 1e-8) / wdata.
dt_modes);
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
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);
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();
349 wdata.
dx0 = wdata.r_rmodes.get_xlen() / (wdata.
n0_sp - 1);
351 wdata.r_rmodes.get_ylen() / amrex::max(wdata.
n1_sp - 1, 1);
353 depth = wdata.r_rmodes.get_depth();
355 wdata.
dimL = wdata.r_rmodes.get_L();
358 (int)((wdata.r_rmodes.get_Tstop() + 1e-8) / wdata.
dt_modes);
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
378 <<
"WARNING (waves2amr_ops): available mode data exceeded, "
379 "resetting to beginning of mode data.\n";
383 if (std::abs(depth - (wdata.
zsl - prob_lo_input[2])) > 1e-3 * depth) {
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";
393 modes_hosgrid::allocate_complex(wdata.
n0, wdata.
n1);
395 modes_hosgrid::allocate_complex(wdata.
n0, wdata.
n1);
397 modes_hosgrid::allocate_complex(wdata.
n0, wdata.
n1);
399 modes_hosgrid::allocate_complex(wdata.
n0, wdata.
n1);
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);
409 wdata.plan_out = modes_hosgrid::allocate_real(wdata.
n0, wdata.
n1);
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")) {
421 <<
"WARNING (waves2amr_ops): invalid fftw_planner_flag "
422 "specified; defaulting to estimate (FFTW_ESTIMATE).\n";
426 wdata.plan_vector.emplace_back(
427 modes_hosgrid::plan_ifftw(
428 wdata.
n0, wdata.
n1, wdata.c_eta_mptr, wdata.plan_out,
431 modes_hosgrid::plan_ifftw_nwt(
432 wdata.
n0, wdata.
n1, wdata.plan_vector, wdata.r_eta_mptr,
433 wdata.plan_out, plan_f);
438 int flag = interp_to_mfab::create_height_vector(
439 wdata.
hvec, nheights, dz0, wdata.
zsl, prob_lo_input[2], nh_above);
443 "Waves2AMR ReadInputsOp: create_height_vector error, failure "
445 std::to_string(flag));
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;
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]
468 << wdata.
hvec[wdata.
hvec.size() - 2] -
472 if (found_z_double) {
474 <<
" interp spacing reaches double that of surface at "
475 << z_double << std::endl;
477 amrex::Print() <<
" interp spacing remains less than double "
485 ? wdata.c_rmodes.get_data(
488 : wdata.r_rmodes.get_data(
496 <<
"WARNING (waves2amr_ops): end of mode data file "
497 "detected, resetting to beginning of mode data.\n";
501 ? wdata.c_rmodes.get_data(
504 : wdata.r_rmodes.get_data(
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 "
519 modes_hosgrid::copy_complex(
520 wdata.
n0, wdata.
n1, wdata.
c_mFS, wdata.c_eta_mptr);
522 static_cast<size_t>(wdata.
n0) *
static_cast<size_t>(wdata.
n1),
524 modes_hosgrid::populate_hos_eta(
525 wdata.c_rmodes, wdata.plan_vector, wdata.c_eta_mptr,
528 const auto n_hts = wdata.
hvec.size();
530 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
532 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
534 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
535 for (
int iht = 0; iht < static_cast<int>(n_hts); ++iht) {
537 amrex::Real ht = wdata.
hvec[iht];
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,
546 modes_hosgrid::copy_real(
547 wdata.
n0, wdata.
n1, wdata.
r_mFS, wdata.r_eta_mptr);
549 static_cast<size_t>(wdata.
n0) *
static_cast<size_t>(wdata.
n1),
551 modes_hosgrid::populate_hos_eta(
552 wdata.r_rmodes, wdata.plan_vector, wdata.r_eta_mptr,
555 const auto n_hts = wdata.
hvec.size();
557 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
559 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
561 static_cast<size_t>(wdata.
n0 * wdata.
n1) * n_hts);
562 for (
int iht = 0; iht < static_cast<int>(n_hts); ++iht) {
564 amrex::Real ht = wdata.
hvec[iht];
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,
579 "w2a_velocity", AMREX_SPACEDIM, 3, 2);
583 w2a_velocity.set_default_fillpatch_bc(data.
sim().
time());
596 const amrex::Geometry& geom,
597 bool multiphase_mode)
600#ifdef AMR_WIND_USE_W2A
601 auto& wdata = data.
meta();
603 auto& sim = data.
sim();
617 Field* levelset{
nullptr};
618 if (multiphase_mode) {
622 const auto& problo = geom.ProbLoArray();
623 const auto& probhi = geom.ProbHiArray();
624 const auto& dx = geom.CellSizeArray();
630 amrex::Vector<int> indvec(wdata.
hvec.size());
631 std::iota(indvec.begin(), indvec.end(), 0);
633 interp_to_mfab::interp_eta_to_levelset_multifab(
636 ow_levelset(level), problo, dx);
637 interp_to_mfab::interp_velocity_to_multifab(
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>();
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();
655 const amrex::Real gen_length = wdata.
gen_length;
657 const amrex::Real zero_sea_level = wdata.
zsl;
659 const bool has_beach = wdata.
has_beach && multiphase_mode;
660 const bool init_wave_field = wdata.
init_wave_field || !multiphase_mode;
662 amrex::MultiFab phi_tmp_fab(
663 velocity(level).boxArray(), velocity(level).DistributionMap(), 1,
665 amrex::MultiFab vel_tmp_fab(
666 velocity(level).boxArray(), velocity(level).DistributionMap(),
669 const auto& phi_tmp = phi_tmp_fab.arrays();
670 const auto& vel_tmp = vel_tmp_fab.arrays();
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];
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};
686 const auto bulk = init_wave_field ? wave_sol : quiescent;
687 const auto outlet = has_beach ? quiescent : wave_sol;
690 x, problo[0], gen_length, probhi[0], beach_length, wave_sol,
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];
698 if (multiphase_mode) {
699 phi_arrs[nbx](i, j, k) = phi_tmp[nbx](i, j, k);
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];
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];
715 amrex::Gpu::streamSynchronize();
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 =
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) {
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) /
748 amrex::Gpu::streamSynchronize();
750 amrex::ignore_unused(data, level, geom, multiphase_mode);
762#ifdef AMR_WIND_USE_W2A
763 auto& wdata = data.
meta();
764 auto& sim = data.
sim();
779 auto nlevels = sim.repo().num_active_levels();
780 auto geom = sim.mesh().Geom();
783 amrex::Real t_last = wdata.
t_last;
786 bool read_flag =
false;
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(
804 <<
"WARNING (waves2amr_ops): available mode data exceeded, "
805 "resetting to beginning of mode data.\n";
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
831 bool flag_xlo =
false;
832 bool flag_xhi =
false;
833 bool flag_ylo =
false;
834 bool flag_yhi =
false;
837 (interp_to_mfab::get_local_height_indices(
838 wdata.
indvec, wdata.
hvec, ow_velocity.vec_ptrs(),
844 (interp_to_mfab::check_lateral_overlap_lo(
845 wdata.
gen_length, 0, ow_velocity.vec_ptrs(), geom) ==
851 (interp_to_mfab::check_lateral_overlap_hi(
858 (interp_to_mfab::check_lateral_overlap_lo(
862 (interp_to_mfab::check_lateral_overlap_hi(
874 if (flag_z && (flag_xlo || flag_xhi || flag_ylo || flag_yhi)) {
879 static_cast<size_t>(wdata.
n0) *
880 static_cast<size_t>(wdata.
n1),
883 static_cast<size_t>(wdata.
n0 * wdata.
n1) *
886 static_cast<size_t>(wdata.
n0 * wdata.
n1) *
889 static_cast<size_t>(wdata.
n0 * wdata.
n1) *
900 amrex::ParallelDescriptor::ReduceIntMax(wdata.
n_offset);
905 if (double_data == 0) {
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) {
916 populate_fields_all_levels(
917 wdata, geom, w2a_lvs_old, w2a_vel_old, -1);
922 for (
int lev = nlevels - 1; lev > 0; --lev) {
924 w2a_vel_old(lev), w2a_vel_old(lev - 1), 0,
925 AMREX_SPACEDIM, sim.mesh().refRatio(lev - 1));
927 w2a_lvs_old(lev), w2a_lvs_old(lev - 1), 0, 1,
928 sim.mesh().refRatio(lev - 1));
931 w2a_vel_old.fillpatch(sim.time().new_time());
932 w2a_lvs_old.fillpatch(sim.time().new_time());
938 populate_fields_all_levels(
939 wdata, geom, w2a_levelset, w2a_velocity);
944 for (
int lev = nlevels - 1; lev > 0; --lev) {
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());
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());
957 time_interpolate_wave_fields(
958 w2a_lvs_old, w2a_vel_old, w2a_levelset, w2a_velocity, t_last, time,
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());
967 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: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 amrex::Real time)
Definition waves2amr_ops.H:759
Definition OceanWavesOps.H:17