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 // Output 2D horizontal mixing coefficients if using scaled_to_grid option
417 ncf.def_var_fill("visc2", ncutils::NCDType::Real, { ny_r_name, nx_r_name }, &netcdf_fill_value);
418 ncf.var("visc2").put_attr("long_name","horizontal harmonic viscosity coefficient at RHO-points");
419 ncf.var("visc2").put_attr("units","meter2 second-1");
420 ncf.var("visc2").put_attr("grid","grid");
421 ncf.var("visc2").put_attr("location","face");
422 ncf.var("visc2").put_attr("coordinates","x_rho y_rho");
423 ncf.var("visc2").put_attr("field","visc2, scalar");
424
425 for (int n = 0; n < Tracer_comp + 1; ++n) {
426 const std::string nm = std::string("diff2_") + cons_names[n];
427 ncf.def_var_fill(nm, ncutils::NCDType::Real, { ny_r_name, nx_r_name }, &netcdf_fill_value);
428 ncf.var(nm).put_attr("long_name", std::string("horizontal harmonic diffusivity coefficient for ") + cons_names[n] + " at RHO-points");
429 ncf.var(nm).put_attr("units","meter2 second-1");
430 ncf.var(nm).put_attr("grid","grid");
431 ncf.var(nm).put_attr("location","face");
432 ncf.var(nm).put_attr("coordinates","x_rho y_rho");
433 ncf.var(nm).put_attr("field", nm + ", scalar");
434 }
435 }
436
437 ncf.def_var_fill("u", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_u_name, nx_u_name }, &netcdf_fill_value);
438 ncf.var("u").put_attr("long_name","u-momentum component");
439 ncf.var("u").put_attr("units","meter second-1");
440 ncf.var("u").put_attr("time","ocean_time");
441 ncf.var("u").put_attr("grid","grid");
442 ncf.var("u").put_attr("location","edge1");
443 ncf.var("u").put_attr("coordinates","x_u y_u s_rho ocean_time");
444 ncf.var("u").put_attr("field","u-velocity, scalar, series");
445
446 ncf.def_var_fill("v", ncutils::NCDType::Real, { nt_name, nz_r_name, ny_v_name, nx_v_name }, &netcdf_fill_value);
447 ncf.var("v").put_attr("long_name","v-momentum component");
448 ncf.var("v").put_attr("units","meter second-1");
449 ncf.var("v").put_attr("time","ocean_time");
450 ncf.var("v").put_attr("grid","grid");
451 ncf.var("v").put_attr("location","edge2");
452 ncf.var("v").put_attr("coordinates","x_v y_v s_rho ocean_time");
453 ncf.var("v").put_attr("field","v-velocity, scalar, series");
454
455 ncf.def_var_fill("ubar", ncutils::NCDType::Real, { nt_name, ny_u_name, nx_u_name }, &netcdf_fill_value);
456 ncf.var("ubar").put_attr("long_name","vertically integrated u-momentum component");
457 ncf.var("ubar").put_attr("units","meter second-1");
458 ncf.var("ubar").put_attr("time","ocean_time");
459 ncf.var("ubar").put_attr("grid","grid");
460 ncf.var("ubar").put_attr("location","edge1");
461 ncf.var("ubar").put_attr("coordinates","x_u y_u ocean_time");
462 ncf.var("ubar").put_attr("field","ubar-velocity, scalar, series");
463
464 ncf.def_var_fill("vbar", ncutils::NCDType::Real, { nt_name, ny_v_name, nx_v_name }, &netcdf_fill_value);
465 ncf.var("vbar").put_attr("long_name","vertically integrated v-momentum component");
466 ncf.var("vbar").put_attr("units","meter second-1");
467 ncf.var("vbar").put_attr("time","ocean_time");
468 ncf.var("vbar").put_attr("grid","grid");
469 ncf.var("vbar").put_attr("location","edge2");
470 ncf.var("vbar").put_attr("coordinates","x_v y_v ocean_time");
471 ncf.var("vbar").put_attr("field","vbar-velocity, scalar, series");
472
473 ncf.def_var("sustr", ncutils::NCDType::Real, { nt_name, ny_u_name, nx_u_name });
474 ncf.var("sustr").put_attr("long_name","surface u-momentum stress");
475 ncf.var("sustr").put_attr("units","newton meter-2");
476 ncf.var("sustr").put_attr("time","ocean_time");
477 ncf.var("sustr").put_attr("grid","grid");
478 ncf.var("sustr").put_attr("location","edge1");
479 ncf.var("sustr").put_attr("coordinates","x_u y_u ocean_time");
480 ncf.var("sustr").put_attr("field","surface u-momentum stress, scalar, series");
481
482 ncf.def_var("svstr", ncutils::NCDType::Real, { nt_name, ny_v_name, nx_v_name });
483 ncf.var("svstr").put_attr("long_name","surface v-momentum stress");
484 ncf.var("svstr").put_attr("units","newton meter-2");
485 ncf.var("svstr").put_attr("time","ocean_time");
486 ncf.var("svstr").put_attr("grid","grid");
487 ncf.var("svstr").put_attr("location","edge2");
488 ncf.var("svstr").put_attr("coordinates","x_v y_v ocean_time");
489 ncf.var("svstr").put_attr("field","surface v-momentum stress, scalar, series");
490
492 // Surface air temperature (Celsius)
493 ncf.def_var("Tair", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
494 ncf.var("Tair").put_attr("long_name","surface air temperature");
495 ncf.var("Tair").put_attr("units","Celsius");
496 ncf.var("Tair").put_attr("time","ocean_time");
497 ncf.var("Tair").put_attr("grid","grid");
498 ncf.var("Tair").put_attr("location","face");
499 ncf.var("Tair").put_attr("coordinates","x_rho y_rho ocean_time");
500 ncf.var("Tair").put_attr("field","Tair, scalar, series");
501
502 // Surface air pressure (Pascal)
503 ncf.def_var("Pair", ncutils::NCDType::Real,{ nt_name, ny_r_name, nx_r_name });
504 ncf.var("Pair").put_attr("long_name","surface air pressure");
505 ncf.var("Pair").put_attr("units","Pascal");
506 ncf.var("Pair").put_attr("time","ocean_time");
507 ncf.var("Pair").put_attr("grid","grid");
508 ncf.var("Pair").put_attr("location","face");
509 ncf.var("Pair").put_attr("coordinates","x_rho y_rho ocean_time");
510 ncf.var("Pair").put_attr("field","Pair, scalar, series");
511
512 // Surface net heat flux (W/m2)
513 ncf.def_var("qnet", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
514 ncf.var("qnet").put_attr("long_name","surface net heat flux");
515 ncf.var("qnet").put_attr("units","watt meter-2");
516 ncf.var("qnet").put_attr("time","ocean_time");
517 ncf.var("qnet").put_attr("grid","grid");
518 ncf.var("qnet").put_attr("location","face");
519 ncf.var("qnet").put_attr("coordinates","x_rho y_rho ocean_time");
520 ncf.var("qnet").put_attr("field","surface heat flux, scalar, series");
521
522 // Surface net salt flux (kinematic)
523 ncf.def_var("ssflux", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
524 ncf.var("ssflux").put_attr("long_name","kinematic surface net salt flux, SALT*(E-P)/rhow");
525 ncf.var("ssflux").put_attr("units","meter second-1");
526 ncf.var("ssflux").put_attr("time","ocean_time");
527 ncf.var("ssflux").put_attr("grid","grid");
528 ncf.var("ssflux").put_attr("location","face");
529 ncf.var("ssflux").put_attr("coordinates","x_rho y_rho ocean_time");
530 ncf.var("ssflux").put_attr("field","surface net salt flux, scalar, series");
531
532 // Latent heat flux (W/m2)
533 ncf.def_var("latent", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
534 ncf.var("latent").put_attr("long_name","net latent heat flux");
535 ncf.var("latent").put_attr("units","watt meter-2");
536 ncf.var("latent").put_attr("time","ocean_time");
537 ncf.var("latent").put_attr("grid","grid");
538 ncf.var("latent").put_attr("location","face");
539 ncf.var("latent").put_attr("coordinates","x_rho y_rho ocean_time");
540 ncf.var("latent").put_attr("field","latent heat flux, scalar, series");
541
542 // Sensible heat flux (W/m2)
543 ncf.def_var("sensible", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
544 ncf.var("sensible").put_attr("long_name","net sensible heat flux");
545 ncf.var("sensible").put_attr("units","watt meter-2");
546 ncf.var("sensible").put_attr("time","ocean_time");
547 ncf.var("sensible").put_attr("grid","grid");
548 ncf.var("sensible").put_attr("location","face");
549 ncf.var("sensible").put_attr("coordinates","x_rho y_rho ocean_time");
550 ncf.var("sensible").put_attr("field","sensible heat flux, scalar, series");
551
552 // Longwave radiation (W/m2)
553 ncf.def_var("lwrad", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
554 ncf.var("lwrad").put_attr("long_name","net longwave radiation flux");
555 ncf.var("lwrad").put_attr("units","watt meter-2");
556 ncf.var("lwrad").put_attr("time","ocean_time");
557 ncf.var("lwrad").put_attr("grid","grid");
558 ncf.var("lwrad").put_attr("location","face");
559 ncf.var("lwrad").put_attr("coordinates","x_rho y_rho ocean_time");
560 ncf.var("lwrad").put_attr("field","longwave radiation, scalar, series");
561
562 // Shortwave radiation (W/m2)
563 ncf.def_var("swrad", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
564 ncf.var("swrad").put_attr("long_name","solar shortwave radiation flux");
565 ncf.var("swrad").put_attr("units","watt meter-2");
566 ncf.var("swrad").put_attr("time","ocean_time");
567 ncf.var("swrad").put_attr("grid","grid");
568 ncf.var("swrad").put_attr("location","face");
569 ncf.var("swrad").put_attr("coordinates","x_rho y_rho ocean_time");
570 ncf.var("swrad").put_attr("field","shortwave radiation, scalar, series");
571
572 // Evaporation rate (kg m-2 s-1)
573 ncf.def_var("evaporation", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
574 ncf.var("evaporation").put_attr("long_name","evaporation rate");
575 ncf.var("evaporation").put_attr("units","kilogram meter-2 second-1");
576 ncf.var("evaporation").put_attr("time","ocean_time");
577 ncf.var("evaporation").put_attr("grid","grid");
578 ncf.var("evaporation").put_attr("location","face");
579 ncf.var("evaporation").put_attr("coordinates","x_rho y_rho ocean_time");
580 ncf.var("evaporation").put_attr("field","evaporation, scalar, series");
581
582 // Rain rate (kg m-2 s-1)
583 ncf.def_var("rain", ncutils::NCDType::Real, {nt_name, ny_r_name, nx_r_name });
584 ncf.var("rain").put_attr("long_name","rain fall rate");
585 ncf.var("rain").put_attr("units","kilogram meter-2 second-1");
586 ncf.var("rain").put_attr("time","ocean_time");
587 ncf.var("rain").put_attr("grid","grid");
588 ncf.var("rain").put_attr("location","face");
589 ncf.var("rain").put_attr("coordinates","x_rho y_rho ocean_time");
590 ncf.var("rain").put_attr("field","rain, scalar, series");
591 }
592 // Right now this is hard-wired to {temp, salt, tracer, u, v}
593 ncf.put_attr("space_dimension", std::vector<int> { AMREX_SPACEDIM });
594// ncf.put_attr("current_time", std::vector<double> { time });
595 ncf.put_attr("start_time", std::vector<double> { start_bdy_time });
596 ncf.put_attr("CurrentLevel", std::vector<int> { flev });
597 ncf.put_attr("DefaultGeometry", std::vector<int> { amrex::DefaultGeometry().Coord() });
598
599 ncf.exit_def_mode();
600
601 // We are doing single-level writes but it doesn't have to be level 0
602 //
603 // Write out the header information.
604 //
605
606 Real dx[AMREX_SPACEDIM];
607 for (int i = 0; i < AMREX_SPACEDIM; i++) {
608 dx[i] = geom[lev].CellSize()[i];
609 }
610 const auto *base = geom[lev].ProbLo();
611 RealBox rb(subdomain, dx, base);
612
613 amrex::Vector<Real> probLo;
614 amrex::Vector<Real> probHi;
615 for (int i = 0; i < AMREX_SPACEDIM; i++) {
616 probLo.push_back(rb.lo(i));
617 probHi.push_back(rb.hi(i));
618 }
619
620 //nc_probLo.par_access(NC_COLLECTIVE);
621 // small variable data written by just the master proc
622 ncmpi_begin_indep_data(ncf.ncid);
623 if (amrex::ParallelDescriptor::IOProcessor()) // only master proc
624 {
625 auto nc_probLo = ncf.var("probLo");
626
627 nc_probLo.put(probLo.data(), { 0 }, { AMREX_SPACEDIM });
628
629 auto nc_probHi = ncf.var("probHi");
630 //nc_probHi.par_access(NC_COLLECTIVE);
631 nc_probHi.put(probHi.data(), { 0 }, { AMREX_SPACEDIM });
632
633 amrex::Vector<int> smallend;
634 amrex::Vector<int> bigend;
635 for (int i = lev; i < flev; i++) {
636 smallend.clear();
637 bigend.clear();
638 for (int j = 0; j < AMREX_SPACEDIM; j++) {
639 smallend.push_back(subdomain.smallEnd(j));
640 bigend.push_back(subdomain.bigEnd(j));
641 }
642 auto nc_Geom_smallend = ncf.var("Geom.smallend");
643 //nc_Geom_smallend.par_access(NC_COLLECTIVE);
644 nc_Geom_smallend.put(smallend.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
645 AMREX_SPACEDIM });
646
647 auto nc_Geom_bigend = ncf.var("Geom.bigend");
648 //nc_Geom_bigend.par_access(NC_COLLECTIVE);
649 nc_Geom_bigend.put(bigend.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
650 AMREX_SPACEDIM });
651 }
652
653 amrex::Vector<Real> CellSize;
654 for (int i = lev; i < flev; i++) {
655 CellSize.clear();
656 for (Real &j : dx) {
657 CellSize.push_back(amrex::Real(j));
658 }
659 auto nc_CellSize = ncf.var("CellSize");
660 //nc_CellSize.par_access(NC_COLLECTIVE);
661 nc_CellSize.put(CellSize.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
662 AMREX_SPACEDIM });
663 }
664 Real hc = solverChoice.tcline;
665 Real theta_s = solverChoice.theta_s;
666 Real theta_b = solverChoice.theta_b;
667 ncf.var("hc").put(&hc);
668 ncf.var("theta_s").put(&theta_s);
669 ncf.var("theta_b").put(&theta_b);
670
671 }
672 ncmpi_end_indep_data(ncf.ncid);
673
674 } // end if write_header
675
676 ncmpi_begin_indep_data(ncf.ncid);
677 //
678 // We compute the offsets based on location of the box within the domain
679 //
680 long long adjusted_history_count = chunk_history_file ? history_count % steps_per_history_file : history_count;
681 long long local_start_nt = (is_history ? static_cast<long long>(adjusted_history_count) : static_cast<long long>(0));
682 long long local_nt = 1; // We write data for only one time
683
684 {
685 auto nc_plot_var = ncf.var("ocean_time");
686 //nc_plot_var.par_access(NC_COLLECTIVE);
687 nc_plot_var.put(&t_new[lev], { local_start_nt }, { local_nt });
688 }
689 // do all independent writes
690 //ncmpi_end_indep_data(ncf.ncid);
691
692 mask_arrays_for_write(lev, (Real) netcdf_fill_value, 0.0_rt);
693
694 // Check whether there are any nans or infs in variables that we will write out
695 if (vec_Zt_avg1[lev]->contains_nan() || vec_Zt_avg1[lev]->contains_inf()) {
696 amrex::Abort("Found while writing output: zeta contains nan or inf");
697 }
698 if (plotMF->contains_nan(Temp_comp,1) || plotMF->contains_inf(Temp_comp,1)) {
699 amrex::Abort("Found while writing output: Temperature contains nan or inf");
700 }
701 if (plotMF->contains_nan(Salt_comp,1) || plotMF->contains_inf(Salt_comp,1)) {
702 amrex::Abort("Found while writing output: Salinity contains nan or inf");
703 }
704 if (plotMF->contains_nan(Tracer_comp,1) || plotMF->contains_inf(Tracer_comp,1)) {
705 amrex::Abort("Found while writing output: Passive tracer contains nan or inf");
706 }
707 if (xvel_new[lev]->contains_nan() || xvel_new[lev]->contains_inf()) {
708 amrex::Abort("Found while writing output: velocity u contains nan or inf");
709 }
710 if (vec_ubar[lev]->contains_nan(0,1) || vec_ubar[lev]->contains_inf(0,1)) {
711 amrex::Abort("Found while writing output: velocity ubar contains nan or inf");
712 }
713 if (yvel_new[lev]->contains_nan() || yvel_new[lev]->contains_inf()) {
714 amrex::Abort("Found while writing output: velocity v contains nan or inf");
715 }
716 if (vec_vbar[lev]->contains_nan(0,1) || vec_vbar[lev]->contains_inf(0,1)) {
717 amrex::Abort("Found while writing output: velocity vbar contains nan or inf");
718 }
719
720 for (MFIter mfi(*plotMF, false); mfi.isValid(); ++mfi) {
721 auto bx = mfi.validbox();
722 if (subdomain.contains(bx)) {
723 //
724 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
725 //
726 Box tmp_bx(bx);
727 if (tmp_bx.smallEnd()[0] == subdomain.smallEnd()[0])
728 tmp_bx.growLo(0, 1);
729 if (tmp_bx.smallEnd()[1] == subdomain.smallEnd()[1])
730 tmp_bx.growLo(1, 1);
731 if (tmp_bx.bigEnd()[0] == subdomain.bigEnd()[0])
732 tmp_bx.growHi(0, 1);
733 if (tmp_bx.bigEnd()[1] == subdomain.bigEnd()[1])
734 tmp_bx.growHi(1, 1);
735 // amrex::Print() << " BX " << bx << std::endl;
736 // amrex::Print() << "TMP_BX " << tmp_bx << std::endl;
737
738 Box tmp_bx_2d(tmp_bx);
739 tmp_bx_2d.makeSlab(2, 0);
740
741 Box tmp_bx_1d(tmp_bx);
742 tmp_bx_1d.makeSlab(0, 0);
743 tmp_bx_1d.makeSlab(1, 0);
744
745 //
746 // These are the dimensions of the data we write for only this box
747 //
748 long long local_nx = tmp_bx.length()[0];
749 long long local_ny = tmp_bx.length()[1];
750 long long local_nz = tmp_bx.length()[2];
751
752 // We do the "+1" because the offset needs to start at 0
753 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0] + 1);
754 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1] + 1);
755 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
756
757 if (write_header) {
758 // Only write out s_rho and s_w at x=0,y=0 to avoid NaNs
759 if (bx.contains(IntVect(0,0,0)))
760 {
761 {
762 amrex::Vector<amrex::Real> tmp_srho(local_nz);
763
764#ifdef AMREX_USE_GPU
765 Gpu::dtoh_memcpy(tmp_srho.data(), s_r.data(), sizeof(amrex::Real)*local_nz);
766#else
767 std::memcpy(tmp_srho.data(), s_r.data(), sizeof(amrex::Real)*local_nz);
768#endif
769 Gpu::streamSynchronize();
770
771 auto nc_plot_var = ncf.var("s_rho");
772 //nc_plot_var.par_access(NC_INDEPENDENT);
773 nc_plot_var.put(tmp_srho.data(), { local_start_z }, { local_nz });
774 }
775 {
776 amrex::Vector<amrex::Real> tmp_sw(local_nz+1);
777
778#ifdef AMREX_USE_GPU
779 Gpu::dtoh_memcpy(tmp_sw.data(), s_w.data(), sizeof(amrex::Real)*(local_nz+1));
780#else
781 std::memcpy(tmp_sw.data(), s_w.data(), sizeof(amrex::Real)*(local_nz+1));
782#endif
783 Gpu::streamSynchronize();
784
785 auto nc_plot_var = ncf.var("s_w");
786 //nc_plot_var.par_access(NC_INDEPENDENT);
787 nc_plot_var.put(tmp_sw.data(), { local_start_z }, { local_nz + 1});
788 }
789 {
790 amrex::Vector<amrex::Real> tmp_Csrho(local_nz);
791
792#ifdef AMREX_USE_GPU
793 Gpu::dtoh_memcpy(tmp_Csrho.data(), Cs_r.data(), sizeof(amrex::Real)*(local_nz));
794#else
795 std::memcpy(tmp_Csrho.data(), Cs_r.data(), sizeof(amrex::Real)*(local_nz));
796#endif
797 Gpu::streamSynchronize();
798
799 auto nc_plot_var = ncf.var("Cs_r");
800 //nc_plot_var.par_access(NC_INDEPENDENT);
801 nc_plot_var.put(tmp_Csrho.data(), { local_start_z }, { local_nz });
802 }
803 {
804 amrex::Vector<amrex::Real> tmp_Csw(local_nz+1);
805
806#ifdef AMREX_USE_GPU
807 Gpu::dtoh_memcpy(tmp_Csw.data(), Cs_w.data(), sizeof(amrex::Real)*(local_nz+1));
808#else
809 std::memcpy(tmp_Csw.data(), Cs_w.data(), sizeof(amrex::Real)*(local_nz+1));
810#endif
811
812 Gpu::streamSynchronize();
813
814 auto nc_plot_var = ncf.var("Cs_w");
815 //nc_plot_var.par_access(NC_INDEPENDENT);
816 nc_plot_var.put(tmp_Csw.data(), { local_start_z }, { local_nz + 1});
817 }
818 }
819
820 {
821 FArrayBox tmp_bathy;
822 tmp_bathy.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
823
824 tmp_bathy.template copy<RunOn::Device>((*vec_h[lev])[mfi.index()], 0, 0, 1);
825 Gpu::streamSynchronize();
826
827 auto nc_plot_var = ncf.var("h");
828 //nc_plot_var.par_access(NC_INDEPENDENT);
829 nc_plot_var.put(tmp_bathy.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
830 }
831
832 {
833 FArrayBox tmp_pm;
834 tmp_pm.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
835
836 tmp_pm.template copy<RunOn::Device>((*vec_pm[lev])[mfi.index()], 0, 0, 1);
837 Gpu::streamSynchronize();
838
839 auto nc_plot_var = ncf.var("pm");
840 //nc_plot_var.par_access(NC_INDEPENDENT);
841
842 nc_plot_var.put(tmp_pm.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
843 }
844
845 {
846 FArrayBox tmp_pn;
847 tmp_pn.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
848
849 tmp_pn.template copy<RunOn::Device>((*vec_pn[lev])[mfi.index()], 0, 0, 1);
850 Gpu::streamSynchronize();
851
852 auto nc_plot_var = ncf.var("pn");
853 //nc_plot_var.par_access(NC_INDEPENDENT);
854 nc_plot_var.put(tmp_pn.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
855 }
856
857 {
858 FArrayBox tmp_f;
859 tmp_f.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
860
861 tmp_f.template copy<RunOn::Device>((*vec_fcor[lev])[mfi.index()], 0, 0, 1);
862 Gpu::streamSynchronize();
863
864 auto nc_plot_var = ncf.var("f");
865 //nc_plot_var.par_access(NC_INDEPENDENT);
866 nc_plot_var.put(tmp_f.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
867 }
868
869 {
870 FArrayBox tmp_xr;
871 tmp_xr.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
872
873 tmp_xr.template copy<RunOn::Device>((*vec_xr[lev])[mfi.index()], 0, 0, 1);
874 Gpu::streamSynchronize();
875
876 auto nc_plot_var = ncf.var("x_rho");
877 //nc_plot_var.par_access(NC_INDEPENDENT);
878
879 nc_plot_var.put(tmp_xr.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
880 }
881
882 {
883 FArrayBox tmp_yr;
884 tmp_yr.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
885
886 tmp_yr.template copy<RunOn::Device>((*vec_yr[lev])[mfi.index()], 0, 0, 1);
887 Gpu::streamSynchronize();
888
889 auto nc_plot_var = ncf.var("y_rho");
890 //nc_plot_var.par_access(NC_INDEPENDENT);
891 nc_plot_var.put(tmp_yr.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
892 }
893 }
894
895 {
896 FArrayBox tmp_zeta;
897 tmp_zeta.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
898 tmp_zeta.template copy<RunOn::Device>((*vec_Zt_avg1[lev])[mfi.index()], 0, 0, 1);
899 Gpu::streamSynchronize();
900
901 auto nc_plot_var = ncf.var("zeta");
902 nc_plot_var.put(tmp_zeta.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny,
903 local_nx });
904 }
905
907 {
908 const Real Hscale = solverChoice.rho0 * Cp;
909 // Tair
910 {
911 FArrayBox tmp_Tair;
912 tmp_Tair.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
913 tmp_Tair.template copy<RunOn::Device>((*vec_Tair[lev])[mfi.index()], 0, 0, 1);
914 Gpu::streamSynchronize();
915
916 auto nc_plot_var = ncf.var("Tair");
917 nc_plot_var.put(tmp_Tair.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
918 }
919 // Pair
920 {
921 FArrayBox tmp_Pair;
922 tmp_Pair.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
923 tmp_Pair.template copy<RunOn::Device>((*vec_Pair[lev])[mfi.index()], 0, 0, 1);
924 Gpu::streamSynchronize();
925
926 auto nc_plot_var = ncf.var("Pair");
927 nc_plot_var.put(tmp_Pair.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
928 }
929 // qnet (stored °C m/s → write W/m²)
930 {
931 FArrayBox tmp;
932 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
933
934 // Copy stflux Temp component
935 tmp.template copy<RunOn::Device>(
936 (*vec_stflux[lev])[mfi.index()],
937 Temp_comp, // source component
938 0, // dest component
939 1 // number of comps
940 );
941
942 Gpu::streamSynchronize();
943
944 // Convert °C·m/s → W/m²
945 tmp.mult<RunOn::Device>(Hscale);
946
947 auto nc_var = ncf.var("qnet");
948 nc_var.put(tmp.dataPtr(),
949 { local_start_nt, local_start_y, local_start_x },
950 { local_nt, local_ny, local_nx });
951 }
952 // ssflux = surface net freshwater flux (kg/m²/s converted to m/s)
953 {
954 FArrayBox tmp;
955 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
956
957 // Copy stflux Salt component
958 tmp.template copy<RunOn::Device>(
959 (*vec_stflux[lev])[mfi.index()],
960 Salt_comp, // source component
961 0, // destination component
962 1 // number of components
963 );
964
965 Gpu::streamSynchronize();
966
967 auto nc_var = ncf.var("ssflux");
968 nc_var.put(tmp.dataPtr(),
969 { local_start_nt, local_start_y, local_start_x },
970 { local_nt, local_ny, local_nx });
971 }
972 // latent (stored °C m/s → write W/m²)
973 {
974 FArrayBox tmp;
975 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
976 tmp.template copy<RunOn::Device>((*vec_lhflx[lev])[mfi.index()], 0, 0, 1);
977 Gpu::streamSynchronize();
978
979 // Convert °C·m/s → W/m²
980 tmp.mult<RunOn::Device>(Hscale);
981
982 auto nc_var = ncf.var("latent");
983 nc_var.put(tmp.dataPtr(),
984 { local_start_nt, local_start_y, local_start_x },
985 { local_nt, local_ny, local_nx });
986 }
987 // sensible (stored °C m/s → write W/m²)
988 {
989 FArrayBox tmp;
990 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
991 tmp.template copy<RunOn::Device>((*vec_shflx[lev])[mfi.index()], 0, 0, 1);
992 Gpu::streamSynchronize();
993
994 // Convert °C·m/s → W/m²
995 tmp.mult<RunOn::Device>(Hscale);
996
997 auto nc_var = ncf.var("sensible");
998 nc_var.put(tmp.dataPtr(),
999 { local_start_nt, local_start_y, local_start_x },
1000 { local_nt, local_ny, local_nx });
1001 }
1002 // lwrad (stored °C m/s → write W/m²)
1003 {
1004 FArrayBox tmp;
1005 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1006 tmp.template copy<RunOn::Device>((*vec_lrflx[lev])[mfi.index()], 0, 0, 1);
1007 Gpu::streamSynchronize();
1008
1009 // Convert °C·m/s → W/m²
1010 tmp.mult<RunOn::Device>(Hscale);
1011
1012 auto nc_var = ncf.var("lwrad");
1013 nc_var.put(tmp.dataPtr(),
1014 { local_start_nt, local_start_y, local_start_x },
1015 { local_nt, local_ny, local_nx });
1016 }
1017 // swrad, note this is stored explicitly as W/m², not degC m/s in REMORA.bulk_flux.cpp
1018 {
1019 FArrayBox tmp;
1020 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1021 tmp.template copy<RunOn::Device>((*vec_srflx[lev])[mfi.index()], 0, 0, 1);
1022 Gpu::streamSynchronize();
1023
1024 auto nc_var = ncf.var("swrad");
1025 nc_var.put(tmp.dataPtr(),
1026 { local_start_nt, local_start_y, local_start_x },
1027 { local_nt, local_ny, local_nx });
1028 }
1029 // evaporation
1030 {
1031 FArrayBox tmp;
1032 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1033 tmp.template copy<RunOn::Device>((*vec_evap[lev])[mfi.index()], 0, 0, 1);
1034 Gpu::streamSynchronize();
1035
1036 auto nc_var = ncf.var("evaporation");
1037 nc_var.put(tmp.dataPtr(),
1038 { local_start_nt, local_start_y, local_start_x },
1039 { local_nt, local_ny, local_nx });
1040 }
1041 // rain
1042 {
1043 FArrayBox tmp;
1044 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1045 tmp.template copy<RunOn::Device>((*vec_rain[lev])[mfi.index()], 0, 0, 1);
1046 Gpu::streamSynchronize();
1047
1048 auto nc_var = ncf.var("rain");
1049 nc_var.put(tmp.dataPtr(),
1050 { local_start_nt, local_start_y, local_start_x },
1051 { local_nt, local_ny, local_nx });
1052 }
1053 } // end output forcing
1054
1055 // **************************************************************************
1056 { // Temp
1057 int comp = -1;
1058 for (int i = 0; i < plot_var_names_3d.size(); i++) {
1059 if (plot_var_names_3d[i] == "temp") comp = i;
1060 }
1061 if (comp >= 0) {
1062 FArrayBox tmp;
1063 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1064 tmp.template copy<RunOn::Device>((*plotMF)[mfi.index()], comp, 0, 1);
1065 Gpu::streamSynchronize();
1066
1067 auto nc_plot_var = ncf.var(plot_var_names_3d[comp]);
1068 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1069 local_nz, local_ny, local_nx });
1070 } // if temp exists in plotMF
1071 } // end temp
1072 // **************************************************************************
1073
1074 // **************************************************************************
1075 { // Salt
1076 int comp = -1;
1077 for (int i = 0; i < plot_var_names_3d.size(); i++) {
1078 if (plot_var_names_3d[i] == "salt") comp = i;
1079 }
1080 if (comp >= 0) {
1081 FArrayBox tmp;
1082 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1083 tmp.template copy<RunOn::Device>((*plotMF)[mfi.index()], comp, 0, 1);
1084 Gpu::streamSynchronize();
1085
1086 auto nc_plot_var = ncf.var(plot_var_names_3d[comp]);
1087 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1088 local_nz, local_ny, local_nx });
1089 } // if salt exists in plotMF
1090 } // end salt
1091 // **************************************************************************
1092
1093 // **************************************************************************
1094 { // Tracer
1095 int comp = -1;
1096 for (int i = 0; i < plot_var_names_3d.size(); i++) {
1097 if (plot_var_names_3d[i] == "tracer") comp = i;
1098 }
1099 if (comp >= 0) {
1100 FArrayBox tmp;
1101 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1102 tmp.template copy<RunOn::Device>((*plotMF)[mfi.index()], comp, 0, 1);
1103 Gpu::streamSynchronize();
1104
1105 auto nc_plot_var = ncf.var(plot_var_names_3d[comp]);
1106 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1107 local_nz, local_ny, local_nx });
1108 } // if tracer exists in plotMF
1109 } // end tracer
1110 // **************************************************************************
1111
1112 // **************************************************************************
1113 { // Vorticity
1114 int comp = -1;
1115 for (int i = 0; i < plot_var_names_3d.size(); i++) {
1116 if (plot_var_names_3d[i] == "vorticity") comp = i;
1117 }
1118 if (comp >= 0) {
1119 FArrayBox tmp;
1120 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1121 tmp.template copy<RunOn::Device>((*plotMF)[mfi.index()], comp, 0, 1);
1122 Gpu::streamSynchronize();
1123
1124 auto nc_plot_var = ncf.var(plot_var_names_3d[comp]);
1125 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1126 local_nz, local_ny, local_nx });
1127 } // if vorticity exists in plotMF
1128 } // end vorticity
1129 // **************************************************************************
1130
1131 // **************************************************************************
1132 // Horizontal mixing coefficients (scaled_to_grid):
1133 // vertically homogeneous and time-invariant -> write static 2D fields once.
1134 // **************************************************************************
1136 { // visc2
1137 FArrayBox tmp;
1138 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1139 tmp.template copy<RunOn::Device>((*vec_visc2_r[lev])[mfi.index()], 0, 0, 1);
1140 Gpu::streamSynchronize();
1141
1142 auto nc_var = ncf.var("visc2");
1143 nc_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1144 }
1145
1146 // diff2_*
1147 for (int n = 0; n < Tracer_comp + 1; ++n) {
1148 const std::string nm = std::string("diff2_") + cons_names[n];
1149 FArrayBox tmp;
1150 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1151 tmp.template copy<RunOn::Device>((*vec_diff2[lev])[mfi.index()], n, 0, 1);
1152 Gpu::streamSynchronize();
1153
1154 auto nc_var = ncf.var(nm);
1155 nc_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1156 }
1157 }
1158
1159 } // subdomain
1160 } // mfi
1161
1162 //ncf.wait_all(irq, &requests[0]);
1163 //requests.resize(0);
1164 //irq = 0;
1165 // Writing u (we loop over cons to get cell-centered box)
1166 for (MFIter mfi(*plotMF, false); mfi.isValid(); ++mfi) {
1167 Box bx = mfi.validbox();
1168
1169 if (subdomain.contains(bx)) {
1170 //
1171 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
1172 //
1173 Box tmp_bx(bx);
1174 tmp_bx.surroundingNodes(0);
1175 if (tmp_bx.smallEnd()[1] == subdomain.smallEnd()[1])
1176 tmp_bx.growLo(1, 1);
1177 if (tmp_bx.bigEnd()[1] == subdomain.bigEnd()[1])
1178 tmp_bx.growHi(1, 1);
1179 Box tmp_bx_2d(tmp_bx);
1180 tmp_bx_2d.makeSlab(2, 0);
1181
1182 //
1183 // These are the dimensions of the data we write for only this box
1184 //
1185 long long local_nx = tmp_bx.length()[0];
1186 long long local_ny = tmp_bx.length()[1];
1187 long long local_nz = tmp_bx.length()[2];
1188
1189 // We do the "+1" because the offset needs to start at 0
1190 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0]);
1191 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1] + 1);
1192 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
1193
1194 if (write_header) {
1195 {
1196 FArrayBox tmp;
1197 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1198 tmp.template copy<RunOn::Device>((*vec_xu[lev])[mfi.index()], 0, 0, 1);
1199 Gpu::streamSynchronize();
1200
1201 auto nc_plot_var = ncf.var("x_u");
1202 //nc_plot_var.par_access(NC_INDEPENDENT);
1203 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1204 }
1205 {
1206 FArrayBox tmp;
1207 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1208 tmp.template copy<RunOn::Device>((*vec_yu[lev])[mfi.index()], 0, 0, 1);
1209 Gpu::streamSynchronize();
1210
1211 auto nc_plot_var = ncf.var("y_u");
1212 //nc_plot_var.par_access(NC_INDEPENDENT);
1213 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1214 }
1215 }
1216
1217 {
1218 FArrayBox tmp;
1219 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1220 tmp.template copy<RunOn::Device>((*xvel_new[lev])[mfi.index()], 0, 0, 1);
1221 Gpu::streamSynchronize();
1222
1223 auto nc_plot_var = ncf.var("u");
1224 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1225 local_nz, local_ny, local_nx });
1226 }
1227
1228 {
1229 FArrayBox tmp;
1230 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1231 tmp.template copy<RunOn::Device>((*vec_ubar[lev])[mfi.index()], 0, 0, 1);
1232 Gpu::streamSynchronize();
1233
1234 auto nc_plot_var = ncf.var("ubar");
1235 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
1236 }
1237 {
1238 FArrayBox tmp;
1239 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1240 tmp.template copy<RunOn::Device>((*vec_sustr[lev])[mfi.index()], 0, 0, 1);
1241 Gpu::streamSynchronize();
1242
1243 auto nc_plot_var = ncf.var("sustr");
1244 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
1245 }
1246 } // in subdomain
1247 } // mfi
1248
1249 // Writing v (we loop over cons to get cell-centered box)
1250 for (MFIter mfi(*plotMF, false); mfi.isValid(); ++mfi) {
1251 Box bx = mfi.validbox();
1252
1253 if (subdomain.contains(bx)) {
1254 //
1255 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
1256 //
1257 Box tmp_bx(bx);
1258 tmp_bx.surroundingNodes(1);
1259 if (tmp_bx.smallEnd()[0] == subdomain.smallEnd()[0])
1260 tmp_bx.growLo(0, 1);
1261 if (tmp_bx.bigEnd()[0] == subdomain.bigEnd()[0])
1262 tmp_bx.growHi(0, 1);
1263 // amrex::Print() << " BX " << bx << std::endl;
1264 // amrex::Print() << "TMP_BX " << tmp_bx << std::endl;
1265
1266 Box tmp_bx_2d(tmp_bx);
1267 tmp_bx_2d.makeSlab(2, 0);
1268
1269 //
1270 // These are the dimensions of the data we write for only this box
1271 //
1272 long long local_nx = tmp_bx.length()[0];
1273 long long local_ny = tmp_bx.length()[1];
1274 long long local_nz = tmp_bx.length()[2];
1275
1276 // We do the "+1" because the offset needs to start at 0
1277 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0] + 1);
1278 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1]);
1279 long long local_start_z = static_cast<long long>(tmp_bx.smallEnd()[2]);
1280
1281 if (write_header) {
1282 {
1283 FArrayBox tmp;
1284 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1285 tmp.template copy<RunOn::Device>((*vec_xv[lev])[mfi.index()], 0, 0, 1);
1286 Gpu::streamSynchronize();
1287
1288 auto nc_plot_var = ncf.var("x_v");
1289 //nc_plot_var.par_access(NC_INDEPENDENT);
1290 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1291 }
1292 {
1293 FArrayBox tmp;
1294 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1295 tmp.template copy<RunOn::Device>((*vec_yv[lev])[mfi.index()], 0, 0, 1);
1296 Gpu::streamSynchronize();
1297
1298 auto nc_plot_var = ncf.var("y_v");
1299 //nc_plot_var.par_access(NC_INDEPENDENT);
1300 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1301 }
1302 }
1303
1304 {
1305 FArrayBox tmp;
1306 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
1307 tmp.template copy<RunOn::Device>((*yvel_new[lev])[mfi.index()], 0, 0, 1);
1308 Gpu::streamSynchronize();
1309
1310 auto nc_plot_var = ncf.var("v");
1311 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
1312 local_nz, local_ny, local_nx });
1313 }
1314
1315 {
1316 FArrayBox tmp;
1317 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1318 tmp.template copy<RunOn::Device>((*vec_vbar[lev])[mfi.index()], 0, 0, 1);
1319 Gpu::streamSynchronize();
1320
1321 auto nc_plot_var = ncf.var("vbar");
1322 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
1323 }
1324
1325 {
1326 FArrayBox tmp;
1327 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1328 tmp.template copy<RunOn::Device>((*vec_svstr[lev])[mfi.index()], 0, 0, 1);
1329 Gpu::streamSynchronize();
1330
1331 auto nc_plot_var = ncf.var("svstr");
1332 nc_plot_var.put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
1333 }
1334
1335 } // in subdomain
1336 } // mfi
1337
1338 for (MFIter mfi(*plotMF, false); mfi.isValid(); ++mfi) {
1339 Box bx = mfi.validbox();
1340
1341 if (subdomain.contains(bx)) {
1342 //
1343 // We only include one grow cell at subdomain boundaries, not internal grid boundaries
1344 //
1345 Box tmp_bx(bx);
1346 tmp_bx.surroundingNodes(0);
1347 tmp_bx.surroundingNodes(1);
1348
1349 Box tmp_bx_2d(tmp_bx);
1350 tmp_bx_2d.makeSlab(2, 0);
1351
1352 //
1353 // These are the dimensions of the data we write for only this box
1354 //
1355 long long local_nx = tmp_bx.length()[0];
1356 long long local_ny = tmp_bx.length()[1];
1357
1358 // We do the "+1" because the offset needs to start at 0
1359 long long local_start_x = static_cast<long long>(tmp_bx.smallEnd()[0]);
1360 long long local_start_y = static_cast<long long>(tmp_bx.smallEnd()[1]);
1361
1362 if (write_header) {
1363 {
1364 FArrayBox tmp;
1365 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1366 tmp.template copy<RunOn::Device>((*vec_xp[lev])[mfi.index()], 0, 0, 1);
1367 Gpu::streamSynchronize();
1368
1369 auto nc_plot_var = ncf.var("x_psi");
1370 //nc_plot_var.par_access(NC_INDEPENDENT);
1371 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1372 }
1373 {
1374 FArrayBox tmp;
1375 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
1376 tmp.template copy<RunOn::Device>((*vec_yp[lev])[mfi.index()], 0, 0, 1);
1377 Gpu::streamSynchronize();
1378
1379 auto nc_plot_var = ncf.var("y_psi");
1380 //nc_plot_var.par_access(NC_INDEPENDENT);
1381 nc_plot_var.put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
1382 }
1383
1384 } // header
1385 } // in subdomain
1386 } // mfi
1387
1389
1390 ncf.close();
1391
1393}
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:353
amrex::Vector< std::string > cons_names
Names of scalars for plotfile output.
Definition REMORA.H:1397
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:409
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:254
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:404
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lrflx
longwave radiation
Definition REMORA.H:333
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yv
y_grid on v-points (2D)
Definition REMORA.H:424
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
Definition REMORA.H:1159
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:296
int history_count
Counter for which time index we are writing to in the netcdf history file.
Definition REMORA.H:1390
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rain
precipitation rate [kg/m^2/s]
Definition REMORA.H:351
bool chunk_history_file
Whether to chunk netcdf history file.
Definition REMORA.H:1383
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:315
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yp
y_grid on psi-points (2D)
Definition REMORA.H:429
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xr
x_grid on rho points (2D)
Definition REMORA.H:412
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xv
x_grid on v-points (2D)
Definition REMORA.H:422
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1266
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:246
amrex::Real start_bdy_time
Start time in the time series of boundary data.
Definition REMORA.H:1154
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:294
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_shflx
sensible heat flux
Definition REMORA.H:339
int steps_per_history_file
Number of time steps per netcdf history file.
Definition REMORA.H:1388
int max_step
maximum number of steps
Definition REMORA.H:1333
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:244
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lhflx
latent heat flux
Definition REMORA.H:337
int plot_int
Plotfile output interval in iterations.
Definition REMORA.H:1372
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:1447
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_r
Harmonic viscosity defined on the rho points (centers)
Definition REMORA.H:284
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:317
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:304
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1273
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1403
static int total_nc_plot_file_step
Definition REMORA.H:1067
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xp
x_grid on psi-points (2D)
Definition REMORA.H:427
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:388
amrex::Gpu::DeviceVector< amrex::Real > Cs_w
Stretching coefficients at w points.
Definition REMORA.H:306
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:386
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yu
y_grid on u-points (2D)
Definition REMORA.H:419
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:1411
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xu
x_grid on u-points (2D)
Definition REMORA.H:417
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:406
amrex::Vector< std::string > plot_var_names_3d
Names of 3D variables to output to AMReX plotfile.
Definition REMORA.H:1393
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflux
Surface tracer flux; input arrays.
Definition REMORA.H:344
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1370
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:312
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_srflx
Shortwave radiation flux [W/m²], defined at rho-points.
Definition REMORA.H:331
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Pair
Air pressure [mb], defined at rho-points.
Definition REMORA.H:328
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Tair
Air temperature [°C], defined at rho-points.
Definition REMORA.H:324
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yr
y_grid on rho points (2D)
Definition REMORA.H:414
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
Definition REMORA.H:286
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)
HorizMixingType horiz_mixing_type
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".