REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_ReadFromInitNetcdf.cpp
Go to the documentation of this file.
1#include "REMORA_NCFile.H"
2#include "AMReX_FArrayBox.H"
3#include "REMORA_DataStruct.H"
4
5using namespace amrex;
6
7#ifdef REMORA_USE_NETCDF
8/**
9 * @param lev level of data to read
10 * @param domain simulation domain
11 * @param fname file name to read from
12 * @param NC_temp_fab container for temperature data
13 * @param NC_salt_fab container for salinity data
14 * @param NC_u_fab container for u velocity data
15 * @param NC_v_fab container for v velocity data
16 * @param NC_ubar_fab container for u_bar velocity data
17 * @param NC_vbar_fab container for v_bar velocity data
18 */
19void
21 const Box& domain,
22 const std::string& fname,
23 FArrayBox& NC_temp_fab, FArrayBox& NC_salt_fab,
24 FArrayBox& NC_xvel_fab, FArrayBox& NC_yvel_fab,
25 FArrayBox& NC_ubar_fab, FArrayBox& NC_vbar_fab)
26{
27 amrex::Print() << "Loading initial solution data from NetCDF file " << fname << std::endl;
28
29 Vector<FArrayBox*> NC_fabs;
30 Vector<std::string> NC_names;
31 Vector<enum NC_Data_Dims_Type> NC_dim_types;
32
33 NC_fabs.push_back(&NC_temp_fab); NC_names.push_back("temp"); NC_dim_types.push_back(NC_Data_Dims_Type::Time_BT_SN_WE); // 0
34 NC_fabs.push_back(&NC_salt_fab); NC_names.push_back("salt"); NC_dim_types.push_back(NC_Data_Dims_Type::Time_BT_SN_WE); // 1
35 NC_fabs.push_back(&NC_xvel_fab); NC_names.push_back("u"); NC_dim_types.push_back(NC_Data_Dims_Type::Time_BT_SN_WE); // 2
36 NC_fabs.push_back(&NC_yvel_fab); NC_names.push_back("v"); NC_dim_types.push_back(NC_Data_Dims_Type::Time_BT_SN_WE); // 3
37 NC_fabs.push_back(&NC_ubar_fab), NC_names.push_back("ubar"); NC_dim_types.push_back(NC_Data_Dims_Type::Time_SN_WE); // 4
38 NC_fabs.push_back(&NC_vbar_fab); NC_names.push_back("vbar"); NC_dim_types.push_back(NC_Data_Dims_Type::Time_SN_WE); // 5
39
40 // Read the netcdf file and fill these FABs
41 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
42}
43
44/**
45 * @param lev level of data to read
46 * @param domain simulation domain
47 * @param fname file name to read from
48 * @param NC_zeta_fab container for sea surface height data
49 */
50void
52 const Box& domain,
53 const std::string& fname,
54 FArrayBox& NC_zeta_fab)
55{
56 amrex::Print() << "Loading initial sea surface height from NetCDF file " << fname << std::endl;
57
58 Vector<FArrayBox*> NC_fabs;
59 Vector<std::string> NC_names;
60 Vector<enum NC_Data_Dims_Type> NC_dim_types;
61
62 NC_fabs.push_back(&NC_zeta_fab ) ; NC_names.push_back("zeta") ; NC_dim_types.push_back(NC_Data_Dims_Type::Time_SN_WE); // 0
63
64 // Read the netcdf file and fill these FABs
65 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
66}
67
68/**
69 * @param lev level of data to read
70 * @param domain simulation domain
71 * @param fname file name to read from
72 * @param NC_h_fab container for bathymetry data
73 * @param NC_pm_fab container for pm data
74 * @param NC_pn_fab container for pn data
75 * @param NC_xr_fab container for x_rho data
76 * @param NC_yr_fab container for y_rho data
77 * @param NC_xu_fab container for x_u data
78 * @param NC_yu_fab container for y_u data
79 * @param NC_xv_fab container for x_v data
80 * @param NC_yv_fab container for y_v data
81 * @param NC_xp_fab container for x_p data
82 * @param NC_yp_fab container for y_p data
83 */
84void
86 const Box& domain,
87 const std::string& fname,
88 FArrayBox& NC_h_fab,
89 FArrayBox& NC_pm_fab, FArrayBox& NC_pn_fab,
90 FArrayBox& NC_xr_fab, FArrayBox& NC_yr_fab,
91 FArrayBox& NC_xu_fab, FArrayBox& NC_yu_fab,
92 FArrayBox& NC_xv_fab, FArrayBox& NC_yv_fab,
93 FArrayBox& NC_xp_fab, FArrayBox& NC_yp_fab)
94{
95 amrex::Print() << "Loading initial bathymetry from NetCDF file " << fname << std::endl;
96
97 Vector<FArrayBox*> NC_fabs;
98 Vector<std::string> NC_names;
99 Vector<enum NC_Data_Dims_Type> NC_dim_types;
100
101 NC_fabs.push_back(&NC_h_fab ) ; NC_names.push_back("h") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 0
102 NC_fabs.push_back(&NC_pm_fab) ; NC_names.push_back("pm") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 1
103 NC_fabs.push_back(&NC_pn_fab) ; NC_names.push_back("pn") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 2
104 NC_fabs.push_back(&NC_xr_fab) ; NC_names.push_back("x_rho") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 3
105 NC_fabs.push_back(&NC_yr_fab) ; NC_names.push_back("y_rho") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 4
106 NC_fabs.push_back(&NC_xu_fab) ; NC_names.push_back("x_u") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 5
107 NC_fabs.push_back(&NC_yu_fab) ; NC_names.push_back("y_u") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 6
108 NC_fabs.push_back(&NC_xv_fab) ; NC_names.push_back("x_v") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 7
109 NC_fabs.push_back(&NC_yv_fab) ; NC_names.push_back("y_v") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 8
110 NC_fabs.push_back(&NC_xp_fab) ; NC_names.push_back("x_psi") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 9
111 NC_fabs.push_back(&NC_yp_fab) ; NC_names.push_back("y_psi") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 10
112
113 // Read the netcdf file and fill these FABs
114 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
115}
116
117/**
118 * @param lev level of data to read
119 * @param domain simulation domain
120 * @param fname file name to read from
121 * @param NC_fcor_fab container for Coriolis parameter data
122 */
123void
125 const Box& domain,
126 const std::string& fname,
127 FArrayBox& NC_fcor_fab)
128{
129 amrex::Print() << "Loading initial coriolis from NetCDF file " << fname << std::endl;
130
131 Vector<FArrayBox*> NC_fabs;
132 Vector<std::string> NC_names;
133 Vector<enum NC_Data_Dims_Type> NC_dim_types;
134
135 NC_fabs.push_back(&NC_fcor_fab ) ; NC_names.push_back("f") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 0
136
137 // Read the netcdf file and fill these FABs
138 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
139}
140
141/**
142 * @param lev level of data to read
143 * @param domain simulation domain
144 * @param fname file name to read from
145 * @param NC_mskr_fab container for rho-point land/sea mask data
146 * @param NC_msku_fab container for u-point land/sea mask data
147 * @param NC_mskv_fab container for v-point land/sea mask data
148 */
149void
151 const Box& domain,
152 const std::string& fname,
153 FArrayBox& NC_mskr_fab,
154 FArrayBox& NC_msku_fab,
155 FArrayBox& NC_mskv_fab)
156{
157 amrex::Print() << "Loading masks from NetCDF file " << fname << std::endl;
158
159 Vector<FArrayBox*> NC_fabs;
160 Vector<std::string> NC_names;
161 Vector<enum NC_Data_Dims_Type> NC_dim_types;
162
163 NC_fabs.push_back(&NC_mskr_fab ) ; NC_names.push_back("mask_rho") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 0
164 NC_fabs.push_back(&NC_msku_fab ) ; NC_names.push_back("mask_u") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 1
165 NC_fabs.push_back(&NC_mskv_fab ) ; NC_names.push_back("mask_v") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 2
166
167 // Read the netcdf file and fill these FABs
168 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
169}
170
171/**
172 * @param lev level of data to read
173 * @param domain simulation domain
174 * @param fname file name to read from
175 * @param do_m2_clim_nudg whether to do 2d momentum climatology nudging
176 * @param do_m3_clim_nudg whether to do 3d momentum climatology nudging
177 * @param do_temp_clim_nudg whether to do temperature climatology nudging
178 * @param do_salt_clim_nudg whether to do salinity climatology nudging
179 * @param NC_M2NC_fab container for 2d momentum climatology data
180 * @param NC_M3NC_fab container for 3d momentum climatology data
181 * @param NC_TempNC_fab container for temperature climatology data
182 * @param NC_SaltNC_fab container for salinity climatology data
183 */
184void
186 const Box& domain,
187 const std::string& fname,
188 bool do_m2_clim_nudg,
189 bool do_m3_clim_nudg,
190 bool do_temp_clim_nudg,
191 bool do_salt_clim_nudg,
192 FArrayBox& NC_M2NC_fab,
193 FArrayBox& NC_M3NC_fab,
194 FArrayBox& NC_TempNC_fab,
195 FArrayBox& NC_SaltNC_fab)
196{
197 amrex::Print() << "Loading nudging coefficients from NetCDF file " << fname << std::endl;
198
199 Vector<FArrayBox*> NC_fabs;
200 Vector<std::string> NC_names;
201 Vector<enum NC_Data_Dims_Type> NC_dim_types;
202
203 if (do_m3_clim_nudg) {
204 NC_fabs.push_back(&NC_M3NC_fab ); NC_names.push_back("M3_NudgeCoef"); NC_dim_types.push_back(NC_Data_Dims_Type::BT_SN_WE);
205 }
206 if (do_m2_clim_nudg) {
207 NC_fabs.push_back(&NC_M2NC_fab ); NC_names.push_back("M2_NudgeCoef"); NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE);
208 }
209 if (do_temp_clim_nudg) {
210 NC_fabs.push_back(&NC_TempNC_fab ); NC_names.push_back("temp_NudgeCoef"); NC_dim_types.push_back(NC_Data_Dims_Type::BT_SN_WE);
211 }
212 if (do_salt_clim_nudg) {
213 NC_fabs.push_back(&NC_SaltNC_fab ); NC_names.push_back("salt_NudgeCoef"); NC_dim_types.push_back(NC_Data_Dims_Type::BT_SN_WE);
214 }
215 // Read the netcdf file and fill these FABs
216 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
217}
218
219/**
220 * @param lev level of data to read
221 * @param fname file name to read from
222 * @param field_name field name to read
223 * @param vec_dat vector to fill data
224 */
225//template <typename DType>
226void read_vec_from_netcdf (int /*lev*/, const std::string& fname, const std::string& field_name, amrex::Vector<int>& vec_dat)
227{
228 amrex::Print() << "Reading " << field_name << " from NetCDF file" << std::endl;
229
230 // get x-positions and put in array
231 using ARRAY = NDArray<int>;
232 amrex::Vector<ARRAY> array_dat(1);
233 ReadNetCDFFile(fname, {field_name}, array_dat); // filled only on proc 0
234 if (amrex::ParallelDescriptor::IOProcessor())
235 {
236 int n = array_dat[0].get_vshape()[0];
237 for (int i(0); i < n; i++)
238 {
239 vec_dat.push_back((*(array_dat[0].get_data() + i)));
240 }
241 }
242 int nvals = vec_dat.size();
243 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
244 amrex::ParallelDescriptor::Bcast(&nvals,1,ioproc);
245 if (!(amrex::ParallelDescriptor::IOProcessor())) {
246 vec_dat.resize(nvals);
247 }
248 amrex::ParallelDescriptor::Bcast(vec_dat.data(), vec_dat.size(), ioproc);
249}
250
251#endif // ROMSX_USE_NETCDF
void ReadNetCDFFile(const std::string &fname, amrex::Vector< std::string > names, amrex::Vector< NDArray< DType > > &arrays, bool one_time=false, int fill_time=0)
Read in data from netcdf file and save to data arrays.
void read_clim_nudg_coeff_from_netcdf(int, const Box &domain, const std::string &fname, bool do_m2_clim_nudg, bool do_m3_clim_nudg, bool do_temp_clim_nudg, bool do_salt_clim_nudg, FArrayBox &NC_M2NC_fab, FArrayBox &NC_M3NC_fab, FArrayBox &NC_TempNC_fab, FArrayBox &NC_SaltNC_fab)
helper function to read climatology nudging from netcdf
void read_data_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_temp_fab, FArrayBox &NC_salt_fab, FArrayBox &NC_xvel_fab, FArrayBox &NC_yvel_fab, FArrayBox &NC_ubar_fab, FArrayBox &NC_vbar_fab)
helper function for reading in initial state data from netcdf
void read_coriolis_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_fcor_fab)
helper function to read coriolis factor from netcdf
void read_masks_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_mskr_fab, FArrayBox &NC_msku_fab, FArrayBox &NC_mskv_fab)
helper function for reading in land-sea masks from netcdf
void read_zeta_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_zeta_fab)
helper function to read sea surface height from netcdf
void read_bathymetry_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_h_fab, FArrayBox &NC_pm_fab, FArrayBox &NC_pn_fab, FArrayBox &NC_xr_fab, FArrayBox &NC_yr_fab, FArrayBox &NC_xu_fab, FArrayBox &NC_yu_fab, FArrayBox &NC_xv_fab, FArrayBox &NC_yv_fab, FArrayBox &NC_xp_fab, FArrayBox &NC_yp_fab)
helper function to read bathymetry from netcdf
void read_vec_from_netcdf(int, const std::string &fname, const std::string &field_name, amrex::Vector< int > &vec_dat)
helper function to read in vector of data from netcdf
NDArray is the datatype designed to hold any data, including scalars, multidimensional arrays,...