REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_Plotfile.cpp
Go to the documentation of this file.
1#include <REMORA.H>
2#include "AMReX_Interp_3D_C.H"
3#include "AMReX_PlotFileUtil.H"
4
5using namespace amrex;
6
7PhysBCFunctNoOp null_bc_for_fill;
8
9template<typename V, typename T>
10bool containerHasElement(const V& iterable, const T& query) {
11 return std::find(iterable.begin(), iterable.end(), query) != iterable.end();
12}
13
14// Write plotfile to disk
15void
16REMORA::WritePlotFile (int istep_for_plot)
17{
18#ifndef REMORA_USE_NETCDF
19 amrex::ignore_unused(istep_for_plot);
20#endif
21 Vector<std::string> varnames_3d;
22 varnames_3d.insert(varnames_3d.end(), plot_var_names_3d.begin(), plot_var_names_3d.end());
23
24 Vector<std::string> varnames_2d;
25 varnames_2d.insert(varnames_2d.end(), plot_var_names_2d.begin(), plot_var_names_2d.end());
26
27 Vector<std::string> varnames_2d_rho;
28 Vector<std::string> varnames_2d_u;
29 Vector<std::string> varnames_2d_v;
30
31 const int ncomp_mf_3d = varnames_3d.size();
32 const auto ngrow_vars = IntVect(NGROW-1,NGROW-1,0);
33
34 // These are the ncomp for the 2D cell-centered, x-face-based, y-face-based MultiFabs respectively
35 int ncomp_mf_2d_rho = 0;
36 int ncomp_mf_2d_u = 0;
37 int ncomp_mf_2d_v = 0;
38
39 // Check to see if we found all the requested variables
40 for (auto plot_name : varnames_2d) {
41 {
42 if (plot_name == "zeta" ) {varnames_2d_rho.push_back(plot_name); ncomp_mf_2d_rho++;}
43 if (plot_name == "h" ) {varnames_2d_rho.push_back(plot_name); ncomp_mf_2d_rho++;}
44 if (plot_name == "ubar" ) {varnames_2d_u.push_back(plot_name); ncomp_mf_2d_u++;}
45 if (plot_name == "sustr") {varnames_2d_u.push_back(plot_name); ncomp_mf_2d_u++;}
46 if (plot_name == "bustr") {varnames_2d_u.push_back(plot_name); ncomp_mf_2d_u++;}
47 if (plot_name == "vbar" ) {varnames_2d_v.push_back(plot_name); ncomp_mf_2d_v++;}
48 if (plot_name == "svstr") {varnames_2d_v.push_back(plot_name); ncomp_mf_2d_v++;}
49 if (plot_name == "bvstr") {varnames_2d_v.push_back(plot_name); ncomp_mf_2d_v++;}
50 }
51 }
52
53 // We fillpatch here because some of the derived quantities require derivatives
54 // which require ghost cells to be filled. Don't fill the boundary, though.
55 for (int lev = 0; lev <= finest_level; ++lev) {
56 FillPatchNoBC(lev, t_new[lev], *cons_new[lev], cons_new, BdyVars::t,0,true,false);
57 FillPatchNoBC(lev, t_new[lev], *xvel_new[lev], xvel_new, BdyVars::u,0,true,false);
58 FillPatchNoBC(lev, t_new[lev], *yvel_new[lev], yvel_new, BdyVars::v,0,true,false);
59 FillPatchNoBC(lev, t_new[lev], *zvel_new[lev], zvel_new, BdyVars::null,0,true,false);
60 }
61
62 for (int lev = 0; lev <= finest_level; ++lev) {
64 }
65
66 // Array of 3D MultiFabs to hold the plotfile data
67 Vector<MultiFab> plotMF(finest_level+1);
68 for (int lev = 0; lev <= finest_level; ++lev) {
69 plotMF[lev].define(grids[lev], dmap[lev], ncomp_mf_3d, ngrow_vars);
70 plotMF[lev].setVal(1.234e20);
71 }
72
73 // Array of 2D MultiFabs to hold the plotfile data
74 Vector<MultiFab> mf_2d_rho(finest_level+1);
75 Vector<MultiFab> mf_2d_u(finest_level+1);
76 Vector<MultiFab> mf_2d_v(finest_level+1);
77 for (int lev = 0; lev <= finest_level; ++lev) {
78 BoxArray ba(grids[lev]);
79 BoxList bl2d = ba.boxList();
80 for (auto& b : bl2d) {
81 b.setRange(2,0);
82 }
83 BoxArray ba2d(std::move(bl2d));
84 mf_2d_rho[lev].define(ba2d, dmap[lev], ncomp_mf_2d_rho, IntVect(0,0,0));
85 mf_2d_u[lev].define(ba2d, dmap[lev], ncomp_mf_2d_u , IntVect(0,0,0));
86 mf_2d_v[lev].define(ba2d, dmap[lev], ncomp_mf_2d_v , IntVect(0,0,0));
87 }
88
89 // Array of MultiFabs for nodal data
90 Vector<MultiFab> mf_nd(finest_level+1);
91 for (int lev = 0; lev <= finest_level; ++lev) {
92 BoxArray nodal_grids(grids[lev]); nodal_grids.surroundingNodes();
93 mf_nd[lev].define(nodal_grids, dmap[lev], AMREX_SPACEDIM, 0);
94 mf_nd[lev].setVal(0.);
95 }
96
97 // Vector of MultiFabs for face-centered velocity
98 Vector<MultiFab> mf_u(finest_level+1);
99 Vector<MultiFab> mf_v(finest_level+1);
100 Vector<MultiFab> mf_w(finest_level+1);
102 for (int lev = 0; lev <= finest_level; ++lev) {
103 BoxArray grid_stag_u(grids[lev]); grid_stag_u.surroundingNodes(0);
104 BoxArray grid_stag_v(grids[lev]); grid_stag_v.surroundingNodes(1);
105 BoxArray grid_stag_w(grids[lev]); grid_stag_w.surroundingNodes(2);
106 mf_u[lev].define(grid_stag_u, dmap[lev], 1, 0);
107 mf_v[lev].define(grid_stag_v, dmap[lev], 1, 0);
108 mf_w[lev].define(grid_stag_w, dmap[lev], 1, 0);
109 MultiFab::Copy(mf_u[lev],*xvel_new[lev],0,0,1,0);
110 MultiFab::Copy(mf_v[lev],*yvel_new[lev],0,0,1,0);
111 MultiFab::Copy(mf_w[lev],*zvel_new[lev],0,0,1,0);
112 }
113 }
114
115 // Array of MultiFabs for cell-centered velocity
116 Vector<MultiFab> mf_cc_vel(finest_level+1);
117
118 if (containerHasElement(plot_var_names_3d, "x_velocity") ||
119 containerHasElement(plot_var_names_3d, "y_velocity") ||
120 containerHasElement(plot_var_names_3d, "z_velocity") ||
121 containerHasElement(plot_var_names_3d, "vorticity") ) {
122
123 for (int lev = 0; lev <= finest_level; ++lev) {
124 mf_cc_vel[lev].define(grids[lev], dmap[lev], AMREX_SPACEDIM, IntVect(1,1,0));
125 mf_cc_vel[lev].setVal(0.0_rt); // zero out velocity in case we have any wall boundaries
126 average_face_to_cellcenter(mf_cc_vel[lev],0,
127 Array<const MultiFab*,3>{xvel_new[lev],yvel_new[lev],zvel_new[lev]},IntVect(1,1,0));
128 mf_cc_vel[lev].FillBoundary(geom[lev].periodicity());
129 } // lev
130
131 // We need ghost cells if computing vorticity
132 amrex::Interpolater* mapper = &cell_cons_interp;
133 if ( containerHasElement(plot_var_names_3d, "vorticity") ) {
134 for (int lev = 1; lev <= finest_level; ++lev) {
135 Vector<MultiFab*> fmf = {&(mf_cc_vel[lev]), &(mf_cc_vel[lev])};
136 Vector<Real> ftime = {t_new[lev], t_new[lev]};
137 Vector<MultiFab*> cmf = {&mf_cc_vel[lev-1], &mf_cc_vel[lev-1]};
138 Vector<Real> ctime = {t_new[lev], t_new[lev]};
139
140 MultiFab mf_to_fill;
141 amrex::FillPatchTwoLevels(mf_cc_vel[lev], mf_cc_vel[lev].nGrowVect(), IntVect(0,0,0),
142 t_new[lev], cmf, ctime, fmf, ftime,
143 0, 0, mf_cc_vel[lev].nComp(), geom[lev-1], geom[lev],
144 refRatio(lev-1), mapper, domain_bcs_type, BCVars::foextrap_bc);
145 } // lev
146 } // if
147 } // if
148
149 int icomp_rho = 0;
150 for (auto plot_name : varnames_2d_rho)
151 {
152 if (plot_name == "zeta" ) {
153 for (int lev = 0; lev <= finest_level; ++lev) { MultiFab::Copy(mf_2d_rho[lev],*vec_Zt_avg1[lev],0,icomp_rho,1,0); }
154 icomp_rho++;
155 }
156 if (plot_name == "h" ) {
157 for (int lev = 0; lev <= finest_level; ++lev) { MultiFab::Copy(mf_2d_rho[lev],*vec_h[lev],0,icomp_rho,1,0); }
158 icomp_rho++;
159 }
160 }
161
162 int icomp_u = 0;
163 for (auto plot_name : varnames_2d_u)
164 {
165 if (plot_name == "ubar" ) {
166 for (int lev = 0; lev <= finest_level; ++lev) {
167 MultiFab::Copy(mf_2d_u[lev],*vec_DU_avg1[lev],0,icomp_u,1,0);
168 }
169 icomp_u++;
170 }
171 if (plot_name == "sustr" ) {
172 for (int lev = 0; lev <= finest_level; ++lev) { MultiFab::Copy(mf_2d_u[lev],*vec_sustr[lev],0,icomp_u,1,0); }
173 icomp_u++;
174 }
175 if (plot_name == "bustr" ) {
176 for (int lev = 0; lev <= finest_level; ++lev) { MultiFab::Copy(mf_2d_u[lev],*vec_bustr[lev],0,icomp_u,1,0); }
177 icomp_u++;
178 }
179 }
180
181 int icomp_v = 0;
182 for (auto plot_name : varnames_2d_v)
183 {
184 if (plot_name == "vbar" ) {
185 for (int lev = 0; lev <= finest_level; ++lev) { MultiFab::Copy(mf_2d_v[lev],*vec_DV_avg1[lev],0,icomp_v,1,0); }
186 icomp_v++;
187 }
188 if (plot_name == "svstr" ) {
189 for (int lev = 0; lev <= finest_level; ++lev) { MultiFab::Copy(mf_2d_v[lev],*vec_svstr[lev],0,icomp_v,1,0); }
190 icomp_v++;
191 }
192 if (plot_name == "bvstr" ) {
193 for (int lev = 0; lev <= finest_level; ++lev) { MultiFab::Copy(mf_2d_v[lev],*vec_bvstr[lev],0,icomp_v,1,0); }
194 icomp_v++;
195 }
196 }
197
198 for (int lev = 0; lev <= finest_level; ++lev)
199 {
200 int mf_comp = 0;
201
202 // First, copy any of the conserved state variables into the output plotfile
203 AMREX_ALWAYS_ASSERT(cons_names.size() == NCONS);
204 for (int i = 0; i < NCONS; ++i) {
206 if (cons_new[lev]->contains_nan() || cons_new[lev]->contains_inf()) {
207 amrex::Abort("Found while writing output: Cons (salt, temp, or tracer, etc) contains nan or inf");
208 }
209 MultiFab::Copy(plotMF[lev],*cons_new[lev],i,mf_comp,1,ngrow_vars);
210 mf_comp++;
211 }
212 } // NCONS
213
214 // Next, check for velocities
215 if (containerHasElement(plot_var_names_3d, "x_velocity")) {
216 if (mf_cc_vel[lev].contains_nan(0,1) || mf_cc_vel[lev].contains_inf(0,1)) {
217 amrex::Abort("Found while writing output: u velocity contains nan or inf");
218 }
219 MultiFab::Copy(plotMF[lev], mf_cc_vel[lev], 0, mf_comp, 1, 0);
220 mf_comp += 1;
221 }
222 if (containerHasElement(plot_var_names_3d, "y_velocity")) {
223 if (mf_cc_vel[lev].contains_nan(1,1) || mf_cc_vel[lev].contains_inf(1,1)) {
224 amrex::Abort("Found while writing output: v velocity contains nan or inf");
225 }
226 MultiFab::Copy(plotMF[lev], mf_cc_vel[lev], 1, mf_comp, 1, 0);
227 mf_comp += 1;
228 }
229 if (containerHasElement(plot_var_names_3d, "z_velocity")) {
230 if (mf_cc_vel[lev].contains_nan(2,1) || mf_cc_vel[lev].contains_inf(2,1)) {
231 amrex::Abort("Found while writing output: z velocity contains nan or inf");
232 }
233 MultiFab::Copy(plotMF[lev], mf_cc_vel[lev], 2, mf_comp, 1, 0);
234 mf_comp += 1;
235 }
236
237 // Define standard process for calling the functions in Derive.cpp
238 auto calculate_derived = [&](const std::string& der_name,
239 decltype(derived::remora_dernull)& der_function)
240 {
241 if (containerHasElement(plot_var_names_3d, der_name)) {
242 MultiFab dmf(plotMF[lev], make_alias, mf_comp, 1);
243#ifdef _OPENMP
244#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
245#endif
246 for (MFIter mfi(dmf, TilingIfNotGPU()); mfi.isValid(); ++mfi)
247 {
248 const Box& bx = mfi.tilebox();
249 auto& dfab = dmf[mfi];
250
251 if (der_name == "vorticity") {
252 auto const& sfab = mf_cc_vel[lev][mfi];
253 der_function(bx, dfab, 0, 1, sfab, vec_pm[lev]->const_array(mfi), vec_pn[lev]->const_array(mfi), vec_mskr[lev]->const_array(mfi), Geom(lev), t_new[0], nullptr, lev);
254 } else {
255 auto const& sfab = (*cons_new[lev])[mfi];
256 der_function(bx, dfab, 0, 1, sfab, vec_pm[lev]->const_array(mfi), vec_pn[lev]->const_array(mfi), vec_mskr[lev]->const_array(mfi), Geom(lev), t_new[0], nullptr, lev);
257 }
258 }
259
260 mf_comp++;
261 }
262 };
263
264 // Note: All derived variables must be computed in order of "derived_names" defined in REMORA.H
265 calculate_derived("vorticity", derived::remora_dervort);
266
267 // Fill cell-centered location
268 Real dx = Geom()[lev].CellSizeArray()[0];
269 Real dy = Geom()[lev].CellSizeArray()[1];
270
271 // Next, check for location names -- if we write one we write all
275 {
276 MultiFab dmf(plotMF[lev], make_alias, mf_comp, AMREX_SPACEDIM);
277#ifdef _OPENMP
278#pragma omp parallel if (Gpu::notInLaunchRegion())
279#endif
280 for (MFIter mfi(dmf, TilingIfNotGPU()); mfi.isValid(); ++mfi) {
281 const Box& bx = mfi.tilebox();
282 const Array4<Real> loc_arr = dmf.array(mfi);
283 const Array4<Real const> zp_arr = vec_z_phys_nd[lev]->const_array(mfi);
284
285 ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
286 loc_arr(i,j,k,0) = (i+0.5_rt) * dx;
287 loc_arr(i,j,k,1) = (j+0.5_rt) * dy;
288 loc_arr(i,j,k,2) = 0.125_rt * (zp_arr(i,j ,k ) + zp_arr(i+1,j ,k ) +
289 zp_arr(i,j+1,k ) + zp_arr(i+1,j+1,k ) +
290 zp_arr(i,j ,k+1) + zp_arr(i+1,j ,k+1) +
291 zp_arr(i,j+1,k+1) + zp_arr(i+1,j+1,k+1) );
292 });
293 } // mfi
294 mf_comp += AMREX_SPACEDIM;
295 } // if containerHasElement
296
297#ifdef REMORA_USE_PARTICLES
298 const auto& particles_namelist( particleData.getNames() );
299 for (ParticlesNamesVector::size_type i = 0; i < particles_namelist.size(); i++) {
300 if (containerHasElement(plot_var_names_3d, std::string(particles_namelist[i]+"_count")))
301 {
302 MultiFab temp_dat(plotMF[lev].boxArray(), plotMF[lev].DistributionMap(), 1, 0);
303 temp_dat.setVal(0);
304 particleData[particles_namelist[i]]->Increment(temp_dat, lev);
305 MultiFab::Copy(plotMF[lev], temp_dat, 0, mf_comp, 1, 0);
306 mf_comp += 1;
307 }
308 }
309
310 Vector<std::string> particle_mesh_plot_names(0);
311 particleData.GetMeshPlotVarNames( particle_mesh_plot_names );
312 for (int i = 0; i < particle_mesh_plot_names.size(); i++) {
313 std::string plot_var_name(particle_mesh_plot_names[i]);
314 if (containerHasElement(plot_var_names_3d, plot_var_name) ) {
315 MultiFab temp_dat(plotMF[lev].boxArray(), plotMF[lev].DistributionMap(), 1, 1);
316 temp_dat.setVal(0);
317 particleData.GetMeshPlotVar(plot_var_name, temp_dat, lev);
318 MultiFab::Copy(plotMF[lev], temp_dat, 0, mf_comp, 1, 0);
319 mf_comp += 1;
320 }
321 }
322#endif
323
324 MultiFab::Copy(mf_nd[lev],*vec_z_phys_nd[lev],0,2,1,0);
325 Real dz = Geom()[lev].CellSizeArray()[2];
326 int N = Geom()[lev].Domain().size()[2];
327
328#ifdef _OPENMP
329#pragma omp parallel if (Gpu::notInLaunchRegion())
330#endif
331 for (MFIter mfi(mf_nd[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
332 {
333 const Box& bx = mfi.tilebox();
334 Array4<Real> mf_arr = mf_nd[lev].array(mfi);
335 ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
336 mf_arr(i,j,k,2) = mf_arr(i,j,k,2) + (N-k) * dz;
337 });
338 } // mfi
339
340 } // lev
341
344 {
345
346 std::string plotfilename = Concatenate(plot_file_name, istep[0], file_min_digits);
347
348 if (finest_level == 0)
349 {
351 amrex::Print() << "Writing plotfile " << plotfilename << "\n";
352 WriteMultiLevelPlotfileWithBathymetry(plotfilename, finest_level+1,
353 GetVecOfConstPtrs(plotMF),
354 GetVecOfConstPtrs(mf_nd),
355 GetVecOfConstPtrs(mf_u),
356 GetVecOfConstPtrs(mf_v),
357 GetVecOfConstPtrs(mf_w),
358 GetVecOfConstPtrs(mf_2d_rho),
359 GetVecOfConstPtrs(mf_2d_u),
360 GetVecOfConstPtrs(mf_2d_v),
361 varnames_3d, varnames_2d_rho,
362 varnames_2d_u, varnames_2d_v,
363 Geom(),
364 t_new[0], istep, refRatio());
365 writeJobInfo(plotfilename);
366
367#ifdef REMORA_USE_PARTICLES
368 particleData.Checkpoint(plotfilename);
369#endif
370
371#ifdef REMORA_USE_HDF5
372 } else if (plotfile_type == PlotfileType::hdf5) {
373 amrex::Print() << "Writing plotfile " << plotfilename+"d01.h5" << "\n";
374 WriteMultiLevelPlotfileHDF5(plotfilename, finest_level+1,
375 GetVecOfConstPtrs(plotMF),
376 varnames_3d,
377 Geom(), t_new[0], istep, refRatio());
378#endif
379 }
380
381 } else { // multilevel
383 amrex::Print() << "Writing plotfile " << plotfilename << "\n";
384 int lev0 = 0;
385 [[maybe_unused]] int desired_ratio = std::max(std::max(ref_ratio[lev0][0],ref_ratio[lev0][1]),ref_ratio[lev0][2]);
386 bool any_ratio_one = ( ( (ref_ratio[lev0][0] == 1) || (ref_ratio[lev0][1] == 1) ) ||
387 (ref_ratio[lev0][2] == 1) );
388 for (int lev = 1; lev < finest_level; lev++) {
389 any_ratio_one = any_ratio_one ||
390 ( ( (ref_ratio[lev][0] == 1) || (ref_ratio[lev][1] == 1) ) ||
391 (ref_ratio[lev][2] == 1) );
392 }
393 if (any_ratio_one && expand_plotvars_to_unif_rr) {
394 Vector<IntVect> r2(finest_level);
395 Vector<Geometry> g2(finest_level+1);
396 Vector<MultiFab> mf2(finest_level+1);
397
398 mf2[0].define(grids[0], dmap[0], ncomp_mf_3d, 0);
399
400 // Copy level 0 as is
401 MultiFab::Copy(mf2[0],plotMF[0],0,0,plotMF[0].nComp(),0);
402
403 // Define a new multi-level array of Geometry's so that we pass the new "domain" at lev > 0
404 Array<int,AMREX_SPACEDIM> periodicity =
405 {Geom()[0].isPeriodic(0),Geom()[0].isPeriodic(1),Geom()[0].isPeriodic(2)};
406 g2[0].define(Geom()[0].Domain(),&(Geom()[0].ProbDomain()),0,periodicity.data());
407
408 r2[0] = IntVect(1,1,ref_ratio[0][0]);
409 for (int lev = 1; lev <= finest_level; ++lev) {
410 if (lev > 1) {
411 r2[lev-1][0] = 1;
412 r2[lev-1][1] = 1;
413 r2[lev-1][2] = r2[lev-2][2] * ref_ratio[lev-1][0];
414 }
415
416 mf2[lev].define(refine(grids[lev],r2[lev-1]), dmap[lev], ncomp_mf_3d, 0);
417
418 // Set the new problem domain
419 Box d2(Geom()[lev].Domain());
420 d2.refine(r2[lev-1]);
421
422 g2[lev].define(d2,&(Geom()[lev].ProbDomain()),0,periodicity.data());
423 }
424
425 // Make a vector of BCRec with default values so we can use it here -- note the values
426 // aren't actually used because we do PCInterp
427 amrex::Vector<amrex::BCRec> null_dom_bcs;
428 null_dom_bcs.resize(mf2[0].nComp());
429 for (int n = 0; n < mf2[0].nComp(); n++) {
430 for (int dir = 0; dir < AMREX_SPACEDIM; dir++) {
431 null_dom_bcs[n].setLo(dir, REMORABCType::int_dir);
432 null_dom_bcs[n].setHi(dir, REMORABCType::int_dir);
433 }
434 }
435
436 // Do piecewise interpolation of mf into mf2
437 for (int lev = 1; lev <= finest_level; ++lev) {
438 Interpolater* mapper_c = &pc_interp;
439 InterpFromCoarseLevel(mf2[lev], t_new[lev], plotMF[lev],
440 0, 0, mf2[lev].nComp(),
441 geom[lev], g2[lev],
443 r2[lev-1], mapper_c, null_dom_bcs, 0);
444 }
445
446 // Define an effective ref_ratio which is isotropic to be passed into WriteMultiLevelPlotfile
447 Vector<IntVect> rr(finest_level);
448 for (int lev = 0; lev < finest_level; ++lev) {
449 rr[lev] = IntVect(ref_ratio[lev][0],ref_ratio[lev][1],ref_ratio[lev][0]);
450 }
451
452 WriteMultiLevelPlotfileWithBathymetry(plotfilename, finest_level+1,
453 GetVecOfConstPtrs(mf2),
454 GetVecOfConstPtrs(mf_nd),
455 GetVecOfConstPtrs(mf_u),
456 GetVecOfConstPtrs(mf_v),
457 GetVecOfConstPtrs(mf_w),
458 GetVecOfConstPtrs(mf_2d_rho),
459 GetVecOfConstPtrs(mf_2d_u),
460 GetVecOfConstPtrs(mf_2d_v),
461 varnames_3d, varnames_2d_rho,
462 varnames_2d_u, varnames_2d_v,
463 g2,
464 t_new[0], istep, rr);
465 writeJobInfo(plotfilename);
466
467#ifdef REMORA_USE_PARTICLES
468 particleData.Checkpoint(plotfilename);
469#endif
470 } else {
471 WriteMultiLevelPlotfileWithBathymetry(plotfilename, finest_level+1,
472 GetVecOfConstPtrs(plotMF),
473 GetVecOfConstPtrs(mf_nd),
474 GetVecOfConstPtrs(mf_u),
475 GetVecOfConstPtrs(mf_v),
476 GetVecOfConstPtrs(mf_w),
477 GetVecOfConstPtrs(mf_2d_rho),
478 GetVecOfConstPtrs(mf_2d_u),
479 GetVecOfConstPtrs(mf_2d_v),
480 varnames_3d, varnames_2d_rho,
481 varnames_2d_u, varnames_2d_v,
482 Geom(),
483 t_new[0], istep, ref_ratio);
484 writeJobInfo(plotfilename);
485#ifdef REMORA_USE_PARTICLES
486 particleData.Checkpoint(plotfilename);
487#endif
488 }
489 }
490 } // end multi-level
491 for (int lev = 0; lev <= finest_level; ++lev) {
493 }
494
495 }
496#ifdef REMORA_USE_NETCDF
498 {
499 // Currently this is hard-coded to plot only level 0
500 AMREX_ASSERT(finest_level == 0);
501 int lev = 0;
502 plotMF[0].FillBoundary(geom[lev].periodicity());
503 WriteNCPlotFile(istep_for_plot,&plotMF[lev]);
504 } // end if plotfile_type == netcdf
505#endif
506}
507
508/**
509 * @param plotfilename name of plotfile to write to
510 * @param nlevels number of levels to write out
511 * @param mf MultiFab of data to write out
512 * @param mf_nd Multifab of nodal data to write out
513 * @param varnames_3d 3D variable names to write out
514 * @param varnames_2d_rho 2D cell-centered variable names to write out
515 * @param varnames_2d_u 2D x-face-based variable names to write out
516 * @param varnames_2d_v 2D y-face-based variable names to write out
517 * @param my_geom geometry to use for writing plotfile
518 * @param time time at which to output
519 * @param level_steps vector over level of iterations
520 * @param rr refinement ratio to use for writing plotfile
521 * @param versionName version string for VisIt
522 * @param levelPrefix string to prepend to level number
523 * @param mfPrefix subdirectory for multifab data
524 * @param extra_dirs additional subdirectories within plotfile
525 */
526 void
527 REMORA::WriteMultiLevelPlotfileWithBathymetry (const std::string& plotfilename, int nlevels,
528 const Vector<const MultiFab*>& mf,
529 const Vector<const MultiFab*>& mf_nd,
530 const Vector<const MultiFab*>& mf_u,
531 const Vector<const MultiFab*>& mf_v,
532 const Vector<const MultiFab*>& mf_w,
533 const Vector<const MultiFab*>& mf_2d_rho,
534 const Vector<const MultiFab*>& mf_2d_u,
535 const Vector<const MultiFab*>& mf_2d_v,
536 const Vector<std::string>& varnames_3d,
537 const Vector<std::string>& varnames_2d_rho,
538 const Vector<std::string>& varnames_2d_u,
539 const Vector<std::string>& varnames_2d_v,
540 const Vector<Geometry>& my_geom,
541 Real time,
542 const Vector<int>& level_steps,
543 const Vector<IntVect>& rr,
544 const std::string &versionName,
545 const std::string &levelPrefix,
546 const std::string &mfPrefix,
547 const Vector<std::string>& extra_dirs) const
548{
549 BL_PROFILE("WriteMultiLevelPlotfileWithBathymetry()");
550
551 AMREX_ASSERT(nlevels <= mf.size());
552 AMREX_ASSERT(nlevels <= ref_ratio.size()+1);
553 AMREX_ASSERT(nlevels <= level_steps.size());
554
555 AMREX_ASSERT(mf[0]->nComp() == varnames_3d.size());
556
557 bool callBarrier(false);
558 PreBuildDirectorHierarchy(plotfilename, levelPrefix, nlevels, callBarrier);
559 if (!extra_dirs.empty()) {
560 for (const auto& d : extra_dirs) {
561 const std::string ed = plotfilename+"/"+d;
562 PreBuildDirectorHierarchy(ed, levelPrefix, nlevels, callBarrier);
563 }
564 }
565 ParallelDescriptor::Barrier();
566
567 if (ParallelDescriptor::MyProc() == ParallelDescriptor::NProcs()-1) {
568 Vector<BoxArray> boxArrays(nlevels);
569 for(int level(0); level < boxArrays.size(); ++level) {
570 boxArrays[level] = mf[level]->boxArray();
571 }
572
573 auto f = [=]() {
574 VisMF::IO_Buffer io_buffer(VisMF::IO_Buffer_Size);
575 std::string HeaderFileName(plotfilename + "/Header");
576 std::ofstream HeaderFile;
577 HeaderFile.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
578 HeaderFile.open(HeaderFileName.c_str(), std::ofstream::out |
579 std::ofstream::trunc |
580 std::ofstream::binary);
581 if( ! HeaderFile.good()) FileOpenFailed(HeaderFileName);
582 WriteGenericPlotfileHeaderWithBathymetry(HeaderFile, nlevels, boxArrays, varnames_3d,
583 varnames_2d_rho, varnames_2d_u, varnames_2d_v,
584 my_geom, time, level_steps, rr, versionName,
585 levelPrefix, mfPrefix);
586 };
587
588 if (AsyncOut::UseAsyncOut()) {
589 AsyncOut::Submit(std::move(f));
590 } else {
591 f();
592 }
593 }
594
595 std::string mf_nodal_prefix = "Nu_nd";
596 std::string mf_uface_prefix = "UFace";
597 std::string mf_vface_prefix = "VFace";
598 std::string mf_wface_prefix = "WFace";
599 std::string mf_2d_rho_prefix = "rho2d";
600 std::string mf_2d_u_prefix = "u2d";
601 std::string mf_2d_v_prefix = "v2d";
602
603 for (int level = 0; level <= finest_level; ++level)
604 {
605 if (AsyncOut::UseAsyncOut()) {
606 VisMF::AsyncWrite(*mf[level],
607 MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mfPrefix),
608 true);
609 VisMF::AsyncWrite(*mf_nd[level],
610 MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_nodal_prefix),
611 true);
613 VisMF::AsyncWrite(*mf_u[level],
614 MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_uface_prefix),
615 true);
616 VisMF::AsyncWrite(*mf_v[level],
617 MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_vface_prefix),
618 true);
619 VisMF::AsyncWrite(*mf_w[level],
620 MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_wface_prefix),
621 true);
622 }
623 if (mf_2d_rho[level]->nComp() > 0) {
624 VisMF::AsyncWrite(*mf_2d_rho[level],
625 MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_2d_rho_prefix),
626 true);
627 }
628 if (mf_2d_u[level]->nComp() > 0) {
629 VisMF::AsyncWrite(*mf_2d_u[level],
630 MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_2d_u_prefix),
631 true);
632 }
633 if (mf_2d_v[level]->nComp() > 0) {
634 VisMF::AsyncWrite(*mf_2d_v[level],
635 MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_2d_v_prefix),
636 true);
637 }
638 } else {
639 const MultiFab* data;
640 std::unique_ptr<MultiFab> mf_tmp;
641 if (mf[level]->nGrowVect() != 0) {
642 mf_tmp = std::make_unique<MultiFab>(mf[level]->boxArray(),
643 mf[level]->DistributionMap(),
644 mf[level]->nComp(), 0, MFInfo(),
645 mf[level]->Factory());
646 MultiFab::Copy(*mf_tmp, *mf[level], 0, 0, mf[level]->nComp(), 0);
647 data = mf_tmp.get();
648 } else {
649 data = mf[level];
650 }
651 VisMF::Write(*data , MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mfPrefix));
652 VisMF::Write(*mf_nd[level], MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_nodal_prefix));
654 VisMF::Write(*mf_u[level], MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_uface_prefix));
655 VisMF::Write(*mf_v[level], MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_vface_prefix));
656 VisMF::Write(*mf_w[level], MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_wface_prefix));
657 }
658 if (mf_2d_rho[level]->nComp() > 0) {
659 VisMF::Write(*mf_2d_rho[level], MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_2d_rho_prefix));
660 }
661 if (mf_2d_u[level]->nComp() > 0) {
662 VisMF::Write(*mf_2d_u[level], MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_2d_u_prefix));
663 }
664 if (mf_2d_v[level]->nComp() > 0) {
665 VisMF::Write(*mf_2d_v[level], MultiFabFileFullPrefix(level, plotfilename, levelPrefix, mf_2d_v_prefix));
666 }
667 }
668 } // level
669}
670
671/**
672 * @param HeaderFile output stream for header
673 * @param nlevels number of levels to write out
674 * @param bArray vector over levels of BoxArrays
675 * @param varnames_3d 3D variable names to write out
676 * @param varnames_2d 2D variable names to write out
677 * @param my_geom geometry to use for writing plotfile
678 * @param time time at which to output
679 * @param level_steps vector over level of iterations
680 * @param my_ref_ratio refinement ratio to use for writing plotfile
681 * @param versionName version string for VisIt
682 * @param levelPrefix string to prepend to level number
683 * @param mfPrefix subdirectory for multifab data
684 */
685void
687 [[maybe_unused]] int nlevels,
688 const Vector<BoxArray> &bArray,
689 const Vector<std::string> &varnames_3d,
690 const Vector<std::string> &varnames_2d_rho,
691 const Vector<std::string> &varnames_2d_u,
692 const Vector<std::string> &varnames_2d_v,
693 const Vector<Geometry>& my_geom,
694 Real time,
695 const Vector<int> &level_steps,
696 const Vector<IntVect>& my_ref_ratio,
697 const std::string &versionName,
698 const std::string &levelPrefix,
699 const std::string &mfPrefix) const
700{
701 AMREX_ASSERT(nlevels <= bArray.size());
702 AMREX_ASSERT(nlevels <= ref_ratio.size()+1);
703 AMREX_ASSERT(nlevels <= level_steps.size());
704
705 int num_extra_mfs = 1; // for nodal, which is always on
707 num_extra_mfs += 3; // for nodal, which is always on
708 }
709
710 HeaderFile.precision(17);
711
712 // ---- this is the generic plot file type name
713 HeaderFile << versionName << '\n';
714
715 HeaderFile << varnames_3d.size() << '\n';
716
717 for (int ivar = 0; ivar < varnames_3d.size(); ++ivar) {
718 HeaderFile << varnames_3d[ivar] << "\n";
719 }
720 HeaderFile << AMREX_SPACEDIM << '\n';
721 HeaderFile << time << '\n';
722 HeaderFile << finest_level << '\n';
723 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
724 HeaderFile << my_geom[0].ProbLo(i) << ' ';
725 }
726 HeaderFile << '\n';
727 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
728 HeaderFile << my_geom[0].ProbHi(i) << ' ';
729 }
730 HeaderFile << '\n';
731 for (int i = 0; i < finest_level; ++i) {
732 HeaderFile << my_ref_ratio[i][0] << ' ';
733 }
734 HeaderFile << '\n';
735 for (int i = 0; i <= finest_level; ++i) {
736 HeaderFile << my_geom[i].Domain() << ' ';
737 }
738 HeaderFile << '\n';
739 for (int i = 0; i <= finest_level; ++i) {
740 HeaderFile << level_steps[i] << ' ';
741 }
742 HeaderFile << '\n';
743 for (int i = 0; i <= finest_level; ++i) {
744 for (int k = 0; k < AMREX_SPACEDIM; ++k) {
745 HeaderFile << my_geom[i].CellSize()[k] << ' ';
746 }
747 HeaderFile << '\n';
748 }
749 HeaderFile << (int) my_geom[0].Coord() << '\n';
750 HeaderFile << "0\n";
751
752 for (int level = 0; level <= finest_level; ++level) {
753 HeaderFile << level << ' ' << bArray[level].size() << ' ' << time << '\n';
754 HeaderFile << level_steps[level] << '\n';
755
756 const IntVect& domain_lo = my_geom[level].Domain().smallEnd();
757 for (int i = 0; i < bArray[level].size(); ++i)
758 {
759 // Need to shift because the RealBox ctor we call takes the
760 // physical location of index (0,0,0). This does not affect
761 // the usual cases where the domain index starts with 0.
762 const Box& b = shift(bArray[level][i], -domain_lo);
763 RealBox loc = RealBox(b, my_geom[level].CellSize(), my_geom[level].ProbLo());
764 for (int n = 0; n < AMREX_SPACEDIM; ++n) {
765 HeaderFile << loc.lo(n) << ' ' << loc.hi(n) << '\n';
766 }
767 }
768
769 HeaderFile << MultiFabHeaderPath(level, levelPrefix, mfPrefix) << '\n';
770 }
771 HeaderFile << num_extra_mfs << "\n";
772 HeaderFile << "3" << "\n";
773 HeaderFile << "amrexvec_nu_x" << "\n";
774 HeaderFile << "amrexvec_nu_y" << "\n";
775 HeaderFile << "amrexvec_nu_z" << "\n";
776 std::string mf_nodal_prefix = "Nu_nd";
777 for (int level = 0; level <= finest_level; ++level) {
778 HeaderFile << MultiFabHeaderPath(level, levelPrefix, mf_nodal_prefix) << '\n';
779 }
781 HeaderFile << "1" << "\n"; // number of components in the multifab
782 HeaderFile << "u_vel" << "\n";
783 std::string mf_uface_prefix = "UFace";
784 for (int level = 0; level <= finest_level; ++level) {
785 HeaderFile << MultiFabHeaderPath(level, levelPrefix, mf_uface_prefix) << '\n';
786 }
787 HeaderFile << "1" << "\n";
788 HeaderFile << "v_vel" << "\n";
789 std::string mf_vface_prefix = "VFace";
790 for (int level = 0; level <= finest_level; ++level) {
791 HeaderFile << MultiFabHeaderPath(level, levelPrefix, mf_vface_prefix) << '\n';
792 }
793 HeaderFile << "1" << "\n";
794 HeaderFile << "w_vel" << "\n";
795 std::string mf_wface_prefix = "WFace";
796 for (int level = 0; level <= finest_level; ++level) {
797 HeaderFile << MultiFabHeaderPath(level, levelPrefix, mf_wface_prefix) << '\n';
798 }
799 }
800
801 if (varnames_2d_rho.size() > 0) {
802 HeaderFile << varnames_2d_rho.size() << "\n"; // number of components in the 2D rho multifab
803 for (int ivar = 0; ivar < varnames_2d_rho.size(); ++ivar) {
804 HeaderFile << varnames_2d_rho[ivar] << "\n";
805 }
806 std::string mf_2d_rho_prefix = "rho2d";
807 for (int level = 0; level <= finest_level; ++level) {
808 HeaderFile << MultiFabHeaderPath(level, levelPrefix, mf_2d_rho_prefix) << "\n";
809 }
810 }
811
812 if (varnames_2d_u.size() > 0) {
813 HeaderFile << varnames_2d_u.size() << "\n"; // number of components in the 2D rho multifab
814 for (int ivar = 0; ivar < varnames_2d_u.size(); ++ivar) {
815 HeaderFile << varnames_2d_u[ivar] << "\n";
816 }
817 std::string mf_2d_u_prefix = "u2d";
818 for (int level = 0; level <= finest_level; ++level) {
819 HeaderFile << MultiFabHeaderPath(level, levelPrefix, mf_2d_u_prefix) << "\n";
820 }
821 }
822
823 if (varnames_2d_v.size() > 0) {
824 HeaderFile << varnames_2d_v.size() << "\n"; // number of components in the 2D v multifab
825 for (int ivar = 0; ivar < varnames_2d_v.size(); ++ivar) {
826 HeaderFile << varnames_2d_v[ivar] << "\n";
827 }
828 std::string mf_2d_v_prefix = "v2d";
829 for (int level = 0; level <= finest_level; ++level) {
830 HeaderFile << MultiFabHeaderPath(level, levelPrefix, mf_2d_v_prefix) << "\n";
831 }
832 }
833}
834
835/**
836 * @param lev level to mask
837 * @param fill_value fill value to mask with
838 * @param fill_where value at cells where we will apply the mask. This is necessary because rivers
839 */
840void
841REMORA::mask_arrays_for_write(int lev, Real fill_value, Real fill_where)
842{
843 for (MFIter mfi(*cons_new[lev],false); mfi.isValid(); ++mfi) {
844 Box gbx1 = mfi.growntilebox(IntVect(NGROW+1,NGROW+1,0));
845 Box ubx = mfi.grownnodaltilebox(0,IntVect(NGROW,NGROW,0));
846 Box vbx = mfi.grownnodaltilebox(1,IntVect(NGROW,NGROW,0));
847
848 Array4<Real> const& Zt_avg1 = vec_Zt_avg1[lev]->array(mfi);
849 Array4<Real> const& ubar = vec_ubar[lev]->array(mfi);
850 Array4<Real> const& vbar = vec_vbar[lev]->array(mfi);
851 Array4<Real> const& xvel = xvel_new[lev]->array(mfi);
852 Array4<Real> const& yvel = yvel_new[lev]->array(mfi);
853 Array4<Real> const& temp = cons_new[lev]->array(mfi,Temp_comp);
854 Array4<Real> const& salt = cons_new[lev]->array(mfi,Salt_comp);
855
856 Array4<Real const> const& mskr = vec_mskr[lev]->array(mfi);
857 Array4<Real const> const& msku = vec_msku[lev]->array(mfi);
858 Array4<Real const> const& mskv = vec_mskv[lev]->array(mfi);
859
860 ParallelFor(makeSlab(gbx1,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int )
861 {
862 if (mskr(i,j,0) == 0.0) { // Explicitly compare to 0.0
863 Zt_avg1(i,j,0) = fill_value;
864 }
865 });
866 ParallelFor(gbx1, [=] AMREX_GPU_DEVICE (int i, int j, int k)
867 {
868 if (mskr(i,j,0) == 0.0) { // Explicitly compare to 0.0
869 temp(i,j,k) = fill_value;
870 salt(i,j,k) = fill_value;
871 }
872 });
873 ParallelFor(makeSlab(ubx,2,0), 3, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
874 {
875 if (msku(i,j,0) == 0.0 && ubar(i,j,0)==fill_where) { // Explicitly compare to 0.0
876 ubar(i,j,0,n) = fill_value;
877 }
878 });
879 ParallelFor(makeSlab(vbx,2,0), 3, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
880 {
881 if (mskv(i,j,0) == 0.0 && vbar(i,j,0)==fill_where) { // Explicitly compare to 0.0
882 vbar(i,j,0,n) = fill_value;
883 }
884 });
885 ParallelFor(ubx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
886 {
887 if (msku(i,j,0) == 0.0 && xvel(i,j,k)==fill_where) { // Explicitly compare to 0.0
888 xvel(i,j,k) = fill_value;
889 }
890 });
891 ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
892 {
893 if (mskv(i,j,0) == 0.0 && yvel(i,j,k)==fill_where) { // Explicitly compare to 0.0
894 yvel(i,j,k) = fill_value;
895 }
896 });
897 } // mfi
898 Gpu::streamSynchronize();
899}
Coord
Coordinates.
#define NGROW
#define Temp_comp
#define Salt_comp
#define NCONS
bool containerHasElement(const V &iterable, const T &query)
PhysBCFunctNoOp null_bc_for_fill
static PlotfileType plotfile_type
Native or NetCDF plotfile output.
Definition REMORA.H:1386
const amrex::Vector< std::string > cons_names
Names of scalars for plotfile output.
Definition REMORA.H:1330
amrex::Vector< amrex::BCRec > domain_bcs_type
vector (over BCVars) of BCRecs
Definition REMORA.H:1233
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:253
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:400
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:241
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:389
amrex::Real plotfile_fill_value
fill value for masked arrays in amrex plotfiles
Definition REMORA.H:1342
void writeJobInfo(const std::string &dir) const
Write job info to stdout.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:311
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:247
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_msku
land/sea mask at x-faces (2D)
Definition REMORA.H:391
void FillPatchNoBC(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true)
Fill a new MultiFab by copying in phi from valid region and filling ghost cells without applying boun...
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:245
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:243
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskv
land/sea mask at y-faces (2D)
Definition REMORA.H:393
amrex::Vector< int > istep
which step?
Definition REMORA.H:1207
void WriteMultiLevelPlotfileWithBathymetry(const std::string &plotfilename, int nlevels, const amrex::Vector< const amrex::MultiFab * > &mf, const amrex::Vector< const amrex::MultiFab * > &mf_nd, const amrex::Vector< const amrex::MultiFab * > &mf_u, const amrex::Vector< const amrex::MultiFab * > &mf_v, const amrex::Vector< const amrex::MultiFab * > &mf_w, const amrex::Vector< const amrex::MultiFab * > &mf_2d_rho, const amrex::Vector< const amrex::MultiFab * > &mf_2d_u, const amrex::Vector< const amrex::MultiFab * > &mf_2d_v, const amrex::Vector< std::string > &varnames_3d, const amrex::Vector< std::string > &varnames_2d_rho, const amrex::Vector< std::string > &varnames_2d_u, const amrex::Vector< std::string > &varnames_2d_v, const amrex::Vector< amrex::Geometry > &my_geom, amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &rr, const std::string &versionName="HyperCLaw-V1.1", const std::string &levelPrefix="Level_", const std::string &mfPrefix="Cell", const amrex::Vector< std::string > &extra_dirs=amrex::Vector< std::string >()) const
write out particular data to an AMReX plotfile
void mask_arrays_for_write(int lev, amrex::Real fill_value, amrex::Real fill_where)
Mask data arrays before writing output.
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
Definition REMORA.H:1380
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:313
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg1
time average of barotropic y velocity flux
Definition REMORA.H:372
void WriteNCPlotFile(int istep, amrex::MultiFab const *plotMF)
Write plotfile using NetCDF (wrapper)
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1211
void WriteGenericPlotfileHeaderWithBathymetry(std::ostream &HeaderFile, int nlevels, const amrex::Vector< amrex::BoxArray > &bArray, const amrex::Vector< std::string > &varnames_3d, const amrex::Vector< std::string > &varnames_2d_rho, const amrex::Vector< std::string > &varnames_2d_u, const amrex::Vector< std::string > &varnames_2d_v, const amrex::Vector< amrex::Geometry > &my_geom, amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &rr, const std::string &versionName, const std::string &levelPrefix, const std::string &mfPrefix) const
write out header data for an AMReX plotfile
static int plot_staggered_vels
Whether to write the staggered velocities (not averaged to cell centers)
Definition REMORA.H:1383
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:384
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg1
time average of barotropic x velocity flux (2D)
Definition REMORA.H:368
bool expand_plotvars_to_unif_rr
whether plotfile variables should be expanded to a uniform refinement ratio
Definition REMORA.H:1339
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:382
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bustr
Bottom stress in the u direction.
Definition REMORA.H:363
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvstr
Bottom stress in the v direction.
Definition REMORA.H:365
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_phys_nd
z coordinates at psi points (cell nodes)
Definition REMORA.H:305
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:402
amrex::Vector< std::string > plot_var_names_3d
Names of 3D variables to output to AMReX plotfile.
Definition REMORA.H:1326
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1303
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:308
void WritePlotFile(int istep)
main driver for writing AMReX plotfiles
amrex::Vector< std::string > plot_var_names_2d
Names of 2D variables to output to AMReX plotfile.
Definition REMORA.H:1328
void remora_dervort(const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Array4< const amrex::Real > &pm, const amrex::Array4< const amrex::Real > &pn, const amrex::Array4< const amrex::Real > &, const amrex::Geometry &, amrex::Real, const int *, const int)
void remora_dernull(const amrex::Box &, amrex::FArrayBox &, int, int, const amrex::FArrayBox &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Array4< const amrex::Real > &, const amrex::Geometry &, amrex::Real, const int *, const int)