Skip to content

Commit

Permalink
expand area of valid target wave velocity for all current wave types
Browse files Browse the repository at this point in the history
  • Loading branch information
mbkuhn committed Oct 9, 2024
1 parent 3ad9250 commit 33f9aed
Show file tree
Hide file tree
Showing 3 changed files with 22 additions and 7 deletions.
12 changes: 9 additions & 3 deletions amr-wind/ocean_waves/relaxation_zones/linear_waves_ops.H
Original file line number Diff line number Diff line change
Expand Up @@ -73,7 +73,9 @@ struct InitDataOp<LinearWaves>

phi(i, j, k) = eta - zc;

if (phi(i, j, k) + 0.5 * dx[2] >= 0) {
const amrex::Real cell_length_2D =
std::sqrt(std::pow(dx[0], 2) + std::pow(dx[2], 2));
if (phi(i, j, k) + cell_length_2D >= 0) {
vel(i, j, k, 0) =
omega * wave_height / 2.0 *
std::cosh(
Expand All @@ -95,7 +97,9 @@ struct InitDataOp<LinearWaves>
gbx3, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
const amrex::Real z = problo[2] + (k + 0.5) * dx[2];
phi(i, j, k) = zsl - z;
if (phi(i, j, k) + 0.5 * dx[2] >= 0) {
const amrex::Real cell_length_2D =
std::sqrt(std::pow(dx[0], 2) + std::pow(dx[2], 2));
if (phi(i, j, k) + cell_length_2D >= 0) {
vel(i, j, k, 0) = 0.0;
vel(i, j, k, 1) = 0.0;
vel(i, j, k, 2) = 0.0;
Expand Down Expand Up @@ -157,7 +161,9 @@ struct UpdateRelaxZonesOp<LinearWaves>

phi(i, j, k) = eta - zc;

if (phi(i, j, k) + 0.5 * dx[2] >= 0) {
const amrex::Real cell_length_2D =
std::sqrt(std::pow(dx[0], 2) + std::pow(dx[2], 2));
if (phi(i, j, k) + cell_length_2D >= 0) {
vel(i, j, k, 0) =
omega * wave_height / 2.0 *
std::cosh(
Expand Down
13 changes: 10 additions & 3 deletions amr-wind/ocean_waves/relaxation_zones/stokes_waves_ops.H
Original file line number Diff line number Diff line change
Expand Up @@ -109,7 +109,9 @@ struct InitDataOp<StokesWaves>
zero_sea_level, g, x, z, 0.0, eta, u_w, v_w, w_w);

phi(i, j, k) = eta - z;
if (phi(i, j, k) + 0.5 * dx[2] >= 0) {
const amrex::Real cell_length_2D =
std::sqrt(std::pow(dx[0], 2) + std::pow(dx[2], 2));
if (phi(i, j, k) + cell_length_2D >= 0) {
vel(i, j, k, 0) = u_w;
vel(i, j, k, 1) = v_w;
vel(i, j, k, 2) = w_w;
Expand All @@ -121,7 +123,9 @@ struct InitDataOp<StokesWaves>
gbx3, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
const amrex::Real z = problo[2] + (k + 0.5) * dx[2];
phi(i, j, k) = zero_sea_level - z;
if (phi(i, j, k) + 0.5 * dx[2] >= 0) {
const amrex::Real cell_length_2D =
std::sqrt(std::pow(dx[0], 2) + std::pow(dx[2], 2));
if (phi(i, j, k) + cell_length_2D >= 0) {
vel(i, j, k, 0) = 0.0;
vel(i, j, k, 1) = 0.0;
vel(i, j, k, 2) = 0.0;
Expand Down Expand Up @@ -178,7 +182,10 @@ struct UpdateRelaxZonesOp<StokesWaves>
zero_sea_level, g, x, z, time, eta, u_w, v_w, w_w);

phi(i, j, k) = eta - z;
if (phi(i, j, k) + 0.5 * dx[2] >= 0) {
const amrex::Real cell_length_2D =
std::sqrt(std::pow(dx[0], 2) + std::pow(dx[2], 2));
if (phi(i, j, k) + cell_length_2D >= 0) {
// Wave velocity within a cell of interface
vel(i, j, k, 0) = u_w;
vel(i, j, k, 1) = v_w;
vel(i, j, k, 2) = w_w;
Expand Down
4 changes: 3 additions & 1 deletion amr-wind/ocean_waves/relaxation_zones/waves2amr_ops.H
Original file line number Diff line number Diff line change
Expand Up @@ -90,7 +90,9 @@ void postprocess_velocity_mfab_liquid(
amrex::ParallelFor(
gbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
// Set velocity to zero if no liquid present
if (phi(i, j, k) + 0.5 * dx[2] < 0.0) {
const amrex::Real cell_length_2D =
std::sqrt(std::pow(dx[0], 2) + std::pow(dx[2], 2));
if (phi(i, j, k) + cell_length_2D >= 0) {
vel(i, j, k, 0) = 0.0;
vel(i, j, k, 1) = 0.0;
vel(i, j, k, 2) = 0.0;
Expand Down

0 comments on commit 33f9aed

Please sign in to comment.