REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_Checkpoint.cpp
Go to the documentation of this file.
1#include <REMORA.H>
2#include "AMReX_PlotFileUtil.H"
3
4using namespace amrex;
5
6// utility to skip to next line in Header
7void
8REMORA::GotoNextLine (std::istream& is)
9{
10 constexpr std::streamsize bl_ignore_max { 100000 };
11 is.ignore(bl_ignore_max, '\n');
12}
13
14void
16{
17
18 // chk00010 write a checkpoint file with this root directory
19 // chk00010/Header this contains information you need to save (e.g., finest_level, t_new, etc.) and also
20 // the BoxArrays at each level
21 // chk00010/Level_0/
22 // chk00010/Level_1/
23 // etc. these subdirectories will hold the MultiFab data at each level of refinement
24
25 // checkpoint file name, e.g., chk00010
26 const std::string& checkpointname = amrex::Concatenate(check_file,istep[0],file_min_digits);
27
28 amrex::Print() << "Writing checkpoint " << checkpointname << "\n";
29
30 const int nlevels = finest_level+1;
31
32 // ---- prebuild a hierarchy of directories
33 // ---- dirName is built first. if dirName exists, it is renamed. then build
34 // ---- dirName/subDirPrefix_0 .. dirName/subDirPrefix_nlevels-1
35 // ---- if callBarrier is true, call ParallelDescriptor::Barrier()
36 // ---- after all directories are built
37 // ---- ParallelDescriptor::IOProcessor() creates the directories
38 amrex::PreBuildDirectorHierarchy(checkpointname, "Level_", nlevels, true);
39
40 // write Header file
41 if (ParallelDescriptor::IOProcessor()) {
42
43 std::string HeaderFileName(checkpointname + "/Header");
44 VisMF::IO_Buffer io_buffer(VisMF::IO_Buffer_Size);
45 std::ofstream HeaderFile;
46 HeaderFile.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
47 HeaderFile.open(HeaderFileName.c_str(), std::ofstream::out |
48 std::ofstream::trunc |
49 std::ofstream::binary);
50 if( ! HeaderFile.good()) {
51 amrex::FileOpenFailed(HeaderFileName);
52 }
53
54 HeaderFile.precision(17);
55
56 // write out title line
57 HeaderFile << "Checkpoint file for REMORA\n";
58
59 // write out finest_level
60 HeaderFile << finest_level << "\n";
61
62 // write the number of components
63 // for each variable we store
64
65 // conservative, cell-centered vars
66 HeaderFile << NCONS << "\n";
67
68 // x-velocity on faces
69 HeaderFile << 1 << "\n";
70
71 // y-velocity on faces
72 HeaderFile << 1 << "\n";
73
74 // z-velocity on faces
75 HeaderFile << 1 << "\n";
76
77 HeaderFile << 2 << "\n";
78
79 HeaderFile << 2 << "\n";
80
81 HeaderFile << 3 << "\n";
82
83 HeaderFile << 3 << "\n";
84
85 // write out array of istep
86 for (int i = 0; i < istep.size(); ++i) {
87 HeaderFile << istep[i] << " ";
88 }
89 HeaderFile << "\n";
90
91 // write out array of dt
92 for (int i = 0; i < dt.size(); ++i) {
93 HeaderFile << dt[i] << " ";
94 }
95 HeaderFile << "\n";
96
97 // write out array of t_new
98 for (int i = 0; i < t_new.size(); ++i) {
99 HeaderFile << t_new[i] << " ";
100 }
101 HeaderFile << "\n";
102
103 // write the BoxArray at each level
104 for (int lev = 0; lev <= finest_level; ++lev) {
105 boxArray(lev).writeOn(HeaderFile);
106 HeaderFile << '\n';
107 }
108 }
109
110 // write the MultiFab data to, e.g., chk00010/Level_0/
111 // Here we make copies of the MultiFab with no ghost cells
112 for (int lev = 0; lev <= finest_level; ++lev)
113 {
114 BoxList bl2d = grids[lev].boxList();
115 for (auto& b : bl2d) {
116 b.setRange(2,0);
117 }
118 BoxArray ba2d(std::move(bl2d));
119
120 MultiFab cons(grids[lev],dmap[lev],NCONS,cons_new[lev]->nGrowVect());
121 MultiFab::Copy(cons,*cons_new[lev],0,0,NCONS,cons_new[lev]->nGrowVect());
122 VisMF::Write(cons, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Cell"));
123
124 MultiFab::Copy(cons,*cons_old[lev],0,0,NCONS,cons_old[lev]->nGrowVect());
125 VisMF::Write(cons, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Cell_old"));
126
127 MultiFab xvel(convert(grids[lev],IntVect(1,0,0)),dmap[lev],1,xvel_new[lev]->nGrowVect());
128 MultiFab::Copy(xvel,*xvel_new[lev],0,0,1,xvel_new[lev]->nGrowVect());
129 VisMF::Write(xvel, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "XFace"));
130
131 MultiFab::Copy(xvel,*xvel_old[lev],0,0,1,xvel_old[lev]->nGrowVect());
132 VisMF::Write(xvel, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "XFace_old"));
133
134 MultiFab yvel(convert(grids[lev],IntVect(0,1,0)),dmap[lev],1,yvel_new[lev]->nGrowVect());
135 MultiFab::Copy(yvel,*yvel_new[lev],0,0,1,yvel_new[lev]->nGrowVect());
136 VisMF::Write(yvel, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "YFace"));
137
138 MultiFab::Copy(yvel,*yvel_old[lev],0,0,1,yvel_old[lev]->nGrowVect());
139 VisMF::Write(yvel, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "YFace_old"));
140
141 MultiFab zvel(convert(grids[lev],IntVect(0,0,1)),dmap[lev],1,zvel_new[lev]->nGrowVect());
142 MultiFab::Copy(zvel,*zvel_new[lev],0,0,1,zvel_new[lev]->nGrowVect());
143 VisMF::Write(zvel, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "ZFace"));
144
145 MultiFab::Copy(zvel,*zvel_old[lev],0,0,1,zvel_old[lev]->nGrowVect());
146 VisMF::Write(zvel, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "ZFace_old"));
147
148 MultiFab mf_ru(convert(grids[lev],IntVect(1,0,0)),dmap[lev],2,(vec_ru[lev])->nGrowVect());
149 MultiFab::Copy(mf_ru,*vec_ru[lev],0,0,2,(vec_ru[lev])->nGrowVect());
150 VisMF::Write(mf_ru, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "XRHS"));
151
152 MultiFab mf_rv(convert(grids[lev],IntVect(0,1,0)),dmap[lev],2,(vec_rv[lev])->nGrowVect());
153 MultiFab::Copy(mf_rv,*vec_rv[lev],0,0,2,(vec_rv[lev])->nGrowVect());
154 VisMF::Write(mf_rv, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "YRHS"));
155
156 MultiFab mf_ubar(convert(ba2d,IntVect(1,0,0)),dmap[lev],3,(vec_ubar[lev])->nGrowVect());
157 MultiFab::Copy(mf_ubar,*(vec_ubar[lev]),0,0,3,(vec_ubar[lev])->nGrowVect());
158 VisMF::Write(mf_ubar, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "XBar"));
159
160 MultiFab mf_vbar(convert(ba2d,IntVect(0,1,0)),dmap[lev],3,(vec_vbar[lev])->nGrowVect());
161 MultiFab::Copy(mf_vbar,*(vec_vbar[lev]),0,0,3,(vec_vbar[lev])->nGrowVect());
162 VisMF::Write(mf_vbar, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "YBar"));
163
164 MultiFab mf_ru2d(convert(ba2d,IntVect(1,0,0)),dmap[lev],2,(vec_ru2d[lev])->nGrowVect());
165 MultiFab::Copy(mf_ru2d,*(vec_ru2d[lev]),0,0,2,(vec_ru2d[lev])->nGrowVect());
166 VisMF::Write(mf_ru2d, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "XRHS2d"));
167
168 MultiFab mf_rv2d(convert(ba2d,IntVect(0,1,0)),dmap[lev],2,(vec_rv2d[lev])->nGrowVect());
169 MultiFab::Copy(mf_rv2d,*(vec_rv2d[lev]),0,0,2,(vec_rv2d[lev])->nGrowVect());
170 VisMF::Write(mf_rv2d, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "YRHS2d"));
171
172 MultiFab mf_mskr(ba2d,dmap[lev],1,(vec_mskr[lev])->nGrowVect());
173 MultiFab::Copy(mf_mskr,*(vec_mskr[lev]),0,0,1,(vec_mskr[lev])->nGrowVect());
174 VisMF::Write(mf_mskr, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Mskr"));
175
176 MultiFab mf_msku(convert(ba2d,IntVect(1,0,0)),dmap[lev],1,(vec_msku[lev])->nGrowVect());
177 MultiFab::Copy(mf_msku,*(vec_msku[lev]),0,0,1,(vec_msku[lev])->nGrowVect());
178 VisMF::Write(mf_msku, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Msku"));
179
180 MultiFab mf_mskv(convert(ba2d,IntVect(0,1,0)),dmap[lev],1,(vec_mskv[lev])->nGrowVect());
181 MultiFab::Copy(mf_mskv,*(vec_mskv[lev]),0,0,1,(vec_mskv[lev])->nGrowVect());
182 VisMF::Write(mf_mskv, amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Mskv"));
183
185 VisMF::Write(*(vec_rdrag[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "rdrag"));
187 VisMF::Write(*(vec_rdrag2[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "rdrag2"));
188 }
189
192 VisMF::Write(*(vec_ZoBot[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "ZoBot"));
193 }
194
195 VisMF::Write(*(vec_DU_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "DU_avg1"));
196 VisMF::Write(*(vec_DU_avg2[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "DU_avg2"));
197 VisMF::Write(*(vec_DV_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "DV_avg1"));
198 VisMF::Write(*(vec_DV_avg2[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "DV_avg2"));
199
200 VisMF::Write(*(vec_zeta[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "zeta"));
201 VisMF::Write(*(vec_Zt_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Zt_avg1"));
202
203 VisMF::Write(*(vec_hOfTheConfusingName[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "h"));
204
205 VisMF::Write(*(vec_tke[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "tke"));
206 VisMF::Write(*(vec_gls[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "gls"));
207 VisMF::Write(*(vec_Lscale[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Lscale"));
208 VisMF::Write(*(vec_Akk[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Akk"));
209 VisMF::Write(*(vec_Akp[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Akp"));
210 VisMF::Write(*(vec_Akv[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Akv"));
211 VisMF::Write(*(vec_Akt[lev]), amrex::MultiFabFileFullPrefix(lev, checkpointname, "Level_", "Akt"));
212 }
213
214#ifdef REMORA_USE_PARTICLES
215 particleData.Checkpoint(checkpointname);
216#endif
217
218#ifdef REMORA_USE_NETCDF
219 // Write bdy_data files
220 if ( ParallelDescriptor::IOProcessor() && (solverChoice.ic_bc_type == IC_BC_Type::netcdf) )
221 {
222
223 // Vector dimensions
224 int num_time = bdy_data_xlo.size();
225 int num_var = bdy_data_xlo[0].size();
226
227 // Open header file and write to it
228 std::ofstream bdy_h_file(amrex::MultiFabFileFullPrefix(0, checkpointname, "Level_", "bdy_H"));
229 bdy_h_file << std::setprecision(1) << std::fixed;
230 bdy_h_file << num_time << "\n";
231 bdy_h_file << num_var << "\n";
232 bdy_h_file << start_bdy_time << "\n";
233 bdy_h_file << bdy_time_interval << "\n";
234 for (int ivar(0); ivar<num_var; ++ivar) {
235 bdy_h_file << bdy_data_xlo[0][ivar].box() << "\n";
236 bdy_h_file << bdy_data_xhi[0][ivar].box() << "\n";
237 bdy_h_file << bdy_data_ylo[0][ivar].box() << "\n";
238 bdy_h_file << bdy_data_yhi[0][ivar].box() << "\n";
239 }
240
241 // Open data file and write to it
242 std::ofstream bdy_d_file(amrex::MultiFabFileFullPrefix(0, checkpointname, "Level_", "bdy_D"));
243 for (int itime(0); itime<num_time; ++itime) {
244 for (int ivar(0); ivar<num_var; ++ivar) {
245 bdy_data_xlo[itime][ivar].writeOn(bdy_d_file,0,1);
246 bdy_data_xhi[itime][ivar].writeOn(bdy_d_file,0,1);
247 bdy_data_ylo[itime][ivar].writeOn(bdy_d_file,0,1);
248 bdy_data_yhi[itime][ivar].writeOn(bdy_d_file,0,1);
249 }
250 }
251 }
252#endif
253}
254
255void
257{
258 amrex::Print() << "Restart from checkpoint " << restart_chkfile << "\n";
259
260 // Header
261 std::string File(restart_chkfile + "/Header");
262
263 VisMF::IO_Buffer io_buffer(VisMF::GetIOBufferSize());
264
265 Vector<char> fileCharPtr;
266 ParallelDescriptor::ReadAndBcastFile(File, fileCharPtr);
267 std::string fileCharPtrString(fileCharPtr.dataPtr());
268 std::istringstream is(fileCharPtrString, std::istringstream::in);
269
270 std::string line, word;
271
272 int chk_ncomp;
273
274 // read in title line
275 std::getline(is, line);
276
277 // read in finest_level
278 is >> finest_level;
279 GotoNextLine(is);
280
281 // read the number of components
282 // for each variable we store
283
284 // conservative, cell-centered vars
285 is >> chk_ncomp;
286 GotoNextLine(is);
287 AMREX_ASSERT(chk_ncomp == NCONS);
288
289 // x-velocity on faces
290 is >> chk_ncomp;
291 GotoNextLine(is);
292 AMREX_ASSERT(chk_ncomp == 1);
293
294 // y-velocity on faces
295 is >> chk_ncomp;
296 GotoNextLine(is);
297 AMREX_ASSERT(chk_ncomp == 1);
298
299 // z-velocity on faces
300 is >> chk_ncomp;
301 GotoNextLine(is);
302 AMREX_ASSERT(chk_ncomp == 1);
303
304 is >> chk_ncomp;
305 GotoNextLine(is);
306 AMREX_ASSERT(chk_ncomp == 2);
307
308 is >> chk_ncomp;
309 GotoNextLine(is);
310 AMREX_ASSERT(chk_ncomp == 2);
311
312 is >> chk_ncomp;
313 GotoNextLine(is);
314 AMREX_ASSERT(chk_ncomp == 3);
315
316 is >> chk_ncomp;
317 GotoNextLine(is);
318 AMREX_ASSERT(chk_ncomp == 3);
319
320 // read in array of istep
321 std::getline(is, line);
322 {
323 std::istringstream lis(line);
324 int i = 0;
325 while (lis >> word) {
326 istep[i++] = std::stoi(word);
327 }
328 }
329
330 // read in array of dt
331 std::getline(is, line);
332 {
333 std::istringstream lis(line);
334 int i = 0;
335 while (lis >> word) {
336#ifdef AMREX_USE_FLOAT
337 dt[i++] = std::stof(word);
338#else
339 dt[i++] = std::stod(word);
340#endif
341 }
342 }
343
344 // read in array of t_new
345 std::getline(is, line);
346 {
347 std::istringstream lis(line);
348 int i = 0;
349 while (lis >> word) {
350#ifdef AMREX_USE_FLOAT
351 t_new[i++] = std::stof(word);
352#else
353 t_new[i++] = std::stod(word);
354#endif
355 }
356 }
357
358 for (int lev = 0; lev <= finest_level; ++lev) {
359
360 // read in level 'lev' BoxArray from Header
361 BoxArray ba;
362 ba.readFrom(is);
363 GotoNextLine(is);
364
365 // create a distribution mapping
366 DistributionMapping dm { ba, ParallelDescriptor::NProcs() };
367
368 MakeNewLevelFromScratch (lev, t_new[lev], ba, dm);
369 }
370
371 // read in the MultiFab data
372 for (int lev = 0; lev <= finest_level; ++lev)
373 {
374 BoxList bl2d = grids[lev].boxList();
375 for (auto& b : bl2d) {
376 b.setRange(2,0);
377 }
378 BoxArray ba2d(std::move(bl2d));
379
380 MultiFab cons(grids[lev],dmap[lev],NCONS,cons_new[lev]->nGrowVect());
381 VisMF::Read(cons, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Cell"));
382 MultiFab::Copy(*cons_new[lev],cons,0,0,NCONS,cons_new[lev]->nGrowVect());
383
384 VisMF::Read(cons, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Cell_old"));
385 MultiFab::Copy(*cons_old[lev],cons,0,0,NCONS,cons_old[lev]->nGrowVect());
386
387 MultiFab xvel(convert(grids[lev],IntVect(1,0,0)),dmap[lev],1,xvel_new[lev]->nGrowVect());
388 VisMF::Read(xvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XFace"));
389 MultiFab::Copy(*xvel_new[lev],xvel,0,0,1,xvel_new[lev]->nGrowVect());
390
391 VisMF::Read(xvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XFace_old"));
392 MultiFab::Copy(*xvel_old[lev],xvel,0,0,1,xvel_old[lev]->nGrowVect());
393
394 MultiFab yvel(convert(grids[lev],IntVect(0,1,0)),dmap[lev],1,yvel_new[lev]->nGrowVect());
395 VisMF::Read(yvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YFace"));
396 MultiFab::Copy(*yvel_new[lev],yvel,0,0,1,yvel_new[lev]->nGrowVect());
397
398 VisMF::Read(yvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YFace_old"));
399 MultiFab::Copy(*yvel_old[lev],yvel,0,0,1,yvel_old[lev]->nGrowVect());
400
401 MultiFab zvel(convert(grids[lev],IntVect(0,0,1)),dmap[lev],1,zvel_new[lev]->nGrowVect());
402 VisMF::Read(zvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "ZFace"));
403 MultiFab::Copy(*zvel_new[lev],zvel,0,0,1,zvel_new[lev]->nGrowVect());
404
405 VisMF::Read(zvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "ZFace_old"));
406 MultiFab::Copy(*zvel_old[lev],zvel,0,0,1,zvel_old[lev]->nGrowVect());
407
408 MultiFab mf_ru(convert(grids[lev],IntVect(1,0,0)),dmap[lev],2,(vec_ru[lev])->nGrowVect());
409 VisMF::Read(mf_ru, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XRHS"));
410 MultiFab::Copy(*(vec_ru[lev]),mf_ru,0,0,2,(vec_ru[lev])->nGrowVect());
411
412 MultiFab mf_rv(convert(grids[lev],IntVect(0,1,0)),dmap[lev],2,(vec_rv[lev])->nGrowVect());
413 VisMF::Read(mf_rv, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YRHS"));
414 MultiFab::Copy(*(vec_rv[lev]),mf_rv,0,0,2,(vec_rv[lev])->nGrowVect());
415
416 MultiFab mf_ubar(convert(ba2d,IntVect(1,0,0)),dmap[lev],3,(vec_ubar[lev])->nGrowVect());
417 VisMF::Read(mf_ubar, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XBar"));
418 MultiFab::Copy(*(vec_ubar[lev]),mf_ubar,0,0,3,(vec_ubar[lev])->nGrowVect());
419
420 MultiFab mf_vbar(convert(ba2d,IntVect(0,1,0)),dmap[lev],3,(vec_vbar[lev])->nGrowVect());
421 VisMF::Read(mf_vbar, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YBar"));
422 MultiFab::Copy(*(vec_vbar[lev]),mf_vbar,0,0,3,(vec_vbar[lev])->nGrowVect());
423
424 MultiFab mf_ru2d(convert(ba2d,IntVect(1,0,0)),dmap[lev],2,(vec_ru2d[lev])->nGrowVect());
425 VisMF::Read(mf_ru2d, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XRHS2d"));
426 MultiFab::Copy(*(vec_ru2d[lev]),mf_ru2d,0,0,2,(vec_ru2d[lev])->nGrowVect());
427
428 MultiFab mf_rv2d(convert(ba2d,IntVect(0,1,0)),dmap[lev],2,(vec_rv2d[lev])->nGrowVect());
429 VisMF::Read(mf_rv2d, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YRHS2d"));
430 MultiFab::Copy(*(vec_rv2d[lev]),mf_rv2d,0,0,2,(vec_rv2d[lev])->nGrowVect());
431
432 MultiFab mf_mskr(ba2d,dmap[lev],1,(vec_mskr[lev])->nGrowVect());
433 VisMF::Read(mf_mskr, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Mskr"));
434 MultiFab::Copy(*(vec_mskr[lev]),mf_mskr,0,0,1,(vec_mskr[lev])->nGrowVect());
435
436 MultiFab mf_msku(convert(ba2d,IntVect(1,0,0)),dmap[lev],1,(vec_msku[lev])->nGrowVect());
437 VisMF::Read(mf_msku, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Msku"));
438 MultiFab::Copy(*(vec_msku[lev]),mf_msku,0,0,1,(vec_msku[lev])->nGrowVect());
439
440 MultiFab mf_mskv(convert(ba2d,IntVect(0,1,0)),dmap[lev],1,(vec_mskv[lev])->nGrowVect());
441 VisMF::Read(mf_mskv, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Mskv"));
442 MultiFab::Copy(*(vec_mskv[lev]),mf_mskv,0,0,1,(vec_mskv[lev])->nGrowVect());
443
445 VisMF::Read(*(vec_rdrag[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "rdrag"));
447 VisMF::Read(*(vec_rdrag2[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "rdrag2"));
448 }
449
452 VisMF::Read(*(vec_ZoBot[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "ZoBot"));
453 }
454
455 VisMF::Read(*(vec_DU_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "DU_avg1"));
456 VisMF::Read(*(vec_DU_avg2[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "DU_avg2"));
457 VisMF::Read(*(vec_DV_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "DV_avg1"));
458 VisMF::Read(*(vec_DV_avg2[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "DV_avg2"));
459
460 VisMF::Read(*(vec_zeta[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "zeta"));
461 VisMF::Read(*(vec_Zt_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Zt_avg1"));
462
463 VisMF::Read(*(vec_hOfTheConfusingName[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "h"));
464
465 VisMF::Read(*(vec_tke[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "tke"));
466 VisMF::Read(*(vec_gls[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "gls"));
467 VisMF::Read(*(vec_Lscale[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Lscale"));
468 VisMF::Read(*(vec_Akk[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Akk"));
469 VisMF::Read(*(vec_Akp[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Akp"));
470 VisMF::Read(*(vec_Akt[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Akt"));
471 VisMF::Read(*(vec_Akv[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Akv"));
472
474
475 }
476
477#ifdef REMORA_USE_PARTICLES
478 particleData.Restart((amrex::ParGDBBase*)GetParGDB(),restart_chkfile);
479#endif
480
481#ifdef REMORA_USE_NETCDF
482 // Read bdy_data files
484 {
485 int ioproc = ParallelDescriptor::IOProcessorNumber(); // I/O rank
486 int num_time;
487 int num_var;
488 Vector<Box> bx_v;
489 if (ParallelDescriptor::IOProcessor()) {
490 // Open header file and read from it
491 std::ifstream bdy_h_file(amrex::MultiFabFileFullPrefix(0, restart_chkfile, "Level_", "bdy_H"));
492 bdy_h_file >> num_time;
493 bdy_h_file >> num_var;
494 bdy_h_file >> start_bdy_time;
495 bdy_h_file >> bdy_time_interval;
496 bx_v.resize(4*num_var);
497 for (int ivar(0); ivar<num_var; ++ivar) {
498 bdy_h_file >> bx_v[4*ivar ];
499 bdy_h_file >> bx_v[4*ivar+1];
500 bdy_h_file >> bx_v[4*ivar+2];
501 bdy_h_file >> bx_v[4*ivar+3];
502 }
503
504 // IO size the FABs
505 bdy_data_xlo.resize(num_time);
506 bdy_data_xhi.resize(num_time);
507 bdy_data_ylo.resize(num_time);
508 bdy_data_yhi.resize(num_time);
509 for (int itime(0); itime<num_time; ++itime) {
510 bdy_data_xlo[itime].resize(num_var);
511 bdy_data_xhi[itime].resize(num_var);
512 bdy_data_ylo[itime].resize(num_var);
513 bdy_data_yhi[itime].resize(num_var);
514 for (int ivar(0); ivar<num_var; ++ivar) {
515 bdy_data_xlo[itime][ivar].resize(bx_v[4*ivar ]);
516 bdy_data_xhi[itime][ivar].resize(bx_v[4*ivar+1]);
517 bdy_data_ylo[itime][ivar].resize(bx_v[4*ivar+2]);
518 bdy_data_yhi[itime][ivar].resize(bx_v[4*ivar+3]);
519 }
520 }
521
522 // Open data file and read from it
523 std::ifstream bdy_d_file(amrex::MultiFabFileFullPrefix(0, restart_chkfile, "Level_", "bdy_D"));
524 for (int itime(0); itime<num_time; ++itime) {
525 for (int ivar(0); ivar<num_var; ++ivar) {
526 bdy_data_xlo[itime][ivar].readFrom(bdy_d_file);
527 bdy_data_xhi[itime][ivar].readFrom(bdy_d_file);
528 bdy_data_ylo[itime][ivar].readFrom(bdy_d_file);
529 bdy_data_yhi[itime][ivar].readFrom(bdy_d_file);
530 }
531 }
532 } // IO
533
534 // Broadcast the data
535 ParallelDescriptor::Barrier();
536 ParallelDescriptor::Bcast(&start_bdy_time,1,ioproc);
537 ParallelDescriptor::Bcast(&bdy_time_interval,1,ioproc);
538 ParallelDescriptor::Bcast(&num_time,1,ioproc);
539 ParallelDescriptor::Bcast(&num_var,1,ioproc);
540
541 // Everyone size their boxes
542 bx_v.resize(4*num_var);
543
544 ParallelDescriptor::Bcast(bx_v.dataPtr(),bx_v.size(),ioproc);
545
546 // Everyone but IO size their FABs
547 if (!ParallelDescriptor::IOProcessor()) {
548 bdy_data_xlo.resize(num_time);
549 bdy_data_xhi.resize(num_time);
550 bdy_data_ylo.resize(num_time);
551 bdy_data_yhi.resize(num_time);
552 for (int itime(0); itime<num_time; ++itime) {
553 bdy_data_xlo[itime].resize(num_var);
554 bdy_data_xhi[itime].resize(num_var);
555 bdy_data_ylo[itime].resize(num_var);
556 bdy_data_yhi[itime].resize(num_var);
557 for (int ivar(0); ivar<num_var; ++ivar) {
558 bdy_data_xlo[itime][ivar].resize(bx_v[4*ivar ]);
559 bdy_data_xhi[itime][ivar].resize(bx_v[4*ivar+1]);
560 bdy_data_ylo[itime][ivar].resize(bx_v[4*ivar+2]);
561 bdy_data_yhi[itime][ivar].resize(bx_v[4*ivar+3]);
562 }
563 }
564 }
565
566 for (int itime(0); itime<num_time; ++itime) {
567 for (int ivar(0); ivar<num_var; ++ivar) {
568 ParallelDescriptor::Bcast(bdy_data_xlo[itime][ivar].dataPtr(),bdy_data_xlo[itime][ivar].box().numPts(),ioproc);
569 ParallelDescriptor::Bcast(bdy_data_xhi[itime][ivar].dataPtr(),bdy_data_xhi[itime][ivar].box().numPts(),ioproc);
570 ParallelDescriptor::Bcast(bdy_data_ylo[itime][ivar].dataPtr(),bdy_data_ylo[itime][ivar].box().numPts(),ioproc);
571 ParallelDescriptor::Bcast(bdy_data_yhi[itime][ivar].dataPtr(),bdy_data_yhi[itime][ivar].box().numPts(),ioproc);
572 }
573 }
574 } // init real
575#endif
576}
#define NCONS
static void GotoNextLine(std::istream &is)
utility to skip to next line in Header
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:244
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_hOfTheConfusingName
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:229
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ZoBot
Bottom roughness length [m], defined at rho points.
Definition REMORA.H:317
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg2
correct time average of barotropic x velocity flux for coupling (2D)
Definition REMORA.H:327
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:217
void stretch_transform(int lev)
Calculate vertical stretched coordinates.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:346
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_tke
Turbulent kinetic energy.
Definition REMORA.H:403
void ReadCheckpointFile()
read checkpoint file from disk
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_gls
Turbulent generic length scale.
Definition REMORA.H:405
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:242
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:223
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Lscale
Vertical mixing turbulent length scale.
Definition REMORA.H:407
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
Definition REMORA.H:252
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_msku
land/sea mask at x-faces (2D)
Definition REMORA.H:348
std::string check_file
Checkpoint file prefix.
Definition REMORA.H:1196
amrex::Vector< amrex::Vector< amrex::FArrayBox > > bdy_data_yhi
Vectors (over time) of Vector (over variables) of FArrayBoxs for holding Northern boundary data from ...
Definition REMORA.H:996
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:210
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:221
amrex::Real bdy_time_interval
Interval between boundary data times.
Definition REMORA.H:1001
amrex::Real start_bdy_time
Start time in the time series of boundary data.
Definition REMORA.H:999
amrex::Vector< amrex::MultiFab * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
Definition REMORA.H:214
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:219
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskv
land/sea mask at y-faces (2D)
Definition REMORA.H:350
amrex::Vector< amrex::Vector< amrex::FArrayBox > > bdy_data_ylo
Vectors (over time) of Vector (over variables) of FArrayBoxs for holding Southern boundary data from ...
Definition REMORA.H:994
amrex::Vector< int > istep
which step?
Definition REMORA.H:1094
void WriteCheckpointFile()
write checkpoint file to disk
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
Definition REMORA.H:1256
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:212
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg1
time average of barotropic y velocity flux
Definition REMORA.H:329
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1098
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1221
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akk
Turbulent kinetic energy vertical diffusion coefficient.
Definition REMORA.H:409
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag2
Quadratic drag coefficient [unitless], defined at rho points.
Definition REMORA.H:315
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:238
amrex::Vector< amrex::Vector< amrex::FArrayBox > > bdy_data_xhi
Vectors (over time) of Vector (over variables) of FArrayBoxs for holding Eastern boundary data from f...
Definition REMORA.H:992
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
Definition REMORA.H:343
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:341
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg1
time average of barotropic x velocity flux (2D)
Definition REMORA.H:325
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:339
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:208
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg2
correct time average of barotropic y velocity flux for coupling (2D)
Definition REMORA.H:331
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
Definition REMORA.H:250
virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Make a new level from scratch using provided BoxArray and DistributionMapping. Only used during initi...
amrex::Vector< amrex::Vector< amrex::FArrayBox > > bdy_data_xlo
Vectors (over time) of Vector (over variables) of FArrayBoxs for holding Western boundary data from f...
Definition REMORA.H:990
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag
Linear drag coefficient [m/s], defined at rho points.
Definition REMORA.H:313
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:279
std::string restart_chkfile
If set, restart from this checkpoint file.
Definition REMORA.H:1166
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:240
amrex::Vector< amrex::Real > dt
time step at each level
Definition REMORA.H:1102
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akp
Turbulent length scale vertical diffusion coefficient.
Definition REMORA.H:411
IC_BC_Type ic_bc_type
BottomStressType bottom_stress_type
VertMixingType vert_mixing_type