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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/ocean_waves/relaxation_zones/linear_waves_ops.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
linear_waves_ops.H
Go to the documentation of this file.
1#ifndef LINEAR_WAVES_OPS_H
2#define LINEAR_WAVES_OPS_H
3
9
11
12template <>
14{
16 LinearWaves::DataType& data, const ::amr_wind::utils::MultiParser& pp)
17 {
18 auto& wdata = data.meta();
19 auto& info = data.info();
20 relaxation_zones::read_inputs(wdata, info, pp);
21 // Get gravity, assume negative z
22 amrex::Vector<amrex::Real> gravity{0.0, 0.0, -9.81};
23 amrex::ParmParse pp_incflo("incflo");
24 pp_incflo.queryarr("gravity", gravity);
25 wdata.g = -gravity[2];
26
27 pp.get("wave_length", wdata.wave_length);
28 pp.get("wave_height", wdata.wave_height);
29 pp.query("wave_phase_offset_radians", wdata.wave_phase_offset);
30 if (!pp.contains("wave_phase_offset_radians")) {
31 pp.query("wave_phase_offset_degrees", wdata.wave_phase_offset);
32 wdata.wave_phase_offset *= M_PI / 180.;
33 } else if (pp.contains("wave_phase_offset_degrees")) {
34 amrex::Abort(
35 "ReadInputsOp<LinearWaves> : wave phase offset is specified in "
36 "both radians and degrees. Please use only one.");
37 }
38 }
39};
40
41template <>
43{
46 int level,
47 const amrex::Geometry& geom,
48 bool multiphase_mode)
49 {
50 const auto& wdata = data.meta();
51
52 auto& sim = data.sim();
53
54 Field* levelset{nullptr};
55 if (multiphase_mode) {
56 levelset = &sim.repo().get_field("levelset");
57 }
58 // cppcheck-suppress constVariableReference
59 auto& velocity = sim.repo().get_field("velocity");
60
61 const auto& problo = geom.ProbLoArray();
62 const auto& probhi = geom.ProbHiArray();
63 const auto& dx = geom.CellSizeArray();
64
65 const auto& vel = velocity(level).arrays();
66 const auto& phi_arrs = multiphase_mode
67 ? (*levelset)(level).arrays()
68 : amrex::MultiArray4<amrex::Real>();
69
70 const amrex::Real zero_sea_level = wdata.zsl;
71 const amrex::Real gen_length = wdata.gen_length;
72 const amrex::Real beach_length = wdata.beach_length;
73 const amrex::Real g = wdata.g;
74 const bool has_beach = wdata.has_beach && multiphase_mode;
75 const bool init_wave_field = wdata.init_wave_field || !multiphase_mode;
76
77 const amrex::Real wave_height = wdata.wave_height;
78 const amrex::Real wave_length = wdata.wave_length;
79 const amrex::Real phase_offset = wdata.wave_phase_offset;
80 const amrex::Real water_depth = wdata.water_depth;
81 const amrex::Real current = wdata.current;
82
83 amrex::ParallelFor(
84 velocity(level), amrex::IntVect(3),
85 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
86 const amrex::Real x = problo[0] + (i + 0.5) * dx[0];
87 const amrex::Real z = problo[2] + (k + 0.5) * dx[2];
88
89 const amrex::Real wave_number = 2. * M_PI / wave_length;
90 const amrex::Real omega = std::pow(
91 wave_number * g * std::tanh(wave_number * water_depth),
92 0.5);
93 const amrex::Real phase = wave_number * x - phase_offset;
94
95 // Wave profile
96 const amrex::Real eta_w =
97 wave_height / 2.0 * std::cos(phase) + zero_sea_level;
98 const amrex::Real u_w =
99 omega * wave_height / 2.0 *
100 std::cosh(
101 wave_number * (z - zero_sea_level + water_depth)) /
102 std::sinh(wave_number * water_depth) * std::cos(phase);
103 const amrex::Real v_w = 0.0;
104 const amrex::Real w_w =
105 omega * wave_height / 2.0 *
106 std::sinh(
107 wave_number * (z - zero_sea_level + water_depth)) /
108 std::sinh(wave_number * water_depth) * std::sin(phase);
109 const utils::WaveVec wave_sol{u_w, v_w, w_w, eta_w};
110
111 // Quiescent profile
112 const utils::WaveVec quiescent{0.0, 0.0, 0.0, zero_sea_level};
113
114 // Specify initial state for each region of domain
115 const auto bulk = init_wave_field ? wave_sol : quiescent;
116 const auto outlet = has_beach ? quiescent : wave_sol;
117
118 const auto local_profile = utils::harmonize_profiles_1d(
119 x, problo[0], gen_length, probhi[0], beach_length, wave_sol,
120 bulk, outlet);
121
122 const amrex::Real phi = local_profile[3] - z;
123 const amrex::Real cell_length_2D =
124 std::sqrt(dx[0] * dx[0] + dx[2] * dx[2]);
125
126 if (phi + cell_length_2D >= 0) {
127 vel[nbx](i, j, k, 0) = current + local_profile[0];
128 vel[nbx](i, j, k, 1) = local_profile[1];
129 vel[nbx](i, j, k, 2) = local_profile[2];
130 }
131 if (multiphase_mode) {
132 phi_arrs[nbx](i, j, k) = phi;
133 }
134 });
135 amrex::Gpu::streamSynchronize();
136 }
137};
138
139template <>
141{
142 void operator()(LinearWaves::DataType& data, const amrex::Real time)
143 {
144 const auto& wdata = data.meta();
145
146 auto& sim = data.sim();
147
148 // cppcheck-suppress constVariableReference
149 auto& ow_levelset = sim.repo().get_field("ow_levelset");
150 // cppcheck-suppress constVariableReference
151 auto& ow_velocity = sim.repo().get_field("ow_velocity");
152
153 auto nlevels = sim.repo().num_active_levels();
154 auto geom = sim.mesh().Geom();
155
156 for (int lev = 0; lev < nlevels; ++lev) {
157 const auto& problo = geom[lev].ProbLoArray();
158 const auto& dx = geom[lev].CellSizeArray();
159
160 const auto& phi = ow_levelset(lev).arrays();
161 const auto& vel = ow_velocity(lev).arrays();
162
163 const amrex::Real wave_height = wdata.wave_height;
164 const amrex::Real wave_length = wdata.wave_length;
165 const amrex::Real phase_offset = wdata.wave_phase_offset;
166 const amrex::Real water_depth = wdata.water_depth;
167 const amrex::Real zero_sea_level = wdata.zsl;
168 const amrex::Real g = wdata.g;
169 const amrex::Real current = wdata.current;
170
171 amrex::ParallelFor(
172 ow_velocity(lev), amrex::IntVect(3),
173 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
174 const amrex::Real x =
175 amrex::max(problo[0], problo[0] + (i + 0.5) * dx[0]);
176 const amrex::Real z = problo[2] + (k + 0.5) * dx[2];
177
178 const amrex::Real wave_number = 2. * M_PI / wave_length;
179 const amrex::Real omega = std::pow(
180 wave_number * g * std::tanh(wave_number * water_depth),
181 0.5);
182 const amrex::Real phase =
183 wave_number * (x - current * time) - omega * time -
184 phase_offset;
185
186 const amrex::Real eta =
187 wave_height / 2.0 * std::cos(phase) + zero_sea_level;
188
189 phi[nbx](i, j, k) = eta - z;
190
191 const amrex::Real cell_length_2D =
192 std::sqrt(dx[0] * dx[0] + dx[2] * dx[2]);
193 if (phi[nbx](i, j, k) + cell_length_2D >= 0) {
194 vel[nbx](i, j, k, 0) =
195 omega * wave_height / 2.0 *
196 std::cosh(
197 wave_number *
198 (z - zero_sea_level + water_depth)) /
199 std::sinh(wave_number * water_depth) *
200 std::cos(phase) +
201 current;
202 vel[nbx](i, j, k, 1) = 0.0;
203 vel[nbx](i, j, k, 2) =
204 omega * wave_height / 2.0 *
205 std::sinh(
206 wave_number *
207 (z - zero_sea_level + water_depth)) /
208 std::sinh(wave_number * water_depth) *
209 std::sin(phase);
210 } else {
211 vel[nbx](i, j, k, 0) = 0.;
212 vel[nbx](i, j, k, 1) = 0.;
213 vel[nbx](i, j, k, 2) = 0.;
214 }
215 });
216 }
217 amrex::Gpu::streamSynchronize();
218 }
219};
220
221} // namespace amr_wind::ocean_waves::ops
222
223#endif /* LINEAR_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
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
Definition LinearWaves.H:16
void operator()(LinearWaves::DataType &data, int level, const amrex::Geometry &geom, bool multiphase_mode)
Definition linear_waves_ops.H:44
Definition OceanWavesOps.H:14
void operator()(LinearWaves::DataType &data, const ::amr_wind::utils::MultiParser &pp)
Definition linear_waves_ops.H:15
Definition OceanWavesOps.H:11
void operator()(LinearWaves::DataType &data, const amrex::Real time)
Definition linear_waves_ops.H:142