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