REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA.cpp
Go to the documentation of this file.
1/**
2 * \file REMORA.cpp
3 */
4
6#include <REMORA.H>
7
8#include <AMReX_buildInfo.H>
9
10using namespace amrex;
11
12amrex::Real REMORA::startCPUTime = 0.0_rt;
13amrex::Real REMORA::previousCPUTimeUsed = 0.0_rt;
14
15Vector<AMRErrorTag> REMORA::ref_tags;
16
18
19// Time step control
20amrex::Real REMORA::cfl = 0.8_rt;
21amrex::Real REMORA::fixed_dt = -1.0_rt;
22amrex::Real REMORA::fixed_fast_dt = -1.0_rt;
23amrex::Real REMORA::change_max = 1.1_rt;
24
26
27// Dictate verbosity in screen output
28int REMORA::verbose = 0;
29
30// Frequency of diagnostic output
32amrex::Real REMORA::sum_per = -1.0_rt;
33
34// Minimum number of digits in plotfile name
36
37// Native AMReX vs NetCDF
39
40#ifdef REMORA_USE_NETCDF
41
43
44// Do we write one file per timestep (false) or one file for all timesteps (true)
46
47// NetCDF initialization file
48amrex::Vector<std::string> REMORA::nc_bdry_file = {""}; // Must provide via input
49amrex::Vector<amrex::Vector<std::string>> REMORA::nc_init_file = {{""}}; // Must provide via input
50amrex::Vector<amrex::Vector<std::string>> REMORA::nc_grid_file = {{""}}; // Must provide via input
51#endif
52
53amrex::Vector<std::string> BCNames = {"xlo", "ylo", "zlo", "xhi", "yhi", "zhi"};
54
55/**
56 * constructor:
57 * - reads in parameters from inputs file
58 * - sizes multilevel arrays and data structures
59 * - initializes BCRec boundary condition object
60 */
62{
63 if (ParallelDescriptor::IOProcessor()) {
64 const char* remora_hash = amrex::buildInfoGetGitHash(1);
65 const char* amrex_hash = amrex::buildInfoGetGitHash(2);
66 const char* buildgithash = amrex::buildInfoGetBuildGitHash();
67 const char* buildgitname = amrex::buildInfoGetBuildGitName();
68
69 if (strlen(remora_hash) > 0) {
70 amrex::Print() << "\n"
71 << "REMORA git hash: " << remora_hash << "\n";
72 }
73 if (strlen(amrex_hash) > 0) {
74 amrex::Print() << "AMReX git hash: " << amrex_hash << "\n";
75 }
76 if (strlen(buildgithash) > 0) {
77 amrex::Print() << buildgitname << " git hash: " << buildgithash << "\n";
78 }
79
80 amrex::Print() << "\n";
81 }
82
84 const std::string& pv1 = "plot_vars"; setPlotVariables(pv1);
85
86 prob = amrex_probinit(geom[0].ProbLo(),geom[0].ProbHi());
87
88 // Geometry on all levels has been defined already.
89
90 // No valid BoxArray and DistributionMapping have been defined.
91 // But the arrays for them have been resized.
92
93 int nlevs_max = max_level + 1;
94
95 istep.resize(nlevs_max, 0);
96 nsubsteps.resize(nlevs_max, 1);
97 for (int lev = 1; lev <= max_level; ++lev) {
98 nsubsteps[lev] = do_substep ? MaxRefRatio(lev-1) : 1;
99 }
100
101 physbcs.resize(nlevs_max);
102
103 t_new.resize(nlevs_max, 0.0_rt);
104 t_old.resize(nlevs_max, -1.e100_rt);
105 dt.resize(nlevs_max, 1.e100_rt);
106
107 cons_new.resize(nlevs_max);
108 cons_old.resize(nlevs_max);
109 xvel_new.resize(nlevs_max);
110 xvel_old.resize(nlevs_max);
111 yvel_new.resize(nlevs_max);
112 yvel_old.resize(nlevs_max);
113 zvel_new.resize(nlevs_max);
114 zvel_old.resize(nlevs_max);
115
116 advflux_reg.resize(nlevs_max);
117
118 // Initialize tagging criteria for mesh refinement
120
121 // We have already read in the ref_Ratio (via amr.ref_ratio =) but we need to enforce
122 // that there is no refinement in the vertical so we test on that here.
123 for (int lev = 0; lev < max_level; ++lev)
124 {
125 amrex::Print() << "Refinement ratio at level " << lev << " set to be " <<
126 ref_ratio[lev][0] << " " << ref_ratio[lev][1] << " " << ref_ratio[lev][2] << std::endl;
127
128 if (ref_ratio[lev][2] != 1)
129 {
130 amrex::Error("We don't allow refinement in the vertical -- make sure to set ref_ratio = 1 in z");
131 }
132 }
133}
134
136{
137}
138
139void
141{
142 Real cur_time = t_new[0];
143
144 // Take one coarse timestep by calling timeStep -- which recursively calls timeStep
145 // for finer levels (with or without subcycling)
146 for (int step = istep[0]; step < max_step && cur_time < stop_time; ++step)
147 {
148 amrex::Print() << "\nCoarse STEP " << step+1 << " starts ..." << std::endl;
149
150 ComputeDt();
151
152 int lev = 0;
153 int iteration = 1;
154
155 if (max_level == 0) {
156 timeStep(lev, cur_time, iteration);
157 }
158 else {
159 timeStepML(cur_time, iteration);
160 }
161
162 cur_time += dt[0];
163
164 amrex::Print() << "Coarse STEP " << step+1 << " ends." << " TIME = " << cur_time
165 << " DT = " << dt[0] << std::endl;
166
167 if ((plot_int > 0 && (step+1 - last_plot_file_step) == plot_int)
168 || (plot_int_time > 0 && cur_time >= (last_plot_file_time + plot_int_time))) {
169 last_plot_file_step = step+1;
170 last_plot_file_time = cur_time;
172
174 }
175#ifdef REMORA_USE_NETCDF
177 WriteNCPlotFile(step+1);
179 }
180#endif
181 }
182
183 if ((check_int > 0 && (step+1 - last_check_file_step) == check_int)
184 || (check_int_time > 0 && cur_time >= (last_check_file_time + check_int_time))) {
185 last_check_file_step = step+1;
186 last_check_file_time = cur_time;
188 }
189
190 post_timestep(step, cur_time, dt[0]);
191
192#ifdef AMREX_MEM_PROFILING
193 {
194 std::ostringstream ss;
195 ss << "[STEP " << step+1 << "]";
196 MemProfiler::report(ss.str());
197 }
198#endif
199
200 if (cur_time >= stop_time - 1.e-6*dt[0]) break;
201 }
202
203 if ((plot_int > 0 || plot_int_time > 0.0) && istep[0] > last_plot_file_step) {
206 }
207#ifdef REMORA_USE_NETCDF
211 }
212#endif
213 }
214
215 if ((check_int > 0 || check_int_time > 0.0) && istep[0] > last_check_file_step) {
217 }
218}
219
220/**
221 * @param[in ] nstep which step we're on
222 * @param[in ] time current time
223 * @param[in ] dt_lev0 time step on level 0
224 */
225void
226REMORA::post_timestep (int nstep, Real time, Real dt_lev0)
227{
228 BL_PROFILE("REMORA::post_timestep()");
229
230#ifdef REMORA_USE_PARTICLES
231 particleData.Redistribute();
232#endif
233
235 {
236 for (int lev = finest_level-1; lev >= 0; lev--)
237 {
238 // This call refluxes from the lev/lev+1 interface onto lev
239 //getAdvFluxReg(lev+1)->Reflux(*cons_new[lev], 0, 0, NCONS);
240
241 // We need to do this before anything else because refluxing changes the
242 // values of coarse cells underneath fine grids with the assumption they'll
243 // be over-written by averaging down
244 //
245 AverageDownTo(lev);
246 }
247 }
248
249 if (is_it_time_for_action(nstep, time, dt_lev0, sum_interval, sum_per)) {
251 }
252}
253
254/**
255 * This is called from main.cpp and handles all initialization, whether from start or restart
256 */
257void
259{
260 // Initialize the start time for our CPU-time tracker
261 startCPUTime = Real(ParallelDescriptor::second());
262
263 // Map the words in the inputs file to BC types, then translate
264 // those types into what they mean for each variable
265 init_bcs();
266
269 last_plot_file_time = -1.0_rt;
270 last_check_file_time = -1.0_rt;
271
272 if (restart_chkfile == "") {
273 // start simulation from the beginning
274
275 InitFromScratch(start_time);
276
278 AverageDown();
279 }
280
281 } else { // Restart from a checkpoint
282
283 restart();
284
285 }
286#ifdef REMORA_USE_MOAB
287 InitMOABMesh();
288#endif
289 // Initialize flux registers (whether we start from scratch or restart)
291 advflux_reg[0] = nullptr;
292 for (int lev = 1; lev <= finest_level; lev++)
293 {
294 advflux_reg[lev] = new YAFluxRegister(grids[lev], grids[lev-1],
295 dmap[lev], dmap[lev-1],
296 geom[lev], geom[lev-1],
297 ref_ratio[lev-1], lev, NCONS);
298 }
299 }
300
301 // Fill ghost cells/faces
302 for (int lev = 0; lev <= finest_level; ++lev)
303 {
304 if (lev > 0 && cf_width >= 0) {
306 }
307
308 if (restart_chkfile == "") {
309 FillPatch(lev, t_new[lev], *cons_new[lev], cons_new, BCVars::cons_bc, BdyVars::t, 0, true, false,0,0,0.0,*cons_new[lev]);
310 FillPatch(lev, t_new[lev], *xvel_new[lev], xvel_new, BCVars::xvel_bc, BdyVars::u, 0, true, false,0,0,0.0,*xvel_new[lev]);
311 FillPatch(lev, t_new[lev], *yvel_new[lev], yvel_new, BCVars::yvel_bc, BdyVars::v, 0, true, false,0,0,0.0,*yvel_new[lev]);
312 FillPatch(lev, t_new[lev], *zvel_new[lev], zvel_new, BCVars::zvel_bc, BdyVars::null, 0, true, false);
313
314 // Copy from new into old just in case when initializing from scratch
315 int ngs = cons_new[lev]->nGrow();
316 int ngvel = xvel_new[lev]->nGrow();
317 MultiFab::Copy(*cons_old[lev],*cons_new[lev],0,0,NCONS,ngs);
318 MultiFab::Copy(*xvel_old[lev],*xvel_new[lev],0,0,1,ngvel);
319 MultiFab::Copy(*yvel_old[lev],*yvel_new[lev],0,0,1,ngvel);
320 MultiFab::Copy(*zvel_old[lev],*zvel_new[lev],0,0,1,IntVect(ngvel,ngvel,0));
321 }
322 } // lev
323
324 // Check for additional plotting variables that are available after
325 // particle containers are setup.
326 const std::string& pv1 = "plot_vars"; appendPlotVariables(pv1);
327
328 if (restart_chkfile == "" && (check_int > 0 || check_int_time > 0.0_rt))
329 {
332 }
333
334 if ( (restart_chkfile == "") ||
336 {
337 if (plot_int > 0 || plot_int_time > 0.0)
338 {
341#ifdef REMORA_USE_NETCDF
343 int step0 = 0;
344 WriteNCPlotFile(step0);
346 }
347#endif
349 }
350 }
351
354 }
355
356 ComputeDt();
357
358}
359
360/**
361 * @param[in ] lev level to operate on
362 */
363void
365{
366 amrex::Print() << ":::Construct_REMORAFillPatchers " << lev << std::endl;
367
368 auto& ba_fine = cons_new[lev ]->boxArray();
369 auto& ba_crse = cons_new[lev-1]->boxArray();
370 auto& dm_fine = cons_new[lev ]->DistributionMap();
371 auto& dm_crse = cons_new[lev-1]->DistributionMap();
372
373 BoxList bl2d_fine = ba_fine.boxList();
374 for (auto& b : bl2d_fine) {
375 b.setRange(2,0);
376 }
377 BoxArray ba2d_fine(std::move(bl2d_fine));
378
379 BoxList bl2d_crse = ba_crse.boxList();
380 for (auto& b : bl2d_crse) {
381 b.setRange(2,0);
382 }
383 BoxArray ba2d_crse(std::move(bl2d_crse));
384
385 int ncomp = cons_new[lev]->nComp();
386
387 FPr_c.emplace_back(ba_fine, dm_fine, geom[lev] ,
388 ba_crse, dm_crse, geom[lev-1],
389 -cf_width, -cf_set_width, ncomp, &cell_cons_interp);
390 FPr_u.emplace_back(convert(ba_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
391 convert(ba_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
392 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
393 FPr_v.emplace_back(convert(ba_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
394 convert(ba_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
395 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
396 FPr_w.emplace_back(convert(ba_fine, IntVect(0,0,1)), dm_fine, geom[lev] ,
397 convert(ba_crse, IntVect(0,0,1)), dm_crse, geom[lev-1],
398 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
399
400 FPr_ubar.emplace_back(convert(ba2d_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
401 convert(ba2d_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
402 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
403 FPr_vbar.emplace_back(convert(ba2d_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
404 convert(ba2d_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
405 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
406}
407
408/**
409 * @param[in ] lev level to operate on
410 */
411void
413{
414 amrex::Print() << ":::Define_REMORAFillPatchers " << lev << std::endl;
415
416 auto& ba_fine = cons_new[lev ]->boxArray();
417 auto& ba_crse = cons_new[lev-1]->boxArray();
418 auto& dm_fine = cons_new[lev ]->DistributionMap();
419 auto& dm_crse = cons_new[lev-1]->DistributionMap();
420
421 BoxList bl2d_fine = ba_fine.boxList();
422 for (auto& b : bl2d_fine) {
423 b.setRange(2,0);
424 }
425 BoxArray ba2d_fine(std::move(bl2d_fine));
426
427 BoxList bl2d_crse = ba_crse.boxList();
428 for (auto& b : bl2d_crse) {
429 b.setRange(2,0);
430 }
431 BoxArray ba2d_crse(std::move(bl2d_crse));
432
433
434 int ncomp = cons_new[lev]->nComp();
435
436 FPr_c[lev-1].Define(ba_fine, dm_fine, geom[lev] ,
437 ba_crse, dm_crse, geom[lev-1],
438 -cf_width, -cf_set_width, ncomp, &cell_cons_interp);
439 FPr_u[lev-1].Define(convert(ba_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
440 convert(ba_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
441 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
442 FPr_v[lev-1].Define(convert(ba_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
443 convert(ba_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
444 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
445 FPr_w[lev-1].Define(convert(ba_fine, IntVect(0,0,1)), dm_fine, geom[lev] ,
446 convert(ba_crse, IntVect(0,0,1)), dm_crse, geom[lev-1],
447 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
448
449 FPr_ubar[lev-1].Define(convert(ba2d_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
450 convert(ba2d_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
451 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
452 FPr_vbar[lev-1].Define(convert(ba2d_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
453 convert(ba2d_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
454 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
455}
456
457void
459{
461
462 // We set this here so that we don't over-write the checkpoint file we just started from
464}
465
466/**
467 * @param[in ] lev level to operate on
468 */
469void
471{
473 prob->init_analytic_zeta(lev, geom[lev], solverChoice, *this, *vec_zeta[lev]);
474
475#ifdef REMORA_USE_NETCDF
477 amrex::Print() << "Calling init_zeta_from_netcdf " << std::endl;
479 amrex::Print() << "Sea surface height loaded from netcdf file \n " << std::endl;
480#endif
481 } else {
482 Abort("Don't know this ic_bc_type!");
483 }
484 vec_zeta[lev]->FillBoundary(geom[lev].periodicity());
485 set_zeta_average(lev);
486}
487
488/**
489 * @param[in ] lev level to operate on
490 */
491void
493{
494 // Only set bathymetry on level 0, and interpolate for finer levels
496 if (lev==0) {
498 prob->init_analytic_bathymetry(lev, geom[lev], solverChoice, *this, *vec_hOfTheConfusingName[lev]);
499
500#ifdef REMORA_USE_NETCDF
502 amrex::Print() << "Calling init_bathymetry_from_netcdf " << std::endl;
504 amrex::Print() << "Bathymetry loaded from netcdf file \n " << std::endl;
505#endif
506 } else {
507 Abort("Don't know this ic_bc_type!");
508 }
509 // Need FillBoundary to fill at grid-grid boundaries, and EnforcePeriodicity
510 // to make sure ghost cells in the domain corners are consistent.
511 vec_hOfTheConfusingName[lev]->FillBoundary(geom[lev].periodicity());
512 vec_hOfTheConfusingName[lev]->EnforcePeriodicity(geom[lev].periodicity());
513 } else {
514 Real dummy_time = 0.0_rt;
515 FillCoarsePatch(lev,dummy_time,vec_hOfTheConfusingName[lev].get(), vec_hOfTheConfusingName[lev-1].get());
516 }
517 } else if (solverChoice.init_ana_h) {
519 prob->init_analytic_bathymetry(lev, geom[lev], solverChoice, *this,*vec_hOfTheConfusingName[lev]);
520
521#ifdef REMORA_USE_NETCDF
523 amrex::Print() << "Calling init_bathymetry_from_netcdf " << std::endl;
525 amrex::Print() << "Bathymetry loaded from netcdf file \n " << std::endl;
526#endif
527 } else {
528 Abort("Don't know this ic_bc_type!");
529 }
530 // Need FillBoundary to fill at grid-grid boundaries, and EnforcePeriodicity
531 // to make sure ghost cells in the domain corners are consistent.
532 vec_hOfTheConfusingName[lev]->FillBoundary(geom[lev].periodicity());
533 vec_hOfTheConfusingName[lev]->EnforcePeriodicity(geom[lev].periodicity());
534 } else if (solverChoice.init_l1ad_h) {
536 prob->init_analytic_bathymetry(lev, geom[lev], solverChoice, *this, *vec_hOfTheConfusingName[lev]);
537
538#ifdef REMORA_USE_NETCDF
540 amrex::Print() << "Calling init_bathymetry_from_netcdf " << std::endl;
542 amrex::Print() << "Bathymetry loaded from netcdf file \n " << std::endl;
543#endif
544 } else {
545 Abort("Don't know this ic_bc_type!");
546 }
547 // Need FillBoundary to fill at grid-grid boundaries, and EnforcePeriodicity
548 // to make sure ghost cells in the domain corners are consistent.
549 vec_hOfTheConfusingName[lev]->FillBoundary(geom[lev].periodicity());
550 vec_hOfTheConfusingName[lev]->EnforcePeriodicity(geom[lev].periodicity());
551 } else {
552 amrex::Abort("Don't know this h init type");
553 }
554}
555
556/**
557 * @param[in ] lev level to operate on
558 */
559void
563 prob->init_analytic_coriolis(lev, geom[lev], solverChoice, *this, *vec_fcor[lev]);
566#ifdef REMORA_USE_NETCDF
568 amrex::Print() << "Calling init_coriolis_from_netcdf " << std::endl;
570 amrex::Print() << "Coriolis loaded from netcdf file \n" << std::endl;
571#endif
572 } else {
573 Abort("Don't know this coriolis_type!");
574 }
575
576 Real time = 0.0_rt;
577 FillPatch(lev, time, *vec_fcor[lev], GetVecOfPtrs(vec_fcor),BCVars::foextrap_bc);
578 vec_fcor[lev]->EnforcePeriodicity(geom[lev].periodicity());
579 }
580}
581
582void
588 // The GLS initialization just sets the multifab to a value, so there's
589 // no need to call FillPatch here
590 } else {
591 Abort("Don't know this vertical mixing type");
592 }
593}
594
595/**
596 * @param[in ] lev level to operate on
597 */
598void
600 Real time = 0.0_rt;
601 prob->init_analytic_vmix(lev, geom[lev], solverChoice, *this,*vec_Akv[lev], *vec_Akt[lev]);
602 FillPatch(lev, time, *vec_Akv[lev], GetVecOfPtrs(vec_Akv),BCVars::zvel_bc,BdyVars::null,0,true,false);
603 for (int n=0; n<NCONS;n++) {
604 FillPatch(lev, time, *vec_Akt[lev], GetVecOfPtrs(vec_Akt),BCVars::zvel_bc,BdyVars::null,0,false,false);
605 }
606}
607
608/**
609 * @param[in ] lev level to operate on
610 */
611void
613{
615 prob->init_analytic_hmix(lev, geom[lev], solverChoice, *this, *vec_visc2_p[lev], *vec_visc2_r[lev], *vec_diff2[lev]);
617 vec_visc2_p[lev]->setVal(solverChoice.visc2);
618 vec_visc2_r[lev]->setVal(solverChoice.visc2);
619 for (int n=0; n<NCONS; n++) {
620 vec_diff2[lev]->setVal(solverChoice.tnu2[n],n,1);
621 }
622 } else {
623 Abort("Don't know this horizontal mixing type");
624 }
625 Real time = 0.0_rt;
626 FillPatch(lev, time, *vec_visc2_p[lev], GetVecOfPtrs(vec_visc2_p),BCVars::foextrap_periodic_bc);
627 FillPatch(lev, time, *vec_visc2_r[lev], GetVecOfPtrs(vec_visc2_r),BCVars::foextrap_periodic_bc);
628 for (int n=0; n<NCONS; n++) {
629 FillPatch(lev, time, *vec_diff2[lev] , GetVecOfPtrs(vec_diff2),BCVars::foextrap_periodic_bc,BdyVars::null,n,false);
630 }
631}
632
633/**
634 * @param[in ] lev level to operate on
635 */
636void
638{
640 prob->init_analytic_smflux(lev, geom[lev], solverChoice, *this,*vec_sustr[lev], *vec_svstr[lev]);
642#ifdef REMORA_USE_NETCDF
645 FillPatch(lev, t_old[lev], *vec_sustr[lev], GetVecOfPtrs(vec_sustr),BCVars::foextrap_periodic_bc,BdyVars::null,0,false);
646 FillPatch(lev, t_old[lev], *vec_svstr[lev], GetVecOfPtrs(vec_svstr),BCVars::foextrap_periodic_bc,BdyVars::null,0,false);
647#endif
648 }
649}
650
651/**
652 * @param[in ] lev level to operate on
653 */
654void
656{
658 prob->init_analytic_wind(lev,geom[lev], solverChoice, *this, *vec_uwind[lev], *vec_vwind[lev]);
660#ifdef REMORA_USE_NETCDF
663 FillPatch(lev, t_old[lev], *vec_uwind[lev], GetVecOfPtrs(vec_uwind),BCVars::foextrap_periodic_bc,BdyVars::null,0,false);
664 FillPatch(lev, t_old[lev], *vec_vwind[lev], GetVecOfPtrs(vec_vwind),BCVars::foextrap_periodic_bc,BdyVars::null,0,false);
665#endif
666 }
667}
668
669/**
670 * @param[in ] lev level to operate on
671 * @param[in ] time current time for initialization
672 */
673void
674REMORA::init_only (int lev, Real time)
675{
676 t_new[lev] = time;
677 t_old[lev] = time - 1.e200_rt;
678
679 cons_new[lev]->setVal(0.0_rt);
680 xvel_new[lev]->setVal(0.0_rt);
681 yvel_new[lev]->setVal(0.0_rt);
682 zvel_new[lev]->setVal(0.0_rt);
683
684 xvel_old[lev]->setVal(0.0_rt);
685 yvel_old[lev]->setVal(0.0_rt);
686 zvel_old[lev]->setVal(0.0_rt);
687
688 vec_ru[lev]->setVal(0.0_rt);
689 vec_rv[lev]->setVal(0.0_rt);
690
691 vec_ru2d[lev]->setVal(0.0_rt);
692 vec_rv2d[lev]->setVal(0.0_rt);
693
694#ifdef REMORA_USE_NETCDF
696 amrex::Print() << "Calling init_masks_from_netcdf " << std::endl;
698 amrex::Print() << "Masks loaded from netcdf file \n " << std::endl;
699
700 amrex::Print() << "Calling init_bdry_from_netcdf " << std::endl;
702 amrex::Print() << "Boundary data loaded from netcdf file \n " << std::endl;
703
705
707 if (nc_clim_his_file.empty()) {
708 amrex::Error("NetCDF climatology file name must be provided via input");
709 }
711 ubar_clim_data_from_file = new NCTimeSeries(nc_clim_his_file, "ubar", clim_ubar_time_varname, geom[lev].Domain(),vec_ubar[lev].get(),true,true);
712 vbar_clim_data_from_file = new NCTimeSeries(nc_clim_his_file, "vbar", clim_ubar_time_varname, geom[lev].Domain(),vec_vbar[lev].get(),true,true);
715 }
717 u_clim_data_from_file = new NCTimeSeries(nc_clim_his_file, "u", clim_u_time_varname, geom[lev].Domain(),xvel_new[lev],false,true);
718 v_clim_data_from_file = new NCTimeSeries(nc_clim_his_file, "v", clim_v_time_varname, geom[lev].Domain(),yvel_new[lev],false,true);
721 }
722 // Since the NCTimeSeries object isn't filling the cons_new MultiFab directly, we don't have to specify a component.
723 // It just needs to know the shape of the MultiFab
725 temp_clim_data_from_file = new NCTimeSeries(nc_clim_his_file, "temp", clim_temp_time_varname,geom[lev].Domain(),cons_new[lev],false,true);
727 }
729 salt_clim_data_from_file = new NCTimeSeries(nc_clim_his_file, "salt", clim_salt_time_varname,geom[lev].Domain(),cons_new[lev],false,true);
731 }
732 }
733 }
734 // This will be a non-op if forcings specified analytically
736 if (nc_frc_file.empty()) {
737 amrex::Error("NetCDF forcing file name must be provided via input for winds");
738 }
739 Uwind_data_from_file = new NCTimeSeries(nc_frc_file, "Uwind", frc_time_varname, geom[lev].Domain(),vec_uwind[lev].get(), true, false);
740 Vwind_data_from_file = new NCTimeSeries(nc_frc_file, "Vwind", frc_time_varname, geom[lev].Domain(),vec_vwind[lev].get(), true, false);
744 if (nc_frc_file.empty()) {
745 amrex::Error("NetCDF forcing file name must be provided via input for surface momentum fluxes");
746 }
747 sustr_data_from_file = new NCTimeSeries(nc_frc_file, "sustr", frc_time_varname, geom[lev].Domain(),vec_sustr[lev].get(), true, false);
748 svstr_data_from_file = new NCTimeSeries(nc_frc_file, "svstr", frc_time_varname, geom[lev].Domain(),vec_svstr[lev].get(), true, false);
751 }
752
754 auto dom = geom[0].Domain();
755 int nz = dom.length(2);
756 river_source_cons.resize(NCONS);
757 Print() << solverChoice.do_rivers_cons[0] << std::endl;
760 river_source_cons[Salt_comp]->Initialize();
761 }
764 river_source_cons[Temp_comp]->Initialize();
765 }
768 river_source_cons[Scalar_comp]->Initialize();
769 }
775 }
776#endif
777
779 set_grid_scale(lev);
780 }
781 set_bathymetry(lev);
782 set_zeta(lev);
784
786 if (lev==0) {
788 {
789 init_analytic(lev);
790#ifdef REMORA_USE_NETCDF
792 amrex::Print() << "Calling init_data_from_netcdf " << std::endl;
795 amrex::Print() << "Initial data loaded from netcdf file \n " << std::endl;
796
797#endif
798 } else {
799 Abort("Need to specify ic_bc_type");
800 }
801 } else {
802 FillCoarsePatch(lev, time, cons_new[lev], cons_new[lev-1]);
803 FillCoarsePatch(lev, time, xvel_new[lev], xvel_new[lev-1]);
804 FillCoarsePatch(lev, time, yvel_new[lev], yvel_new[lev-1]);
805 FillCoarsePatch(lev, time, zvel_new[lev], zvel_new[lev-1]);
806 }
809 {
810 init_analytic(lev);
811#ifdef REMORA_USE_NETCDF
813 amrex::Print() << "Calling init_data_from_netcdf " << std::endl;
816 amrex::Print() << "Initial data loaded from netcdf file \n " << std::endl;
817
818#endif
819 } else {
820 Abort("Need to specify ic_bc_type");
821 }
822 } else {
823 amrex::Abort("Need to specify T init procedure");
824 }
825
826 // Ensure that the face-based data are the same on both sides of a periodic domain.
827 // The data associated with the lower grid ID is considered the correct value.
828 xvel_new[lev]->OverrideSync(geom[lev].periodicity());
829 yvel_new[lev]->OverrideSync(geom[lev].periodicity());
830 zvel_new[lev]->OverrideSync(geom[lev].periodicity());
831
832 set_2darrays(lev);
833
834 init_set_vmix(lev);
835 set_hmixcoef(lev);
836 set_coriolis(lev);
837
838 // Previously set smflux here with OverrideSync:
839// set_smflux(lev);
840// prob->init_analytic_smflux(lev, geom[lev], solverChoice, *this, *vec_sustr[lev], *vec_svstr[lev]);
841// vec_sustr[lev]->OverrideSync(geom[lev].periodicity());
842// vec_svstr[lev]->OverrideSync(geom[lev].periodicity());
843
844}
845
846void
848{
849 {
850 ParmParse pp; // Traditionally, max_step and stop_time do not have prefix, so allow it for now.
851 bool noprefix_max_step = pp.query("max_step", max_step);
852 bool noprefix_stop_time = pp.query("stop_time", stop_time);
853 bool remora_max_step = pp.query("remora.max_step", max_step);
854 bool remora_stop_time = pp.query("remora.stop_time", stop_time);
855 if (remora_max_step and noprefix_max_step) {
856 Abort("remora.max_step and max_step are both specified. Please use only one!");
857 }
858 if (remora_stop_time and noprefix_stop_time) {
859 Abort("remora.stop_time and stop_time are both specified. Please use only one!");
860 }
861 }
862
863 ParmParse pp(pp_prefix);
864 ParmParse pp_amr("amr");
865 {
866 pp_amr.query("regrid_int", regrid_int);
867 pp.query("check_file", check_file);
868 pp.query("check_int", check_int);
869 pp_amr.query("check_int", check_int);
870 pp.query("check_int_time", check_int_time);
871 pp_amr.query("check_int_time", check_int_time);
872
873 pp.query("restart", restart_chkfile);
874 pp_amr.query("restart", restart_chkfile);
875 pp.query("start_time",start_time);
876
877 if (pp.contains("data_log"))
878 {
879 int num_datalogs = pp.countval("data_log");
880 datalog.resize(num_datalogs);
881 datalogname.resize(num_datalogs);
882 pp.queryarr("data_log",datalogname,0,num_datalogs);
883 for (int i = 0; i < num_datalogs; i++)
885 }
886
887 // Verbosity
888 pp.query("v", verbose);
889
890 // Frequency of diagnostic output
891 pp.query("sum_interval", sum_interval);
892 pp.query("sum_period" , sum_per);
893 pp.query("file_min_digits", file_min_digits);
894
895 if (file_min_digits < 0) {
896 amrex::Abort("remora.file_min_digits must be non-negative");
897 }
898
899 // Time step controls
900 pp.query("cfl", cfl);
901 pp.query("change_max", change_max);
902
903 pp.query("fixed_dt", fixed_dt);
904 pp.query("fixed_fast_dt", fixed_fast_dt);
905
906 pp.query("fixed_ndtfast_ratio", fixed_ndtfast_ratio);
907
908 // If all three are specified, they must be consistent
909 if (fixed_dt > 0. && fixed_fast_dt > 0. && fixed_ndtfast_ratio > 0)
910 {
912 {
913 amrex::Abort("Dt is over-specfied");
914 }
915 }
916 // If two are specified, initialize fixed_ndtfast_ratio
917 else if (fixed_dt > 0. && fixed_fast_dt > 0. && fixed_ndtfast_ratio <= 0)
918 {
919 fixed_ndtfast_ratio = static_cast<int>(fixed_dt / fixed_fast_dt);
920 }
921
922 AMREX_ASSERT(cfl > 0. || fixed_dt > 0.);
923
924 pp_amr.query("do_substep", do_substep);
925 if (do_substep) {
926 amrex::Abort("Time substepping is not yet implemented. amr.do_substep must be 0");
927 }
928
929 // We use this to keep track of how many boxes we read in from WRF initialization
930 num_files_at_level.resize(max_level+1,0);
931
932 // We use this to keep track of how many boxes are specified thru the refinement indicators
933 num_boxes_at_level.resize(max_level+1,0);
934 boxes_at_level.resize(max_level+1);
935
936 // We always have exactly one file at level 0
937 num_boxes_at_level[0] = 1;
938 boxes_at_level[0].resize(1);
939 boxes_at_level[0][0] = geom[0].Domain();
940
941 // Plotfile name and frequency
942 pp.query("plot_file", plot_file_name);
943 pp.query("plot_int", plot_int);
944 pp.query("plot_int_time", plot_int_time);
945 // Output format
946 std::string plotfile_type_str = "amrex";
947 pp.query("plotfile_type", plotfile_type_str);
948 if (plotfile_type_str == "amrex") {
950 } else if (plotfile_type_str == "netcdf" || plotfile_type_str == "NetCDF") {
952#ifdef REMORA_USE_NETCDF
953 pp.query("write_history_file",write_history_file);
954 pp.query("chunk_history_file",chunk_history_file);
955 pp.query("steps_per_history_file",steps_per_history_file);
956 // Estimate size of domain for one timestep of netcdf
957 auto dom = geom[0].Domain();
958 int nx = dom.length(0) + 2;
959 int ny = dom.length(1) + 2;
960 int nz = dom.length(2);
962 // Estimate number of steps that will fit into a 2GB file.
963 steps_per_history_file = int((1.6e10 - NCH2D * nx * ny * 64.0_rt)
964 / (nx * ny * 64.0_rt * (NC3D*nz + NC2D)));
965 // If we calculate that a single step will exceed 2GB and the user has
966 // requested automatic history file sizing, warn about a possible impending
967 // error, and set steps_per_history_file = 1 to attempt output anyway.
968 if (steps_per_history_file == 0) {
969 amrex::Warning("NetCDF output for a single timestep appears to exceed 2GB. NetCDF output may not work. See Documentation for information about tested MPICH versions.");
971 }
972 } else if (write_history_file and !chunk_history_file) {
973 // Estimate number of output steps we'll need
974 int nt_out = int((max_step) / plot_int) + 1;
975 Real est_hist_file_size = NCH2D * nx * ny * 64.0_rt + nt_out * nx * ny * 64.0_rt * (NC3D*nz + NC2D);
976 if (est_hist_file_size > 1.6e10) {
977 amrex::Warning("WARNING: NetCDF history file may be larger than 2GB limit. Consider setting remora.chunk_history_file=true");
978 }
979 }
981 Print() << "NetCDF history files will have " << steps_per_history_file << " steps per file." << std::endl;
982 }
983#endif
984 } else {
985 amrex::Print() << "User selected plotfile_type = " << plotfile_type_str << std::endl;
986 amrex::Abort("Dont know this plotfile_type");
987 }
988#ifndef REMORA_USE_NETCDF
990 {
991 amrex::Abort("Please compile with NetCDF in order to enable NetCDF plotfiles");
992 }
993
994#endif
995
996#ifdef REMORA_USE_NETCDF
997 nc_init_file.resize(max_level+1);
998 nc_grid_file.resize(max_level+1);
999
1000 // NetCDF initialization files -- possibly multiple files at each of multiple levels
1001 // but we always have exactly one file at level 0
1002 for (int lev = 0; lev <= max_level; lev++)
1003 {
1004 const std::string nc_file_names = amrex::Concatenate("nc_init_file_",lev,1);
1005 const std::string nc_bathy_file_names = amrex::Concatenate("nc_grid_file_",lev,1);
1006
1007 if (pp.contains(nc_file_names.c_str()))
1008 {
1009 int num_files = pp.countval(nc_file_names.c_str());
1010 int num_bathy_files = pp.countval(nc_bathy_file_names.c_str());
1011 if (num_files != num_bathy_files) {
1012 amrex::Error("Must have same number of netcdf files for grid info as for solution");
1013 }
1014
1015 num_files_at_level[lev] = num_files;
1016 nc_init_file[lev].resize(num_files);
1017 nc_grid_file[lev].resize(num_files);
1018
1019 pp.queryarr(nc_file_names.c_str() , nc_init_file[lev] ,0,num_files);
1020 pp.queryarr(nc_bathy_file_names.c_str(), nc_grid_file[lev],0,num_files);
1021 }
1022 }
1023 // We only read boundary data at level 0
1024 pp.queryarr("nc_bdry_file", nc_bdry_file);
1025
1026 // Also only read forcings at level 0 (for now)
1027 pp.query("nc_frc_file", nc_frc_file);
1028
1029 // Get river file
1030 pp.query("nc_river_file", nc_riv_file);
1031
1032 // Read in file names for climatology history and nudging weights
1033 pp.query("nc_clim_his_file", nc_clim_his_file);
1034 pp.query("nc_clim_coeff_file", nc_clim_coeff_file);
1035
1036 pp.query("bdy_time_varname",bdry_time_varname);
1037 pp.query("frc_time_varname",frc_time_varname);
1038 pp.query("riv_time_varname",riv_time_varname);
1039 pp.query("clim_ubar_time_varname",clim_ubar_time_varname);
1040 pp.query("clim_vbar_time_varname",clim_vbar_time_varname);
1041 pp.query("clim_u_time_varname",clim_u_time_varname);
1042 pp.query("clim_v_time_varname",clim_v_time_varname);
1043 pp.query("clim_salt_time_varname",clim_salt_time_varname);
1044 pp.query("clim_temp_time_varname",clim_temp_time_varname);
1045
1046#endif
1047
1048#ifdef REMORA_USE_PARTICLES
1049 readTracersParams();
1050#endif
1051 }
1052
1054}
1055
1056void
1058{
1059 for (int lev = finest_level-1; lev >= 0; --lev)
1060 {
1061 AverageDownTo(lev);
1062 }
1063}
1064
1065/**
1066 * @param[in ] crse_lev level to average down to
1067 */
1068void
1070{
1071 average_down(*cons_new[crse_lev+1], *cons_new[crse_lev],
1072 0, cons_new[crse_lev]->nComp(), refRatio(crse_lev));
1073 average_down(*vec_zeta[crse_lev+1].get(), *vec_zeta[crse_lev].get(),
1074 0, vec_zeta[crse_lev]->nComp(), refRatio(crse_lev));
1075
1076 Array<MultiFab*,AMREX_SPACEDIM> faces_crse;
1077 Array<MultiFab*,AMREX_SPACEDIM> faces_fine;
1078 faces_crse[0] = xvel_new[crse_lev];
1079 faces_crse[1] = yvel_new[crse_lev];
1080 faces_crse[2] = zvel_new[crse_lev];
1081
1082 faces_fine[0] = xvel_new[crse_lev+1];
1083 faces_fine[1] = yvel_new[crse_lev+1];
1084 faces_fine[2] = zvel_new[crse_lev+1];
1085
1086 average_down_faces(GetArrOfConstPtrs(faces_fine), faces_crse,
1087 refRatio(crse_lev),geom[crse_lev]);
1088}
1089
1090/**
1091 * @param[in ] lev level at which to get time
1092 */
1093amrex::Real REMORA::get_t_old(int lev) const
1094{
1095 return t_old[lev];
1096}
amrex::Vector< std::string > BCNames
Definition REMORA.cpp:53
PlotfileType
plotfile format
#define Temp_comp
#define Scalar_comp
#define NC3D
#define Salt_comp
#define NC2D
#define NCH2D
#define NCONS
std::unique_ptr< ProblemBase > amrex_probinit(const amrex_real *problo, const amrex_real *probhi) AMREX_ATTRIBUTE_WEAK
Function to init the physical bounds of the domain and instantiate a Problem derived from ProblemBase...
A class to hold and interpolate time series data read from a NetCDF file.
void Initialize()
Read in time array from file and allocate data arrays.
void update_interpolated_to_time(amrex::Real time)
Calculate interpolated values at time, reading in data as necessary.
static PlotfileType plotfile_type
Native or NetCDF plotfile output.
Definition REMORA.H:1259
NCTimeSeriesRiver * river_source_transportbar
Data container for vertically integrated momentum transport in rivers.
Definition REMORA.H:1033
std::string riv_time_varname
Name of time field for river time.
Definition REMORA.H:1292
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:244
static amrex::Real fixed_dt
User specified fixed baroclinic time step.
Definition REMORA.H:1174
amrex::Real last_plot_file_time
Simulation time when we last output a plotfile.
Definition REMORA.H:1145
amrex::Vector< NCTimeSeriesRiver * > river_source_cons
Vector of data containers for scalar data in rivers.
Definition REMORA.H:1029
void init_bathymetry_from_netcdf(int lev)
Bathymetry data initialization from NetCDF file.
void init_bcs()
Read in boundary parameters from input file and set up data structures.
static amrex::Real previousCPUTimeUsed
Accumulator variable for CPU time used thusfar.
Definition REMORA.H:1359
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:360
void init_gls_vmix(int lev, SolverChoice solver_choice)
Initialize GLS variables.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_hOfTheConfusingName
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:229
amrex::Vector< REMORAFillPatcher > FPr_v
Vector over levels of FillPatchers for v (3D)
Definition REMORA.H:1055
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:217
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
Definition REMORA.H:1004
void stretch_transform(int lev)
Calculate vertical stretched coordinates.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vwind
Wind in the v direction, defined at rho-points.
Definition REMORA.H:289
std::unique_ptr< ProblemBase > prob
Pointer to container of analytical functions for problem definition.
Definition REMORA.H:1084
void Construct_REMORAFillPatchers(int lev)
Construct FillPatchers.
Definition REMORA.cpp:364
static int sum_interval
Diagnostic sum output interval in number of steps.
Definition REMORA.H:1251
int history_count
Counter for which time index we are writing to in the netcdf history file.
Definition REMORA.H:1210
amrex::Real stop_time
Time to stop.
Definition REMORA.H:1160
int do_substep
Whether to substep fine levels in time.
Definition REMORA.H:1183
NCTimeSeries * vbar_clim_data_from_file
Data container for vbar climatology data read from file.
Definition REMORA.H:1018
void Evolve()
Advance solution to final time.
Definition REMORA.cpp:140
NCTimeSeries * Uwind_data_from_file
Data container for u-direction wind read from file.
Definition REMORA.H:1011
std::string bdry_time_varname
Name of time field for boundary data.
Definition REMORA.H:1278
void ReadCheckpointFile()
read checkpoint file from disk
bool chunk_history_file
Whether to chunk netcdf history file.
Definition REMORA.H:1203
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:282
amrex::Real get_t_old(int lev) const
Accessor method for t_old to expose to outside classes.
Definition REMORA.cpp:1093
NCTimeSeries * v_clim_data_from_file
Data container for v-velocity climatology data read from file.
Definition REMORA.H:1022
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:242
amrex::Vector< std::string > datalogname
Definition REMORA.H:1391
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:223
void init_only(int lev, amrex::Real time)
Init (NOT restart or regrid)
Definition REMORA.cpp:674
void init_set_vmix(int lev)
Initialize vertical mixing coefficients from file or analytic.
Definition REMORA.cpp:583
std::string clim_u_time_varname
Name of time field for u climatology data.
Definition REMORA.H:1284
void set_grid_scale(int lev)
Set pm and pn arrays on level lev. Only works if using analytic initialization.
void set_coriolis(int lev)
Initialize Coriolis factor from file or analytic.
Definition REMORA.cpp:560
amrex::Vector< REMORAFillPatcher > FPr_u
Vector over levels of FillPatchers for u (3D)
Definition REMORA.H:1053
std::string clim_temp_time_varname
Name of time field for temperature climatology data.
Definition REMORA.H:1290
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1091
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Holds info for dynamically generated tagging criteria.
Definition REMORA.H:1301
std::string nc_frc_file
NetCDF forcing file.
Definition REMORA.H:1269
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
Definition REMORA.H:252
std::string clim_ubar_time_varname
Name of time field for ubar climatology data.
Definition REMORA.H:1280
void WriteNCPlotFile(int istep)
Write plotfile using NetCDF (wrapper)
std::string check_file
Checkpoint file prefix.
Definition REMORA.H:1196
static amrex::Real startCPUTime
Variable for CPU timing.
Definition REMORA.H:1357
void setPlotVariables(const std::string &pp_plot_var_names)
NCTimeSeries * u_clim_data_from_file
Data container for u-velocity climatology data read from file.
Definition REMORA.H:1020
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:210
amrex::Real start_time
Time of the start of the simulation, in seconds.
Definition REMORA.H:1163
NCTimeSeries * ubar_clim_data_from_file
Data container for ubar climatology data read from file.
Definition REMORA.H:1016
void init_data_from_netcdf(int lev)
Problem initialization from NetCDF file.
void init_masks_from_netcdf(int lev)
Mask data initialization from NetCDF file.
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:221
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_uwind
Wind in the u direction, defined at rho-points.
Definition REMORA.H:287
static amrex::Real fixed_fast_dt
User specified fixed barotropic time step.
Definition REMORA.H:1176
int regrid_int
how often each level regrids the higher levels of refinement (after a level advances that many time s...
Definition REMORA.H:1186
amrex::Real check_int_time
Checkpoint output interval in seconds.
Definition REMORA.H:1200
void init_bdry_from_netcdf()
Boundary data initialization from NetCDF file.
NCTimeSeries * sustr_data_from_file
Data container for u-component surface momentum flux read from file.
Definition REMORA.H:1007
void Define_REMORAFillPatchers(int lev)
Define FillPatchers.
Definition REMORA.cpp:412
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_p
Harmonic viscosity defined on the psi points (corners of horizontal grid cells)
Definition REMORA.H:254
amrex::Real plot_int_time
Plotfile output interval in seconds.
Definition REMORA.H:1194
amrex::Vector< int > num_files_at_level
how many netcdf input files specified at each level
Definition REMORA.H:1089
amrex::Vector< REMORAFillPatcher > FPr_vbar
Vector over levels of FillPatchers for vbar (2D)
Definition REMORA.H:1063
void AverageDownTo(int crse_lev)
more flexible version of AverageDown() that lets you average down across multiple levels
Definition REMORA.cpp:1069
int steps_per_history_file
Number of time steps per netcdf history file.
Definition REMORA.H:1208
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Called after every level 0 timestep.
Definition REMORA.cpp:226
int max_step
maximum number of steps
Definition REMORA.H:1158
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< int > num_boxes_at_level
how many boxes specified at each level by tagging criteria
Definition REMORA.H:1087
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:219
void refinement_criteria_setup()
Set refinement criteria.
int last_check_file_step
Step when we last output a checkpoint file.
Definition REMORA.H:1148
void init_beta_plane_coriolis(int lev)
Calculate Coriolis parameters from beta plane parametrization.
std::string clim_vbar_time_varname
Name of time field for vbar climatology data.
Definition REMORA.H:1282
amrex::Vector< int > nsubsteps
How many substeps on each level?
Definition REMORA.H:1096
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
Definition REMORA.H:1108
void ComputeDt()
a wrapper for estTimeStep()
int plot_int
Plotfile output interval in iterations.
Definition REMORA.H:1192
void FillPatch(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true, const int n_not_fill=0, const int icomp_calc=0, const amrex::Real dt=amrex::Real(0.0), const amrex::MultiFab &mf_calc=amrex::MultiFab())
Fill a new MultiFab by copying in phi from valid region and filling ghost cells.
void InitData()
Initialize multilevel data.
Definition REMORA.cpp:258
amrex::Vector< int > istep
which step?
Definition REMORA.H:1094
void WriteCheckpointFile()
write checkpoint file to disk
std::string nc_clim_coeff_file
NetCDF climatology coefficient file.
Definition REMORA.H:1275
void setRecordDataInfo(int i, const std::string &filename)
Definition REMORA.H:1377
void set_analytic_vmix(int lev)
Set vertical mixing coefficients from analytic.
Definition REMORA.cpp:599
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
Definition REMORA.H:1256
void init_riv_pos_from_netcdf(int lev)
static amrex::Vector< std::string > nc_bdry_file
NetCDF boundary data.
Definition REMORA.H:48
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_r
Harmonic viscosity defined on the rho points (centers)
Definition REMORA.H:256
std::string nc_riv_file
Definition REMORA.H:1270
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:284
REMORA()
Definition REMORA.cpp:61
void set_zeta(int lev)
Initialize zeta from file or analytic.
Definition REMORA.cpp:470
static amrex::Real change_max
Fraction maximum change in subsequent time steps.
Definition REMORA.H:1172
void init_zeta_from_netcdf(int lev)
Sea-surface height data initialization from NetCDF file.
void set_zeta_average(int lev)
Set Zt_avg1 to zeta.
void init_coriolis_from_netcdf(int lev)
Coriolis parameter data initialization from NetCDF file.
std::string pp_prefix
default prefix for input file parameters
Definition REMORA.H:205
void set_bathymetry(int lev)
Initialize bathymetry from file or analytic.
Definition REMORA.cpp:492
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:212
amrex::Vector< REMORAFillPatcher > FPr_c
Vector over levels of FillPatchers for scalars.
Definition REMORA.H:1051
std::string clim_v_time_varname
Name of time field for v climatology data.
Definition REMORA.H:1286
amrex::Vector< REMORAFillPatcher > FPr_w
Vector over levels of FillPatchers for w.
Definition REMORA.H:1057
amrex::Vector< amrex::YAFluxRegister * > advflux_reg
array of flux registers for refluxing in multilevel
Definition REMORA.H:1111
std::string clim_salt_time_varname
Name of time field for salinity climatology data.
Definition REMORA.H:1288
NCTimeSeries * Vwind_data_from_file
Data container for v-direction wind read from file.
Definition REMORA.H:1013
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1098
NCTimeSeriesRiver * river_source_transport
Data container for momentum transport in rivers.
Definition REMORA.H:1031
NCTimeSeries * temp_clim_data_from_file
Data container for temperature climatology data read from file.
Definition REMORA.H:1024
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1221
int cf_set_width
Width for fixing values at coarse-fine interface.
Definition REMORA.H:1048
NCTimeSeries * salt_clim_data_from_file
Data container for salinity climatology data read from file.
Definition REMORA.H:1026
void ReadParameters()
read in some parameters from inputs file
Definition REMORA.cpp:847
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:238
NCTimeSeries * svstr_data_from_file
Data container for v-component surface momentum flux read from file.
Definition REMORA.H:1009
void sum_integrated_quantities(amrex::Real time)
Integrate conserved quantities for diagnostics.
static int total_nc_plot_file_step
Definition REMORA.H:905
static amrex::Vector< amrex::Vector< std::string > > nc_grid_file
NetCDF grid file.
Definition REMORA.H:50
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
void set_zeta_to_Ztavg(int lev)
Set zeta components to be equal to time-averaged Zt_avg1.
int plot_file_on_restart
Whether to output a plotfile on restart from checkpoint.
Definition REMORA.H:1152
void set_2darrays(int lev)
Set 2D momentum arrays from 3D momentum.
void init_analytic(int lev)
Initialize initial problem data from analytic functions.
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
std::string frc_time_varname
Name of time field for forcing data.
Definition REMORA.H:1294
void set_wind(int lev)
Initialize or calculate wind speed from file or analytic.
Definition REMORA.cpp:655
amrex::Vector< REMORAFillPatcher > FPr_ubar
Vector over levels of FillPatchers for ubar (2D)
Definition REMORA.H:1061
bool is_it_time_for_action(int nstep, amrex::Real time, amrex::Real dt, int action_interval, amrex::Real action_per)
Decide if it is time to take an action.
void appendPlotVariables(const std::string &pp_plot_var_names)
void set_hmixcoef(int lev)
Initialize horizontal mixing coefficients.
Definition REMORA.cpp:612
std::string nc_clim_his_file
NetCDF climatology history file.
Definition REMORA.H:1273
void timeStep(int lev, amrex::Real time, int iteration)
advance a level by dt, includes a recursive call for finer levels
void FillCoarsePatch(int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse, const int icomp=0, const bool fill_all=true)
fill an entire multifab by interpolating from the coarser level
void AverageDown()
set covered coarse cells to be the average of overlying fine cells
Definition REMORA.cpp:1057
static int fixed_ndtfast_ratio
User specified, number of barotropic steps per baroclinic step.
Definition REMORA.H:1178
void timeStepML(amrex::Real time, int iteration)
advance all levels by dt, loops over finer levels
amrex::Vector< std::unique_ptr< std::fstream > > datalog
Definition REMORA.H:1390
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
Definition REMORA.H:250
static amrex::Real cfl
CFL condition.
Definition REMORA.H:1170
static int verbose
Verbosity level of output.
Definition REMORA.H:1248
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1190
int check_int
Checkpoint output interval in iterations.
Definition REMORA.H:1198
void WritePlotFile()
main driver for writing AMReX plotfiles
void set_smflux(int lev)
Initialize or calculate surface momentum flux from file or analytic.
Definition REMORA.cpp:637
std::string restart_chkfile
If set, restart from this checkpoint file.
Definition REMORA.H:1166
void init_clim_nudg_coeff(int lev)
Wrapper to initialize climatology nudging coefficient.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:240
int cf_width
Nudging width at coarse-fine interface.
Definition REMORA.H:1046
static amrex::Vector< amrex::Vector< std::string > > nc_init_file
NetCDF initialization file.
Definition REMORA.H:49
int last_plot_file_step
Step when we last output a plotfile.
Definition REMORA.H:1143
amrex::Vector< amrex::Real > t_old
old time at each level
Definition REMORA.H:1100
amrex::Real last_check_file_time
Simulation time when we last output a checkpoint file.
Definition REMORA.H:1150
amrex::Vector< amrex::Real > dt
time step at each level
Definition REMORA.H:1102
static amrex::Real sum_per
Diagnostic sum output interval in time.
Definition REMORA.H:1253
virtual ~REMORA()
Definition REMORA.cpp:135
void restart()
Definition REMORA.cpp:458
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
Definition REMORA.H:258
const char * buildInfoGetGitHash(int i)
HorizMixingType horiz_mixing_type
IC_BC_Type ic_bc_type
amrex::Vector< amrex::Real > tnu2
amrex::Vector< int > do_rivers_cons
void init_params()
read in and initialize parameters
amrex::Real visc2
SMFluxType smflux_type
VertMixingType vert_mixing_type
CouplingType coupling_type