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_h[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
219void
221{
222 amrex::Print() << "Restart from checkpoint " << restart_chkfile << "\n";
223
224 // Header
225 std::string File(restart_chkfile + "/Header");
226
227 VisMF::IO_Buffer io_buffer(VisMF::GetIOBufferSize());
228
229 Vector<char> fileCharPtr;
230 ParallelDescriptor::ReadAndBcastFile(File, fileCharPtr);
231 std::string fileCharPtrString(fileCharPtr.dataPtr());
232 std::istringstream is(fileCharPtrString, std::istringstream::in);
233
234 std::string line, word;
235
236 int chk_ncomp;
237
238 // read in title line
239 std::getline(is, line);
240
241 // read in finest_level
242 is >> finest_level;
243 GotoNextLine(is);
244
245 // read the number of components
246 // for each variable we store
247
248 // conservative, cell-centered vars
249 is >> chk_ncomp;
250 GotoNextLine(is);
251 if (chk_ncomp != ncons) {
252 amrex::Abort("Checkpoint scalar component count does not match remora.nscalar");
253 }
254
255 // x-velocity on faces
256 is >> chk_ncomp;
257 GotoNextLine(is);
258 AMREX_ASSERT(chk_ncomp == 1);
259
260 // y-velocity on faces
261 is >> chk_ncomp;
262 GotoNextLine(is);
263 AMREX_ASSERT(chk_ncomp == 1);
264
265 // z-velocity on faces
266 is >> chk_ncomp;
267 GotoNextLine(is);
268 AMREX_ASSERT(chk_ncomp == 1);
269
270 is >> chk_ncomp;
271 GotoNextLine(is);
272 AMREX_ASSERT(chk_ncomp == 2);
273
274 is >> chk_ncomp;
275 GotoNextLine(is);
276 AMREX_ASSERT(chk_ncomp == 2);
277
278 is >> chk_ncomp;
279 GotoNextLine(is);
280 AMREX_ASSERT(chk_ncomp == 3);
281
282 is >> chk_ncomp;
283 GotoNextLine(is);
284 AMREX_ASSERT(chk_ncomp == 3);
285
286 // read in array of istep
287 std::getline(is, line);
288 {
289 std::istringstream lis(line);
290 int i = 0;
291 while (lis >> word) {
292 istep[i++] = std::stoi(word);
293 }
294 }
295
296 // read in array of dt
297 std::getline(is, line);
298 {
299 std::istringstream lis(line);
300 int i = 0;
301 while (lis >> word) {
302#ifdef AMREX_USE_FLOAT
303 dt[i++] = std::stof(word);
304#else
305 dt[i++] = std::stod(word);
306#endif
307 }
308 }
309
310 // read in array of t_new
311 std::getline(is, line);
312 {
313 std::istringstream lis(line);
314 int i = 0;
315 while (lis >> word) {
316#ifdef AMREX_USE_FLOAT
317 t_new[i++] = std::stof(word);
318#else
319 t_new[i++] = std::stod(word);
320#endif
321 }
322 }
323
324 for (int lev = 0; lev <= finest_level; ++lev) {
325
326 // read in level 'lev' BoxArray from Header
327 BoxArray ba;
328 ba.readFrom(is);
329 GotoNextLine(is);
330
331 // create a distribution mapping
332 DistributionMapping dm { ba, ParallelDescriptor::NProcs() };
333
334 MakeNewLevelFromScratch (lev, t_new[lev], ba, dm);
335 }
336
337 // read in the MultiFab data
338 for (int lev = 0; lev <= finest_level; ++lev)
339 {
340 BoxList bl2d = grids[lev].boxList();
341 for (auto& b : bl2d) {
342 b.setRange(2,0);
343 }
344 BoxArray ba2d(std::move(bl2d));
345
346 MultiFab cons(grids[lev],dmap[lev],ncons,cons_new[lev]->nGrowVect());
347 VisMF::Read(cons, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Cell"));
348 MultiFab::Copy(*cons_new[lev],cons,0,0,ncons,cons_new[lev]->nGrowVect());
349
350 VisMF::Read(cons, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Cell_old"));
351 MultiFab::Copy(*cons_old[lev],cons,0,0,ncons,cons_old[lev]->nGrowVect());
352
353 MultiFab xvel(convert(grids[lev],IntVect(1,0,0)),dmap[lev],1,xvel_new[lev]->nGrowVect());
354 VisMF::Read(xvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XFace"));
355 MultiFab::Copy(*xvel_new[lev],xvel,0,0,1,xvel_new[lev]->nGrowVect());
356
357 VisMF::Read(xvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XFace_old"));
358 MultiFab::Copy(*xvel_old[lev],xvel,0,0,1,xvel_old[lev]->nGrowVect());
359
360 MultiFab yvel(convert(grids[lev],IntVect(0,1,0)),dmap[lev],1,yvel_new[lev]->nGrowVect());
361 VisMF::Read(yvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YFace"));
362 MultiFab::Copy(*yvel_new[lev],yvel,0,0,1,yvel_new[lev]->nGrowVect());
363
364 VisMF::Read(yvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YFace_old"));
365 MultiFab::Copy(*yvel_old[lev],yvel,0,0,1,yvel_old[lev]->nGrowVect());
366
367 MultiFab zvel(convert(grids[lev],IntVect(0,0,1)),dmap[lev],1,zvel_new[lev]->nGrowVect());
368 VisMF::Read(zvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "ZFace"));
369 MultiFab::Copy(*zvel_new[lev],zvel,0,0,1,zvel_new[lev]->nGrowVect());
370
371 VisMF::Read(zvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "ZFace_old"));
372 MultiFab::Copy(*zvel_old[lev],zvel,0,0,1,zvel_old[lev]->nGrowVect());
373
374 MultiFab mf_ru(convert(grids[lev],IntVect(1,0,0)),dmap[lev],2,(vec_ru[lev])->nGrowVect());
375 VisMF::Read(mf_ru, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XRHS"));
376 MultiFab::Copy(*(vec_ru[lev]),mf_ru,0,0,2,(vec_ru[lev])->nGrowVect());
377
378 MultiFab mf_rv(convert(grids[lev],IntVect(0,1,0)),dmap[lev],2,(vec_rv[lev])->nGrowVect());
379 VisMF::Read(mf_rv, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YRHS"));
380 MultiFab::Copy(*(vec_rv[lev]),mf_rv,0,0,2,(vec_rv[lev])->nGrowVect());
381
382 MultiFab mf_ubar(convert(ba2d,IntVect(1,0,0)),dmap[lev],3,(vec_ubar[lev])->nGrowVect());
383 VisMF::Read(mf_ubar, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XBar"));
384 MultiFab::Copy(*(vec_ubar[lev]),mf_ubar,0,0,3,(vec_ubar[lev])->nGrowVect());
385
386 MultiFab mf_vbar(convert(ba2d,IntVect(0,1,0)),dmap[lev],3,(vec_vbar[lev])->nGrowVect());
387 VisMF::Read(mf_vbar, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YBar"));
388 MultiFab::Copy(*(vec_vbar[lev]),mf_vbar,0,0,3,(vec_vbar[lev])->nGrowVect());
389
390 MultiFab mf_ru2d(convert(ba2d,IntVect(1,0,0)),dmap[lev],2,(vec_ru2d[lev])->nGrowVect());
391 VisMF::Read(mf_ru2d, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XRHS2d"));
392 MultiFab::Copy(*(vec_ru2d[lev]),mf_ru2d,0,0,2,(vec_ru2d[lev])->nGrowVect());
393
394 MultiFab mf_rv2d(convert(ba2d,IntVect(0,1,0)),dmap[lev],2,(vec_rv2d[lev])->nGrowVect());
395 VisMF::Read(mf_rv2d, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "YRHS2d"));
396 MultiFab::Copy(*(vec_rv2d[lev]),mf_rv2d,0,0,2,(vec_rv2d[lev])->nGrowVect());
397
398 MultiFab mf_mskr(ba2d,dmap[lev],1,(vec_mskr[lev])->nGrowVect());
399 VisMF::Read(mf_mskr, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Mskr"));
400 MultiFab::Copy(*(vec_mskr[lev]),mf_mskr,0,0,1,(vec_mskr[lev])->nGrowVect());
401
402 MultiFab mf_msku(convert(ba2d,IntVect(1,0,0)),dmap[lev],1,(vec_msku[lev])->nGrowVect());
403 VisMF::Read(mf_msku, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Msku"));
404 MultiFab::Copy(*(vec_msku[lev]),mf_msku,0,0,1,(vec_msku[lev])->nGrowVect());
405
406 MultiFab mf_mskv(convert(ba2d,IntVect(0,1,0)),dmap[lev],1,(vec_mskv[lev])->nGrowVect());
407 VisMF::Read(mf_mskv, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Mskv"));
408 MultiFab::Copy(*(vec_mskv[lev]),mf_mskv,0,0,1,(vec_mskv[lev])->nGrowVect());
409
411 VisMF::Read(*(vec_rdrag[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "rdrag"));
413 VisMF::Read(*(vec_rdrag2[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "rdrag2"));
414 }
415
418 VisMF::Read(*(vec_ZoBot[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "ZoBot"));
419 }
420
421 VisMF::Read(*(vec_DU_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "DU_avg1"));
422 VisMF::Read(*(vec_DU_avg2[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "DU_avg2"));
423 VisMF::Read(*(vec_DV_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "DV_avg1"));
424 VisMF::Read(*(vec_DV_avg2[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "DV_avg2"));
425
426 VisMF::Read(*(vec_zeta[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "zeta"));
427 VisMF::Read(*(vec_Zt_avg1[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Zt_avg1"));
428
429 VisMF::Read(*(vec_h[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "h"));
430
431 VisMF::Read(*(vec_tke[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "tke"));
432 VisMF::Read(*(vec_gls[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "gls"));
433 VisMF::Read(*(vec_Lscale[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Lscale"));
434 VisMF::Read(*(vec_Akk[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Akk"));
435 VisMF::Read(*(vec_Akp[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Akp"));
436 VisMF::Read(*(vec_Akt[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Akt"));
437 VisMF::Read(*(vec_Akv[lev]), amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Akv"));
438
440
441 }
442
443#ifdef REMORA_USE_PARTICLES
444 particleData.Restart((amrex::ParGDBBase*)GetParGDB(),restart_chkfile);
445#endif
446}
static void GotoNextLine(std::istream &is)
utility to skip to next line in Header
int ncons
Number of conserved scalars in the state (temperature + salt + passive scalars)
Definition REMORA.H:1346
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:272
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_ZoBot
Bottom roughness length [m], defined at rho points.
Definition REMORA.H:364
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg2
correct time average of barotropic x velocity flux for coupling (2D)
Definition REMORA.H:374
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:242
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:393
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_tke
Turbulent kinetic energy.
Definition REMORA.H:452
void ReadCheckpointFile()
read checkpoint file from disk
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_gls
Turbulent generic length scale.
Definition REMORA.H:454
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:270
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:248
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Lscale
Vertical mixing turbulent length scale.
Definition REMORA.H:456
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
Definition REMORA.H:280
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_msku
land/sea mask at x-faces (2D)
Definition REMORA.H:395
std::string check_file
Checkpoint file prefix.
Definition REMORA.H:1376
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:235
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 * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
Definition REMORA.H:239
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 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:1447
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:237
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg1
time average of barotropic y velocity flux
Definition REMORA.H:376
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1273
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1403
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akk
Turbulent kinetic energy vertical diffusion coefficient.
Definition REMORA.H:458
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag2
Quadratic drag coefficient [unitless], defined at rho points.
Definition REMORA.H:362
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:266
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
Definition REMORA.H:390
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
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:386
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:233
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg2
correct time average of barotropic y velocity flux for coupling (2D)
Definition REMORA.H:378
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
Definition REMORA.H:278
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< std::unique_ptr< amrex::MultiFab > > vec_rdrag
Linear drag coefficient [m/s], defined at rho points.
Definition REMORA.H:360
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:312
std::string restart_chkfile
If set, restart from this checkpoint file.
Definition REMORA.H:1341
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:268
amrex::Vector< amrex::Real > dt
time step at each level
Definition REMORA.H:1277
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akp
Turbulent length scale vertical diffusion coefficient.
Definition REMORA.H:460
BottomStressType bottom_stress_type
VertMixingType vert_mixing_type