/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
10
12
13template <>
15{
17 LinearWaves::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 // Get gravity from MultiPhase physics, assume negative z
23 wdata.g = -data.sim().physics_manager().get<MultiPhase>().gravity()[2];
24
25 pp.get("wave_length", wdata.wave_length);
26 pp.get("wave_height", wdata.wave_height);
27 pp.query("wave_phase_offset_radians", wdata.wave_phase_offset);
28 if (!pp.contains("wave_phase_offset_radians")) {
29 pp.query("wave_phase_offset_degrees", wdata.wave_phase_offset);
30 wdata.wave_phase_offset *= M_PI / 180.;
31 } else if (pp.contains("wave_phase_offset_degrees")) {
32 amrex::Abort(
33 "ReadInputsOp<LinearWaves> : wave phase offset is specified in "
34 "both radians and degrees. Please use only one.");
35 }
36 }
37};
38
39template <>
41{
43 LinearWaves::DataType& data, int level, const amrex::Geometry& geom)
44 {
45 const auto& wdata = data.meta();
46
47 auto& sim = data.sim();
48
49 // cppcheck-suppress constVariableReference
50 auto& m_levelset = sim.repo().get_field("levelset");
51 // cppcheck-suppress constVariableReference
52 auto& m_velocity = sim.repo().get_field("velocity");
53
54 const auto& problo = geom.ProbLoArray();
55 const auto& probhi = geom.ProbHiArray();
56 const auto& dx = geom.CellSizeArray();
57
58 const auto& phi = m_levelset(level).arrays();
59 const auto& vel = m_velocity(level).arrays();
60
61 const amrex::Real zero_sea_level = wdata.zsl;
62 const amrex::Real gen_length = wdata.gen_length;
63 const amrex::Real beach_length = wdata.beach_length;
64 const amrex::Real g = wdata.g;
65 const bool has_beach = wdata.has_beach;
66 const bool init_wave_field = wdata.init_wave_field;
67
68 const amrex::Real wave_height = wdata.wave_height;
69 const amrex::Real wave_length = wdata.wave_length;
70 const amrex::Real phase_offset = wdata.wave_phase_offset;
71 const amrex::Real water_depth = wdata.water_depth;
72
73 amrex::ParallelFor(
74 m_velocity(level), amrex::IntVect(3),
75 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
76 const amrex::Real x = problo[0] + (i + 0.5) * dx[0];
77 const amrex::Real z = problo[2] + (k + 0.5) * dx[2];
78
79 const amrex::Real wave_number = 2. * M_PI / wave_length;
80 const amrex::Real omega = std::pow(
81 wave_number * g * std::tanh(wave_number * water_depth),
82 0.5);
83 const amrex::Real phase = wave_number * x - phase_offset;
84
85 // Wave profile
86 const amrex::Real eta_w =
87 wave_height / 2.0 * std::cos(phase) + zero_sea_level;
88 const amrex::Real u_w =
89 omega * wave_height / 2.0 *
90 std::cosh(
91 wave_number * (z - zero_sea_level + water_depth)) /
92 std::sinh(wave_number * water_depth) * std::cos(phase);
93 const amrex::Real v_w = 0.0;
94 const amrex::Real w_w =
95 omega * wave_height / 2.0 *
96 std::sinh(
97 wave_number * (z - zero_sea_level + water_depth)) /
98 std::sinh(wave_number * water_depth) * std::sin(phase);
99 const utils::WaveVec wave_sol{u_w, v_w, w_w, eta_w};
100
101 // Quiescent profile
102 const utils::WaveVec quiescent{0.0, 0.0, 0.0, zero_sea_level};
103
104 // Specify initial state for each region of domain
105 const auto bulk = init_wave_field ? wave_sol : quiescent;
106 const auto outlet = has_beach ? quiescent : wave_sol;
107
108 const auto local_profile = utils::harmonize_profiles_1d(
109 x, problo[0], gen_length, probhi[0], beach_length, wave_sol,
110 bulk, outlet);
111
112 phi[nbx](i, j, k) = local_profile[3] - z;
113 const amrex::Real cell_length_2D =
114 std::sqrt(dx[0] * dx[0] + dx[2] * dx[2]);
115 if (phi[nbx](i, j, k) + cell_length_2D >= 0) {
116 vel[nbx](i, j, k, 0) = local_profile[0];
117 vel[nbx](i, j, k, 1) = local_profile[1];
118 vel[nbx](i, j, k, 2) = local_profile[2];
119 }
120 });
121 amrex::Gpu::synchronize();
122 }
123};
124
125template <>
127{
129 {
130 const auto& wdata = data.meta();
131
132 auto& sim = data.sim();
133
134 const auto& time = sim.time().new_time();
135 // cppcheck-suppress constVariableReference
136 auto& m_ow_levelset = sim.repo().get_field("ow_levelset");
137 // cppcheck-suppress constVariableReference
138 auto& m_ow_velocity = sim.repo().get_field("ow_velocity");
139
140 auto nlevels = sim.repo().num_active_levels();
141 auto geom = sim.mesh().Geom();
142
143 for (int lev = 0; lev < nlevels; ++lev) {
144 const auto& problo = geom[lev].ProbLoArray();
145 const auto& dx = geom[lev].CellSizeArray();
146
147 const auto& phi = m_ow_levelset(lev).arrays();
148 const auto& vel = m_ow_velocity(lev).arrays();
149
150 const amrex::Real wave_height = wdata.wave_height;
151 const amrex::Real wave_length = wdata.wave_length;
152 const amrex::Real phase_offset = wdata.wave_phase_offset;
153 const amrex::Real water_depth = wdata.water_depth;
154 const amrex::Real zero_sea_level = wdata.zsl;
155 const amrex::Real g = wdata.g;
156
157 amrex::ParallelFor(
158 m_ow_velocity(lev), amrex::IntVect(3),
159 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
160 const amrex::Real x = problo[0] + (i + 0.5) * dx[0];
161 const amrex::Real z = problo[2] + (k + 0.5) * dx[2];
162
163 const amrex::Real wave_number = 2. * M_PI / wave_length;
164 const amrex::Real omega = std::pow(
165 wave_number * g * std::tanh(wave_number * water_depth),
166 0.5);
167 const amrex::Real phase =
168 wave_number * x - omega * time - phase_offset;
169
170 const amrex::Real eta =
171 wave_height / 2.0 * std::cos(phase) + zero_sea_level;
172
173 phi[nbx](i, j, k) = eta - z;
174
175 const amrex::Real cell_length_2D =
176 std::sqrt(dx[0] * dx[0] + dx[2] * dx[2]);
177 if (phi[nbx](i, j, k) + cell_length_2D >= 0) {
178 vel[nbx](i, j, k, 0) =
179 omega * wave_height / 2.0 *
180 std::cosh(
181 wave_number *
182 (z - zero_sea_level + water_depth)) /
183 std::sinh(wave_number * water_depth) *
184 std::cos(phase);
185 vel[nbx](i, j, k, 1) = 0.0;
186 vel[nbx](i, j, k, 2) =
187 omega * wave_height / 2.0 *
188 std::sinh(
189 wave_number *
190 (z - zero_sea_level + water_depth)) /
191 std::sinh(wave_number * water_depth) *
192 std::sin(phase);
193 } else {
194 vel[nbx](i, j, k, 0) = 0.;
195 vel[nbx](i, j, k, 1) = 0.;
196 vel[nbx](i, j, k, 2) = 0.;
197 }
198 });
199 }
200 amrex::Gpu::synchronize();
201 }
202};
203
204} // namespace amr_wind::ocean_waves::ops
205
206#endif /* LINEAR_WAVES_OPS_H */
FieldRepo & repo()
Return the field repository.
Definition CFDSim.H:62
PhysicsMgr & physics_manager()
Definition CFDSim.H:68
SimTime & time()
Return simulation time control.
Definition CFDSim.H:58
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
Definition MultiPhase.H:27
T & get()
Return a concrete physics instance.
Definition Physics.H:104
AMREX_FORCE_INLINE amrex::Real new_time() const
Definition SimTime.H:108
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)
Definition linear_waves_ops.H:42
Definition OceanWavesOps.H:32
void operator()(LinearWaves::DataType &data, const ::amr_wind::utils::MultiParser &pp)
Definition linear_waves_ops.H:16
Definition OceanWavesOps.H:19
void operator()(LinearWaves::DataType &data)
Definition linear_waves_ops.H:128