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) {
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, 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, 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, bool write_header, ncutils::NCFile &ncf, bool is_history) {
110 // Number of cells in this "domain" at this level
111 std::vector<int> n_cells;
112
113 // We only do single-level writes when using NetCDF format
114 int flev = 1; //max_level;
115
116 Box subdomain;
117 if (lev == 0) {
118 subdomain = geom[lev].Domain();
119 } else {
120 subdomain = boxes_at_level[lev][which_subdomain];
121 }
122
123 int nx = subdomain.length(0);
124 int ny = subdomain.length(1);
125 int nz = subdomain.length(2);
126
127 // unsigned long int nt= NC_UNLIMITED;
128 if (is_history && max_step < 0)
129 amrex::Abort("Need to know max_step if writing history file");
130 long long int nt = is_history ? static_cast<long long int>(max_step / std::min(plot_int, max_step)) + 1 : 1;
131 if (chunk_history_file) {
132 // First index of the last history file
133 int last_file_index = REMORA::steps_per_history_file * int(nt / REMORA::steps_per_history_file);
134 if (history_count >= last_file_index) {
135 nt = nt - last_file_index;
136 } else {
138 }
139 }
140
141 n_cells.push_back(nx);
142 n_cells.push_back(ny);
143 n_cells.push_back(nz);
144
145 const std::string nt_name = "ocean_time";
146 const std::string ndim_name = "num_geo_dimensions";
147
148 const std::string flev_name = "FINEST_LEVEL";
149
150 const std::string nx_name = "NX";
151 const std::string ny_name = "NY";
152 const std::string nz_name = "NZ";
153
154 const std::string nx_r_name = "xi_rho";
155 const std::string ny_r_name = "eta_rho";
156 const std::string nz_r_name = "s_rho";
157
158 const std::string nx_u_name = "xi_u";
159 const std::string ny_u_name = "eta_u";
160
161 const std::string nx_v_name = "xi_v";
162 const std::string ny_v_name = "eta_v";
163
164 const std::string nx_p_name = "xi_psi";
165 const std::string ny_p_name = "eta_psi";
166 const std::string nz_w_name = "s_w";
167
168 const Real fill_value = 1.0e37_rt;
169
170 if (write_header) {
171 ncf.enter_def_mode();
172 ncf.put_attr("title", "REMORA data ");
173 ncf.def_dim(nt_name, nt);
174 ncf.def_dim(ndim_name, AMREX_SPACEDIM);
175
176 ncf.def_dim(nx_r_name, nx + 2);
177 ncf.def_dim(ny_r_name, ny + 2);
178 ncf.def_dim(nz_r_name, nz);
179
180 ncf.def_dim(nx_u_name, nx + 1);
181 ncf.def_dim(ny_u_name, ny + 2);
182
183 ncf.def_dim(nx_v_name, nx + 2);
184 ncf.def_dim(ny_v_name, ny + 1);
185
186 ncf.def_dim(nx_p_name, nx + 1);
187 ncf.def_dim(ny_p_name, ny + 1);
188
189 ncf.def_dim(nz_w_name, nz + 1);
190
191 ncf.def_dim(flev_name, flev);
192
193 ncf.def_dim(nx_name, n_cells[0]);
194 ncf.def_dim(ny_name, n_cells[1]);
195 ncf.def_dim(nz_name, n_cells[2]);
196
197 ncf.def_var("probLo", ncutils::NCDType::Real, { ndim_name });
198 ncf.var("probLo").put_attr("long_name","Low side of problem domain in internal AMReX grid");
199 ncf.var("probLo").put_attr("units","meter");
200 ncf.def_var("probHi", ncutils::NCDType::Real, { ndim_name });
201 ncf.var("probHi").put_attr("long_name","High side of problem domain in internal AMReX grid");
202 ncf.var("probHi").put_attr("units","meter");
203
204 ncf.def_var("Geom.smallend", NC_INT, { flev_name, ndim_name });
205 ncf.var("Geom.smallend").put_attr("long_name","Low side of problem domain in index space");
206 ncf.def_var("Geom.bigend", NC_INT, { flev_name, ndim_name });
207 ncf.var("Geom.bigend").put_attr("long_name","High side of problem domain in index space");
208 ncf.def_var("CellSize", ncutils::NCDType::Real, { flev_name, ndim_name });
209 ncf.var("CellSize").put_attr("long_name","Cell size on internal AMReX grid");
210 ncf.var("CellSize").put_attr("units","meter");
211
212 ncf.def_var("theta_s",ncutils::NCDType::Real,{});
213 ncf.var("theta_s").put_attr("long_name","S-coordinate surface control parameter");
214 ncf.def_var("theta_b",ncutils::NCDType::Real,{});
215 ncf.var("theta_b").put_attr("long_name","S-coordinate bottom control parameter");
216 ncf.def_var("hc",ncutils::NCDType::Real,{});
217 ncf.var("hc").put_attr("long_name","S-coordinate parameter, critical depth");
218 ncf.var("hc").put_attr("units","meter");
219
220 ncf.def_var("grid",NC_INT, {});
221 ncf.var("grid").put_attr("cf_role","grid_topology");
222 ncf.var("grid").put_attr("topology_dimension",std::vector({2}));
223 ncf.var("grid").put_attr("node_dimensions", "xi_psi eta_psi");
224 ncf.var("grid").put_attr("face_dimensions", "xi_rho: xi_psi (padding: both) eta_rho: eta_psi (padding: both)");
225 ncf.var("grid").put_attr("edge1_dimensions", "xi_u: xi_psi eta_u: eta_psi (padding: both)");
226 ncf.var("grid").put_attr("edge2_dimensions", "xi_v: xi_psi (padding: both) eta_v: eta_psi");
227 ncf.var("grid").put_attr("node_coordinates", "x_psi y_psi");
228 ncf.var("grid").put_attr("face_coordinates", "x_rho y_rho");
229 ncf.var("grid").put_attr("edge1_coordinates", "x_u y_u");
230 ncf.var("grid").put_attr("edge2_coordinates", "x_v y_v");
231 ncf.var("grid").put_attr("vertical_dimensions", "s_rho: s_w (padding: none)");
232
233 ncf.def_var("s_rho",ncutils::NCDType::Real, {nz_r_name});
234 ncf.var("s_rho").put_attr("long_name","S-coordinate at RHO-points");
235 ncf.var("s_rho").put_attr("field","s_rho, scalar");
236
237 ncf.def_var("s_w",ncutils::NCDType::Real, {nz_w_name});
238 ncf.var("s_w").put_attr("long_name","S-coordinate at W-points");
239 ncf.var("s_w").put_attr("field","s_w, scalar");
240
241 ncf.def_var("pm",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
242 ncf.var("pm").put_attr("long_name","curvilinear coordinate metric in XI");
243 ncf.var("pm").put_attr("units","meter-1");
244 ncf.var("pm").put_attr("grid","grid");
245 ncf.var("pm").put_attr("location","face");
246 ncf.var("pm").put_attr("coordinates","x_rho y_rho");
247 ncf.var("pm").put_attr("field","pm, scalar");
248
249 ncf.def_var("pn",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
250 ncf.var("pn").put_attr("long_name","curvilinear coordinate metric in ETA");
251 ncf.var("pn").put_attr("units","meter-1");
252 ncf.var("pn").put_attr("grid","grid");
253 ncf.var("pn").put_attr("location","face");
254 ncf.var("pn").put_attr("coordinates","x_rho y_rho");
255 ncf.var("pn").put_attr("field","pn, scalar");
256
257 ncf.def_var("f",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
258 ncf.var("f").put_attr("long_name","Coriolis parameter at RHO-points");
259 ncf.var("f").put_attr("units","second-1");
260 ncf.var("f").put_attr("grid","grid");
261 ncf.var("f").put_attr("location","face");
262 ncf.var("f").put_attr("coordinates","x_rho y_rho");
263 ncf.var("f").put_attr("field","coriolis, scalar");
264
265 ncf.def_var("x_rho",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
266 ncf.var("x_rho").put_attr("long_name","x-locations of RHO-points");
267 ncf.var("x_rho").put_attr("units","meter");
268 ncf.var("x_rho").put_attr("field","x_rho, scalar");
269
270 ncf.def_var("y_rho",ncutils::NCDType::Real, {ny_r_name, nx_r_name});
271 ncf.var("y_rho").put_attr("long_name","y-locations of RHO-points");
272 ncf.var("y_rho").put_attr("units","meter");
273 ncf.var("y_rho").put_attr("field","y_rho, scalar");
274
275 ncf.def_var("x_u",ncutils::NCDType::Real, {ny_u_name, nx_u_name});
276 ncf.var("x_u").put_attr("long_name","x-locations of U-points");
277 ncf.var("x_u").put_attr("units","meter");
278 ncf.var("x_u").put_attr("field","x_u, scalar");
279
280 ncf.def_var("y_u",ncutils::NCDType::Real, {ny_u_name, nx_u_name});
281 ncf.var("y_u").put_attr("long_name","y-locations of U-points");
282 ncf.var("y_u").put_attr("units","meter");
283 ncf.var("y_u").put_attr("field","y_u, scalar");
284
285 ncf.def_var("x_v",ncutils::NCDType::Real, {ny_v_name, nx_v_name});
286 ncf.var("x_v").put_attr("long_name","x-locations of V-points");
287 ncf.var("x_v").put_attr("units","meter");
288 ncf.var("x_v").put_attr("field","x_v, scalar");
289
290 ncf.def_var("y_v",ncutils::NCDType::Real, {ny_v_name, nx_v_name});
291 ncf.var("y_v").put_attr("long_name","y-locations of V-points");
292 ncf.var("y_v").put_attr("units","meter");
293 ncf.var("y_v").put_attr("field","y_v, scalar");
294
295 ncf.def_var("x_psi",ncutils::NCDType::Real, {ny_p_name, nx_p_name});
296 ncf.var("x_psi").put_attr("long_name","x-locations of PSI-points");
297 ncf.var("x_psi").put_attr("units","meter");
298 ncf.var("x_psi").put_attr("field","x_psi, scalar");
299
300 ncf.def_var("y_psi",ncutils::NCDType::Real, {ny_p_name, nx_p_name});
301 ncf.var("y_psi").put_attr("long_name","y-locations of PSI-points");
302 ncf.var("y_psi").put_attr("units","meter");
303 ncf.var("y_psi").put_attr("field","y_psi, scalar");
304
305 ncf.def_var("ocean_time", ncutils::NCDType::Real, { nt_name });
306 ncf.var("ocean_time").put_attr("long_name","time since initialization");
307 ncf.var("ocean_time").put_attr("units","seconds since 0001-01-01 00:00:00");
308 ncf.var("ocean_time").put_attr("field","time, scalar, series");
309
310 ncf.def_var("Cs_r", ncutils::NCDType::Real, {nz_r_name});
311 ncf.var("Cs_r").put_attr("long_name", "S-coordinate stretching curves at RHO points");
312 ncf.var("Cs_r").put_attr("valid_min",std::vector({-1.}));
313 ncf.var("Cs_r").put_attr("valid_max",std::vector({0.}));
314 ncf.var("Cs_r").put_attr("field","Cs_r, scalar");
315
316 ncf.def_var("Cs_w", ncutils::NCDType::Real, {nz_w_name});
317 ncf.var("Cs_w").put_attr("long_name", "S-coordinate stretching curves at W points");
318 ncf.var("Cs_w").put_attr("valid_min",std::vector({-1.}));
319 ncf.var("Cs_w").put_attr("valid_max",std::vector({0.}));
320 ncf.var("Cs_w").put_attr("field","Cs_w, scalar");
321
322 ncf.def_var("h", ncutils::NCDType::Real, { ny_r_name, nx_r_name });
323 ncf.var("h").put_attr("long_name","bathymetry at RHO-points");
324 ncf.var("h").put_attr("units","meter");
325 ncf.var("h").put_attr("grid","grid");
326 ncf.var("h").put_attr("location","face");
327 ncf.var("h").put_attr("coordinates","x_rho y_rho");
328 ncf.var("h").put_attr("field","bath, scalar");
329
330 ncf.def_var_fill("zeta", ncutils::NCDType::Real, { nt_name, ny_r_name, nx_r_name }, &fill_value);
331 ncf.var("zeta").put_attr("long_name","free-surface");
332 ncf.var("zeta").put_attr("units","meter");
333 ncf.var("zeta").put_attr("time","ocean_time");
334 ncf.var("zeta").put_attr("grid","grid");
335 ncf.var("zeta").put_attr("location","face");
336 ncf.var("zeta").put_attr("coordinates","x_rho y_rho ocean_time");
337 ncf.var("zeta").put_attr("field","free-surface, scalar, series");
338
339 ncf.def_var_fill("temp", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_r_name, nx_r_name }, &fill_value);
340 ncf.var("temp").put_attr("long_name","potential temperature");
341 ncf.var("temp").put_attr("units","Celsius");
342 ncf.var("temp").put_attr("time","ocean_time");
343 ncf.var("temp").put_attr("grid","grid");
344 ncf.var("temp").put_attr("location","face");
345 ncf.var("temp").put_attr("coordinates","x_rho y_rho s_rho ocean_time");
346 ncf.var("temp").put_attr("field","temperature, scalar, series");
347
348 ncf.def_var_fill("salt", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_r_name, nx_r_name }, &fill_value);
349 ncf.var("salt").put_attr("long_name","salinity");
350 ncf.var("salt").put_attr("time","ocean_time");
351 ncf.var("salt").put_attr("grid","grid");
352 ncf.var("salt").put_attr("location","face");
353 ncf.var("salt").put_attr("coordinates","x_rho y_rho s_rho ocean_time");
354 ncf.var("salt").put_attr("field","salinity, scalar, series");
355
356 ncf.def_var_fill("tracer", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_r_name, nx_r_name }, &fill_value);
357 ncf.var("tracer").put_attr("long_name","passive tracer");
358 ncf.var("tracer").put_attr("time","ocean_time");
359 ncf.var("tracer").put_attr("grid","grid");
360 ncf.var("tracer").put_attr("location","face");
361 ncf.var("tracer").put_attr("coordinates","x_rho y_rho s_rho ocean_time");
362 ncf.var("tracer").put_attr("field","tracer, scalar, series");
363
364 ncf.def_var_fill("u", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_u_name, nx_u_name }, &fill_value);
365 ncf.var("u").put_attr("long_name","u-momentum component");
366 ncf.var("u").put_attr("units","meter second-1");
367 ncf.var("u").put_attr("time","ocean_time");
368 ncf.var("u").put_attr("grid","grid");
369 ncf.var("u").put_attr("location","edge1");
370 ncf.var("u").put_attr("coordinates","x_u y_u s_rho ocean_time");
371 ncf.var("u").put_attr("field","u-velocity, scalar, series");
372
373 ncf.def_var_fill("v", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_v_name, nx_v_name }, &fill_value);
374 ncf.var("v").put_attr("long_name","v-momentum component");
375 ncf.var("v").put_attr("units","meter second-1");
376 ncf.var("v").put_attr("time","ocean_time");
377 ncf.var("v").put_attr("grid","grid");
378 ncf.var("v").put_attr("location","edge2");
379 ncf.var("v").put_attr("coordinates","x_v y_v s_rho ocean_time");
380 ncf.var("v").put_attr("field","v-velocity, scalar, series");
381
382 ncf.def_var_fill("ubar", ncutils::NCDType::Real, { nt_name, ny_u_name, nx_u_name }, &fill_value);
383 ncf.var("ubar").put_attr("long_name","vertically integrated u-momentum component");
384 ncf.var("ubar").put_attr("units","meter second-1");
385 ncf.var("ubar").put_attr("time","ocean_time");
386 ncf.var("ubar").put_attr("grid","grid");
387 ncf.var("ubar").put_attr("location","edge1");
388 ncf.var("ubar").put_attr("coordinates","x_u y_u ocean_time");
389 ncf.var("ubar").put_attr("field","ubar-velocity, scalar, series");
390
391 ncf.def_var_fill("vbar", ncutils::NCDType::Real, { nt_name, ny_v_name, nx_v_name }, &fill_value);
392 ncf.var("vbar").put_attr("long_name","vertically integrated v-momentum component");
393 ncf.var("vbar").put_attr("units","meter second-1");
394 ncf.var("vbar").put_attr("time","ocean_time");
395 ncf.var("vbar").put_attr("grid","grid");
396 ncf.var("vbar").put_attr("location","edge2");
397 ncf.var("vbar").put_attr("coordinates","x_v y_v ocean_time");
398 ncf.var("vbar").put_attr("field","vbar-velocity, scalar, series");
399
400 ncf.def_var("sustr", ncutils::NCDType::Real, { nt_name, ny_u_name, nx_u_name });
401 ncf.var("sustr").put_attr("long_name","surface u-momentum stress");
402 ncf.var("sustr").put_attr("units","newton meter-2");
403 ncf.var("sustr").put_attr("time","ocean_time");
404 ncf.var("sustr").put_attr("grid","grid");
405 ncf.var("sustr").put_attr("location","edge1");
406 ncf.var("sustr").put_attr("coordinates","x_u y_u ocean_time");
407 ncf.var("sustr").put_attr("field","surface u-momentum stress, scalar, series");
408
409 ncf.def_var("svstr", ncutils::NCDType::Real, { nt_name, ny_v_name, nx_v_name });
410 ncf.var("svstr").put_attr("long_name","surface v-momentum stress");
411 ncf.var("svstr").put_attr("units","newton meter-2");
412 ncf.var("svstr").put_attr("time","ocean_time");
413 ncf.var("svstr").put_attr("grid","grid");
414 ncf.var("svstr").put_attr("location","edge2");
415 ncf.var("svstr").put_attr("coordinates","x_v y_v ocean_time");
416 ncf.var("svstr").put_attr("field","surface v-momentum stress, scalar, series");
417
418 // Right now this is hard-wired to {temp, salt, tracer, u, v}
419 ncf.put_attr("space_dimension", std::vector<int> { AMREX_SPACEDIM });
420// ncf.put_attr("current_time", std::vector<double> { time });
421 ncf.put_attr("start_time", std::vector<double> { start_bdy_time });
422 ncf.put_attr("CurrentLevel", std::vector<int> { flev });
423 ncf.put_attr("DefaultGeometry", std::vector<int> { amrex::DefaultGeometry().Coord() });
424
425 ncf.exit_def_mode();
426
427 // We are doing single-level writes but it doesn't have to be level 0
428 //
429 // Write out the header information.
430 //
431
432 Real dx[AMREX_SPACEDIM];
433 for (int i = 0; i < AMREX_SPACEDIM; i++) {
434 dx[i] = geom[lev].CellSize()[i];
435 }
436 const auto *base = geom[lev].ProbLo();
437 RealBox rb(subdomain, dx, base);
438
439 amrex::Vector<Real> probLo;
440 amrex::Vector<Real> probHi;
441 for (int i = 0; i < AMREX_SPACEDIM; i++) {
442 probLo.push_back(rb.lo(i));
443 probHi.push_back(rb.hi(i));
444 }
445
446 //nc_probLo.par_access(NC_COLLECTIVE);
447 // small variable data written by just the master proc
448 ncmpi_begin_indep_data(ncf.ncid);
449 if (amrex::ParallelDescriptor::IOProcessor()) // only master proc
450 {
451 auto nc_probLo = ncf.var("probLo");
452
453 nc_probLo.put(probLo.data(), { 0 }, { AMREX_SPACEDIM });
454
455 auto nc_probHi = ncf.var("probHi");
456 //nc_probHi.par_access(NC_COLLECTIVE);
457 nc_probHi.put(probHi.data(), { 0 }, { AMREX_SPACEDIM });
458
459 amrex::Vector<int> smallend;
460 amrex::Vector<int> bigend;
461 for (int i = lev; i < flev; i++) {
462 smallend.clear();
463 bigend.clear();
464 for (int j = 0; j < AMREX_SPACEDIM; j++) {
465 smallend.push_back(subdomain.smallEnd(j));
466 bigend.push_back(subdomain.bigEnd(j));
467 }
468 auto nc_Geom_smallend = ncf.var("Geom.smallend");
469 //nc_Geom_smallend.par_access(NC_COLLECTIVE);
470 nc_Geom_smallend.put(smallend.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
471 AMREX_SPACEDIM });
472
473 auto nc_Geom_bigend = ncf.var("Geom.bigend");
474 //nc_Geom_bigend.par_access(NC_COLLECTIVE);
475 nc_Geom_bigend.put(bigend.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
476 AMREX_SPACEDIM });
477 }
478
479 amrex::Vector<Real> CellSize;
480 for (int i = lev; i < flev; i++) {
481 CellSize.clear();
482 for (Real &j : dx) {
483 CellSize.push_back(amrex::Real(j));
484 }
485 auto nc_CellSize = ncf.var("CellSize");
486 //nc_CellSize.par_access(NC_COLLECTIVE);
487 nc_CellSize.put(CellSize.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
488 AMREX_SPACEDIM });
489 }
490 Real hc = solverChoice.tcline;
491 Real theta_s = solverChoice.theta_s;
492 Real theta_b = solverChoice.theta_b;
493 ncf.var("hc").put(&hc);
494 ncf.var("theta_s").put(&theta_s);
495 ncf.var("theta_b").put(&theta_b);
496
497 }
498 ncmpi_end_indep_data(ncf.ncid);
499
500 } // end if write_header
501
502 ncmpi_begin_indep_data(ncf.ncid);
503 //
504 // We compute the offsets based on location of the box within the domain
505 //
506 long long adjusted_history_count = chunk_history_file ? history_count % steps_per_history_file : history_count;
507 long long local_start_nt = (is_history ? static_cast<long long>(adjusted_history_count) : static_cast<long long>(0));
508 long long local_nt = 1; // We write data for only one time
509
510 {
511 auto nc_plot_var = ncf.var("ocean_time");
512 //nc_plot_var.par_access(NC_COLLECTIVE);
513 nc_plot_var.put(&t_new[lev], { local_start_nt }, { local_nt });
514 }
515 // do all independent writes
516 //ncmpi_end_indep_data(ncf.ncid);
517
518 cons_new[lev]->FillBoundary(geom[lev].periodicity());
519
520 mask_arrays_for_write(lev, (Real) fill_value, 0.0_rt);
521
522 // Check whether there are any nans or infs in variables that we will write out
523 if (vec_Zt_avg1[lev]->contains_nan() || vec_Zt_avg1[lev]->contains_inf()) {
524 amrex::Abort("Found while writing output: zeta contains nan or inf");
525 }
526 if (cons_new[lev]->contains_nan(Temp_comp,1) || cons_new[lev]->contains_inf(Temp_comp,1)) {
527 amrex::Abort("Found while writing output: Temperature contains nan or inf");
528 }
529 if (cons_new[lev]->contains_nan(Salt_comp,1) || cons_new[lev]->contains_inf(Salt_comp,1)) {
530 amrex::Abort("Found while writing output: Salinity contains nan or inf");
531 }
532 if (cons_new[lev]->contains_nan(Scalar_comp,1) || cons_new[lev]->contains_inf(Scalar_comp,1)) {
533 amrex::Abort("Found while writing output: Passive tracer contains nan or inf");
534 }
535 if (xvel_new[lev]->contains_nan() || xvel_new[lev]->contains_inf()) {
536 amrex::Abort("Found while writing output: velocity u contains nan or inf");
537 }
538 if (vec_ubar[lev]->contains_nan(0,1) || vec_ubar[lev]->contains_inf(0,1)) {
539 amrex::Abort("Found while writing output: velocity ubar contains nan or inf");
540 }
541 if (yvel_new[lev]->contains_nan() || yvel_new[lev]->contains_inf()) {
542 amrex::Abort("Found while writing output: velocity v contains nan or inf");
543 }
544 if (vec_vbar[lev]->contains_nan(0,1) || vec_vbar[lev]->contains_inf(0,1)) {
545 amrex::Abort("Found while writing output: velocity vbar contains nan or inf");
546 }
547
548 for (MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi) {
549 auto bx = mfi.validbox();
550 if (subdomain.contains(bx)) {
551 //
552 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
553 //
554 Box tmp_bx(bx);
555 if (tmp_bx.smallEnd()[0] == subdomain.smallEnd()[0])
556 tmp_bx.growLo(0, 1);
557 if (tmp_bx.smallEnd()[1] == subdomain.smallEnd()[1])
558 tmp_bx.growLo(1, 1);
559 if (tmp_bx.bigEnd()[0] == subdomain.bigEnd()[0])
560 tmp_bx.growHi(0, 1);
561 if (tmp_bx.bigEnd()[1] == subdomain.bigEnd()[1])
562 tmp_bx.growHi(1, 1);
563 // amrex::Print() << " BX " << bx << std::endl;
564 // amrex::Print() << "TMP_BX " << tmp_bx << std::endl;
565
566 Box tmp_bx_2d(tmp_bx);
567 tmp_bx_2d.makeSlab(2, 0);
568
569 Box tmp_bx_1d(tmp_bx);
570 tmp_bx_1d.makeSlab(0, 0);
571 tmp_bx_1d.makeSlab(1, 0);
572
573 //
574 // These are the dimensions of the data we write for only this box
575 //
576 long long local_nx = tmp_bx.length()[0];
577 long long local_ny = tmp_bx.length()[1];
578 long long local_nz = tmp_bx.length()[2];
579
580 // We do the "+1" because the offset needs to start at 0
581 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0] + 1);
582 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1] + 1);
583 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
584
585 if (write_header) {
586 // Only write out s_rho and s_w at x=0,y=0 to avoid NaNs
587 if (bx.contains(IntVect(0,0,0)))
588 {
589 {
590 amrex::Gpu::DeviceVector<amrex::Real> tmp_srho(local_nz);
591
592#ifdef AMREX_USE_GPU
593 Gpu::htod_memcpy(tmp_srho.data(), s_r.data(), sizeof(amrex::Real)*local_nz);
594#else
595 std::memcpy(tmp_srho.data(), s_r.data(), sizeof(amrex::Real)*local_nz);
596#endif
597 Gpu::streamSynchronize();
598
599 auto nc_plot_var = ncf.var("s_rho");
600 //nc_plot_var.par_access(NC_INDEPENDENT);
601 nc_plot_var.put(tmp_srho.data(), { local_start_z }, { local_nz });
602 }
603 {
604 amrex::Gpu::DeviceVector<amrex::Real> tmp_sw(local_nz+1);
605
606#ifdef AMREX_USE_GPU
607 Gpu::htod_memcpy(tmp_sw.data(), s_w.data(), sizeof(amrex::Real)*(local_nz+1));
608#else
609 std::memcpy(tmp_sw.data(), s_w.data(), sizeof(amrex::Real)*(local_nz+1));
610#endif
611 Gpu::streamSynchronize();
612
613 auto nc_plot_var = ncf.var("s_w");
614 //nc_plot_var.par_access(NC_INDEPENDENT);
615 nc_plot_var.put(tmp_sw.data(), { local_start_z }, { local_nz + 1});
616 }
617 {
618 amrex::Gpu::DeviceVector<amrex::Real> tmp_Csrho(local_nz);
619
620#ifdef AMREX_USE_GPU
621 Gpu::htod_memcpy(tmp_Csrho.data(), Cs_r.data(), sizeof(amrex::Real)*(local_nz));
622#else
623 std::memcpy(tmp_Csrho.data(), Cs_r.data(), sizeof(amrex::Real)*(local_nz));
624#endif
625 Gpu::streamSynchronize();
626
627 auto nc_plot_var = ncf.var("Cs_r");
628 //nc_plot_var.par_access(NC_INDEPENDENT);
629 nc_plot_var.put(tmp_Csrho.data(), { local_start_z }, { local_nz });
630 }
631 {
632 amrex::Gpu::DeviceVector<amrex::Real> tmp_Csw(local_nz+1);
633
634#ifdef AMREX_USE_GPU
635 Gpu::htod_memcpy(tmp_Csw.data(), Cs_w.data(), sizeof(amrex::Real)*(local_nz+1));
636#else
637 std::memcpy(tmp_Csw.data(), Cs_w.data(), sizeof(amrex::Real)*(local_nz+1));
638#endif
639
640 Gpu::streamSynchronize();
641
642 auto nc_plot_var = ncf.var("Cs_w");
643 //nc_plot_var.par_access(NC_INDEPENDENT);
644 nc_plot_var.put(tmp_Csw.data(), { local_start_z }, { local_nz + 1});
645 }
646 }
647
648 {
649 FArrayBox tmp_bathy;
650 tmp_bathy.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
651
652 tmp_bathy.template copy<RunOn::Device>((*vec_hOfTheConfusingName[lev])[mfi.index()], 0, 0, 1);
653 Gpu::streamSynchronize();
654
655 auto nc_plot_var = ncf.var("h");
656 //nc_plot_var.par_access(NC_INDEPENDENT);
657 nc_plot_var.put(tmp_bathy.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
658 }
659
660 {
661 FArrayBox tmp_pm;
662 tmp_pm.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
663
664 tmp_pm.template copy<RunOn::Device>((*vec_pm[lev])[mfi.index()], 0, 0, 1);
665 Gpu::streamSynchronize();
666
667 auto nc_plot_var = ncf.var("pm");
668 //nc_plot_var.par_access(NC_INDEPENDENT);
669
670 nc_plot_var.put(tmp_pm.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
671 }
672
673 {
674 FArrayBox tmp_pn;
675 tmp_pn.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
676
677 tmp_pn.template copy<RunOn::Device>((*vec_pn[lev])[mfi.index()], 0, 0, 1);
678 Gpu::streamSynchronize();
679
680 auto nc_plot_var = ncf.var("pn");
681 //nc_plot_var.par_access(NC_INDEPENDENT);
682 nc_plot_var.put(tmp_pn.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
683 }
684
685 {
686 FArrayBox tmp_f;
687 tmp_f.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
688
689 tmp_f.template copy<RunOn::Device>((*vec_fcor[lev])[mfi.index()], 0, 0, 1);
690 Gpu::streamSynchronize();
691
692 auto nc_plot_var = ncf.var("f");
693 //nc_plot_var.par_access(NC_INDEPENDENT);
694 nc_plot_var.put(tmp_f.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
695 }
696
697 {
698 FArrayBox tmp_xr;
699 tmp_xr.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
700
701 tmp_xr.template copy<RunOn::Device>((*vec_xr[lev])[mfi.index()], 0, 0, 1);
702 Gpu::streamSynchronize();
703
704 auto nc_plot_var = ncf.var("x_rho");
705 //nc_plot_var.par_access(NC_INDEPENDENT);
706
707 nc_plot_var.put(tmp_xr.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
708 }
709
710 {
711 FArrayBox tmp_yr;
712 tmp_yr.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
713
714 tmp_yr.template copy<RunOn::Device>((*vec_yr[lev])[mfi.index()], 0, 0, 1);
715 Gpu::streamSynchronize();
716
717 auto nc_plot_var = ncf.var("y_rho");
718 //nc_plot_var.par_access(NC_INDEPENDENT);
719 nc_plot_var.put(tmp_yr.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
720 }
721 }
722
723 {
724 FArrayBox tmp_zeta;
725 tmp_zeta.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
726 tmp_zeta.template copy<RunOn::Device>((*vec_Zt_avg1[lev])[mfi.index()], 0, 0, 1);
727 Gpu::streamSynchronize();
728
729 auto nc_plot_var = ncf.var("zeta");
730 nc_plot_var.put(tmp_zeta.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny,
731 local_nx });
732 }
733
734 {
735 FArrayBox tmp_temp;
736 tmp_temp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
737 tmp_temp.template copy<RunOn::Device>((*cons_new[lev])[mfi.index()], Temp_comp, 0, 1);
738 Gpu::streamSynchronize();
739
740 auto nc_plot_var = ncf.var("temp");
741 nc_plot_var.put(tmp_temp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
742 local_nz, local_ny, local_nx });
743 }
744
745 {
746 FArrayBox tmp_salt;
747 tmp_salt.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
748 tmp_salt.template copy<RunOn::Device>((*cons_new[lev])[mfi.index()], Salt_comp, 0, 1);
749 Gpu::streamSynchronize();
750
751 auto nc_plot_var = ncf.var("salt");
752 nc_plot_var.put(tmp_salt.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
753 local_nz, local_ny, local_nx });
754 }
755
756 {
757 FArrayBox tmp_tracer;
758 tmp_tracer.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
759 tmp_tracer.template copy<RunOn::Device>((*cons_new[lev])[mfi.index()], Scalar_comp, 0, 1);
760 Gpu::streamSynchronize();
761
762 auto nc_plot_var = ncf.var("tracer");
763 nc_plot_var.put(tmp_tracer.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
764 local_nz, local_ny, local_nx });
765 }
766 } // subdomain
767 } // mfi
768
769 //ncf.wait_all(irq, &requests[0]);
770 //requests.resize(0);
771 //irq = 0;
772 // Writing u (we loop over cons to get cell-centered box)
773 for (MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi) {
774 Box bx = mfi.validbox();
775
776 if (subdomain.contains(bx)) {
777 //
778 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
779 //
780 Box tmp_bx(bx);
781 tmp_bx.surroundingNodes(0);
782 if (tmp_bx.smallEnd()[1] == subdomain.smallEnd()[1])
783 tmp_bx.growLo(1, 1);
784 if (tmp_bx.bigEnd()[1] == subdomain.bigEnd()[1])
785 tmp_bx.growHi(1, 1);
786 Box tmp_bx_2d(tmp_bx);
787 tmp_bx_2d.makeSlab(2, 0);
788
789 //
790 // These are the dimensions of the data we write for only this box
791 //
792 long long local_nx = tmp_bx.length()[0];
793 long long local_ny = tmp_bx.length()[1];
794 long long local_nz = tmp_bx.length()[2];
795
796 // We do the "+1" because the offset needs to start at 0
797 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0]);
798 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1] + 1);
799 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
800
801 if (write_header) {
802 {
803 FArrayBox tmp;
804 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
805 tmp.template copy<RunOn::Device>((*vec_xu[lev])[mfi.index()], 0, 0, 1);
806 Gpu::streamSynchronize();
807
808 auto nc_plot_var = ncf.var("x_u");
809 //nc_plot_var.par_access(NC_INDEPENDENT);
810 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
811 }
812 {
813 FArrayBox tmp;
814 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
815 tmp.template copy<RunOn::Device>((*vec_yu[lev])[mfi.index()], 0, 0, 1);
816 Gpu::streamSynchronize();
817
818 auto nc_plot_var = ncf.var("y_u");
819 //nc_plot_var.par_access(NC_INDEPENDENT);
820 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
821 }
822 }
823
824 {
825 FArrayBox tmp;
826 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
827 tmp.template copy<RunOn::Device>((*xvel_new[lev])[mfi.index()], 0, 0, 1);
828 Gpu::streamSynchronize();
829
830 auto nc_plot_var = ncf.var("u");
831 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
832 local_nz, local_ny, local_nx });
833 }
834
835 {
836 FArrayBox tmp;
837 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
838 tmp.template copy<RunOn::Device>((*vec_ubar[lev])[mfi.index()], 0, 0, 1);
839 Gpu::streamSynchronize();
840
841 auto nc_plot_var = ncf.var("ubar");
842 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
843 }
844 {
845 FArrayBox tmp;
846 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
847 tmp.template copy<RunOn::Device>((*vec_sustr[lev])[mfi.index()], 0, 0, 1);
848 Gpu::streamSynchronize();
849
850 auto nc_plot_var = ncf.var("sustr");
851 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
852 }
853 } // in subdomain
854 } // mfi
855
856 // Writing v (we loop over cons to get cell-centered box)
857 for (MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi) {
858 Box bx = mfi.validbox();
859
860 if (subdomain.contains(bx)) {
861 //
862 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
863 //
864 Box tmp_bx(bx);
865 tmp_bx.surroundingNodes(1);
866 if (tmp_bx.smallEnd()[0] == subdomain.smallEnd()[0])
867 tmp_bx.growLo(0, 1);
868 if (tmp_bx.bigEnd()[0] == subdomain.bigEnd()[0])
869 tmp_bx.growHi(0, 1);
870 // amrex::Print() << " BX " << bx << std::endl;
871 // amrex::Print() << "TMP_BX " << tmp_bx << std::endl;
872
873 Box tmp_bx_2d(tmp_bx);
874 tmp_bx_2d.makeSlab(2, 0);
875
876 //
877 // These are the dimensions of the data we write for only this box
878 //
879 long long local_nx = tmp_bx.length()[0];
880 long long local_ny = tmp_bx.length()[1];
881 long long local_nz = tmp_bx.length()[2];
882
883 // We do the "+1" because the offset needs to start at 0
884 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0] + 1);
885 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1]);
886 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
887
888 if (write_header) {
889 {
890 FArrayBox tmp;
891 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
892 tmp.template copy<RunOn::Device>((*vec_xv[lev])[mfi.index()], 0, 0, 1);
893 Gpu::streamSynchronize();
894
895 auto nc_plot_var = ncf.var("x_v");
896 //nc_plot_var.par_access(NC_INDEPENDENT);
897 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
898 }
899 {
900 FArrayBox tmp;
901 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
902 tmp.template copy<RunOn::Device>((*vec_yv[lev])[mfi.index()], 0, 0, 1);
903 Gpu::streamSynchronize();
904
905 auto nc_plot_var = ncf.var("y_v");
906 //nc_plot_var.par_access(NC_INDEPENDENT);
907 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
908 }
909
910 }
911 {
912 FArrayBox tmp;
913 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
914 tmp.template copy<RunOn::Device>((*yvel_new[lev])[mfi.index()], 0, 0, 1);
915 Gpu::streamSynchronize();
916
917 auto nc_plot_var = ncf.var("v");
918 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
919 local_nz, local_ny, local_nx });
920 }
921
922 {
923 FArrayBox tmp;
924 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
925 tmp.template copy<RunOn::Device>((*vec_vbar[lev])[mfi.index()], 0, 0, 1);
926 Gpu::streamSynchronize();
927
928 auto nc_plot_var = ncf.var("vbar");
929 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
930 }
931 {
932 FArrayBox tmp;
933 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
934 tmp.template copy<RunOn::Device>((*vec_svstr[lev])[mfi.index()], 0, 0, 1);
935 Gpu::streamSynchronize();
936
937 auto nc_plot_var = ncf.var("svstr");
938 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
939 }
940
941 } // in subdomain
942 } // mfi
943
944 for (MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi) {
945 Box bx = mfi.validbox();
946
947 if (subdomain.contains(bx)) {
948 //
949 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
950 //
951 Box tmp_bx(bx);
952 tmp_bx.surroundingNodes(0);
953 tmp_bx.surroundingNodes(1);
954
955 Box tmp_bx_2d(tmp_bx);
956 tmp_bx_2d.makeSlab(2, 0);
957
958 //
959 // These are the dimensions of the data we write for only this box
960 //
961 long long local_nx = tmp_bx.length()[0];
962 long long local_ny = tmp_bx.length()[1];
963
964 // We do the "+1" because the offset needs to start at 0
965 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0]);
966 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1]);
967
968 if (write_header) {
969 {
970 FArrayBox tmp;
971 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
972 tmp.template copy<RunOn::Device>((*vec_xp[lev])[mfi.index()], 0, 0, 1);
973 Gpu::streamSynchronize();
974
975 auto nc_plot_var = ncf.var("x_psi");
976 //nc_plot_var.par_access(NC_INDEPENDENT);
977 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
978 }
979 {
980 FArrayBox tmp;
981 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
982 tmp.template copy<RunOn::Device>((*vec_yp[lev])[mfi.index()], 0, 0, 1);
983 Gpu::streamSynchronize();
984
985 auto nc_plot_var = ncf.var("y_psi");
986 //nc_plot_var.par_access(NC_INDEPENDENT);
987 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
988 }
989
990 } // header
991 } // in subdomain
992 } // mfi
993
994 mask_arrays_for_write(lev, 0.0_rt, (Real) fill_value);
995
996 ncf.close();
997
999}
1000
#define Temp_comp
#define Scalar_comp
#define Salt_comp
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:371
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_hOfTheConfusingName
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:235
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:366
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yv
y_grid on v-points (2D)
Definition REMORA.H:386
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:223
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
Definition REMORA.H:1062
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:274
int history_count
Counter for which time index we are writing to in the netcdf history file.
Definition REMORA.H:1268
bool chunk_history_file
Whether to chunk netcdf history file.
Definition REMORA.H:1261
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:293
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yp
y_grid on psi-points (2D)
Definition REMORA.H:391
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xr
x_grid on rho points (2D)
Definition REMORA.H:374
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xv
x_grid on v-points (2D)
Definition REMORA.H:384
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1149
void WriteNCPlotFile(int istep)
Write plotfile using NetCDF (wrapper)
void WriteNCPlotFile_which(int lev, int which_subdomain, bool write_header, ncutils::NCFile &ncf, bool is_history)
Write a particular NetCDF plotfile.
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:227
amrex::Real start_bdy_time
Start time in the time series of boundary data.
Definition REMORA.H:1057
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:272
int steps_per_history_file
Number of time steps per netcdf history file.
Definition REMORA.H:1266
int max_step
maximum number of steps
Definition REMORA.H:1216
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:225
int plot_int
Plotfile output interval in iterations.
Definition REMORA.H:1250
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:1314
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:295
amrex::Gpu::DeviceVector< amrex::Real > Cs_r
Stretching coefficients at rho points.
Definition REMORA.H:282
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1156
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1279
static int total_nc_plot_file_step
Definition REMORA.H:968
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xp
x_grid on psi-points (2D)
Definition REMORA.H:389
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:352
amrex::Gpu::DeviceVector< amrex::Real > Cs_w
Stretching coefficients at w points.
Definition REMORA.H:284
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:350
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yu
y_grid on u-points (2D)
Definition REMORA.H:381
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xu
x_grid on u-points (2D)
Definition REMORA.H:379
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:368
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1248
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:290
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yr
y_grid on rho points (2D)
Definition REMORA.H:376
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".