111 std::vector<int> n_cells;
118 subdomain = geom[lev].Domain();
123 int nx = subdomain.length(0);
124 int ny = subdomain.length(1);
125 int nz = subdomain.length(2);
129 amrex::Abort(
"Need to know max_step if writing history file");
135 nt = nt - last_file_index;
141 n_cells.push_back(nx);
142 n_cells.push_back(ny);
143 n_cells.push_back(nz);
145 const std::string nt_name =
"ocean_time";
146 const std::string ndim_name =
"num_geo_dimensions";
148 const std::string flev_name =
"FINEST_LEVEL";
150 const std::string nx_name =
"NX";
151 const std::string ny_name =
"NY";
152 const std::string nz_name =
"NZ";
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";
158 const std::string nx_u_name =
"xi_u";
159 const std::string ny_u_name =
"eta_u";
161 const std::string nx_v_name =
"xi_v";
162 const std::string ny_v_name =
"eta_v";
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";
168 const Real fill_value = 1.0e37_rt;
172 ncf.
put_attr(
"title",
"REMORA data ");
174 ncf.
def_dim(ndim_name, AMREX_SPACEDIM);
176 ncf.
def_dim(nx_r_name, nx + 2);
177 ncf.
def_dim(ny_r_name, ny + 2);
180 ncf.
def_dim(nx_u_name, nx + 1);
181 ncf.
def_dim(ny_u_name, ny + 2);
183 ncf.
def_dim(nx_v_name, nx + 2);
184 ncf.
def_dim(ny_v_name, ny + 1);
186 ncf.
def_dim(nx_p_name, nx + 1);
187 ncf.
def_dim(ny_p_name, ny + 1);
189 ncf.
def_dim(nz_w_name, nz + 1);
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]);
198 ncf.
var(
"probLo").
put_attr(
"long_name",
"Low side of problem domain in internal AMReX grid");
201 ncf.
var(
"probHi").
put_attr(
"long_name",
"High side of problem domain in internal AMReX grid");
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");
209 ncf.
var(
"CellSize").
put_attr(
"long_name",
"Cell size on internal AMReX grid");
213 ncf.
var(
"theta_s").
put_attr(
"long_name",
"S-coordinate surface control parameter");
215 ncf.
var(
"theta_b").
put_attr(
"long_name",
"S-coordinate bottom control parameter");
217 ncf.
var(
"hc").
put_attr(
"long_name",
"S-coordinate parameter, critical depth");
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)");
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");
238 ncf.
var(
"s_w").
put_attr(
"long_name",
"S-coordinate at W-points");
242 ncf.
var(
"pm").
put_attr(
"long_name",
"curvilinear coordinate metric in XI");
246 ncf.
var(
"pm").
put_attr(
"coordinates",
"x_rho y_rho");
250 ncf.
var(
"pn").
put_attr(
"long_name",
"curvilinear coordinate metric in ETA");
254 ncf.
var(
"pn").
put_attr(
"coordinates",
"x_rho y_rho");
258 ncf.
var(
"f").
put_attr(
"long_name",
"Coriolis parameter at RHO-points");
262 ncf.
var(
"f").
put_attr(
"coordinates",
"x_rho y_rho");
266 ncf.
var(
"x_rho").
put_attr(
"long_name",
"x-locations of RHO-points");
268 ncf.
var(
"x_rho").
put_attr(
"field",
"x_rho, scalar");
271 ncf.
var(
"y_rho").
put_attr(
"long_name",
"y-locations of RHO-points");
273 ncf.
var(
"y_rho").
put_attr(
"field",
"y_rho, scalar");
276 ncf.
var(
"x_u").
put_attr(
"long_name",
"x-locations of U-points");
281 ncf.
var(
"y_u").
put_attr(
"long_name",
"y-locations of U-points");
286 ncf.
var(
"x_v").
put_attr(
"long_name",
"x-locations of V-points");
291 ncf.
var(
"y_v").
put_attr(
"long_name",
"y-locations of V-points");
296 ncf.
var(
"x_psi").
put_attr(
"long_name",
"x-locations of PSI-points");
298 ncf.
var(
"x_psi").
put_attr(
"field",
"x_psi, scalar");
301 ncf.
var(
"y_psi").
put_attr(
"long_name",
"y-locations of PSI-points");
303 ncf.
var(
"y_psi").
put_attr(
"field",
"y_psi, scalar");
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");
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.}));
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.}));
323 ncf.
var(
"h").
put_attr(
"long_name",
"bathymetry at RHO-points");
327 ncf.
var(
"h").
put_attr(
"coordinates",
"x_rho y_rho");
331 ncf.
var(
"zeta").
put_attr(
"long_name",
"free-surface");
336 ncf.
var(
"zeta").
put_attr(
"coordinates",
"x_rho y_rho ocean_time");
337 ncf.
var(
"zeta").
put_attr(
"field",
"free-surface, scalar, series");
340 ncf.
var(
"temp").
put_attr(
"long_name",
"potential temperature");
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");
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");
357 ncf.
var(
"tracer").
put_attr(
"long_name",
"passive tracer");
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");
365 ncf.
var(
"u").
put_attr(
"long_name",
"u-momentum component");
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");
374 ncf.
var(
"v").
put_attr(
"long_name",
"v-momentum component");
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");
383 ncf.
var(
"ubar").
put_attr(
"long_name",
"vertically integrated u-momentum component");
384 ncf.
var(
"ubar").
put_attr(
"units",
"meter second-1");
388 ncf.
var(
"ubar").
put_attr(
"coordinates",
"x_u y_u ocean_time");
389 ncf.
var(
"ubar").
put_attr(
"field",
"ubar-velocity, scalar, series");
392 ncf.
var(
"vbar").
put_attr(
"long_name",
"vertically integrated v-momentum component");
393 ncf.
var(
"vbar").
put_attr(
"units",
"meter second-1");
397 ncf.
var(
"vbar").
put_attr(
"coordinates",
"x_v y_v ocean_time");
398 ncf.
var(
"vbar").
put_attr(
"field",
"vbar-velocity, scalar, series");
401 ncf.
var(
"sustr").
put_attr(
"long_name",
"surface u-momentum stress");
402 ncf.
var(
"sustr").
put_attr(
"units",
"newton meter-2");
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");
410 ncf.
var(
"svstr").
put_attr(
"long_name",
"surface v-momentum stress");
411 ncf.
var(
"svstr").
put_attr(
"units",
"newton meter-2");
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");
419 ncf.
put_attr(
"space_dimension", std::vector<int> { AMREX_SPACEDIM });
422 ncf.
put_attr(
"CurrentLevel", std::vector<int> { flev });
423 ncf.
put_attr(
"DefaultGeometry", std::vector<int> { amrex::DefaultGeometry().Coord() });
432 Real dx[AMREX_SPACEDIM];
433 for (
int i = 0; i < AMREX_SPACEDIM; i++) {
434 dx[i] = geom[lev].CellSize()[i];
436 const auto *base = geom[lev].ProbLo();
437 RealBox rb(subdomain, dx, base);
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));
448 ncmpi_begin_indep_data(ncf.
ncid);
449 if (amrex::ParallelDescriptor::IOProcessor())
451 auto nc_probLo = ncf.
var(
"probLo");
453 nc_probLo.
put(probLo.data(), { 0 }, { AMREX_SPACEDIM });
455 auto nc_probHi = ncf.
var(
"probHi");
457 nc_probHi.
put(probHi.data(), { 0 }, { AMREX_SPACEDIM });
459 amrex::Vector<int> smallend;
460 amrex::Vector<int> bigend;
461 for (
int i = lev; i < flev; i++) {
464 for (
int j = 0; j < AMREX_SPACEDIM; j++) {
465 smallend.push_back(subdomain.smallEnd(j));
466 bigend.push_back(subdomain.bigEnd(j));
468 auto nc_Geom_smallend = ncf.
var(
"Geom.smallend");
470 nc_Geom_smallend.
put(smallend.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
473 auto nc_Geom_bigend = ncf.
var(
"Geom.bigend");
475 nc_Geom_bigend.
put(bigend.data(), { static_cast<long long int>(i - lev), 0 }, { 1,
479 amrex::Vector<Real> CellSize;
480 for (
int i = lev; i < flev; i++) {
483 CellSize.push_back(amrex::Real(j));
485 auto nc_CellSize = ncf.
var(
"CellSize");
487 nc_CellSize.
put(CellSize.data(), {
static_cast<long long int>(i - lev), 0 }, { 1,
494 ncf.
var(
"theta_s").
put(&theta_s);
495 ncf.
var(
"theta_b").
put(&theta_b);
498 ncmpi_end_indep_data(ncf.
ncid);
502 ncmpi_begin_indep_data(ncf.
ncid);
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;
511 auto nc_plot_var = ncf.
var(
"ocean_time");
513 nc_plot_var.
put(&
t_new[lev], { local_start_nt }, { local_nt });
518 cons_new[lev]->FillBoundary(geom[lev].periodicity());
524 amrex::Abort(
"Found while writing output: zeta contains nan or inf");
527 amrex::Abort(
"Found while writing output: Temperature contains nan or inf");
530 amrex::Abort(
"Found while writing output: Salinity contains nan or inf");
533 amrex::Abort(
"Found while writing output: Passive tracer contains nan or inf");
536 amrex::Abort(
"Found while writing output: velocity u contains nan or inf");
539 amrex::Abort(
"Found while writing output: velocity ubar contains nan or inf");
542 amrex::Abort(
"Found while writing output: velocity v contains nan or inf");
545 amrex::Abort(
"Found while writing output: velocity vbar contains nan or inf");
548 for (MFIter mfi(*
cons_new[lev],
false); mfi.isValid(); ++mfi) {
549 auto bx = mfi.validbox();
550 if (subdomain.contains(bx)) {
555 if (tmp_bx.smallEnd()[0] == subdomain.smallEnd()[0])
557 if (tmp_bx.smallEnd()[1] == subdomain.smallEnd()[1])
559 if (tmp_bx.bigEnd()[0] == subdomain.bigEnd()[0])
561 if (tmp_bx.bigEnd()[1] == subdomain.bigEnd()[1])
566 Box tmp_bx_2d(tmp_bx);
567 tmp_bx_2d.makeSlab(2, 0);
569 Box tmp_bx_1d(tmp_bx);
570 tmp_bx_1d.makeSlab(0, 0);
571 tmp_bx_1d.makeSlab(1, 0);
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];
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]);
587 if (bx.contains(IntVect(0,0,0)))
591 tmp_srho.resize(tmp_bx_1d, 1, amrex::The_Pinned_Arena());
593 tmp_srho.template copy<RunOn::Device>((*
vec_s_r[lev])[mfi.index()], 0, 0, 1);
594 Gpu::streamSynchronize();
596 auto nc_plot_var = ncf.
var(
"s_rho");
598 nc_plot_var.
put(tmp_srho.dataPtr(), { local_start_z }, { local_nz });
602 tmp_sw.resize(convert(tmp_bx_1d,IntVect(0,0,1)), 1, amrex::The_Pinned_Arena());
604 tmp_sw.template copy<RunOn::Device>((*
vec_s_w[lev])[mfi.index()], 0, 0, 1);
605 Gpu::streamSynchronize();
607 auto nc_plot_var = ncf.
var(
"s_w");
609 nc_plot_var.
put(tmp_sw.dataPtr(), { local_start_z }, { local_nz + 1});
613 tmp_Csrho.resize(tmp_bx_1d, 1, amrex::The_Pinned_Arena());
615 tmp_Csrho.template copy<RunOn::Device>((*
vec_Cs_r[lev])[mfi.index()], 0, 0, 1);
616 Gpu::streamSynchronize();
618 auto nc_plot_var = ncf.
var(
"Cs_r");
620 nc_plot_var.
put(tmp_Csrho.dataPtr(), { local_start_z }, { local_nz });
624 tmp_Csw.resize(convert(tmp_bx_1d,IntVect(0,0,1)), 1, amrex::The_Pinned_Arena());
626 tmp_Csw.template copy<RunOn::Device>((*
vec_Cs_w[lev])[mfi.index()], 0, 0, 1);
627 Gpu::streamSynchronize();
629 auto nc_plot_var = ncf.
var(
"Cs_w");
631 nc_plot_var.
put(tmp_Csw.dataPtr(), { local_start_z }, { local_nz + 1});
637 tmp_bathy.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
640 Gpu::streamSynchronize();
642 auto nc_plot_var = ncf.
var(
"h");
644 nc_plot_var.
put(tmp_bathy.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
649 tmp_pm.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
651 tmp_pm.template copy<RunOn::Device>((*
vec_pm[lev])[mfi.index()], 0, 0, 1);
652 Gpu::streamSynchronize();
654 auto nc_plot_var = ncf.
var(
"pm");
657 nc_plot_var.
put(tmp_pm.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
662 tmp_pn.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
664 tmp_pn.template copy<RunOn::Device>((*
vec_pn[lev])[mfi.index()], 0, 0, 1);
665 Gpu::streamSynchronize();
667 auto nc_plot_var = ncf.
var(
"pn");
669 nc_plot_var.
put(tmp_pn.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
674 tmp_f.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
676 tmp_f.template copy<RunOn::Device>((*
vec_fcor[lev])[mfi.index()], 0, 0, 1);
677 Gpu::streamSynchronize();
679 auto nc_plot_var = ncf.
var(
"f");
681 nc_plot_var.
put(tmp_f.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
686 tmp_xr.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
688 tmp_xr.template copy<RunOn::Device>((*
vec_xr[lev])[mfi.index()], 0, 0, 1);
689 Gpu::streamSynchronize();
691 auto nc_plot_var = ncf.
var(
"x_rho");
694 nc_plot_var.
put(tmp_xr.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
699 tmp_yr.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
701 tmp_yr.template copy<RunOn::Device>((*
vec_yr[lev])[mfi.index()], 0, 0, 1);
702 Gpu::streamSynchronize();
704 auto nc_plot_var = ncf.
var(
"y_rho");
706 nc_plot_var.
put(tmp_yr.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
712 tmp_zeta.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
713 tmp_zeta.template copy<RunOn::Device>((*
vec_Zt_avg1[lev])[mfi.index()], 0, 0, 1);
714 Gpu::streamSynchronize();
716 auto nc_plot_var = ncf.
var(
"zeta");
717 nc_plot_var.
put(tmp_zeta.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny,
723 tmp_temp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
724 tmp_temp.template copy<RunOn::Device>((*
cons_new[lev])[mfi.index()],
Temp_comp, 0, 1);
725 Gpu::streamSynchronize();
727 auto nc_plot_var = ncf.
var(
"temp");
728 nc_plot_var.
put(tmp_temp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
729 local_nz, local_ny, local_nx });
734 tmp_salt.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
735 tmp_salt.template copy<RunOn::Device>((*
cons_new[lev])[mfi.index()],
Salt_comp, 0, 1);
736 Gpu::streamSynchronize();
738 auto nc_plot_var = ncf.
var(
"salt");
739 nc_plot_var.
put(tmp_salt.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
740 local_nz, local_ny, local_nx });
744 FArrayBox tmp_tracer;
745 tmp_tracer.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
747 Gpu::streamSynchronize();
749 auto nc_plot_var = ncf.
var(
"tracer");
750 nc_plot_var.
put(tmp_tracer.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
751 local_nz, local_ny, local_nx });
760 for (MFIter mfi(*
cons_new[lev],
false); mfi.isValid(); ++mfi) {
761 Box bx = mfi.validbox();
763 if (subdomain.contains(bx)) {
768 tmp_bx.surroundingNodes(0);
769 if (tmp_bx.smallEnd()[1] == subdomain.smallEnd()[1])
771 if (tmp_bx.bigEnd()[1] == subdomain.bigEnd()[1])
773 Box tmp_bx_2d(tmp_bx);
774 tmp_bx_2d.makeSlab(2, 0);
779 long long local_nx = tmp_bx.length()[0];
780 long long local_ny = tmp_bx.length()[1];
781 long long local_nz = tmp_bx.length()[2];
784 long long local_start_x =
static_cast<long long>(tmp_bx.smallEnd()[0]);
785 long long local_start_y =
static_cast<long long>(tmp_bx.smallEnd()[1] + 1);
786 long long local_start_z =
static_cast<long long>(tmp_bx.smallEnd()[2]);
791 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
792 tmp.template copy<RunOn::Device>((*
vec_xu[lev])[mfi.index()], 0, 0, 1);
793 Gpu::streamSynchronize();
795 auto nc_plot_var = ncf.
var(
"x_u");
797 nc_plot_var.
put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
801 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
802 tmp.template copy<RunOn::Device>((*
vec_yu[lev])[mfi.index()], 0, 0, 1);
803 Gpu::streamSynchronize();
805 auto nc_plot_var = ncf.
var(
"y_u");
807 nc_plot_var.
put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
813 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
814 tmp.template copy<RunOn::Device>((*
xvel_new[lev])[mfi.index()], 0, 0, 1);
815 Gpu::streamSynchronize();
817 auto nc_plot_var = ncf.
var(
"u");
818 nc_plot_var.
put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
819 local_nz, local_ny, local_nx });
824 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
825 tmp.template copy<RunOn::Device>((*
vec_ubar[lev])[mfi.index()], 0, 0, 1);
826 Gpu::streamSynchronize();
828 auto nc_plot_var = ncf.
var(
"ubar");
829 nc_plot_var.
put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
833 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
834 tmp.template copy<RunOn::Device>((*
vec_sustr[lev])[mfi.index()], 0, 0, 1);
835 Gpu::streamSynchronize();
837 auto nc_plot_var = ncf.
var(
"sustr");
838 nc_plot_var.
put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
844 for (MFIter mfi(*
cons_new[lev],
false); mfi.isValid(); ++mfi) {
845 Box bx = mfi.validbox();
847 if (subdomain.contains(bx)) {
852 tmp_bx.surroundingNodes(1);
853 if (tmp_bx.smallEnd()[0] == subdomain.smallEnd()[0])
855 if (tmp_bx.bigEnd()[0] == subdomain.bigEnd()[0])
860 Box tmp_bx_2d(tmp_bx);
861 tmp_bx_2d.makeSlab(2, 0);
866 long long local_nx = tmp_bx.length()[0];
867 long long local_ny = tmp_bx.length()[1];
868 long long local_nz = tmp_bx.length()[2];
871 long long local_start_x =
static_cast<long long>(tmp_bx.smallEnd()[0] + 1);
872 long long local_start_y =
static_cast<long long>(tmp_bx.smallEnd()[1]);
873 long long local_start_z =
static_cast<long long>(tmp_bx.smallEnd()[2]);
878 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
879 tmp.template copy<RunOn::Device>((*
vec_xv[lev])[mfi.index()], 0, 0, 1);
880 Gpu::streamSynchronize();
882 auto nc_plot_var = ncf.
var(
"x_v");
884 nc_plot_var.
put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
888 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
889 tmp.template copy<RunOn::Device>((*
vec_yv[lev])[mfi.index()], 0, 0, 1);
890 Gpu::streamSynchronize();
892 auto nc_plot_var = ncf.
var(
"y_v");
894 nc_plot_var.
put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
900 tmp.resize(tmp_bx, 1, amrex::The_Pinned_Arena());
901 tmp.template copy<RunOn::Device>((*
yvel_new[lev])[mfi.index()], 0, 0, 1);
902 Gpu::streamSynchronize();
904 auto nc_plot_var = ncf.
var(
"v");
905 nc_plot_var.
put(tmp.dataPtr(), { local_start_nt, local_start_z, local_start_y, local_start_x }, { local_nt,
906 local_nz, local_ny, local_nx });
911 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
912 tmp.template copy<RunOn::Device>((*
vec_vbar[lev])[mfi.index()], 0, 0, 1);
913 Gpu::streamSynchronize();
915 auto nc_plot_var = ncf.
var(
"vbar");
916 nc_plot_var.
put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
920 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
921 tmp.template copy<RunOn::Device>((*
vec_svstr[lev])[mfi.index()], 0, 0, 1);
922 Gpu::streamSynchronize();
924 auto nc_plot_var = ncf.
var(
"svstr");
925 nc_plot_var.
put(tmp.dataPtr(), { local_start_nt, local_start_y, local_start_x }, { local_nt, local_ny, local_nx });
931 for (MFIter mfi(*
cons_new[lev],
false); mfi.isValid(); ++mfi) {
932 Box bx = mfi.validbox();
934 if (subdomain.contains(bx)) {
939 tmp_bx.surroundingNodes(0);
940 tmp_bx.surroundingNodes(1);
942 Box tmp_bx_2d(tmp_bx);
943 tmp_bx_2d.makeSlab(2, 0);
948 long long local_nx = tmp_bx.length()[0];
949 long long local_ny = tmp_bx.length()[1];
952 long long local_start_x =
static_cast<long long>(tmp_bx.smallEnd()[0]);
953 long long local_start_y =
static_cast<long long>(tmp_bx.smallEnd()[1]);
958 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
959 tmp.template copy<RunOn::Device>((*
vec_xp[lev])[mfi.index()], 0, 0, 1);
960 Gpu::streamSynchronize();
962 auto nc_plot_var = ncf.
var(
"x_psi");
964 nc_plot_var.
put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });
968 tmp.resize(tmp_bx_2d, 1, amrex::The_Pinned_Arena());
969 tmp.template copy<RunOn::Device>((*
vec_yp[lev])[mfi.index()], 0, 0, 1);
970 Gpu::streamSynchronize();
972 auto nc_plot_var = ncf.
var(
"y_psi");
974 nc_plot_var.
put(tmp.dataPtr(), { local_start_y, local_start_x }, { local_ny, local_nx });