REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_NCPlotFile.cpp
Go to the documentation of this file.
1#include <iomanip>
2#include <iostream>
3#include <string>
4#include <ctime>
5
6#ifdef _OPENMP
7#include <omp.h>
8#endif
9
10#include <AMReX_Utility.H>
11#include <AMReX_buildInfo.H>
12#include <AMReX_ParmParse.H>
13
14#include "REMORA.H"
15#include "REMORA_NCInterface.H"
16#include "REMORA_NCPlotFile.H"
17#include "REMORA_IndexDefines.H"
18
19using namespace amrex;
20
21/**
22 * @param which_step current step for output
23 */
24void REMORA::WriteNCPlotFile(int which_step, MultiFab const* plotMF) {
25 AMREX_ASSERT(max_level == 0);
26 // For right now we assume single level -- we will generalize this later to multilevel
27 int lev = 0;
28 int which_subdomain = 0;
29 int which_step_in_chunk = -1;
30
31 // Create filename
32 std::string plt_string;
33 std::string plotfilename;
35 plotfilename = plot_file_name + "_his";
36 } else {
37 plotfilename = Concatenate(plot_file_name, which_step, file_min_digits);
38 }
39 // If chunking, concatenate with which file we're in
42 plotfilename = Concatenate(plotfilename, which_chunk, file_min_digits);
43 which_step_in_chunk = history_count - which_chunk * REMORA::steps_per_history_file;
44 }
45
46 // Set the full IO path for NetCDF output
47 std::string FullPath = plotfilename;
48 if (lev == 0) {
49 const std::string &extension = amrex::Concatenate("_d", lev + 1, 2);
50 FullPath += extension + ".nc";
51 } else {
52 const std::string &extension = amrex::Concatenate("_d", lev + 1 + which_subdomain, 2);
53 FullPath += extension + ".nc";
54 }
55
56 //
57 // Check if this file/directory already exists and if so,
58 // have the IOProcessor move the existing
59 // file/directory to filename.old
60 //
61 if ((!REMORA::write_history_file) || (which_step == 0) || (which_step_in_chunk == 0)) {
62 if (amrex::ParallelDescriptor::IOProcessor()) {
63 if (amrex::FileExists(FullPath)) {
64 std::string newoldname(FullPath + ".old." + amrex::UniqueString());
65 amrex::Print() << "WriteNCPlotFile: " << FullPath << " exists. Renaming to: " << newoldname << std::endl;
66 if (std::rename(FullPath.c_str(), newoldname.c_str())) {
67 amrex::Abort("WriteNCPlotfile:: std::rename failed");
68 }
69 }
70 }
71 ParallelDescriptor::Barrier();
72 }
73
74 bool is_history;
75
77 is_history = true;
78 bool write_header = !(amrex::FileExists(FullPath));
79
80 auto ncf =
81 write_header ?
82 ncutils::NCFile::create(FullPath, NC_CLOBBER|NC_64BIT_DATA, amrex::ParallelContext::CommunicatorSub(), MPI_INFO_NULL) :
83 ncutils::NCFile::open(FullPath, NC_WRITE, amrex::ParallelContext::CommunicatorSub(), MPI_INFO_NULL);
84
85 amrex::Print() << "Writing into level " << lev << " NetCDF history file " << FullPath << std::endl;
86
87 WriteNCPlotFile_which(lev, which_subdomain, plotMF, write_header, ncf, is_history);
88
89 } else {
90
91 is_history = false;
92 bool write_header = true;
93
94 // Open new netcdf file to write data
95 auto ncf = ncutils::NCFile::create(FullPath, NC_CLOBBER|NC_64BIT_DATA, amrex::ParallelContext::CommunicatorSub(), MPI_INFO_NULL);
96 amrex::Print() << "Writing level " << lev << " NetCDF plot file " << FullPath << std::endl;
97
98 WriteNCPlotFile_which(lev, which_subdomain, plotMF, write_header, ncf, is_history);
99 }
100}
101
102/**
103 * @param lev level of data to output
104 * @param which_subdomain index of subdomain if lev != 0
105 * @param write_header whether to write a header
106 * @param ncf netcdf file object
107 * @param is_history whether the file being written is a history file
108 */
109void REMORA::WriteNCPlotFile_which(int lev, int which_subdomain, MultiFab const* plotMF,
110 bool write_header, ncutils::NCFile &ncf, bool is_history)
111{
112 // Number of cells in this "domain" at this level
113 std::vector<int> n_cells;
114
115 // We only do single-level writes when using NetCDF format
116 int flev = 1; //max_level;
117
118 Box subdomain;
119 if (lev == 0) {
120 subdomain = geom[lev].Domain();
121 } else {
122 subdomain = boxes_at_level[lev][which_subdomain];
123 }
124
125 int nx = subdomain.length(0);
126 int ny = subdomain.length(1);
127 int nz = subdomain.length(2);
128
129 if (is_history && max_step < 0) {
130 amrex::Abort("Need to know max_step if writing history file");
131 }
132
133 long long int nt;
134 if (is_history) {
135 if (max_step > 0) {
136 nt = static_cast<long long int>(max_step / std::min(plot_int, max_step)) + 1;
137 } else {
138 nt = 1;
139 }
140 } else {
141 nt = 1;
142 }
143
144 if (chunk_history_file) {
145 // First index of the last history file
146 int last_file_index = REMORA::steps_per_history_file * int(nt / REMORA::steps_per_history_file);
147 if (history_count >= last_file_index) {
148 nt = nt - last_file_index;
149 } else {
151 }
152 }
153
154 n_cells.push_back(nx);
155 n_cells.push_back(ny);
156 n_cells.push_back(nz);
157
158 const std::string nt_name = "ocean_time";
159 const std::string ndim_name = "num_geo_dimensions";
160
161 const std::string flev_name = "FINEST_LEVEL";
162
163 const std::string nx_name = "NX";
164 const std::string ny_name = "NY";
165 const std::string nz_name = "NZ";
166
167 const std::string nx_r_name = "xi_rho";
168 const std::string ny_r_name = "eta_rho";
169 const std::string nz_r_name = "s_rho";
170
171 const std::string nx_u_name = "xi_u";
172 const std::string ny_u_name = "eta_u";
173
174 const std::string nx_v_name = "xi_v";
175 const std::string ny_v_name = "eta_v";
176
177 const std::string nx_p_name = "xi_psi";
178 const std::string ny_p_name = "eta_psi";
179 const std::string nz_w_name = "s_w";
180
181 if (write_header) {
182 ncf.enter_def_mode();
183 ncf.put_attr("title", "REMORA data ");
184 ncf.def_dim(nt_name, nt);
185 ncf.def_dim(ndim_name, AMREX_SPACEDIM);
186
187 ncf.def_dim(nx_r_name, nx + 2);
188 ncf.def_dim(ny_r_name, ny + 2);
189 ncf.def_dim(nz_r_name, nz);
190
191 ncf.def_dim(nx_u_name, nx + 1);
192 ncf.def_dim(ny_u_name, ny + 2);
193
194 ncf.def_dim(nx_v_name, nx + 2);
195 ncf.def_dim(ny_v_name, ny + 1);
196
197 ncf.def_dim(nx_p_name, nx + 1);
198 ncf.def_dim(ny_p_name, ny + 1);
199
200 ncf.def_dim(nz_w_name, nz + 1);
201
202 ncf.def_dim(flev_name, flev);
203
204 ncf.def_dim(nx_name, n_cells[0]);
205 ncf.def_dim(ny_name, n_cells[1]);
206 ncf.def_dim(nz_name, n_cells[2]);
207
208 ncf.def_var("probLo", ncutils::NCDType::Real, { ndim_name });
209 ncf.var("probLo").put_attr("long_name","Low side of problem domain in internal AMReX grid");
210 ncf.var("probLo").put_attr("units","meter");
211 ncf.def_var("probHi", ncutils::NCDType::Real, { ndim_name });
212 ncf.var("probHi").put_attr("long_name","High side of problem domain in internal AMReX grid");
213 ncf.var("probHi").put_attr("units","meter");
214
215 ncf.def_var("Geom.smallend", NC_INT, { flev_name, ndim_name });
216 ncf.var("Geom.smallend").put_attr("long_name","Low side of problem domain in index space");
217 ncf.def_var("Geom.bigend", NC_INT, { flev_name, ndim_name });
218 ncf.var("Geom.bigend").put_attr("long_name","High side of problem domain in index space");
219 ncf.def_var("CellSize", ncutils::NCDType::Real, { flev_name, ndim_name });
220 ncf.var("CellSize").put_attr("long_name","Cell size on internal AMReX grid");
221 ncf.var("CellSize").put_attr("units","meter");
222
223 ncf.def_var("theta_s",ncutils::NCDType::Real,{});
224 ncf.var("theta_s").put_attr("long_name","S-coordinate surface control parameter");
225 ncf.def_var("theta_b",ncutils::NCDType::Real,{});
226 ncf.var("theta_b").put_attr("long_name","S-coordinate bottom control parameter");
227 ncf.def_var("hc",ncutils::NCDType::Real,{});
228 ncf.var("hc").put_attr("long_name","S-coordinate parameter, critical depth");
229 ncf.var("hc").put_attr("units","meter");
230
231 ncf.def_var("grid",NC_INT, {});
232 ncf.var("grid").put_attr("cf_role","grid_topology");
233 ncf.var("grid").put_attr("topology_dimension",std::vector({2}));
234 ncf.var("grid").put_attr("node_dimensions", "xi_psi eta_psi");
235 ncf.var("grid").put_attr("face_dimensions", "xi_rho: xi_psi (padding: both) eta_rho: eta_psi (padding: both)");
236 ncf.var("grid").put_attr("edge1_dimensions", "xi_u: xi_psi eta_u: eta_psi (padding: both)");
237 ncf.var("grid").put_attr("edge2_dimensions", "xi_v: xi_psi (padding: both) eta_v: eta_psi");
238 ncf.var("grid").put_attr("node_coordinates", "x_psi y_psi");
239 ncf.var("grid").put_attr("face_coordinates", "x_rho y_rho");
240 ncf.var("grid").put_attr("edge1_coordinates", "x_u y_u");
241 ncf.var("grid").put_attr("edge2_coordinates", "x_v y_v");
242 ncf.var("grid").put_attr("vertical_dimensions", "s_rho: s_w (padding: none)");
243
244 ncf.def_var("s_rho",ncutils::NCDType::Real, {nz_r_name});
245 ncf.var("s_rho").put_attr("long_name","S-coordinate at RHO-points");
246 ncf.var("s_rho").put_attr("field","s_rho, scalar");
247
248 ncf.def_var("s_w",ncutils::NCDType::Real, {nz_w_name});
249 ncf.var("s_w").put_attr("long_name","S-coordinate at W-points");
250 ncf.var("s_w").put_attr("field","s_w, scalar");
251
252 ncf.def_var("pm",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
253 ncf.var("pm").put_attr("long_name","curvilinear coordinate metric in XI");
254 ncf.var("pm").put_attr("units","meter-1");
255 ncf.var("pm").put_attr("grid","grid");
256 ncf.var("pm").put_attr("location","face");
257 ncf.var("pm").put_attr("coordinates","x_rho y_rho");
258 ncf.var("pm").put_attr("field","pm, scalar");
259
260 ncf.def_var("pn",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
261 ncf.var("pn").put_attr("long_name","curvilinear coordinate metric in ETA");
262 ncf.var("pn").put_attr("units","meter-1");
263 ncf.var("pn").put_attr("grid","grid");
264 ncf.var("pn").put_attr("location","face");
265 ncf.var("pn").put_attr("coordinates","x_rho y_rho");
266 ncf.var("pn").put_attr("field","pn, scalar");
267
268 ncf.def_var("f",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
269 ncf.var("f").put_attr("long_name","Coriolis parameter at RHO-points");
270 ncf.var("f").put_attr("units","second-1");
271 ncf.var("f").put_attr("grid","grid");
272 ncf.var("f").put_attr("location","face");
273 ncf.var("f").put_attr("coordinates","x_rho y_rho");
274 ncf.var("f").put_attr("field","coriolis, scalar");
275
276 ncf.def_var("x_rho",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
277 ncf.var("x_rho").put_attr("long_name","x-locations of RHO-points");
278 ncf.var("x_rho").put_attr("units","meter");
279 ncf.var("x_rho").put_attr("field","x_rho, scalar");
280
281 ncf.def_var("y_rho",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
282 ncf.var("y_rho").put_attr("long_name","y-locations of RHO-points");
283 ncf.var("y_rho").put_attr("units","meter");
284 ncf.var("y_rho").put_attr("field","y_rho, scalar");
285
286 ncf.def_var("x_u",ncutils::NCDType::Real, {ny_u_name, nx_u_name});
287 ncf.var("x_u").put_attr("long_name","x-locations of U-points");
288 ncf.var("x_u").put_attr("units","meter");
289 ncf.var("x_u").put_attr("field","x_u, scalar");
290
291 ncf.def_var("y_u",ncutils::NCDType::Real, {ny_u_name, nx_u_name});
292 ncf.var("y_u").put_attr("long_name","y-locations of U-points");
293 ncf.var("y_u").put_attr("units","meter");
294 ncf.var("y_u").put_attr("field","y_u, scalar");
295
296 ncf.def_var("x_v",ncutils::NCDType::Real, {ny_v_name, nx_v_name});
297 ncf.var("x_v").put_attr("long_name","x-locations of V-points");
298 ncf.var("x_v").put_attr("units","meter");
299 ncf.var("x_v").put_attr("field","x_v, scalar");
300
301 ncf.def_var("y_v",ncutils::NCDType::Real, {ny_v_name, nx_v_name});
302 ncf.var("y_v").put_attr("long_name","y-locations of V-points");
303 ncf.var("y_v").put_attr("units","meter");
304 ncf.var("y_v").put_attr("field","y_v, scalar");
305
306 ncf.def_var("x_psi",ncutils::NCDType::Real, {ny_p_name, nx_p_name});
307 ncf.var("x_psi").put_attr("long_name","x-locations of PSI-points");
308 ncf.var("x_psi").put_attr("units","meter");
309 ncf.var("x_psi").put_attr("field","x_psi, scalar");
310
311 ncf.def_var("y_psi",ncutils::NCDType::Real, {ny_p_name, nx_p_name});
312 ncf.var("y_psi").put_attr("long_name","y-locations of PSI-points");
313 ncf.var("y_psi").put_attr("units","meter");
314 ncf.var("y_psi").put_attr("field","y_psi, scalar");
315
316 ncf.def_var("ocean_time", ncutils::NCDType::Real, { nt_name });
317 ncf.var("ocean_time").put_attr("long_name","time since initialization");
318 ncf.var("ocean_time").put_attr("units","seconds since 0001-01-01 00:00:00");
319 ncf.var("ocean_time").put_attr("field","time, scalar, series");
320
321 ncf.def_var("Cs_r", ncutils::NCDType::Real, {nz_r_name});
322 ncf.var("Cs_r").put_attr("long_name", "S-coordinate stretching curves at RHO points");
323 ncf.var("Cs_r").put_attr("valid_min",std::vector({-1.}));
324 ncf.var("Cs_r").put_attr("valid_max",std::vector({0.}));
325 ncf.var("Cs_r").put_attr("field","Cs_r, scalar");
326
327 ncf.def_var("Cs_w", ncutils::NCDType::Real, {nz_w_name});
328 ncf.var("Cs_w").put_attr("long_name", "S-coordinate stretching curves at W points");
329 ncf.var("Cs_w").put_attr("valid_min",std::vector({-1.}));
330 ncf.var("Cs_w").put_attr("valid_max",std::vector({0.}));
331 ncf.var("Cs_w").put_attr("field","Cs_w, scalar");
332
333 ncf.def_var("h", ncutils::NCDType::Real, { ny_r_name, nx_r_name });
334 ncf.var("h").put_attr("long_name","bathymetry at RHO-points");
335 ncf.var("h").put_attr("units","meter");
336 ncf.var("h").put_attr("grid","grid");
337 ncf.var("h").put_attr("location","face");
338 ncf.var("h").put_attr("coordinates","x_rho y_rho");
339 ncf.var("h").put_attr("field","bath, scalar");
340
341 ncf.def_var_fill("zeta", ncutils::NCDType::Real, { nt_name, ny_r_name, nx_r_name }, &netcdf_fill_value);
342 ncf.var("zeta").put_attr("long_name","free-surface");
343 ncf.var("zeta").put_attr("units","meter");
344 ncf.var("zeta").put_attr("time","ocean_time");
345 ncf.var("zeta").put_attr("grid","grid");
346 ncf.var("zeta").put_attr("location","face");
347 ncf.var("zeta").put_attr("coordinates","x_rho y_rho ocean_time");
348 ncf.var("zeta").put_attr("field","free-surface, scalar, series");
349
350 {
351 int comp = -1;
352 for (int i = 0; i < plot_var_names_3d.size(); i++) {
353 if (plot_var_names_3d[i] == "temp") comp = i;
354 }
355 if (comp >= 0) {
356 ncf.def_var_fill("temp", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_r_name, nx_r_name }, &netcdf_fill_value);
357 ncf.var("temp").put_attr("long_name","potential temperature");
358 ncf.var("temp").put_attr("units","Celsius");
359 ncf.var("temp").put_attr("time","ocean_time");
360 ncf.var("temp").put_attr("grid","grid");
361 ncf.var("temp").put_attr("location","face");
362 ncf.var("temp").put_attr("coordinates","x_rho y_rho s_rho ocean_time");
363 ncf.var("temp").put_attr("field","temperature, scalar, series");
364 }
365 } // end temp
366
367 {
368 int comp = -1;
369 for (int i = 0; i < plot_var_names_3d.size(); i++) {
370 if (plot_var_names_3d[i] == "salt") comp = i;
371 }
372 if (comp >= 0) {
373 ncf.def_var_fill("salt", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_r_name, nx_r_name }, &netcdf_fill_value);
374 ncf.var("salt").put_attr("long_name","salinity");
375 ncf.var("salt").put_attr("time","ocean_time");
376 ncf.var("salt").put_attr("grid","grid");
377 ncf.var("salt").put_attr("location","face");
378 ncf.var("salt").put_attr("coordinates","x_rho y_rho s_rho ocean_time");
379 ncf.var("salt").put_attr("field","salinity, scalar, series");
380 }
381 } // end salt
382
383 {
384 int comp = -1;
385 for (int i = 0; i < plot_var_names_3d.size(); i++) {
386 if (plot_var_names_3d[i] == "tracer") comp = i;
387 }
388 if (comp >= 0) {
389 ncf.def_var_fill("tracer", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_r_name, nx_r_name }, &netcdf_fill_value);
390 ncf.var("tracer").put_attr("long_name","passive tracer");
391 ncf.var("tracer").put_attr("time","ocean_time");
392 ncf.var("tracer").put_attr("grid","grid");
393 ncf.var("tracer").put_attr("location","face");
394 ncf.var("tracer").put_attr("coordinates","x_rho y_rho s_rho ocean_time");
395 ncf.var("tracer").put_attr("field","tracer, scalar, series");
396 }
397 } // end tracer
398
399 {
400 int comp = -1;
401 for (int i = 0; i < plot_var_names_3d.size(); i++) {
402 if (plot_var_names_3d[i] == "vorticity") comp = i;
403 }
404 if (comp >= 0) {
405 ncf.def_var_fill("vorticity", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_r_name, nx_r_name }, &netcdf_fill_value);
406 ncf.var("vorticity").put_attr("long_name","vorticity");
407 ncf.var("vorticity").put_attr("time","ocean_time");
408 ncf.var("vorticity").put_attr("grid","grid");
409 ncf.var("vorticity").put_attr("location","face");
410 ncf.var("vorticity").put_attr("coordinates","x_rho y_rho s_rho ocean_time");
411 ncf.var("vorticity").put_attr("field","vorticity, scalar, series");
412 }
413 } // end vorticity
414
415 ncf.def_var_fill("u", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_u_name, nx_u_name }, &netcdf_fill_value);
416 ncf.var("u").put_attr("long_name","u-momentum component");
417 ncf.var("u").put_attr("units","meter second-1");
418 ncf.var("u").put_attr("time","ocean_time");
419 ncf.var("u").put_attr("grid","grid");
420 ncf.var("u").put_attr("location","edge1");
421 ncf.var("u").put_attr("coordinates","x_u y_u s_rho ocean_time");
422 ncf.var("u").put_attr("field","u-velocity, scalar, series");
423
424 ncf.def_var_fill("v", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_v_name, nx_v_name }, &netcdf_fill_value);
425 ncf.var("v").put_attr("long_name","v-momentum component");
426 ncf.var("v").put_attr("units","meter second-1");
427 ncf.var("v").put_attr("time","ocean_time");
428 ncf.var("v").put_attr("grid","grid");
429 ncf.var("v").put_attr("location","edge2");
430 ncf.var("v").put_attr("coordinates","x_v y_v s_rho ocean_time");
431 ncf.var("v").put_attr("field","v-velocity, scalar, series");
432
433 ncf.def_var_fill("ubar", ncutils::NCDType::Real, { nt_name, ny_u_name, nx_u_name }, &netcdf_fill_value);
434 ncf.var("ubar").put_attr("long_name","vertically integrated u-momentum component");
435 ncf.var("ubar").put_attr("units","meter second-1");
436 ncf.var("ubar").put_attr("time","ocean_time");
437 ncf.var("ubar").put_attr("grid","grid");
438 ncf.var("ubar").put_attr("location","edge1");
439 ncf.var("ubar").put_attr("coordinates","x_u y_u ocean_time");
440 ncf.var("ubar").put_attr("field","ubar-velocity, scalar, series");
441
442 ncf.def_var_fill("vbar", ncutils::NCDType::Real, { nt_name, ny_v_name, nx_v_name }, &netcdf_fill_value);
443 ncf.var("vbar").put_attr("long_name","vertically integrated v-momentum component");
444 ncf.var("vbar").put_attr("units","meter second-1");
445 ncf.var("vbar").put_attr("time","ocean_time");
446 ncf.var("vbar").put_attr("grid","grid");
447 ncf.var("vbar").put_attr("location","edge2");
448 ncf.var("vbar").put_attr("coordinates","x_v y_v ocean_time");
449 ncf.var("vbar").put_attr("field","vbar-velocity, scalar, series");
450
451 ncf.def_var("sustr", ncutils::NCDType::Real, { nt_name, ny_u_name, nx_u_name });
452 ncf.var("sustr").put_attr("long_name","surface u-momentum stress");
453 ncf.var("sustr").put_attr("units","newton meter-2");
454 ncf.var("sustr").put_attr("time","ocean_time");
455 ncf.var("sustr").put_attr("grid","grid");
456 ncf.var("sustr").put_attr("location","edge1");
457 ncf.var("sustr").put_attr("coordinates","x_u y_u ocean_time");
458 ncf.var("sustr").put_attr("field","surface u-momentum stress, scalar, series");
459
460 ncf.def_var("svstr", ncutils::NCDType::Real, { nt_name, ny_v_name, nx_v_name });
461 ncf.var("svstr").put_attr("long_name","surface v-momentum stress");
462 ncf.var("svstr").put_attr("units","newton meter-2");
463 ncf.var("svstr").put_attr("time","ocean_time");
464 ncf.var("svstr").put_attr("grid","grid");
465 ncf.var("svstr").put_attr("location","edge2");
466 ncf.var("svstr").put_attr("coordinates","x_v y_v ocean_time");
467 ncf.var("svstr").put_attr("field","surface v-momentum stress, scalar, series");
468
470 // Surface air temperature (Celsius)
471 ncf.def_var("Tair", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
472 ncf.var("Tair").put_attr("long_name","surface air temperature");
473 ncf.var("Tair").put_attr("units","Celsius");
474 ncf.var("Tair").put_attr("time","ocean_time");
475 ncf.var("Tair").put_attr("grid","grid");
476 ncf.var("Tair").put_attr("location","face");
477 ncf.var("Tair").put_attr("coordinates","x_rho y_rho ocean_time");
478 ncf.var("Tair").put_attr("field","Tair, scalar, series");
479
480 // Surface air pressure (Pascal)
481 ncf.def_var("Pair", ncutils::NCDType::Real,{ nt_name, ny_r_name, nx_r_name });
482 ncf.var("Pair").put_attr("long_name","surface air pressure");
483 ncf.var("Pair").put_attr("units","Pascal");
484 ncf.var("Pair").put_attr("time","ocean_time");
485 ncf.var("Pair").put_attr("grid","grid");
486 ncf.var("Pair").put_attr("location","face");
487 ncf.var("Pair").put_attr("coordinates","x_rho y_rho ocean_time");
488 ncf.var("Pair").put_attr("field","Pair, scalar, series");
489
490 // Surface net heat flux (W/m2)
491 ncf.def_var("qnet", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
492 ncf.var("qnet").put_attr("long_name","surface net heat flux");
493 ncf.var("qnet").put_attr("units","watt meter-2");
494 ncf.var("qnet").put_attr("time","ocean_time");
495 ncf.var("qnet").put_attr("grid","grid");
496 ncf.var("qnet").put_attr("location","face");
497 ncf.var("qnet").put_attr("coordinates","x_rho y_rho ocean_time");
498 ncf.var("qnet").put_attr("field","surface heat flux, scalar, series");
499
500 // Surface net salt flux (kinematic)
501 ncf.def_var("ssflux", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
502 ncf.var("ssflux").put_attr("long_name","kinematic surface net salt flux, SALT*(E-P)/rhow");
503 ncf.var("ssflux").put_attr("units","meter second-1");
504 ncf.var("ssflux").put_attr("time","ocean_time");
505 ncf.var("ssflux").put_attr("grid","grid");
506 ncf.var("ssflux").put_attr("location","face");
507 ncf.var("ssflux").put_attr("coordinates","x_rho y_rho ocean_time");
508 ncf.var("ssflux").put_attr("field","surface net salt flux, scalar, series");
509
510 // Latent heat flux (W/m2)
511 ncf.def_var("latent", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
512 ncf.var("latent").put_attr("long_name","net latent heat flux");
513 ncf.var("latent").put_attr("units","watt meter-2");
514 ncf.var("latent").put_attr("time","ocean_time");
515 ncf.var("latent").put_attr("grid","grid");
516 ncf.var("latent").put_attr("location","face");
517 ncf.var("latent").put_attr("coordinates","x_rho y_rho ocean_time");
518 ncf.var("latent").put_attr("field","latent heat flux, scalar, series");
519
520 // Sensible heat flux (W/m2)
521 ncf.def_var("sensible", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
522 ncf.var("sensible").put_attr("long_name","net sensible heat flux");
523 ncf.var("sensible").put_attr("units","watt meter-2");
524 ncf.var("sensible").put_attr("time","ocean_time");
525 ncf.var("sensible").put_attr("grid","grid");
526 ncf.var("sensible").put_attr("location","face");
527 ncf.var("sensible").put_attr("coordinates","x_rho y_rho ocean_time");
528 ncf.var("sensible").put_attr("field","sensible heat flux, scalar, series");
529
530 // Longwave radiation (W/m2)
531 ncf.def_var("lwrad", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
532 ncf.var("lwrad").put_attr("long_name","net longwave radiation flux");
533 ncf.var("lwrad").put_attr("units","watt meter-2");
534 ncf.var("lwrad").put_attr("time","ocean_time");
535 ncf.var("lwrad").put_attr("grid","grid");
536 ncf.var("lwrad").put_attr("location","face");
537 ncf.var("lwrad").put_attr("coordinates","x_rho y_rho ocean_time");
538 ncf.var("lwrad").put_attr("field","longwave radiation, scalar, series");
539
540 // Shortwave radiation (W/m2)
541 ncf.def_var("swrad", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
542 ncf.var("swrad").put_attr("long_name","solar shortwave radiation flux");
543 ncf.var("swrad").put_attr("units","watt meter-2");
544 ncf.var("swrad").put_attr("time","ocean_time");
545 ncf.var("swrad").put_attr("grid","grid");
546 ncf.var("swrad").put_attr("location","face");
547 ncf.var("swrad").put_attr("coordinates","x_rho y_rho ocean_time");
548 ncf.var("swrad").put_attr("field","shortwave radiation, scalar, series");
549
550 // Evaporation rate (kg m-2 s-1)
551 ncf.def_var("evaporation", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
552 ncf.var("evaporation").put_attr("long_name","evaporation rate");
553 ncf.var("evaporation").put_attr("units","kilogram meter-2 second-1");
554 ncf.var("evaporation").put_attr("time","ocean_time");
555 ncf.var("evaporation").put_attr("grid","grid");
556 ncf.var("evaporation").put_attr("location","face");
557 ncf.var("evaporation").put_attr("coordinates","x_rho y_rho ocean_time");
558 ncf.var("evaporation").put_attr("field","evaporation, scalar, series");
559
560 // Rain rate (kg m-2 s-1)
561 ncf.def_var("rain", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
562 ncf.var("rain").put_attr("long_name","rain fall rate");
563 ncf.var("rain").put_attr("units","kilogram meter-2 second-1");
564 ncf.var("rain").put_attr("time","ocean_time");
565 ncf.var("rain").put_attr("grid","grid");
566 ncf.var("rain").put_attr("location","face");
567 ncf.var("rain").put_attr("coordinates","x_rho y_rho ocean_time");
568 ncf.var("rain").put_attr("field","rain, scalar, series");
569 }
570 // Right now this is hard-wired to {temp, salt, tracer, u, v}
571 ncf.put_attr("space_dimension", std::vector<int> { AMREX_SPACEDIM });
572// ncf.put_attr("current_time", std::vector<double> { time });
573 ncf.put_attr("start_time", std::vector<double> { start_bdy_time });
574 ncf.put_attr("CurrentLevel", std::vector<int> { flev });
575 ncf.put_attr("DefaultGeometry", std::vector<int> { amrex::DefaultGeometry().Coord() });
576
577 ncf.exit_def_mode();
578
579 // We are doing single-level writes but it doesn't have to be level 0
580 //
581 // Write out the header information.
582 //
583
584 Real dx[AMREX_SPACEDIM];
585 for (int i = 0; i < AMREX_SPACEDIM; i++) {
586 dx[i] = geom[lev].CellSize()[i];
587 }
588 const auto *base = geom[lev].ProbLo();
589 RealBox rb(subdomain, dx, base);
590
591 amrex::Vector<Real> probLo;
592 amrex::Vector<Real> probHi;
593 for (int i = 0; i < AMREX_SPACEDIM; i++) {
594 probLo.push_back(rb.lo(i));
595 probHi.push_back(rb.hi(i));
596 }
597
598 //nc_probLo.par_access(NC_COLLECTIVE);
599 // small variable data written by just the master proc
600 ncmpi_begin_indep_data(ncf.ncid);
601 if (amrex::ParallelDescriptor::IOProcessor()) // only master proc
602 {
603 auto nc_probLo = ncf.var("probLo");
604
605 nc_probLo.put(probLo.data(), { 0 }, { AMREX_SPACEDIM });
606
607 auto nc_probHi = ncf.var("probHi");
608 //nc_probHi.par_access(NC_COLLECTIVE);
609 nc_probHi.put(probHi.data(), { 0 }, { AMREX_SPACEDIM });
610
611 amrex::Vector<int> smallend;
612 amrex::Vector<int> bigend;
613 for (int i = lev; i < flev; i++) {
614 smallend.clear();
615 bigend.clear();
616 for (int j = 0; j < AMREX_SPACEDIM; j++) {
617 smallend.push_back(subdomain.smallEnd(j));
618 bigend.push_back(subdomain.bigEnd(j));
619 }
620 auto nc_Geom_smallend = ncf.var("Geom.smallend");
621 //nc_Geom_smallend.par_access(NC_COLLECTIVE);
622 nc_Geom_smallend.put(smallend.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
623 AMREX_SPACEDIM });
624
625 auto nc_Geom_bigend = ncf.var("Geom.bigend");
626 //nc_Geom_bigend.par_access(NC_COLLECTIVE);
627 nc_Geom_bigend.put(bigend.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
628 AMREX_SPACEDIM });
629 }
630
631 amrex::Vector<Real> CellSize;
632 for (int i = lev; i < flev; i++) {
633 CellSize.clear();
634 for (Real &j : dx) {
635 CellSize.push_back(amrex::Real(j));
636 }
637 auto nc_CellSize = ncf.var("CellSize");
638 //nc_CellSize.par_access(NC_COLLECTIVE);
639 nc_CellSize.put(CellSize.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
640 AMREX_SPACEDIM });
641 }
642 Real hc = solverChoice.tcline;
643 Real theta_s = solverChoice.theta_s;
644 Real theta_b = solverChoice.theta_b;
645 ncf.var("hc").put(&hc);
646 ncf.var("theta_s").put(&theta_s);
647 ncf.var("theta_b").put(&theta_b);
648
649 }
650 ncmpi_end_indep_data(ncf.ncid);
651
652 } // end if write_header
653
654 ncmpi_begin_indep_data(ncf.ncid);
655 //
656 // We compute the offsets based on location of the box within the domain
657 //
658 long long adjusted_history_count = chunk_history_file ? history_count % steps_per_history_file : history_count;
659 long long local_start_nt = (is_history ? static_cast<long long>(adjusted_history_count) : static_cast<long long>(0));
660 long long local_nt = 1; // We write data for only one time
661
662 {
663 auto nc_plot_var = ncf.var("ocean_time");
664 //nc_plot_var.par_access(NC_COLLECTIVE);
665 nc_plot_var.put(&t_new[lev], { local_start_nt }, { local_nt });
666 }
667 // do all independent writes
668 //ncmpi_end_indep_data(ncf.ncid);
669
670 mask_arrays_for_write(lev, (Real) netcdf_fill_value, 0.0_rt);
671
672 // Check whether there are any nans or infs in variables that we will write out
673 if (vec_Zt_avg1[lev]->contains_nan() || vec_Zt_avg1[lev]->contains_inf()) {
674 amrex::Abort("Found while writing output: zeta contains nan or inf");
675 }
676 if (plotMF->contains_nan(Temp_comp,1) || plotMF->contains_inf(Temp_comp,1)) {
677 amrex::Abort("Found while writing output: Temperature contains nan or inf");
678 }
679 if (plotMF->contains_nan(Salt_comp,1) || plotMF->contains_inf(Salt_comp,1)) {
680 amrex::Abort("Found while writing output: Salinity contains nan or inf");
681 }
682 if (plotMF->contains_nan(Tracer_comp,1) || plotMF->contains_inf(Tracer_comp,1)) {
683 amrex::Abort("Found while writing output: Passive tracer contains nan or inf");
684 }
685 if (xvel_new[lev]->contains_nan() || xvel_new[lev]->contains_inf()) {
686 amrex::Abort("Found while writing output: velocity u contains nan or inf");
687 }
688 if (vec_ubar[lev]->contains_nan(0,1) || vec_ubar[lev]->contains_inf(0,1)) {
689 amrex::Abort("Found while writing output: velocity ubar contains nan or inf");
690 }
691 if (yvel_new[lev]->contains_nan() || yvel_new[lev]->contains_inf()) {
692 amrex::Abort("Found while writing output: velocity v contains nan or inf");
693 }
694 if (vec_vbar[lev]->contains_nan(0,1) || vec_vbar[lev]->contains_inf(0,1)) {
695 amrex::Abort("Found while writing output: velocity vbar contains nan or inf");
696 }
697
698 for (MFIter mfi(*plotMF, false); mfi.isValid(); ++mfi) {
699 auto bx = mfi.validbox();
700 if (subdomain.contains(bx)) {
701 //
702 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
703 //
704 Box tmp_bx(bx);
705 if (tmp_bx.smallEnd()[0] == subdomain.smallEnd()[0])
706 tmp_bx.growLo(0, 1);
707 if (tmp_bx.smallEnd()[1] == subdomain.smallEnd()[1])
708 tmp_bx.growLo(1, 1);
709 if (tmp_bx.bigEnd()[0] == subdomain.bigEnd()[0])
710 tmp_bx.growHi(0, 1);
711 if (tmp_bx.bigEnd()[1] == subdomain.bigEnd()[1])
712 tmp_bx.growHi(1, 1);
713 // amrex::Print() << " BX " << bx << std::endl;
714 // amrex::Print() << "TMP_BX " << tmp_bx << std::endl;
715
716 Box tmp_bx_2d(tmp_bx);
717 tmp_bx_2d.makeSlab(2, 0);
718
719 Box tmp_bx_1d(tmp_bx);
720 tmp_bx_1d.makeSlab(0, 0);
721 tmp_bx_1d.makeSlab(1, 0);
722
723 //
724 // These are the dimensions of the data we write for only this box
725 //
726 long long local_nx = tmp_bx.length()[0];
727 long long local_ny = tmp_bx.length()[1];
728 long long local_nz = tmp_bx.length()[2];
729
730 // We do the "+1" because the offset needs to start at 0
731 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0] + 1);
732 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1] + 1);
733 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
734
735 if (write_header) {
736 // Only write out s_rho and s_w at x=0,y=0 to avoid NaNs
737 if (bx.contains(IntVect(0,0,0)))
738 {
739 {
740 amrex::Vector<amrex::Real> tmp_srho(local_nz);
741
742#ifdef AMREX_USE_GPU
743 Gpu::dtoh_memcpy(tmp_srho.data(), s_r.data(), sizeof(amrex::Real)*local_nz);
744#else
745 std::memcpy(tmp_srho.data(), s_r.data(), sizeof(amrex::Real)*local_nz);
746#endif
747 Gpu::streamSynchronize();
748
749 auto nc_plot_var = ncf.var("s_rho");
750 //nc_plot_var.par_access(NC_INDEPENDENT);
751 nc_plot_var.put(tmp_srho.data(), { local_start_z }, { local_nz });
752 }
753 {
754 amrex::Vector<amrex::Real> tmp_sw(local_nz+1);
755
756#ifdef AMREX_USE_GPU
757 Gpu::dtoh_memcpy(tmp_sw.data(), s_w.data(), sizeof(amrex::Real)*(local_nz+1));
758#else
759 std::memcpy(tmp_sw.data(), s_w.data(), sizeof(amrex::Real)*(local_nz+1));
760#endif
761 Gpu::streamSynchronize();
762
763 auto nc_plot_var = ncf.var("s_w");
764 //nc_plot_var.par_access(NC_INDEPENDENT);
765 nc_plot_var.put(tmp_sw.data(), { local_start_z }, { local_nz + 1});
766 }
767 {
768 amrex::Vector<amrex::Real> tmp_Csrho(local_nz);
769
770#ifdef AMREX_USE_GPU
771 Gpu::dtoh_memcpy(tmp_Csrho.data(), Cs_r.data(), sizeof(amrex::Real)*(local_nz));
772#else
773 std::memcpy(tmp_Csrho.data(), Cs_r.data(), sizeof(amrex::Real)*(local_nz));
774#endif
775 Gpu::streamSynchronize();
776
777 auto nc_plot_var = ncf.var("Cs_r");
778 //nc_plot_var.par_access(NC_INDEPENDENT);
779 nc_plot_var.put(tmp_Csrho.data(), { local_start_z }, { local_nz });
780 }
781 {
782 amrex::Vector<amrex::Real> tmp_Csw(local_nz+1);
783
784#ifdef AMREX_USE_GPU
785 Gpu::dtoh_memcpy(tmp_Csw.data(), Cs_w.data(), sizeof(amrex::Real)*(local_nz+1));
786#else
787 std::memcpy(tmp_Csw.data(), Cs_w.data(), sizeof(amrex::Real)*(local_nz+1));
788#endif
789
790 Gpu::streamSynchronize();
791
792 auto nc_plot_var = ncf.var("Cs_w");
793 //nc_plot_var.par_access(NC_INDEPENDENT);
794 nc_plot_var.put(tmp_Csw.data(), { local_start_z }, { local_nz + 1});
795 }
796 }
797
798 {
799 FArrayBox tmp_bathy;
800 tmp_bathy.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
801
802 tmp_bathy.template copy<RunOn::Device>((*vec_h[lev])[mfi.index()], 0, 0, 1);
803 Gpu::streamSynchronize();
804
805 auto nc_plot_var = ncf.var("h");
806 //nc_plot_var.par_access(NC_INDEPENDENT);
807 nc_plot_var.put(tmp_bathy.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
808 }
809
810 {
811 FArrayBox tmp_pm;
812 tmp_pm.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
813
814 tmp_pm.template copy<RunOn::Device>((*vec_pm[lev])[mfi.index()], 0, 0, 1);
815 Gpu::streamSynchronize();
816
817 auto nc_plot_var = ncf.var("pm");
818 //nc_plot_var.par_access(NC_INDEPENDENT);
819
820 nc_plot_var.put(tmp_pm.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
821 }
822
823 {
824 FArrayBox tmp_pn;
825 tmp_pn.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
826
827 tmp_pn.template copy<RunOn::Device>((*vec_pn[lev])[mfi.index()], 0, 0, 1);
828 Gpu::streamSynchronize();
829
830 auto nc_plot_var = ncf.var("pn");
831 //nc_plot_var.par_access(NC_INDEPENDENT);
832 nc_plot_var.put(tmp_pn.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
833 }
834
835 {
836 FArrayBox tmp_f;
837 tmp_f.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
838
839 tmp_f.template copy<RunOn::Device>((*vec_fcor[lev])[mfi.index()], 0, 0, 1);
840 Gpu::streamSynchronize();
841
842 auto nc_plot_var = ncf.var("f");
843 //nc_plot_var.par_access(NC_INDEPENDENT);
844 nc_plot_var.put(tmp_f.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
845 }
846
847 {
848 FArrayBox tmp_xr;
849 tmp_xr.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
850
851 tmp_xr.template copy<RunOn::Device>((*vec_xr[lev])[mfi.index()], 0, 0, 1);
852 Gpu::streamSynchronize();
853
854 auto nc_plot_var = ncf.var("x_rho");
855 //nc_plot_var.par_access(NC_INDEPENDENT);
856
857 nc_plot_var.put(tmp_xr.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
858 }
859
860 {
861 FArrayBox tmp_yr;
862 tmp_yr.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
863
864 tmp_yr.template copy<RunOn::Device>((*vec_yr[lev])[mfi.index()], 0, 0, 1);
865 Gpu::streamSynchronize();
866
867 auto nc_plot_var = ncf.var("y_rho");
868 //nc_plot_var.par_access(NC_INDEPENDENT);
869 nc_plot_var.put(tmp_yr.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
870 }
871 }
872
873 {
874 FArrayBox tmp_zeta;
875 tmp_zeta.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
876 tmp_zeta.template copy<RunOn::Device>((*vec_Zt_avg1[lev])[mfi.index()], 0, 0, 1);
877 Gpu::streamSynchronize();
878
879 auto nc_plot_var = ncf.var("zeta");
880 nc_plot_var.put(tmp_zeta.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny,
881 local_nx });
882 }
883
885 {
886 const Real Hscale = solverChoice.rho0 * Cp;
887 // Tair
888 {
889 FArrayBox tmp_Tair;
890 tmp_Tair.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
891 tmp_Tair.template copy<RunOn::Device>((*vec_Tair[lev])[mfi.index()], 0, 0, 1);
892 Gpu::streamSynchronize();
893
894 auto nc_plot_var = ncf.var("Tair");
895 nc_plot_var.put(tmp_Tair.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
896 }
897 // Pair
898 {
899 FArrayBox tmp_Pair;
900 tmp_Pair.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
901 tmp_Pair.template copy<RunOn::Device>((*vec_Pair[lev])[mfi.index()], 0, 0, 1);
902 Gpu::streamSynchronize();
903
904 auto nc_plot_var = ncf.var("Pair");
905 nc_plot_var.put(tmp_Pair.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
906 }
907 // qnet (stored °C m/s → write W/m²)
908 {
909 FArrayBox tmp;
910 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
911
912 // Copy stflux Temp component
913 tmp.template copy<RunOn::Device>(
914 (*vec_stflux[lev])[mfi.index()],
915 Temp_comp, // source component
916 0, // dest component
917 1 // number of comps
918 );
919
920 Gpu::streamSynchronize();
921
922 // Convert °C·m/s → W/m²
923 tmp.mult<RunOn::Device>(Hscale);
924
925 auto nc_var = ncf.var("qnet");
926 nc_var.put(tmp.dataPtr(),
927 { local_start_nt, local_start_y, local_start_x },
928 { local_nt, local_ny, local_nx });
929 }
930 // ssflux = surface net freshwater flux (kg/m²/s converted to m/s)
931 {
932 FArrayBox tmp;
933 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
934
935 // Copy stflux Salt component
936 tmp.template copy<RunOn::Device>(
937 (*vec_stflux[lev])[mfi.index()],
938 Salt_comp, // source component
939 0, // destination component
940 1 // number of components
941 );
942
943 Gpu::streamSynchronize();
944
945 auto nc_var = ncf.var("ssflux");
946 nc_var.put(tmp.dataPtr(),
947 { local_start_nt, local_start_y, local_start_x },
948 { local_nt, local_ny, local_nx });
949 }
950 // latent (stored °C m/s → write W/m²)
951 {
952 FArrayBox tmp;
953 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
954 tmp.template copy<RunOn::Device>((*vec_lhflx[lev])[mfi.index()], 0, 0, 1);
955 Gpu::streamSynchronize();
956
957 // Convert °C·m/s → W/m²
958 tmp.mult<RunOn::Device>(Hscale);
959
960 auto nc_var = ncf.var("latent");
961 nc_var.put(tmp.dataPtr(),
962 { local_start_nt, local_start_y, local_start_x },
963 { local_nt, local_ny, local_nx });
964 }
965 // sensible (stored °C m/s → write W/m²)
966 {
967 FArrayBox tmp;
968 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
969 tmp.template copy<RunOn::Device>((*vec_shflx[lev])[mfi.index()], 0, 0, 1);
970 Gpu::streamSynchronize();
971
972 // Convert °C·m/s → W/m²
973 tmp.mult<RunOn::Device>(Hscale);
974
975 auto nc_var = ncf.var("sensible");
976 nc_var.put(tmp.dataPtr(),
977 { local_start_nt, local_start_y, local_start_x },
978 { local_nt, local_ny, local_nx });
979 }
980 // lwrad (stored °C m/s → write W/m²)
981 {
982 FArrayBox tmp;
983 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
984 tmp.template copy<RunOn::Device>((*vec_lrflx[lev])[mfi.index()], 0, 0, 1);
985 Gpu::streamSynchronize();
986
987 // Convert °C·m/s → W/m²
988 tmp.mult<RunOn::Device>(Hscale);
989
990 auto nc_var = ncf.var("lwrad");
991 nc_var.put(tmp.dataPtr(),
992 { local_start_nt, local_start_y, local_start_x },
993 { local_nt, local_ny, local_nx });
994 }
995 // swrad, note this is stored explicitly as W/m², not degC m/s in REMORA.bulk_flux.cpp
996 {
997 FArrayBox tmp;
998 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
999 tmp.template copy<RunOn::Device>((*vec_srflx[lev])[mfi.index()], 0, 0, 1);
1000 Gpu::streamSynchronize();
1001
1002 auto nc_var = ncf.var("swrad");
1003 nc_var.put(tmp.dataPtr(),
1004 { local_start_nt, local_start_y, local_start_x },
1005 { local_nt, local_ny, local_nx });
1006 }
1007 // evaporation
1008 {
1009 FArrayBox tmp;
1010 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1011 tmp.template copy<RunOn::Device>((*vec_evap[lev])[mfi.index()], 0, 0, 1);
1012 Gpu::streamSynchronize();
1013
1014 auto nc_var = ncf.var("evaporation");
1015 nc_var.put(tmp.dataPtr(),
1016 { local_start_nt, local_start_y, local_start_x },
1017 { local_nt, local_ny, local_nx });
1018 }
1019 // rain
1020 {
1021 FArrayBox tmp;
1022 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1023 tmp.template copy<RunOn::Device>((*vec_rain[lev])[mfi.index()], 0, 0, 1);
1024 Gpu::streamSynchronize();
1025
1026 auto nc_var = ncf.var("rain");
1027 nc_var.put(tmp.dataPtr(),
1028 { local_start_nt, local_start_y, local_start_x },
1029 { local_nt, local_ny, local_nx });
1030 }
1031 } // end output forcing
1032
1033 // **************************************************************************
1034 { // Temp
1035 int comp = -1;
1036 for (int i = 0; i < plot_var_names_3d.size(); i++) {
1037 if (plot_var_names_3d[i] == "temp") comp = i;
1038 }
1039 if (comp >= 0) {
1040 FArrayBox tmp;
1041 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1042 tmp.template copy<RunOn::Device>((*plotMF)[mfi.index()], comp, 0, 1);
1043 Gpu::streamSynchronize();
1044
1045 auto nc_plot_var = ncf.var(plot_var_names_3d[comp]);
1046 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1047 local_nz, local_ny, local_nx });
1048 } // if temp exists in plotMF
1049 } // end temp
1050 // **************************************************************************
1051
1052 // **************************************************************************
1053 { // Salt
1054 int comp = -1;
1055 for (int i = 0; i < plot_var_names_3d.size(); i++) {
1056 if (plot_var_names_3d[i] == "salt") comp = i;
1057 }
1058 if (comp >= 0) {
1059 FArrayBox tmp;
1060 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1061 tmp.template copy<RunOn::Device>((*plotMF)[mfi.index()], comp, 0, 1);
1062 Gpu::streamSynchronize();
1063
1064 auto nc_plot_var = ncf.var(plot_var_names_3d[comp]);
1065 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1066 local_nz, local_ny, local_nx });
1067 } // if salt exists in plotMF
1068 } // end salt
1069 // **************************************************************************
1070
1071 // **************************************************************************
1072 { // Tracer
1073 int comp = -1;
1074 for (int i = 0; i < plot_var_names_3d.size(); i++) {
1075 if (plot_var_names_3d[i] == "tracer") comp = i;
1076 }
1077 if (comp >= 0) {
1078 FArrayBox tmp;
1079 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1080 tmp.template copy<RunOn::Device>((*plotMF)[mfi.index()], comp, 0, 1);
1081 Gpu::streamSynchronize();
1082
1083 auto nc_plot_var = ncf.var(plot_var_names_3d[comp]);
1084 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1085 local_nz, local_ny, local_nx });
1086 } // if tracer exists in plotMF
1087 } // end tracer
1088 // **************************************************************************
1089
1090 // **************************************************************************
1091 { // Vorticity
1092 int comp = -1;
1093 for (int i = 0; i < plot_var_names_3d.size(); i++) {
1094 if (plot_var_names_3d[i] == "vorticity") comp = i;
1095 }
1096 if (comp >= 0) {
1097 FArrayBox tmp;
1098 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1099 tmp.template copy<RunOn::Device>((*plotMF)[mfi.index()], comp, 0, 1);
1100 Gpu::streamSynchronize();
1101
1102 auto nc_plot_var = ncf.var(plot_var_names_3d[comp]);
1103 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1104 local_nz, local_ny, local_nx });
1105 } // if vorticity exists in plotMF
1106 } // end vorticity
1107 // **************************************************************************
1108
1109 } // subdomain
1110 } // mfi
1111
1112 //ncf.wait_all(irq, &requests[0]);
1113 //requests.resize(0);
1114 //irq = 0;
1115 // Writing u (we loop over cons to get cell-centered box)
1116 for (MFIter mfi(*plotMF, false); mfi.isValid(); ++mfi) {
1117 Box bx = mfi.validbox();
1118
1119 if (subdomain.contains(bx)) {
1120 //
1121 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
1122 //
1123 Box tmp_bx(bx);
1124 tmp_bx.surroundingNodes(0);
1125 if (tmp_bx.smallEnd()[1] == subdomain.smallEnd()[1])
1126 tmp_bx.growLo(1, 1);
1127 if (tmp_bx.bigEnd()[1] == subdomain.bigEnd()[1])
1128 tmp_bx.growHi(1, 1);
1129 Box tmp_bx_2d(tmp_bx);
1130 tmp_bx_2d.makeSlab(2, 0);
1131
1132 //
1133 // These are the dimensions of the data we write for only this box
1134 //
1135 long long local_nx = tmp_bx.length()[0];
1136 long long local_ny = tmp_bx.length()[1];
1137 long long local_nz = tmp_bx.length()[2];
1138
1139 // We do the "+1" because the offset needs to start at 0
1140 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0]);
1141 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1] + 1);
1142 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
1143
1144 if (write_header) {
1145 {
1146 FArrayBox tmp;
1147 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1148 tmp.template copy<RunOn::Device>((*vec_xu[lev])[mfi.index()], 0, 0, 1);
1149 Gpu::streamSynchronize();
1150
1151 auto nc_plot_var = ncf.var("x_u");
1152 //nc_plot_var.par_access(NC_INDEPENDENT);
1153 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1154 }
1155 {
1156 FArrayBox tmp;
1157 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1158 tmp.template copy<RunOn::Device>((*vec_yu[lev])[mfi.index()], 0, 0, 1);
1159 Gpu::streamSynchronize();
1160
1161 auto nc_plot_var = ncf.var("y_u");
1162 //nc_plot_var.par_access(NC_INDEPENDENT);
1163 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1164 }
1165 }
1166
1167 {
1168 FArrayBox tmp;
1169 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1170 tmp.template copy<RunOn::Device>((*xvel_new[lev])[mfi.index()], 0, 0, 1);
1171 Gpu::streamSynchronize();
1172
1173 auto nc_plot_var = ncf.var("u");
1174 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1175 local_nz, local_ny, local_nx });
1176 }
1177
1178 {
1179 FArrayBox tmp;
1180 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1181 tmp.template copy<RunOn::Device>((*vec_ubar[lev])[mfi.index()], 0, 0, 1);
1182 Gpu::streamSynchronize();
1183
1184 auto nc_plot_var = ncf.var("ubar");
1185 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
1186 }
1187 {
1188 FArrayBox tmp;
1189 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1190 tmp.template copy<RunOn::Device>((*vec_sustr[lev])[mfi.index()], 0, 0, 1);
1191 Gpu::streamSynchronize();
1192
1193 auto nc_plot_var = ncf.var("sustr");
1194 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
1195 }
1196 } // in subdomain
1197 } // mfi
1198
1199 // Writing v (we loop over cons to get cell-centered box)
1200 for (MFIter mfi(*plotMF, false); mfi.isValid(); ++mfi) {
1201 Box bx = mfi.validbox();
1202
1203 if (subdomain.contains(bx)) {
1204 //
1205 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
1206 //
1207 Box tmp_bx(bx);
1208 tmp_bx.surroundingNodes(1);
1209 if (tmp_bx.smallEnd()[0] == subdomain.smallEnd()[0])
1210 tmp_bx.growLo(0, 1);
1211 if (tmp_bx.bigEnd()[0] == subdomain.bigEnd()[0])
1212 tmp_bx.growHi(0, 1);
1213 // amrex::Print() << " BX " << bx << std::endl;
1214 // amrex::Print() << "TMP_BX " << tmp_bx << std::endl;
1215
1216 Box tmp_bx_2d(tmp_bx);
1217 tmp_bx_2d.makeSlab(2, 0);
1218
1219 //
1220 // These are the dimensions of the data we write for only this box
1221 //
1222 long long local_nx = tmp_bx.length()[0];
1223 long long local_ny = tmp_bx.length()[1];
1224 long long local_nz = tmp_bx.length()[2];
1225
1226 // We do the "+1" because the offset needs to start at 0
1227 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0] + 1);
1228 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1]);
1229 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
1230
1231 if (write_header) {
1232 {
1233 FArrayBox tmp;
1234 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1235 tmp.template copy<RunOn::Device>((*vec_xv[lev])[mfi.index()], 0, 0, 1);
1236 Gpu::streamSynchronize();
1237
1238 auto nc_plot_var = ncf.var("x_v");
1239 //nc_plot_var.par_access(NC_INDEPENDENT);
1240 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1241 }
1242 {
1243 FArrayBox tmp;
1244 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1245 tmp.template copy<RunOn::Device>((*vec_yv[lev])[mfi.index()], 0, 0, 1);
1246 Gpu::streamSynchronize();
1247
1248 auto nc_plot_var = ncf.var("y_v");
1249 //nc_plot_var.par_access(NC_INDEPENDENT);
1250 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1251 }
1252 }
1253
1254 {
1255 FArrayBox tmp;
1256 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1257 tmp.template copy<RunOn::Device>((*yvel_new[lev])[mfi.index()], 0, 0, 1);
1258 Gpu::streamSynchronize();
1259
1260 auto nc_plot_var = ncf.var("v");
1261 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1262 local_nz, local_ny, local_nx });
1263 }
1264
1265 {
1266 FArrayBox tmp;
1267 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1268 tmp.template copy<RunOn::Device>((*vec_vbar[lev])[mfi.index()], 0, 0, 1);
1269 Gpu::streamSynchronize();
1270
1271 auto nc_plot_var = ncf.var("vbar");
1272 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
1273 }
1274
1275 {
1276 FArrayBox tmp;
1277 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1278 tmp.template copy<RunOn::Device>((*vec_svstr[lev])[mfi.index()], 0, 0, 1);
1279 Gpu::streamSynchronize();
1280
1281 auto nc_plot_var = ncf.var("svstr");
1282 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
1283 }
1284
1285 } // in subdomain
1286 } // mfi
1287
1288 for (MFIter mfi(*plotMF, false); mfi.isValid(); ++mfi) {
1289 Box bx = mfi.validbox();
1290
1291 if (subdomain.contains(bx)) {
1292 //
1293 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
1294 //
1295 Box tmp_bx(bx);
1296 tmp_bx.surroundingNodes(0);
1297 tmp_bx.surroundingNodes(1);
1298
1299 Box tmp_bx_2d(tmp_bx);
1300 tmp_bx_2d.makeSlab(2, 0);
1301
1302 //
1303 // These are the dimensions of the data we write for only this box
1304 //
1305 long long local_nx = tmp_bx.length()[0];
1306 long long local_ny = tmp_bx.length()[1];
1307
1308 // We do the "+1" because the offset needs to start at 0
1309 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0]);
1310 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1]);
1311
1312 if (write_header) {
1313 {
1314 FArrayBox tmp;
1315 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1316 tmp.template copy<RunOn::Device>((*vec_xp[lev])[mfi.index()], 0, 0, 1);
1317 Gpu::streamSynchronize();
1318
1319 auto nc_plot_var = ncf.var("x_psi");
1320 //nc_plot_var.par_access(NC_INDEPENDENT);
1321 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1322 }
1323 {
1324 FArrayBox tmp;
1325 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1326 tmp.template copy<RunOn::Device>((*vec_yp[lev])[mfi.index()], 0, 0, 1);
1327 Gpu::streamSynchronize();
1328
1329 auto nc_plot_var = ncf.var("y_psi");
1330 //nc_plot_var.par_access(NC_INDEPENDENT);
1331 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1332 }
1333
1334 } // header
1335 } // in subdomain
1336 } // mfi
1337
1339
1340 ncf.close();
1341
1343}
1344
constexpr amrex::Real Cp
#define Temp_comp
#define Tracer_comp
#define Salt_comp
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_evap
evaporation rate [kg/m^2/s]
Definition REMORA.H:349
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:405
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:253
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:400
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lrflx
longwave radiation
Definition REMORA.H:329
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yv
y_grid on v-points (2D)
Definition REMORA.H:420
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
Definition REMORA.H:1097
amrex::Gpu::DeviceVector< amrex::Real > s_w
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at w-points (cell faces)
Definition REMORA.H:292
int history_count
Counter for which time index we are writing to in the netcdf history file.
Definition REMORA.H:1323
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rain
precipitation rate [kg/m^2/s]
Definition REMORA.H:347
bool chunk_history_file
Whether to chunk netcdf history file.
Definition REMORA.H:1316
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:311
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yp
y_grid on psi-points (2D)
Definition REMORA.H:425
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xr
x_grid on rho points (2D)
Definition REMORA.H:408
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xv
x_grid on v-points (2D)
Definition REMORA.H:418
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1204
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:245
amrex::Real start_bdy_time
Start time in the time series of boundary data.
Definition REMORA.H:1092
amrex::Gpu::DeviceVector< amrex::Real > s_r
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at rho points (cell centers)
Definition REMORA.H:290
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_shflx
sensible heat flux
Definition REMORA.H:335
int steps_per_history_file
Number of time steps per netcdf history file.
Definition REMORA.H:1321
int max_step
maximum number of steps
Definition REMORA.H:1271
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:243
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lhflx
latent heat flux
Definition REMORA.H:333
int plot_int
Plotfile output interval in iterations.
Definition REMORA.H:1305
void mask_arrays_for_write(int lev, amrex::Real fill_value, amrex::Real fill_where)
Mask data arrays before writing output.
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
Definition REMORA.H:1380
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:313
void WriteNCPlotFile(int istep, amrex::MultiFab const *plotMF)
Write plotfile using NetCDF (wrapper)
amrex::Gpu::DeviceVector< amrex::Real > Cs_r
Stretching coefficients at rho points.
Definition REMORA.H:300
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1211
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1336
static int total_nc_plot_file_step
Definition REMORA.H:1008
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xp
x_grid on psi-points (2D)
Definition REMORA.H:423
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:384
amrex::Gpu::DeviceVector< amrex::Real > Cs_w
Stretching coefficients at w points.
Definition REMORA.H:302
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:382
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yu
y_grid on u-points (2D)
Definition REMORA.H:415
void WriteNCPlotFile_which(int lev, int which_subdomain, amrex::MultiFab const *plotMF, bool write_header, ncutils::NCFile &ncf, bool is_history)
Write a particular NetCDF plotfile.
amrex::Real netcdf_fill_value
fill value for masked arrays in netcdf output
Definition REMORA.H:1344
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xu
x_grid on u-points (2D)
Definition REMORA.H:413
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:402
amrex::Vector< std::string > plot_var_names_3d
Names of 3D variables to output to AMReX plotfile.
Definition REMORA.H:1326
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflux
Surface tracer flux; input arrays.
Definition REMORA.H:340
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1303
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:308
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_srflx
Shortwave radiation flux [W/m²], defined at rho-points.
Definition REMORA.H:327
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Pair
Air pressure [mb], defined at rho-points.
Definition REMORA.H:324
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Tair
Air temperature [°C], defined at rho-points.
Definition REMORA.H:320
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yr
y_grid on rho points (2D)
Definition REMORA.H:410
void put_attr(const std::string &name, const std::string &value) const
Set file attribute to value.
void enter_def_mode() const
Enter definition mode (not needed for NetCDF4 format)
NCVar def_var(const std::string &name, const nc_type dtype, const std::vector< std::string > &dnames) const
Define a variable (wrapper for def_array)
NCDim def_dim(const std::string &, const size_t len) const
Define new dimension.
static NCFile create(const std::string &name, const int cmode=NC_CLOBBER|NC_MPIIO, MPI_Comm comm=MPI_COMM_WORLD, MPI_Info info=MPI_INFO_NULL)
Create a file.
void exit_def_mode() const
Exit definition mode.
static NCFile open(const std::string &name, const int cmode=NC_NOWRITE, MPI_Comm comm=MPI_COMM_WORLD, MPI_Info info=MPI_INFO_NULL)
Open an existing file.
void close()
Close file object.
NCVar var(const std::string &) const
Get the variable instance by name.
NCVar def_var_fill(const std::string &name, const nc_type dtype, const std::vector< std::string > &dnames, const void *fill_val) const
Define a variable (wrapper for def_array)
amrex::Real theta_b
amrex::Real theta_s
amrex::Real tcline
static constexpr nc_type Real
void put(const double *ptr) const
Write out the entire double variable.
void put_attr(const std::string &name, const std::string &value) const
Set attribute "name" to "value".