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