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 */
17void
19 const Box& domain,
20 const std::string& fname,
21 FArrayBox& NC_temp_fab, FArrayBox& NC_salt_fab,
22 FArrayBox& NC_xvel_fab, FArrayBox& NC_yvel_fab)
23{
24 amrex::Print() << "Loading initial solution data from NetCDF file " << fname << std::endl;
25
26 Vector<FArrayBox*> NC_fabs;
27 Vector<std::string> NC_names;
28 Vector<enum NC_Data_Dims_Type> NC_dim_types;
29
30 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
31 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
32 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
33 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
34
35 // Read the netcdf file and fill these FABs
36 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
37}
38
39/**
40 * @param lev level of data to read
41 * @param domain simulation domain
42 * @param fname file name to read from
43 * @param NC_temp_fab container for temperature data
44 * @param NC_salt_fab container for salinity data
45 * @param NC_u_fab container for u velocity data
46 * @param NC_v_fab container for v velocity data
47 * @param ngrow number of grow cells to read
48 */
49void
51 const Box& domain,
52 const std::string& fname,
53 FArrayBox& NC_temp_fab, FArrayBox& NC_salt_fab,
54 FArrayBox& NC_xvel_fab, FArrayBox& NC_yvel_fab,
55 IntVect ngrow)
56{
57 amrex::Print() << "Loading initial solution data from NetCDF file " << fname << std::endl;
58
59 Vector<FArrayBox*> NC_fabs;
60 Vector<std::string> NC_names;
61 Vector<enum NC_Data_Dims_Type> NC_dim_types;
62
63 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
64 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
65 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
66 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
67
68 // Read the netcdf file and fill these FABs
69 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs, false, 0, ngrow);
70}
71
72/**
73 * @param lev level of data to read
74 * @param domain simulation domain
75 * @param fname file name to read from
76 * @param NC_zeta_fab container for sea surface height data
77 */
78void
80 const Box& domain,
81 const std::string& fname,
82 FArrayBox& NC_zeta_fab)
83{
84 amrex::Print() << "Loading initial sea surface height from NetCDF file " << fname << std::endl;
85
86 Vector<FArrayBox*> NC_fabs;
87 Vector<std::string> NC_names;
88 Vector<enum NC_Data_Dims_Type> NC_dim_types;
89
90 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
91
92 // Read the netcdf file and fill these FABs
93 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
94}
95
96/**
97 * @param lev level of data to read
98 * @param domain simulation domain
99 * @param fname file name to read from
100 * @param NC_zeta_fab container for sea surface height data
101 * @param ngrow number of grow cells to read in, if not default
102 */
103void
105 const Box& domain,
106 const std::string& fname,
107 FArrayBox& NC_zeta_fab,
108 IntVect ngrow)
109{
110 amrex::Print() << "Loading initial sea surface height from NetCDF file " << fname << std::endl;
111
112 Vector<FArrayBox*> NC_fabs;
113 Vector<std::string> NC_names;
114 Vector<enum NC_Data_Dims_Type> NC_dim_types;
115
116 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
117
118 // Read the netcdf file and fill these FABs
119 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs, false, 0, ngrow);
120}
121
122/**
123 * @param lev level of data to read
124 * @param domain simulation domain
125 * @param fname file name to read from
126 * @param NC_h_fab container for bathymetry data
127 */
128void
130 const Box& domain,
131 const std::string& fname,
132 FArrayBox& NC_h_fab)
133{
134 amrex::Print() << "Loading initial bathymetry from NetCDF file " << fname << std::endl;
135
136 Vector<FArrayBox*> NC_fabs;
137 Vector<std::string> NC_names;
138 Vector<enum NC_Data_Dims_Type> NC_dim_types;
139
140 NC_fabs.push_back(&NC_h_fab ) ; NC_names.push_back("h") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 0
141
142 // Read the netcdf file and fill these FABs
143 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
144}
145
146/**
147 * @param lev level of data to read
148 * @param domain simulation domain
149 * @param fname file name to read from
150 * @param NC_pm_fab container for pm data
151 * @param NC_pn_fab container for pn data
152 * @param NC_xr_fab container for x_rho data
153 * @param NC_yr_fab container for y_rho data
154 * @param NC_xu_fab container for x_u data
155 * @param NC_yu_fab container for y_u data
156 * @param NC_xv_fab container for x_v data
157 * @param NC_yv_fab container for y_v data
158 * @param NC_xp_fab container for x_p data
159 * @param NC_yp_fab container for y_p data
160 */
161void
163 const Box& domain,
164 const std::string& fname,
165 FArrayBox& NC_pm_fab, FArrayBox& NC_pn_fab,
166 FArrayBox& NC_xr_fab, FArrayBox& NC_yr_fab,
167 FArrayBox& NC_xu_fab, FArrayBox& NC_yu_fab,
168 FArrayBox& NC_xv_fab, FArrayBox& NC_yv_fab,
169 FArrayBox& NC_xp_fab, FArrayBox& NC_yp_fab)
170{
171 amrex::Print() << "Loading grid variables from NetCDF file " << fname << std::endl;
172
173 Vector<FArrayBox*> NC_fabs;
174 Vector<std::string> NC_names;
175 Vector<enum NC_Data_Dims_Type> NC_dim_types;
176
177 NC_fabs.push_back(&NC_pm_fab) ; NC_names.push_back("pm") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 1
178 NC_fabs.push_back(&NC_pn_fab) ; NC_names.push_back("pn") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 2
179 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
180 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
181 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
182 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
183 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
184 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
185 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
186 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
187
188 // Read the netcdf file and fill these FABs
189 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
190}
191
192/**
193 * @param domain simulation domain at nc_hires_grid_level
194 * @param fname file name to read from
195 * @param NC_h_fab container for bathymetry data
196 * @param ngrow number of grow cells to read in, if not default
197 */
198void
200 const std::string& fname,
201 FArrayBox& NC_h_fab, IntVect ngrow)
202{
203 amrex::Print() << "Loading high resolution bathymetry from NetCDF file " << fname << std::endl;
204
205 Vector<FArrayBox*> NC_fabs;
206 Vector<std::string> NC_names;
207 Vector<enum NC_Data_Dims_Type> NC_dim_types;
208
209 NC_fabs.push_back(&NC_h_fab ) ; NC_names.push_back("h") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 0
210
211 // Read the netcdf file and fill these FABs
212 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs, false, 0, ngrow);
213}
214
215/**
216 * @param lev level of data to read
217 * @param domain simulation domain
218 * @param fname file name to read from
219 * @param NC_pm_fab container for pm data
220 * @param NC_pn_fab container for pn data
221 * @param ngrow number of grow cells to read in, if not default
222 */
223void
225 const std::string& fname,
226 FArrayBox& NC_pm_fab, FArrayBox& NC_pn_fab,
227 IntVect ngrow)
228{
229 amrex::Print() << "Loading high resolution grid variables from NetCDF file " << fname << std::endl;
230
231 Vector<FArrayBox*> NC_fabs;
232 Vector<std::string> NC_names;
233 Vector<enum NC_Data_Dims_Type> NC_dim_types;
234
235 NC_fabs.push_back(&NC_pm_fab) ; NC_names.push_back("pm") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 0
236 NC_fabs.push_back(&NC_pn_fab) ; NC_names.push_back("pn") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 1
237 // Read the netcdf file and fill these FABs
238 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs, false, 0, ngrow);
239}
240
241/**
242 * @param lev level of data to read
243 * @param domain simulation domain
244 * @param fname file name to read from
245 * @param NC_fcor_fab container for Coriolis parameter data
246 */
247void
249 const Box& domain,
250 const std::string& fname,
251 FArrayBox& NC_fcor_fab)
252{
253 amrex::Print() << "Loading initial coriolis from NetCDF file " << fname << std::endl;
254
255 Vector<FArrayBox*> NC_fabs;
256 Vector<std::string> NC_names;
257 Vector<enum NC_Data_Dims_Type> NC_dim_types;
258
259 NC_fabs.push_back(&NC_fcor_fab ) ; NC_names.push_back("f") ; NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE); // 0
260
261 // Read the netcdf file and fill these FABs
262 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
263}
264
265/**
266 * @param lev level of data to read
267 * @param domain simulation domain
268 * @param fname file name to read from
269 * @param NC_mskr_fab container for rho-point land/sea mask data
270 * @param NC_msku_fab container for u-point land/sea mask data
271 * @param NC_mskv_fab container for v-point land/sea mask data
272 */
273void
275 const Box& domain,
276 const std::string& fname,
277 FArrayBox& NC_mskr_fab,
278 FArrayBox& NC_msku_fab,
279 FArrayBox& NC_mskv_fab)
280{
281 amrex::Print() << "Loading masks from NetCDF file " << fname << std::endl;
282
283 Vector<FArrayBox*> NC_fabs;
284 Vector<std::string> NC_names;
285 Vector<enum NC_Data_Dims_Type> NC_dim_types;
286
287 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
288 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
289 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
290
291 // Read the netcdf file and fill these FABs
292 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
293}
294
295/**
296 * @param lev level of data to read
297 * @param domain simulation domain
298 * @param fname file name to read from
299 * @param do_m2_clim_nudg whether to do 2d momentum climatology nudging
300 * @param do_m3_clim_nudg whether to do 3d momentum climatology nudging
301 * @param do_temp_clim_nudg whether to do temperature climatology nudging
302 * @param do_salt_clim_nudg whether to do salinity climatology nudging
303 * @param NC_M2NC_fab container for 2d momentum climatology data
304 * @param NC_M3NC_fab container for 3d momentum climatology data
305 * @param NC_TempNC_fab container for temperature climatology data
306 * @param NC_SaltNC_fab container for salinity climatology data
307 */
308void
310 const Box& domain,
311 const std::string& fname,
312 bool do_m2_clim_nudg,
313 bool do_m3_clim_nudg,
314 bool do_temp_clim_nudg,
315 bool do_salt_clim_nudg,
316 FArrayBox& NC_M2NC_fab,
317 FArrayBox& NC_M3NC_fab,
318 FArrayBox& NC_TempNC_fab,
319 FArrayBox& NC_SaltNC_fab)
320{
321 amrex::Print() << "Loading nudging coefficients from NetCDF file " << fname << std::endl;
322
323 Vector<FArrayBox*> NC_fabs;
324 Vector<std::string> NC_names;
325 Vector<enum NC_Data_Dims_Type> NC_dim_types;
326
327 if (do_m3_clim_nudg) {
328 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);
329 }
330 if (do_m2_clim_nudg) {
331 NC_fabs.push_back(&NC_M2NC_fab ); NC_names.push_back("M2_NudgeCoef"); NC_dim_types.push_back(NC_Data_Dims_Type::SN_WE);
332 }
333 if (do_temp_clim_nudg) {
334 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);
335 }
336 if (do_salt_clim_nudg) {
337 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);
338 }
339 // Read the netcdf file and fill these FABs
340 BuildFABsFromNetCDFFile<FArrayBox,Real>(domain, fname, NC_names, NC_dim_types, NC_fabs);
341}
342
343/**
344 * @param lev level of data to read
345 * @param fnames file name(s) to read from
346 * @param field_name field name to read
347 * @param vec_dat vector to fill data
348 */
349//template <typename DType>
350void read_vec_from_netcdf (int /*lev*/, const amrex::Vector<std::string>& fnames, const std::string& field_name, amrex::Vector<int>& vec_dat)
351{
352 AMREX_ALWAYS_ASSERT(!fnames.empty());
353 const std::string& fname = fnames[0];
354
355 amrex::Print() << "Reading " << field_name << " from NetCDF file" << std::endl;
356
357 // get x-positions and put in array
358 using ARRAY = NDArray<int>;
359 amrex::Vector<ARRAY> array_dat(1);
360 ReadNetCDFFile(fname, {field_name}, array_dat); // filled only on proc 0
361 if (amrex::ParallelDescriptor::IOProcessor())
362 {
363 int n = array_dat[0].get_vshape()[0];
364 for (int i(0); i < n; i++)
365 {
366 vec_dat.push_back((*(array_dat[0].get_data() + i)));
367 }
368 }
369 int nvals = vec_dat.size();
370 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
371 amrex::ParallelDescriptor::Bcast(&nvals,1,ioproc);
372 if (!(amrex::ParallelDescriptor::IOProcessor())) {
373 vec_dat.resize(nvals);
374 }
375 amrex::ParallelDescriptor::Bcast(vec_dat.data(), vec_dat.size(), ioproc);
376}
377
378#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_zeta_full_domain_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_zeta_fab, IntVect ngrow)
helper function to read high-resolution full-domain sea surface height from netcdf
void read_bathymetry_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_h_fab)
helper function to read bathymetry from netcdf
void read_vec_from_netcdf(int, const amrex::Vector< std::string > &fnames, const std::string &field_name, amrex::Vector< int > &vec_dat)
helper function to read in vector of data from netcdf
void read_data_full_domain_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, IntVect ngrow)
helper function for reading in full domain high-resolution initial state data from netcdf
void read_grid_vars_from_netcdf(int, const Box &domain, const std::string &fname, 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 grid variables from netcdf
void read_bathymetry_full_domain_from_netcdf(const Box &domain, const std::string &fname, FArrayBox &NC_h_fab, IntVect ngrow)
helper function to read full-domain high resolution bathymetry 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_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)
helper function for reading in initial state data from netcdf
void read_grid_vars_full_domain_from_netcdf(const Box &domain, const std::string &fname, FArrayBox &NC_pm_fab, FArrayBox &NC_pn_fab, IntVect ngrow)
helper function to read full-domain high resolution grid variables from netcdf
NDArray is the datatype designed to hold any data, including scalars, multidimensional arrays,...