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