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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/ocean_waves/relaxation_zones/W2AWaves.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
W2AWaves.H
Go to the documentation of this file.
1#ifndef W2AWAVES_H
2#define W2AWAVES_H
3
5#include <vector>
6#include <complex>
7#ifdef AMR_WIND_USE_W2A
8#include "Waves2AMR.h"
9#include "AMReX_REAL.H"
10
11using namespace amrex::literals;
12#endif
13
14namespace amr_wind::ocean_waves {
15
17{
18 // Prefix for HOS files
19 std::string modes_file{"modes_HOS_SWENSE.dat"};
20 // Number of modes in each direction (related to FFT)
21 int n0{0};
22 int n1{0};
23 int n2{0};
24 // Number of points in each direction (spatial output)
25 int n0_sp{0};
26 int n1_sp{0};
27 // Spatial resolution in each direction
28 amrex::Real dx0{0.};
29 amrex::Real dx1{0.};
30 // Offsets for spatial data (between AMR-Wind domain and HOS domain)
31 amrex::Real xlo0{0.};
32 amrex::Real xlo1{0.};
33 // Length to dimensionalize nondim data
34 amrex::Real dimL{0.};
35 // Timestep size of reference data
36 amrex::Real dt_modes{0.0_rt};
37
38 // Time (in HOS sim) to begin reading from
39 // (this corresponds to t = 0 in AMR-Wind sim, functions as offset)
40 amrex::Real t_winit{-1.0_rt};
41 // Timestep (in HOS sim) of latest imported wave data
42 int ntime{0};
43 // Timestep (in HOS sim) of first imported wave data
44 int n_winit{0};
45 // Timestep (in HOS sim) of last available wave data
46 int n_wstop{0};
47 // Timestep offset (in HOS sim) for resetting to earlier data
48 int n_offset{0};
49
50 // Time (in AMR-Wind sim) when ow arrays were last changed
51 // (this is from the last timestep, goes with time-interpolated data)
52 amrex::Real t_last{-1.0_rt};
53 // Time (in AMR-Wind sim) of latest imported wave data
54 amrex::Real t{0.0_rt};
55
56 // Vectors to store complex or real modes
57 std::vector<std::complex<amrex::Real>> c_mX, c_mY, c_mZ, c_mFS;
58 std::vector<amrex::Real> r_mX, r_mY, r_mZ, r_mFS, r_mAdd;
59
60// Struct variables that have special types
61#ifdef AMR_WIND_USE_W2A
62 // FFTW plan vector
63 std::vector<fftw_plan> plan_vector;
64 // FFTW pointers to store modes ready for ifft
65 fftw_complex *c_eta_mptr = nullptr, *c_u_mptr = nullptr,
66 *c_v_mptr = nullptr, *c_w_mptr = nullptr;
67 // FFTW pointers to store modes ready for ifft
68 amrex::Real *r_eta_mptr = nullptr, *r_u_mptr = nullptr, *r_v_mptr = nullptr,
69 *r_w_mptr = nullptr, *au_mptr = nullptr, *av_mptr = nullptr,
70 *aw_mptr = nullptr, *plan_out = nullptr;
71 // ReadModes objects
72 ReadModes<std::complex<amrex::Real>> c_rmodes;
73 ReadModes<amrex::Real> r_rmodes;
74#endif
75
76 // Height vector (where velocity is sampled) stuff
77 amrex::Vector<amrex::Real> hvec;
78 // Index vector (where hvec overlaps with local boxes)
79 amrex::Vector<int> indvec;
80 // Flag indicating whether the source simulation used HOS-Ocean or HOS-NWT
81 bool is_ocean{true};
82 // Flag indicating interpolation should take place on this processor
83 bool do_interp{true};
84 // Flag indicating regrid has occurred since the last resizing of vecs
85 bool resize_flag{false};
86 // Vectors for spatial data from transformed modes
87 amrex::Gpu::DeviceVector<amrex::Real> sp_eta_vec, sp_u_vec, sp_v_vec,
89
91 {
92#ifdef AMR_WIND_USE_W2A
93 if (c_eta_mptr) {
94 delete[] c_eta_mptr;
95 }
96 if (c_u_mptr) {
97 delete[] c_u_mptr;
98 }
99 if (c_v_mptr) {
100 delete[] c_v_mptr;
101 }
102 if (c_w_mptr) {
103 delete[] c_w_mptr;
104 }
105 if (r_eta_mptr) {
106 delete[] r_eta_mptr;
107 }
108 if (r_u_mptr) {
109 delete[] r_u_mptr;
110 }
111 if (r_v_mptr) {
112 delete[] r_v_mptr;
113 }
114 if (r_w_mptr) {
115 delete[] r_w_mptr;
116 }
117 if (au_mptr) {
118 delete[] au_mptr;
119 }
120 if (av_mptr) {
121 delete[] av_mptr;
122 }
123 if (aw_mptr) {
124 delete[] aw_mptr;
125 }
126 if (plan_out) {
127 delete[] plan_out;
128 }
129 for (int n = 0; n < (int)plan_vector.size(); ++n) {
130 fftw_destroy_plan(plan_vector[n]);
131 }
132 plan_vector.clear();
133#endif
134 }
135};
136
138{
142
143 static std::string identifier() { return "W2AWaves"; }
144};
145
146} // namespace amr_wind::ocean_waves
147
148#endif /* W2AWAVES_H */
Definition OceanWavesTypes.H:59
Definition OceanWaves.cpp:14
Definition OceanWavesTypes.H:35
Definition RelaxationZones.H:19
Definition RelaxationZones.H:49
Definition W2AWaves.H:17
amrex::Real xlo1
Definition W2AWaves.H:32
amrex::Real dimL
Definition W2AWaves.H:34
amrex::Real dt_modes
Definition W2AWaves.H:36
amrex::Gpu::DeviceVector< amrex::Real > sp_v_vec
Definition W2AWaves.H:87
std::vector< amrex::Real > r_mFS
Definition W2AWaves.H:58
std::vector< std::complex< amrex::Real > > c_mFS
Definition W2AWaves.H:57
amrex::Vector< int > indvec
Definition W2AWaves.H:79
std::vector< amrex::Real > r_mAdd
Definition W2AWaves.H:58
int n1
Definition W2AWaves.H:22
amrex::Gpu::DeviceVector< amrex::Real > sp_u_vec
Definition W2AWaves.H:87
int n_winit
Definition W2AWaves.H:44
amrex::Vector< amrex::Real > hvec
Definition W2AWaves.H:77
int n_offset
Definition W2AWaves.H:48
amrex::Real dx0
Definition W2AWaves.H:28
int ntime
Definition W2AWaves.H:42
~W2AWavesData()
Definition W2AWaves.H:90
int n1_sp
Definition W2AWaves.H:26
std::vector< std::complex< amrex::Real > > c_mZ
Definition W2AWaves.H:57
amrex::Real t_last
Definition W2AWaves.H:52
amrex::Real xlo0
Definition W2AWaves.H:31
amrex::Gpu::DeviceVector< amrex::Real > sp_w_vec
Definition W2AWaves.H:88
std::vector< std::complex< amrex::Real > > c_mX
Definition W2AWaves.H:57
int n0
Definition W2AWaves.H:21
int n0_sp
Definition W2AWaves.H:25
std::vector< amrex::Real > r_mY
Definition W2AWaves.H:58
amrex::Real dx1
Definition W2AWaves.H:29
amrex::Real t_winit
Definition W2AWaves.H:40
int n2
Definition W2AWaves.H:23
int n_wstop
Definition W2AWaves.H:46
std::vector< amrex::Real > r_mZ
Definition W2AWaves.H:58
std::string modes_file
Definition W2AWaves.H:19
std::vector< std::complex< amrex::Real > > c_mY
Definition W2AWaves.H:57
bool is_ocean
Definition W2AWaves.H:81
amrex::Gpu::DeviceVector< amrex::Real > sp_eta_vec
Definition W2AWaves.H:87
bool do_interp
Definition W2AWaves.H:83
bool resize_flag
Definition W2AWaves.H:85
std::vector< amrex::Real > r_mX
Definition W2AWaves.H:58
amrex::Real t
Definition W2AWaves.H:54
Definition W2AWaves.H:138
static std::string identifier()
Definition W2AWaves.H:143
OceanWavesDataHolder< W2AWaves > DataType
Definition W2AWaves.H:141
OceanWavesInfo InfoType
Definition W2AWaves.H:139
W2AWavesData MetaType
Definition W2AWaves.H:140