/home/runner/work/amr-wind/amr-wind/amr-wind/ocean_waves/relaxation_zones/stokes_waves_ops.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/ocean_waves/relaxation_zones/stokes_waves_ops.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
stokes_waves_ops.H
Go to the documentation of this file.
1#ifndef STOKES_WAVES_OPS_H
2#define STOKES_WAVES_OPS_H
3
10#include "AMReX_REAL.H"
11
12using namespace amrex::literals;
13
15
16template <>
18{
20 StokesWaves::DataType& data, const ::amr_wind::utils::MultiParser& pp)
21 {
22 auto& wdata = data.meta();
23 auto& info = data.info();
24 relaxation_zones::read_inputs(wdata, info, pp);
25
26 // Get gravity, assume negative z
27 amrex::Vector<amrex::Real> gravity{0.0_rt, 0.0_rt, -9.81_rt};
28 amrex::ParmParse pp_incflo("incflo");
29 pp_incflo.queryarr("gravity", gravity);
30 wdata.g = -gravity[2];
31
32 if (wdata.current > constants::TIGHT_TOL) {
33 amrex::Abort(
34 "Current is specified as nonzero, but current is not yet "
35 "implemented for Stokes Waves.");
36 }
37
38 // Get wave attributes
39 pp.get("wave_height", wdata.wave_height);
40 pp.get("order", wdata.order);
41 if (pp.contains("wave_length")) {
42 pp.get("wave_length", wdata.wave_length);
43 } else {
44 // Wave length will be calculated using wave period
45 amrex::Real wave_period = 1.0_rt;
46 pp.get("wave_period", wave_period);
47 // User can specify a different order for wavelength calculation,
48 // e.g., an order of 1 will use the linear dispersion relation.
49 // Default is to use the same order as the waves themselves.
50 int order_lambda = wdata.order;
51 pp.query("stokes_wavelength_order", order_lambda);
52 // Get user-specified convergence criteria, if supplied
53 amrex::Real tol_lambda = 1.0e-10_rt;
54 int itmax_lambda = 40;
55 pp.query("stokes_wavelength_tolerance", tol_lambda);
56 pp.query("stokes_wavelength_iter_max", itmax_lambda);
57 wdata.wave_length = relaxation_zones::stokes_wave_length(
58 wave_period, wdata.water_depth, wdata.wave_height, order_lambda,
59 wdata.g, tol_lambda, itmax_lambda);
60 // Abort if wave length is negative
61 if (wdata.wave_length < 0.0_rt) {
62 if (wdata.wave_length == -1) {
63 amrex::Print() << "Stokes wave length too close to 0.\n";
64 }
65 if (wdata.wave_length == -2) {
66 amrex::Print() << "Stokes wave length is not a number.\n";
67 }
68 if (wdata.wave_length == -3) {
69 amrex::Print() << "Stokes wave length calculation used "
70 "maximum iterations before converging.\n";
71 }
72 amrex::Abort(
73 "Failed to properly calculate a wave length in "
74 "stokes_waves_ops.H");
75 } else {
76 amrex::Print() << "Stokes wave length calculated to be "
77 << wdata.wave_length << " m.\n";
78 }
79 }
80 pp.query("wave_phase_offset_radians", wdata.wave_phase_offset);
81 if (!pp.contains("wave_phase_offset_radians")) {
82 pp.query("wave_phase_offset_degrees", wdata.wave_phase_offset);
83 wdata.wave_phase_offset *=
84 static_cast<amrex::Real>(M_PI) / 180.0_rt;
85 } else if (pp.contains("wave_phase_offset_degrees")) {
86 amrex::Abort(
87 "ReadInputsOp<StokesWaves> : wave phase offset is specified in "
88 "both radians and degrees. Please use only one.");
89 }
90 }
91};
92
93template <>
95{
98 int level,
99 const amrex::Geometry& geom,
100 bool multiphase_mode)
101 {
102 const auto& wdata = data.meta();
103
104 auto& sim = data.sim();
105 Field* levelset{nullptr};
106 if (multiphase_mode) {
107 levelset = &sim.repo().get_field("levelset");
108 }
109 // cppcheck-suppress constVariableReference
110 auto& velocity = sim.repo().get_field("velocity");
111 const auto& problo = geom.ProbLoArray();
112 const auto& probhi = geom.ProbHiArray();
113 const auto& dx = geom.CellSizeArray();
114
115 const auto& vel = velocity(level).arrays();
116 const auto& phi_arrs = multiphase_mode
117 ? (*levelset)(level).arrays()
118 : amrex::MultiArray4<amrex::Real>();
119
120 const amrex::Real wave_height = wdata.wave_height;
121 const amrex::Real wave_length = wdata.wave_length;
122 const amrex::Real phase_offset = wdata.wave_phase_offset;
123 const amrex::Real water_depth = wdata.water_depth;
124 const amrex::Real zero_sea_level = wdata.zsl;
125 const amrex::Real gen_length = wdata.gen_length;
126 const amrex::Real beach_length = wdata.beach_length;
127 const amrex::Real g = wdata.g;
128 const int order = wdata.order;
129 const bool has_beach = wdata.has_beach && multiphase_mode;
130 const bool init_wave_field = wdata.init_wave_field || !multiphase_mode;
131
132 amrex::ParallelFor(
133 velocity(level), amrex::IntVect(3),
134 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
135 const amrex::Real x = problo[0] + (i + 0.5_rt) * dx[0];
136 const amrex::Real z = problo[2] + (k + 0.5_rt) * dx[2];
137
138 // Wave profile
139 amrex::Real eta_w{0.0_rt}, u_w{0.0_rt}, v_w{0.0_rt},
140 w_w{0.0_rt};
142 order, wave_length, water_depth, wave_height,
143 zero_sea_level, g, x, z, 0.0_rt, phase_offset, eta_w, u_w,
144 v_w, w_w);
145 const utils::WaveVec wave_sol{u_w, v_w, w_w, eta_w};
146
147 // Quiescent profile
148 const utils::WaveVec quiescent{
149 0.0_rt, 0.0_rt, 0.0_rt, zero_sea_level};
150
151 // Specify initial state for each region of domain
152 const auto bulk = init_wave_field ? wave_sol : quiescent;
153 const auto outlet = has_beach ? quiescent : wave_sol;
154
155 const auto local_profile = utils::harmonize_profiles_1d(
156 x, problo[0], gen_length, probhi[0], beach_length, wave_sol,
157 bulk, outlet);
158
159 const amrex::Real phi = local_profile[3] - z;
160 const amrex::Real cell_length_2D =
161 std::sqrt(dx[0] * dx[0] + dx[2] * dx[2]);
162 if (phi + cell_length_2D >= 0) {
163 vel[nbx](i, j, k, 0) = local_profile[0];
164 vel[nbx](i, j, k, 1) = local_profile[1];
165 vel[nbx](i, j, k, 2) = local_profile[2];
166 }
167 if (multiphase_mode) {
168 phi_arrs[nbx](i, j, k) = phi;
169 }
170 });
171 amrex::Gpu::streamSynchronize();
172 }
173};
174
175template <>
177{
178 void operator()(StokesWaves::DataType& data, const amrex::Real time)
179 {
180 const auto& wdata = data.meta();
181
182 auto& sim = data.sim();
183
184 // cppcheck-suppress constVariableReference
185 auto& ow_levelset = sim.repo().get_field("ow_levelset");
186 // cppcheck-suppress constVariableReference
187 auto& ow_velocity = sim.repo().get_field("ow_velocity");
188
189 auto nlevels = sim.repo().num_active_levels();
190 auto geom = sim.mesh().Geom();
191
192 for (int lev = 0; lev < nlevels; ++lev) {
193 const auto& problo = geom[lev].ProbLoArray();
194 const auto& dx = geom[lev].CellSizeArray();
195
196 const auto& phi = ow_levelset(lev).arrays();
197 const auto& vel = ow_velocity(lev).arrays();
198
199 const amrex::Real wave_height = wdata.wave_height;
200 const amrex::Real wave_length = wdata.wave_length;
201 const amrex::Real phase_offset = wdata.wave_phase_offset;
202 const amrex::Real water_depth = wdata.water_depth;
203 const amrex::Real zero_sea_level = wdata.zsl;
204 const amrex::Real g = wdata.g;
205 const int order = wdata.order;
206
207 amrex::ParallelFor(
208 ow_velocity(lev), amrex::IntVect(3),
209 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
210 const amrex::Real x = amrex::max<amrex::Real>(
211 problo[0], problo[0] + (i + 0.5_rt) * dx[0]);
212 const amrex::Real z = problo[2] + (k + 0.5_rt) * dx[2];
213
214 amrex::Real eta{0.0_rt}, u_w{0.0_rt}, v_w{0.0_rt},
215 w_w{0.0_rt};
216
218 order, wave_length, water_depth, wave_height,
219 zero_sea_level, g, x, z, time, phase_offset, eta, u_w,
220 v_w, w_w);
221
222 phi[nbx](i, j, k) = eta - z;
223 const amrex::Real cell_length_2D =
224 std::sqrt(dx[0] * dx[0] + dx[2] * dx[2]);
225 if (phi[nbx](i, j, k) + cell_length_2D >= 0) {
226 // Wave velocity within a cell of interface
227 vel[nbx](i, j, k, 0) = u_w;
228 vel[nbx](i, j, k, 1) = v_w;
229 vel[nbx](i, j, k, 2) = w_w;
230 } else {
231 vel[nbx](i, j, k, 0) = 0.;
232 vel[nbx](i, j, k, 1) = 0.;
233 vel[nbx](i, j, k, 2) = 0.;
234 }
235 });
236 }
237 amrex::Gpu::streamSynchronize();
238 }
239};
240
241} // namespace amr_wind::ocean_waves::ops
242
243#endif /* STOKES_WAVES_OPS_H */
FieldRepo & repo()
Return the field repository.
Definition CFDSim.H:75
Definition Field.H:116
FieldRepo & repo() const
FieldRepo instance that manages this field.
Definition Field.H:159
Field & get_field(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:152
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
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:17
Definition OceanWavesOps.H:8
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void stokes_waves(int stokes_order, amrex::Real wavelength, amrex::Real water_depth, amrex::Real wave_height, amrex::Real zsl, amrex::Real g, amrex::Real x, amrex::Real z, amrex::Real time, amrex::Real phase_offset, amrex::Real &eta, amrex::Real &u_w, amrex::Real &v_w, amrex::Real &w_w)
Definition stokes_waves_K.H:283
void read_inputs(RelaxZonesBaseData &wdata, OceanWavesInfo &, const ::amr_wind::utils::MultiParser &pp)
Definition relaxation_zones_ops.cpp:16
AMREX_FORCE_INLINE amrex::Real stokes_wave_length(const amrex::Real T, const amrex::Real d, const amrex::Real H, const int order, const amrex::Real g, const amrex::Real tol, const int iter_max)
Definition stokes_waves_K.H:26
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:74
amrex::GpuArray< amrex::Real, 4 > WaveVec
Definition wave_utils_K.H:13
Definition StokesWaves.H:20
OceanWavesDataHolder< StokesWaves > DataType
Definition StokesWaves.H:23
void operator()(StokesWaves::DataType &data, int level, const amrex::Geometry &geom, bool multiphase_mode)
Definition stokes_waves_ops.H:96
Definition OceanWavesOps.H:14
void operator()(StokesWaves::DataType &data, const ::amr_wind::utils::MultiParser &pp)
Definition stokes_waves_ops.H:19
Definition OceanWavesOps.H:11
void operator()(StokesWaves::DataType &data, const amrex::Real time)
Definition stokes_waves_ops.H:178