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