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// Do we include staggered velocities in the plotfile?
39
40// Native AMReX vs NetCDF
42
43#ifdef REMORA_USE_NETCDF
44
46
47// Do we write one file per timestep (false) or one file for all timesteps (true)
49
50// NetCDF initialization file
51amrex::Vector<std::string> REMORA::nc_bdry_file = {""}; // Must provide via input
52amrex::Vector<amrex::Vector<std::string>> REMORA::nc_init_file = {{""}}; // Must provide via input
53amrex::Vector<amrex::Vector<std::string>> REMORA::nc_grid_file = {{""}}; // Must provide via input
54#endif
55
56amrex::Vector<std::string> BCNames = {"xlo", "ylo", "zlo", "xhi", "yhi", "zhi"};
57
58/**
59 * constructor:
60 * - reads in parameters from inputs file
61 * - sizes multilevel arrays and data structures
62 * - initializes BCRec boundary condition object
63 */
65{
66 BL_PROFILE("REMORA::REMORA()");
67 if (ParallelDescriptor::IOProcessor()) {
68 const char* remora_hash = amrex::buildInfoGetGitHash(1);
69 const char* amrex_hash = amrex::buildInfoGetGitHash(2);
70 const char* buildgithash = amrex::buildInfoGetBuildGitHash();
71 const char* buildgitname = amrex::buildInfoGetBuildGitName();
72
73 if (strlen(remora_hash) > 0) {
74 amrex::Print() << "\n"
75 << "REMORA git hash: " << remora_hash << "\n";
76 }
77 if (strlen(amrex_hash) > 0) {
78 amrex::Print() << "AMReX git hash: " << amrex_hash << "\n";
79 }
80 if (strlen(buildgithash) > 0) {
81 amrex::Print() << buildgitname << " git hash: " << buildgithash << "\n";
82 }
83
84 amrex::Print() << "\n";
85 }
86
88
89 const std::string& pv3d = "plot_vars_3d"; set3DPlotVariables(pv3d);
90 const std::string& pv2d = "plot_vars_2d"; set2DPlotVariables(pv2d);
91
92 prob = amrex_probinit(geom[0].ProbLo(),geom[0].ProbHi());
93
94 // Geometry on all levels has been defined already.
95
96 // No valid BoxArray and DistributionMapping have been defined.
97 // But the arrays for them have been resized.
98
99 int nlevs_max = max_level + 1;
100
101 istep.resize(nlevs_max, 0);
102 nsubsteps.resize(nlevs_max, 1);
103 for (int lev = 1; lev <= max_level; ++lev) {
104 nsubsteps[lev] = do_substep ? MaxRefRatio(lev-1) : 1;
105 }
106
107 physbcs.resize(nlevs_max);
108
109 t_new.resize(nlevs_max, 0.0_rt);
110 t_old.resize(nlevs_max, -1.e100_rt);
111 dt.resize(nlevs_max, 1.e100_rt);
112
113 cons_new.resize(nlevs_max);
114 cons_old.resize(nlevs_max);
115 xvel_new.resize(nlevs_max);
116 xvel_old.resize(nlevs_max);
117 yvel_new.resize(nlevs_max);
118 yvel_old.resize(nlevs_max);
119 zvel_new.resize(nlevs_max);
120 zvel_old.resize(nlevs_max);
121
122 advflux_reg.resize(nlevs_max);
123
124 // Initialize tagging criteria for mesh refinement
126
127 // We have already read in the ref_Ratio (via amr.ref_ratio =) but we need to enforce
128 // that there is no refinement in the vertical so we test on that here.
129 for (int lev = 0; lev < max_level; ++lev)
130 {
131 amrex::Print() << "Refinement ratio at level " << lev << " set to be " <<
132 ref_ratio[lev][0] << " " << ref_ratio[lev][1] << " " << ref_ratio[lev][2] << std::endl;
133
134 if (ref_ratio[lev][2] != 1)
135 {
136 amrex::Print() << "********************************************************************************" << std::endl;
137 amrex::Print() << "We don't allow refinement in the vertical -- make sure to set ref_ratio = 1 in z" << std::endl;
138 amrex::Print() << "It's possible you set amr.ref_ratio when you meant to set amr.ref_ratio_vect " << std::endl;
139 amrex::Print() << "********************************************************************************" << std::endl;
140 amrex::Abort();
141 }
142 }
143}
144
146{
147}
148
149void
151{
152 BL_PROFILE("REMORA::Evolve()");
153 Real cur_time = t_new[0];
154
155 // Take one coarse timestep by calling timeStep -- which recursively calls timeStep
156 // for finer levels (with or without subcycling)
157 for (int step = istep[0]; step < max_step && cur_time < stop_time; ++step)
158 {
159 amrex::Print() << "\nCoarse STEP " << step+1 << " starts ..." << std::endl;
160
161 ComputeDt();
162
163 int lev = 0;
164 int iteration = 1;
165
166 if (max_level == 0) {
167 timeStep(lev, cur_time, iteration);
168 }
169 else {
170 timeStepML(cur_time, iteration);
171 }
172
173 cur_time += dt[0];
174
175 amrex::Print() << "Coarse STEP " << step+1 << " ends." << " TIME = " << cur_time
176 << " DT = " << dt[0] << std::endl;
177
178 if ( (plot_int > 0 && (step+1 - last_plot_file_step) == plot_int ) ||
179 (plot_int_time > 0 && (cur_time >= (last_plot_file_time + plot_int_time))) )
180 {
181 last_plot_file_step = step+1;
182 last_plot_file_time = cur_time;
183 WritePlotFile(step+1);
185 }
186
187 if ((check_int > 0 && (step+1 - last_check_file_step) == check_int)
188 || (check_int_time > 0 && cur_time >= (last_check_file_time + check_int_time))) {
189 last_check_file_step = step+1;
190 last_check_file_time = cur_time;
192 }
193
194 post_timestep(step, cur_time, dt[0]);
195
196#ifdef AMREX_MEM_PROFILING
197 {
198 std::ostringstream ss;
199 ss << "[STEP " << step+1 << "]";
200 MemProfiler::report(ss.str());
201 }
202#endif
203
204 if (cur_time >= stop_time - 1.e-6*dt[0]) break;
205 }
206
207 if ( (plot_int > 0 || plot_int_time > 0.0) && istep[0] > last_plot_file_step)
208 {
211 }
212
213 if ((check_int > 0 || check_int_time > 0.0) && istep[0] > last_check_file_step) {
215 }
216}
217
218/**
219 * @param[in ] nstep which step we're on
220 * @param[in ] time current time
221 * @param[in ] dt_lev0 time step on level 0
222 */
223void
224REMORA::post_timestep (int nstep, Real time, Real dt_lev0)
225{
226 BL_PROFILE("REMORA::post_timestep()");
227
228#ifdef REMORA_USE_PARTICLES
229 particleData.Redistribute();
230#endif
231
233 {
234 for (int lev = finest_level-1; lev >= 0; lev--)
235 {
236 // This call refluxes from the lev/lev+1 interface onto lev
237 //getAdvFluxReg(lev+1)->Reflux(*cons_new[lev], 0, 0, NCONS);
238
239 // We need to do this before anything else because refluxing changes the
240 // values of coarse cells underneath fine grids with the assumption they'll
241 // be over-written by averaging down
242 //
243 AverageDownTo(lev);
244 }
245 }
246
247 if (is_it_time_for_action(nstep, time, dt_lev0, sum_interval, sum_per)) {
249 }
250}
251
252/**
253 * This is called from main.cpp and handles all initialization, whether from start or restart
254 */
255void
257{
258 BL_PROFILE("REMORA::InitData()");
259 // Initialize the start time for our CPU-time tracker
260 startCPUTime = Real(ParallelDescriptor::second());
261
262 // Map the words in the inputs file to BC types, then translate
263 // those types into what they mean for each variable
264 init_bcs();
265
266 // Init vertical stretching coeffs
268
271 last_plot_file_time = -1.0_rt;
272 last_check_file_time = -1.0_rt;
273
274 if (restart_chkfile == "") {
275 // start simulation from the beginning
276
277 InitFromScratch(start_time);
278
280 AverageDown();
281 }
282
283 } else { // Restart from a checkpoint
284
285 restart();
286
287 }
288#ifdef REMORA_USE_MOAB
289 InitMOABMesh();
290#endif
291 // Initialize flux registers (whether we start from scratch or restart)
293 advflux_reg[0] = nullptr;
294 for (int lev = 1; lev <= finest_level; lev++)
295 {
296 advflux_reg[lev].reset( new YAFluxRegister(grids[lev], grids[lev-1],
297 dmap[lev], dmap[lev-1],
298 geom[lev], geom[lev-1],
299 ref_ratio[lev-1], lev, NCONS));
300 }
301 }
302
303 // Fill ghost cells/faces
304 for (int lev = 0; lev <= finest_level; ++lev)
305 {
306 if (lev > 0 && cf_width >= 0) {
308 }
309
310 if (restart_chkfile == "") {
311 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]);
312 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]);
313 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]);
314 FillPatch(lev, t_new[lev], *zvel_new[lev], zvel_new, BCVars::zvel_bc, BdyVars::null, 0, true, false);
315
316 // Copy from new into old just in case when initializing from scratch
317 int ngs = cons_new[lev]->nGrow();
318 int ngvel = xvel_new[lev]->nGrow();
319 MultiFab::Copy(*cons_old[lev],*cons_new[lev],0,0,NCONS,ngs);
320 MultiFab::Copy(*xvel_old[lev],*xvel_new[lev],0,0,1,ngvel);
321 MultiFab::Copy(*yvel_old[lev],*yvel_new[lev],0,0,1,ngvel);
322 MultiFab::Copy(*zvel_old[lev],*zvel_new[lev],0,0,1,IntVect(ngvel,ngvel,0));
323 }
324 } // lev
325
326 // Check for additional plotting variables that are available after
327 // particle containers are setup.
328 const std::string& pv3d = "plot_vars_3d"; append3DPlotVariables(pv3d);
329 const std::string& pv2d = "plot_vars_2d"; append2DPlotVariables(pv2d);
330
331 if (restart_chkfile == "" && (check_int > 0 || check_int_time > 0.0_rt))
332 {
335 }
336
337 if ( (restart_chkfile == "") ||
339 {
340 if (plot_int > 0 || plot_int_time > 0.0)
341 {
342 int step0 = 0;
343 WritePlotFile(step0);
346 }
347 }
348
351 }
352
353 ComputeDt();
354
355}
356
357/**
358 * @param[in ] lev level to operate on
359 */
360void
362{
363 BL_PROFILE("REMORA::Construct_REMORAFillPatchers()");
364 amrex::Print() << ":::Construct_REMORAFillPatchers " << lev << std::endl;
365
366 auto& ba_fine = cons_new[lev ]->boxArray();
367 auto& ba_crse = cons_new[lev-1]->boxArray();
368 auto& dm_fine = cons_new[lev ]->DistributionMap();
369 auto& dm_crse = cons_new[lev-1]->DistributionMap();
370
371 BoxList bl2d_fine = ba_fine.boxList();
372 for (auto& b : bl2d_fine) {
373 b.setRange(2,0);
374 }
375 BoxArray ba2d_fine(std::move(bl2d_fine));
376
377 BoxList bl2d_crse = ba_crse.boxList();
378 for (auto& b : bl2d_crse) {
379 b.setRange(2,0);
380 }
381 BoxArray ba2d_crse(std::move(bl2d_crse));
382
383 int ncomp = cons_new[lev]->nComp();
384
385 FPr_c.emplace_back(ba_fine, dm_fine, geom[lev] ,
386 ba_crse, dm_crse, geom[lev-1],
387 -cf_width, -cf_set_width, ncomp, &cell_cons_interp);
388 FPr_u.emplace_back(convert(ba_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
389 convert(ba_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
390 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
391 FPr_v.emplace_back(convert(ba_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
392 convert(ba_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
393 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
394 FPr_w.emplace_back(convert(ba_fine, IntVect(0,0,1)), dm_fine, geom[lev] ,
395 convert(ba_crse, IntVect(0,0,1)), dm_crse, geom[lev-1],
396 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
397
398 FPr_ubar.emplace_back(convert(ba2d_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
399 convert(ba2d_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
400 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
401 FPr_vbar.emplace_back(convert(ba2d_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
402 convert(ba2d_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
403 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
404}
405
406/**
407 * @param[in ] lev level to operate on
408 */
409void
411{
412 BL_PROFILE("REMORA::Define_REMORAFillPatchers()");
413 amrex::Print() << ":::Define_REMORAFillPatchers " << lev << std::endl;
414
415 auto& ba_fine = cons_new[lev ]->boxArray();
416 auto& ba_crse = cons_new[lev-1]->boxArray();
417 auto& dm_fine = cons_new[lev ]->DistributionMap();
418 auto& dm_crse = cons_new[lev-1]->DistributionMap();
419
420 BoxList bl2d_fine = ba_fine.boxList();
421 for (auto& b : bl2d_fine) {
422 b.setRange(2,0);
423 }
424 BoxArray ba2d_fine(std::move(bl2d_fine));
425
426 BoxList bl2d_crse = ba_crse.boxList();
427 for (auto& b : bl2d_crse) {
428 b.setRange(2,0);
429 }
430 BoxArray ba2d_crse(std::move(bl2d_crse));
431
432
433 int ncomp = cons_new[lev]->nComp();
434
435 FPr_c[lev-1].Define(ba_fine, dm_fine, geom[lev] ,
436 ba_crse, dm_crse, geom[lev-1],
437 -cf_width, -cf_set_width, ncomp, &cell_cons_interp);
438 FPr_u[lev-1].Define(convert(ba_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
439 convert(ba_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
440 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
441 FPr_v[lev-1].Define(convert(ba_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
442 convert(ba_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
443 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
444 FPr_w[lev-1].Define(convert(ba_fine, IntVect(0,0,1)), dm_fine, geom[lev] ,
445 convert(ba_crse, IntVect(0,0,1)), dm_crse, geom[lev-1],
446 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
447
448 FPr_ubar[lev-1].Define(convert(ba2d_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
449 convert(ba2d_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
450 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
451 FPr_vbar[lev-1].Define(convert(ba2d_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
452 convert(ba2d_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
453 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
454}
455
456void
458{
459 BL_PROFILE("REMORA::restart()");
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{
472 BL_PROFILE("REMORA::set_zeta()");
474 prob->init_analytic_zeta(lev, geom[lev], solverChoice, *this, *vec_zeta[lev]);
475
476#ifdef REMORA_USE_NETCDF
477 } else if (solverChoice.ic_type == IC_Type::netcdf) {
478 if (lev == 0) {
479 amrex::Print() << "Calling init_zeta_from_netcdf on level " << lev << std::endl;
481 amrex::Print() << "Sea surface height loaded from netcdf file \n " << std::endl;
482 } else {
483 Real dummy_time = 0.0_rt;
484 FillCoarsePatch(lev,dummy_time,vec_zeta[lev].get(), vec_zeta[lev-1].get(),BCVars::cons_bc);
485 }
486#endif
487 } else {
488 Abort("Don't know this ic_type!");
489 }
490 vec_zeta[lev]->FillBoundary(geom[lev].periodicity());
491 set_zeta_average(lev);
492}
493
494/**
495 * @param[in ] lev level to operate on
496 */
497void
499{
500 BL_PROFILE("REMORA::bathymetry()");
501 // Only set bathymetry on level 0, and interpolate for finer levels
503 if (lev==0) {
506 prob->init_analytic_bathymetry(lev, geom[lev], solverChoice, *this, *vec_h[lev]);
507 } else {
509 }
510#ifdef REMORA_USE_NETCDF
511 } else if (solverChoice.ic_type == IC_Type::netcdf) {
513 amrex::Print() << "Calling init_bathymetry_from_netcdf " << std::endl;
515 amrex::Print() << "Bathymetry loaded from netcdf file \n " << std::endl;
516 } else {
518 }
519#endif
520 } else {
521 Abort("Don't know this ic_type!");
522 }
523 // Need FillBoundary to fill at grid-grid boundaries, and EnforcePeriodicity
524 // to make sure ghost cells in the domain corners are consistent.
525 vec_h[lev]->FillBoundary(geom[lev].periodicity());
526 vec_h[lev]->EnforcePeriodicity(geom[lev].periodicity());
527 } else {
528 Real dummy_time = 0.0_rt;
529 FillCoarsePatch(lev,dummy_time,vec_h[lev].get(), vec_h[lev-1].get(),BCVars::cons_bc);
531 set_grid_scale(lev);
532 }
533 }
534 } else if (solverChoice.init_ana_h) {
537 prob->init_analytic_bathymetry(lev, geom[lev], solverChoice, *this,*vec_h[lev]);
538 } else {
540 }
541#ifdef REMORA_USE_NETCDF
542 } else if (solverChoice.ic_type == IC_Type::netcdf) {
543 amrex::Print() << "Calling init_bathymetry_from_netcdf level " << lev << std::endl;
545 amrex::Print() << "Bathymetry loaded from netcdf file \n " << std::endl;
546#endif
547 } else {
548 Abort("Don't know this ic_type!");
549 }
550 // Need FillBoundary to fill at grid-grid boundaries, and EnforcePeriodicity
551 // to make sure ghost cells in the domain corners are consistent.
552 vec_h[lev]->FillBoundary(geom[lev].periodicity());
553 vec_h[lev]->EnforcePeriodicity(geom[lev].periodicity());
554 } else if (solverChoice.init_l1ad_h) {
557 prob->init_analytic_bathymetry(lev, geom[lev], solverChoice, *this, *vec_h[lev]);
558 } else {
560 }
561#ifdef REMORA_USE_NETCDF
562 } else if (solverChoice.ic_type == IC_Type::netcdf) {
564 amrex::Print() << "Calling init_bathymetry_from_netcdf " << std::endl;
566 amrex::Print() << "Bathymetry loaded from netcdf file \n " << std::endl;
567 } else {
569 }
570#endif
571 } else {
572 Abort("Don't know this ic_type!");
573 }
574 // Need FillBoundary to fill at grid-grid boundaries, and EnforcePeriodicity
575 // to make sure ghost cells in the domain corners are consistent.
576 vec_h[lev]->FillBoundary(geom[lev].periodicity());
577 vec_h[lev]->EnforcePeriodicity(geom[lev].periodicity());
578 } else {
579 amrex::Abort("Don't know this h init type");
580 }
581}
582
583/**
584 * @param[in ] lev level to operate on
585 */
586void
588 BL_PROFILE("REMORA::set_coriolis()");
591 prob->init_analytic_coriolis(lev, geom[lev], solverChoice, *this, *vec_fcor[lev]);
594#ifdef REMORA_USE_NETCDF
596 if (lev == 0) {
597 amrex::Print() << "Calling init_coriolis_from_netcdf " << std::endl;
599 amrex::Print() << "Coriolis loaded from netcdf file \n" << std::endl;
600 } else {
601 Real dummy_time = 0.0_rt;
602 FillCoarsePatch(lev,dummy_time,vec_fcor[lev].get(), vec_fcor[lev-1].get(),BCVars::cons_bc);
603 }
604#endif
605 } else {
606 Abort("Don't know this coriolis_type!");
607 }
608
609 Real time = 0.0_rt;
610 FillPatch(lev, time, *vec_fcor[lev], GetVecOfPtrs(vec_fcor),BCVars::foextrap_bc);
611 vec_fcor[lev]->EnforcePeriodicity(geom[lev].periodicity());
612 }
613}
614
615void
617 BL_PROFILE("REMORA::init_set_vmix()");
622 // The GLS initialization just sets the multifab to a value, so there's
623 // no need to call FillPatch here
624 } else {
625 Abort("Don't know this vertical mixing type");
626 }
627}
628
629/**
630 * @param[in ] lev level to operate on
631 */
632void
634 BL_PROFILE("REMORA::set_analytic_vmix()");
635 Real time = 0.0_rt;
636 prob->init_analytic_vmix(lev, geom[lev], solverChoice, *this,*vec_Akv[lev], *vec_Akt[lev]);
637 FillPatch(lev, time, *vec_Akv[lev], GetVecOfPtrs(vec_Akv),BCVars::zvel_bc,BdyVars::null,0,true,false);
638 for (int n=0; n<NCONS;n++) {
639 FillPatch(lev, time, *vec_Akt[lev], GetVecOfPtrs(vec_Akt),BCVars::zvel_bc,BdyVars::null,0,false,false);
640 }
641}
642
643/**
644 * @param[in ] lev level to operate on
645 */
646void
648{
649 BL_PROFILE("REMORA::set_hmixcoef()");
651 prob->init_analytic_hmix(lev, geom[lev], solverChoice, *this, *vec_visc2_p[lev], *vec_visc2_r[lev], *vec_diff2[lev]);
653 vec_visc2_p[lev]->setVal(solverChoice.visc2);
654 vec_visc2_r[lev]->setVal(solverChoice.visc2);
655 for (int n=0; n<NCONS; n++) {
656 vec_diff2[lev]->setVal(solverChoice.tnu2[n],n,1);
657 }
658 } else {
659 Abort("Don't know this horizontal mixing type");
660 }
661 Real time = 0.0_rt;
662 FillPatch(lev, time, *vec_visc2_p[lev], GetVecOfPtrs(vec_visc2_p),BCVars::foextrap_periodic_bc);
663 FillPatch(lev, time, *vec_visc2_r[lev], GetVecOfPtrs(vec_visc2_r),BCVars::foextrap_periodic_bc);
664 for (int n=0; n<NCONS; n++) {
665 FillPatch(lev, time, *vec_diff2[lev] , GetVecOfPtrs(vec_diff2),BCVars::foextrap_periodic_bc,BdyVars::null,n,false);
666 }
667}
668
669/**
670 * @param[in ] lev level to operate on
671 */
672void
674{
675 BL_PROFILE("REMORA::init_flat_bathymetry()");
676 vec_h[lev]->setVal(-geom[0].ProbLo()[2]);
677}
678
679/**
680 * @param[in ] lev level to operate on
681 */
682void
684{
685 BL_PROFILE("REMORA::set_smflux()");
687 prob->init_analytic_smflux(lev, geom[lev], solverChoice, *this,*vec_sustr[lev], *vec_svstr[lev]);
689#ifdef REMORA_USE_NETCDF
690 sustr_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_sustr[lev].get(), geom, ref_ratio);
691 svstr_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_svstr[lev].get(), geom, ref_ratio);
692 FillPatch(lev, t_old[lev], *vec_sustr[lev], GetVecOfPtrs(vec_sustr),BCVars::foextrap_periodic_bc,BdyVars::null,0,false);
693 FillPatch(lev, t_old[lev], *vec_svstr[lev], GetVecOfPtrs(vec_svstr),BCVars::foextrap_periodic_bc,BdyVars::null,0,false);
694#endif
695 }
696}
697
698/**
699 * @param[in ] lev level to operate on
700 */
701void
703{
704 BL_PROFILE("REMORA::set_wind()");
706 prob->init_analytic_wind(lev,geom[lev], solverChoice, *this, *vec_uwind[lev], *vec_vwind[lev]);
708#ifdef REMORA_USE_NETCDF
709 Uwind_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_uwind[lev].get(), geom, ref_ratio);
710 Vwind_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_vwind[lev].get(), geom, ref_ratio);
711 FillPatch(lev, t_old[lev], *vec_uwind[lev], GetVecOfPtrs(vec_uwind),
713 FillPatch(lev, t_old[lev], *vec_vwind[lev], GetVecOfPtrs(vec_vwind),
715
716 // Conditionally update atmospheric fields if loaded from NetCDF
718 Tair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_Tair[lev].get(), geom, ref_ratio);
719 FillPatch(lev, t_old[lev], *vec_Tair[lev], GetVecOfPtrs(vec_Tair),
721 }
723 qair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_qair[lev].get(), geom, ref_ratio);
724 FillPatch(lev, t_old[lev], *vec_qair[lev], GetVecOfPtrs(vec_qair),
726
727 // Convert qair from percentage (0-100) to specific humidity (0-1) if needed
729 vec_qair[lev]->mult(0.01);
730
731 // Update ghost cells after modification
732 vec_qair[lev]->FillBoundary(geom[lev].periodicity());
733 }
734 }
736 Pair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_Pair[lev].get(), geom, ref_ratio);
737 FillPatch(lev, t_old[lev], *vec_Pair[lev], GetVecOfPtrs(vec_Pair),
739 }
741 srflx_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_srflx[lev].get(), geom, ref_ratio);
742 FillPatch(lev, t_old[lev], *vec_srflx[lev], GetVecOfPtrs(vec_srflx),
744 }
746 longwave_down_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_longwave_down[lev].get(), geom, ref_ratio);
747 FillPatch(lev, t_old[lev], *vec_longwave_down[lev], GetVecOfPtrs(vec_longwave_down),
749 }
751 rain_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_rain[lev].get(), geom, ref_ratio);
752 FillPatch(lev, t_old[lev], *vec_rain[lev], GetVecOfPtrs(vec_rain),
754 }
756 cloud_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_cloud[lev].get(), geom, ref_ratio);
757 FillPatch(lev, t_old[lev], *vec_cloud[lev], GetVecOfPtrs(vec_cloud),
759 }
761 EminusP_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_EminusP[lev].get(), geom, ref_ratio);
762 FillPatch(lev, t_old[lev], *vec_EminusP[lev], GetVecOfPtrs(vec_EminusP),
764 }
765#endif
766 }
767}
768/**
769 * @param[in ] lev level to operate on
770 */
771void
773{
775 prob->init_analytic_masks(lev,geom[lev], solverChoice, *this, *vec_mskr[lev]);
778#ifdef REMORA_USE_NETCDF
779 if (lev == 0) {
780 amrex::Print() << "Calling init_masks_from_netcdf level " << lev << std::endl;
782 amrex::Print() << "Masks loaded from netcdf file \n " << std::endl;
783 } else {
784 Real dummy_time = 0.0_rt;
785 FillCoarsePatchPC(lev, dummy_time, vec_mskr[lev].get(), vec_mskr[lev-1].get(),
788 }
789#endif
790 }
791 fill_3d_masks(lev);
792}
793
794/**
795 * @param[in ] lev level to operate on
796 * @param[in ] time current time for initialization
797 */
798void
799REMORA::init_only (int lev, Real time)
800{
801 BL_PROFILE("REMORA::init_only()");
802 t_new[lev] = time;
803 t_old[lev] = time - 1.e200_rt;
804
805 cons_new[lev]->setVal(0.0_rt);
806 xvel_new[lev]->setVal(0.0_rt);
807 yvel_new[lev]->setVal(0.0_rt);
808 zvel_new[lev]->setVal(0.0_rt);
809
810 xvel_old[lev]->setVal(0.0_rt);
811 yvel_old[lev]->setVal(0.0_rt);
812 zvel_old[lev]->setVal(0.0_rt);
813
814 vec_ru[lev]->setVal(0.0_rt);
815 vec_rv[lev]->setVal(0.0_rt);
816
817 vec_ru2d[lev]->setVal(0.0_rt);
818 vec_rv2d[lev]->setVal(0.0_rt);
819
821 set_grid_scale(lev);
822 }
823 set_masks(lev);
824
825#ifdef REMORA_USE_NETCDF
828
829 if (solverChoice.do_any_clim_nudg && lev == 0) {
830 if (nc_clim_his_file.empty()) {
831 amrex::Error("NetCDF climatology file name must be provided via input");
832 }
835 clim_ubar_time_varname, geom[lev].Domain(),vec_ubar[lev].get(),true,true));
837 clim_ubar_time_varname, geom[lev].Domain(),vec_vbar[lev].get(),true,true));
838 ubar_clim_data_from_file->Initialize();
839 vbar_clim_data_from_file->Initialize();
840 }
842 u_clim_data_from_file.reset(new NCTimeSeries(nc_clim_his_file, "u", clim_u_time_varname, geom[lev].Domain(),xvel_new[lev],false,true));
843 v_clim_data_from_file.reset(new NCTimeSeries(nc_clim_his_file, "v", clim_v_time_varname, geom[lev].Domain(),yvel_new[lev],false,true));
844 u_clim_data_from_file->Initialize();
845 v_clim_data_from_file->Initialize();
846 }
847 // Since the NCTimeSeries object isn't filling the cons_new MultiFab directly, we don't have to specify a component.
848 // It just needs to know the shape of the MultiFab
850 temp_clim_data_from_file.reset(new NCTimeSeries(nc_clim_his_file, "temp", clim_temp_time_varname,geom[lev].Domain(),cons_new[lev],false,true));
851 temp_clim_data_from_file->Initialize();
852 }
854 salt_clim_data_from_file.reset(new NCTimeSeries(nc_clim_his_file, "salt", clim_salt_time_varname,geom[lev].Domain(),cons_new[lev],false,true));
855 salt_clim_data_from_file->Initialize();
856 }
857 }
858 }
859
861 amrex::Print() << "Calling init_bdry_from_netcdf at level " << lev << std::endl;
863 amrex::Print() << "Boundary data loaded from netcdf file \n " << std::endl;
864 }
865
866 // This will be a non-op if forcings specified analytically
867 if (solverChoice.wind_type == WindType::netcdf && lev == 0) {
868 if (nc_frc_file.empty()) {
869 amrex::Error("NetCDF forcing file name must be provided via input for winds");
870 }
871 Uwind_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Uwind", frc_time_varname, geom[lev].Domain(),vec_uwind[lev].get(), true, false));
872 Vwind_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Vwind", frc_time_varname, geom[lev].Domain(),vec_vwind[lev].get(), true, false));
873 Uwind_data_from_file->Initialize();
874 Vwind_data_from_file->Initialize();
875
876 // Conditionally load atmospheric forcing fields from NetCDF based on user flags
878 Tair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Tair", frc_time_varname, geom[lev].Domain(),vec_Tair[lev].get(), true, false));
879 Tair_data_from_file->Initialize();
880 }
882 qair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "qair", frc_time_varname, geom[lev].Domain(),vec_qair[lev].get(), true, false));
883 qair_data_from_file->Initialize();
884 }
886 Pair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Pair", frc_time_varname, geom[lev].Domain(),vec_Pair[lev].get(), true, false));
887 Pair_data_from_file->Initialize();
888 }
890 srflx_data_from_file.reset(new NCTimeSeries(nc_frc_file, "swrad", frc_time_varname, geom[lev].Domain(),vec_srflx[lev].get(), true, false));
891 srflx_data_from_file->Initialize();
892 }
894 rain_data_from_file.reset(new NCTimeSeries(nc_frc_file, "rain", frc_time_varname, geom[lev].Domain(),vec_rain[lev].get(), true, false));
895 rain_data_from_file->Initialize();
896 }
898 cloud_data_from_file.reset(new NCTimeSeries(nc_frc_file, "cloud", frc_time_varname, geom[lev].Domain(),vec_cloud[lev].get(), true, false));
899 cloud_data_from_file->Initialize();
900 }
902 EminusP_data_from_file.reset(new NCTimeSeries(nc_frc_file, "EminusP", frc_time_varname, geom[lev].Domain(),vec_EminusP[lev].get(), true, false));
903 EminusP_data_from_file->Initialize();
904 }
905 } else if (solverChoice.smflux_type == SMFluxType::netcdf && lev == 0) {
906 if (nc_frc_file.empty()) {
907 amrex::Error("NetCDF forcing file name must be provided via input for surface momentum fluxes");
908 }
909 sustr_data_from_file.reset(new NCTimeSeries(nc_frc_file, "sustr", frc_time_varname, geom[lev].Domain(),vec_sustr[lev].get(), true, false));
910 svstr_data_from_file.reset(new NCTimeSeries(nc_frc_file, "svstr", frc_time_varname, geom[lev].Domain(),vec_svstr[lev].get(), true, false));
911 sustr_data_from_file->Initialize();
912 svstr_data_from_file->Initialize();
913 }
914 if (solverChoice.longwave_down_from_netcdf && lev == 0) {
915 if (nc_frc_file.empty()) {
916 amrex::Error("NetCDF forcing file name must be provided via input for longwave radiation");
917 }
919 geom[lev].Domain(), vec_longwave_down[lev].get(), true, false));
920 longwave_down_data_from_file->Initialize();
921 }
922
924 auto dom = geom[0].Domain();
925 int nz = dom.length(2);
926 river_source_cons.resize(NCONS);
927 Print() << solverChoice.do_rivers_cons[0] << std::endl;
930 river_source_cons[Salt_comp]->Initialize();
931 }
934 river_source_cons[Temp_comp]->Initialize();
935 }
938 river_source_cons[Tracer_comp]->Initialize();
939 }
940 river_source_transport.reset(new NCTimeSeriesRiver(nc_riv_file, "river_transport", riv_time_varname, nz));
941 river_source_transport->Initialize();
942 river_source_transportbar.reset(new NCTimeSeriesRiver(nc_riv_file, "river_transport", riv_time_varname, nz, 1));
943 river_source_transportbar->Initialize();
945 }
946#else
948 Abort("Not compiled with NetCDF, but selected boundary conditions require NetCDF");
949 }
951 Abort("Not compiled with NetCDF, but using river sources requires NetCDF");
952 }
953#endif
954
955 set_bathymetry(lev);
956 set_zeta(lev);
958
960 if (lev==0) {
962 {
963 init_analytic(lev);
964#ifdef REMORA_USE_NETCDF
965 } else if (solverChoice.ic_type == IC_Type::netcdf) {
966 amrex::Print() << "Calling init_data_from_netcdf " << std::endl;
969 amrex::Print() << "Initial data loaded from netcdf file \n " << std::endl;
970
971#endif
972 } else {
973 Abort("Need to specify ic_type");
974 }
975 } else {
980 }
983 {
984 init_analytic(lev);
985#ifdef REMORA_USE_NETCDF
986 } else if (solverChoice.ic_type == IC_Type::netcdf) {
987 amrex::Print() << "Calling init_data_from_netcdf on level " << lev << std::endl;
990 amrex::Print() << "Initial data loaded from netcdf file on level " << lev << std::endl;
991
992#endif
993 } else {
994 Abort("Need to specify ic_type");
995 }
996 } else {
997 amrex::Abort("Need to specify T init procedure");
998 }
999
1000 // Ensure that the face-based data are the same on both sides of a periodic domain.
1001 // The data associated with the lower grid ID is considered the correct value.
1002 xvel_new[lev]->OverrideSync(geom[lev].periodicity());
1003 yvel_new[lev]->OverrideSync(geom[lev].periodicity());
1004 zvel_new[lev]->OverrideSync(geom[lev].periodicity());
1005
1006 set_2darrays(lev);
1007
1008 init_set_vmix(lev);
1009 set_hmixcoef(lev);
1010 set_coriolis(lev);
1011
1012 // Previously set smflux here with OverrideSync:
1013// set_smflux(lev);
1014// prob->init_analytic_smflux(lev, geom[lev], solverChoice, *this, *vec_sustr[lev], *vec_svstr[lev]);
1015// vec_sustr[lev]->OverrideSync(geom[lev].periodicity());
1016// vec_svstr[lev]->OverrideSync(geom[lev].periodicity());
1017
1018}
1019
1020void
1022{
1023 BL_PROFILE("REMORA::ReadParameters()");
1024 {
1025 ParmParse pp; // Traditionally, max_step and stop_time do not have prefix, so allow it for now.
1026 bool noprefix_max_step = pp.queryAdd("max_step", max_step);
1027 bool noprefix_stop_time = pp.queryAdd("stop_time", stop_time);
1028 bool remora_max_step = pp.queryAdd("remora.max_step", max_step);
1029 bool remora_stop_time = pp.queryAdd("remora.stop_time", stop_time);
1030 if (remora_max_step and noprefix_max_step) {
1031 Abort("remora.max_step and max_step are both specified. Please use only one!");
1032 }
1033 if (remora_stop_time and noprefix_stop_time) {
1034 Abort("remora.stop_time and stop_time are both specified. Please use only one!");
1035 }
1036 }
1037
1038 ParmParse pp(pp_prefix);
1039 ParmParse pp_amr("amr");
1040 {
1041 pp_amr.queryAdd("regrid_int", regrid_int);
1042 pp.queryAdd("check_file", check_file);
1043 pp.queryAdd("check_int", check_int);
1044 pp_amr.queryAdd("check_int", check_int);
1045 pp.queryAdd("check_int_time", check_int_time);
1046 pp_amr.queryAdd("check_int_time", check_int_time);
1047
1048 pp.queryAdd("expand_plotvars_to_unif_rr", expand_plotvars_to_unif_rr);
1049
1050 pp.query("plotfile_fill_value", plotfile_fill_value);
1051 pp.query("netcdf_fill_value", netcdf_fill_value);
1052
1053 pp.queryAdd("restart", restart_chkfile);
1054 pp_amr.queryAdd("restart", restart_chkfile);
1055 pp.queryAdd("start_time",start_time);
1056
1057 if (pp.contains("data_log"))
1058 {
1059 int num_datalogs = pp.countval("data_log");
1060 datalog.resize(num_datalogs);
1061 datalogname.resize(num_datalogs);
1062 pp.queryarr("data_log",datalogname,0,num_datalogs);
1063 for (int i = 0; i < num_datalogs; i++)
1065 }
1066
1067 // Verbosity
1068 pp.queryAdd("v", verbose);
1069
1070 // Frequency of diagnostic output
1071 pp.queryAdd("sum_interval", sum_interval);
1072 pp.queryAdd("sum_period" , sum_per);
1073 pp.queryAdd("file_min_digits", file_min_digits);
1074
1075 if (file_min_digits < 0) {
1076 amrex::Abort("remora.file_min_digits must be non-negative");
1077 }
1078
1079 // Time step controls
1080 pp.queryAdd("cfl", cfl);
1081 pp.queryAdd("change_max", change_max);
1082
1083 pp.queryAdd("fixed_dt", fixed_dt);
1084 pp.queryAdd("fixed_fast_dt", fixed_fast_dt);
1085
1086 pp.queryAdd("fixed_ndtfast_ratio", fixed_ndtfast_ratio);
1087
1088 // If all three are specified, they must be consistent
1089 if (fixed_dt > 0. && fixed_fast_dt > 0. && fixed_ndtfast_ratio > 0)
1090 {
1092 {
1093 amrex::Abort("Dt is over-specfied");
1094 }
1095 }
1096 // If two are specified, initialize fixed_ndtfast_ratio
1097 else if (fixed_dt > 0. && fixed_fast_dt > 0. && fixed_ndtfast_ratio <= 0)
1098 {
1099 fixed_ndtfast_ratio = static_cast<int>(fixed_dt / fixed_fast_dt);
1100 }
1101
1102 AMREX_ASSERT(cfl > 0. || fixed_dt > 0.);
1103
1104 pp_amr.queryAdd("do_substep", do_substep);
1105 if (do_substep) {
1106 amrex::Abort("Time substepping is not yet implemented. amr.do_substep must be 0");
1107 }
1108
1109 // We use this to keep track of how many boxes we read in from WRF initialization
1110 num_files_at_level.resize(max_level+1,0);
1111
1112 // We use this to keep track of how many boxes are specified thru the refinement indicators
1113 num_boxes_at_level.resize(max_level+1,0);
1114 boxes_at_level.resize(max_level+1);
1115
1116 // We always have exactly one file at level 0
1117 num_boxes_at_level[0] = 1;
1118 boxes_at_level[0].resize(1);
1119 boxes_at_level[0][0] = geom[0].Domain();
1120
1121 // Plotfile name and frequency
1122 pp.queryAdd("plot_file", plot_file_name);
1123 pp.queryAdd("plot_int", plot_int);
1124 pp.queryAdd("plot_int_time", plot_int_time);
1125
1126 // Should we plot the staggered face velocities (without averaging to cell centers)
1127 pp.queryAdd("plot_staggered_vels", plot_staggered_vels);
1128
1129 // Output format
1130 std::string plotfile_type_str = "amrex";
1131 pp.queryAdd("plotfile_type", plotfile_type_str);
1132 if (plotfile_type_str == "amrex") {
1134 } else if (plotfile_type_str == "netcdf" || plotfile_type_str == "NetCDF") {
1136#ifdef REMORA_USE_NETCDF
1137 pp.queryAdd("write_history_file",write_history_file);
1138 pp.queryAdd("chunk_history_file",chunk_history_file);
1139 pp.queryAdd("steps_per_history_file",steps_per_history_file);
1140 // Estimate size of domain for one timestep of netcdf
1141 auto dom = geom[0].Domain();
1142 int nx = dom.length(0) + 2;
1143 int ny = dom.length(1) + 2;
1144 int nz = dom.length(2);
1146 // Estimate number of steps that will fit into a 2GB file.
1147 steps_per_history_file = int((1.6e10 - NCH2D * nx * ny * 64.0_rt)
1148 / (nx * ny * 64.0_rt * (NC3D*nz + NC2D)));
1149 // If we calculate that a single step will exceed 2GB and the user has
1150 // requested automatic history file sizing, warn about a possible impending
1151 // error, and set steps_per_history_file = 1 to attempt output anyway.
1152 if (steps_per_history_file == 0) {
1153 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.");
1155 }
1156 } else if (write_history_file and !chunk_history_file) {
1157 // Estimate number of output steps we'll need
1158 int nt_out = int((max_step) / plot_int) + 1;
1159 Real est_hist_file_size = NCH2D * nx * ny * 64.0_rt + nt_out * nx * ny * 64.0_rt * (NC3D*nz + NC2D);
1160 if (est_hist_file_size > 1.6e10) {
1161 amrex::Warning("WARNING: NetCDF history file may be larger than 2GB limit. Consider setting remora.chunk_history_file=true");
1162 }
1163 }
1165 Print() << "NetCDF history files will have " << steps_per_history_file << " steps per file." << std::endl;
1166 }
1167#endif
1168 } else {
1169 amrex::Print() << "User selected plotfile_type = " << plotfile_type_str << std::endl;
1170 amrex::Abort("Dont know this plotfile_type");
1171 }
1172#ifndef REMORA_USE_NETCDF
1174 {
1175 amrex::Abort("Please compile with NetCDF in order to enable NetCDF plotfiles");
1176 }
1177
1178#endif
1179
1180#ifdef REMORA_USE_NETCDF
1181 nc_init_file.resize(max_level+1);
1182 nc_grid_file.resize(max_level+1);
1183
1184 boundary_series.resize(max_level+1);
1185
1186 // NetCDF initialization files -- possibly multiple files at each of multiple levels
1187 // but we always have exactly one file at level 0
1188 for (int lev = 0; lev <= max_level; lev++)
1189 {
1190 const std::string nc_file_names = amrex::Concatenate("nc_init_file_",lev,1);
1191 const std::string nc_bathy_file_names = amrex::Concatenate("nc_grid_file_",lev,1);
1192
1193 if (pp.contains(nc_file_names.c_str()))
1194 {
1195 int num_files = pp.countval(nc_file_names.c_str());
1196 int num_bathy_files = pp.countval(nc_bathy_file_names.c_str());
1197 if (num_files != num_bathy_files) {
1198 amrex::Error("Must have same number of netcdf files for grid info as for solution");
1199 }
1200
1201 num_files_at_level[lev] = num_files;
1202 nc_init_file[lev].resize(num_files);
1203 nc_grid_file[lev].resize(num_files);
1204
1205 pp.queryarr(nc_file_names.c_str() , nc_init_file[lev] ,0,num_files);
1206 pp.queryarr(nc_bathy_file_names.c_str(), nc_grid_file[lev],0,num_files);
1207 }
1208 }
1209 // We only read boundary data at level 0
1210 pp.queryarr("nc_bdry_file", nc_bdry_file);
1211
1212 // Also only read forcings at level 0 (for now)
1213 pp.queryAdd("nc_frc_file", nc_frc_file);
1214
1215 // Get river file
1216 pp.queryAdd("nc_river_file", nc_riv_file);
1217
1218 // Read in file names for climatology history and nudging weights
1219 pp.queryAdd("nc_clim_his_file", nc_clim_his_file);
1220 pp.queryAdd("nc_clim_coeff_file", nc_clim_coeff_file);
1221
1222 for (int i=0; i<BdyVars::NumTypes; i++) {
1223 bdry_time_name_byvar.push_back("");
1224 }
1225 pp.queryAdd("bdy_time_varname",bdry_time_varname);
1226 pp.queryAdd("bdy_temp_time_varname",bdry_time_name_byvar[BdyVars::t]);
1227 pp.queryAdd("bdy_salt_time_varname",bdry_time_name_byvar[BdyVars::s]);
1228 pp.queryAdd("bdy_u_time_varname",bdry_time_name_byvar[BdyVars::u]);
1229 pp.queryAdd("bdy_v_time_varname",bdry_time_name_byvar[BdyVars::v]);
1230 pp.queryAdd("bdy_ubar_time_varname",bdry_time_name_byvar[BdyVars::ubar]);
1231 pp.queryAdd("bdy_vbar_time_varname",bdry_time_name_byvar[BdyVars::vbar]);
1232 pp.queryAdd("bdy_zeta_time_varname",bdry_time_name_byvar[BdyVars::zeta]);
1233
1234 // If not specified per variable, populate with the default
1235 for (int i=0; i<BdyVars::NumTypes; i++) {
1236 if (bdry_time_name_byvar[i] == "") {
1238 }
1239 }
1240
1241 pp.queryAdd("frc_time_varname",frc_time_varname);
1242
1243 pp.queryAdd("riv_time_varname",riv_time_varname);
1244
1245 pp.queryAdd("clim_ubar_time_varname",clim_ubar_time_varname);
1246 pp.queryAdd("clim_vbar_time_varname",clim_vbar_time_varname);
1247 pp.queryAdd("clim_u_time_varname",clim_u_time_varname);
1248 pp.queryAdd("clim_v_time_varname",clim_v_time_varname);
1249 pp.queryAdd("clim_salt_time_varname",clim_salt_time_varname);
1250 pp.queryAdd("clim_temp_time_varname",clim_temp_time_varname);
1251
1252#endif
1253
1254#ifdef REMORA_USE_PARTICLES
1255 readTracersParams();
1256#endif
1257 }
1258
1260}
1261
1262void
1264{
1265 BL_PROFILE("REMORA::AverageDown()");
1266 for (int lev = finest_level-1; lev >= 0; --lev)
1267 {
1268 AverageDownTo(lev);
1269 }
1270}
1271
1272/**
1273 * @param[in ] crse_lev level to average down to
1274 */
1275void
1277{
1278 BL_PROFILE("REMORA::AverageDownTo()");
1279 average_down(*cons_new[crse_lev+1], *cons_new[crse_lev],
1280 0, cons_new[crse_lev]->nComp(), refRatio(crse_lev));
1281 average_down(*vec_Zt_avg1[crse_lev+1].get(), *vec_Zt_avg1[crse_lev].get(),
1282 0, vec_Zt_avg1[crse_lev]->nComp(), refRatio(crse_lev));
1283
1284 Array<MultiFab*,AMREX_SPACEDIM> faces_crse;
1285 Array<MultiFab*,AMREX_SPACEDIM> faces_fine;
1286 faces_crse[0] = xvel_new[crse_lev];
1287 faces_crse[1] = yvel_new[crse_lev];
1288 faces_crse[2] = zvel_new[crse_lev];
1289
1290 faces_fine[0] = xvel_new[crse_lev+1];
1291 faces_fine[1] = yvel_new[crse_lev+1];
1292 faces_fine[2] = zvel_new[crse_lev+1];
1293
1294 average_down_faces(GetArrOfConstPtrs(faces_fine), faces_crse,
1295 refRatio(crse_lev),geom[crse_lev]);
1296 stretch_transform(crse_lev);
1297}
1298
1299/**
1300 * @param[in ] lev level at which to get time
1301 */
1302amrex::Real REMORA::get_t_old(int lev) const
1303{
1304 return t_old[lev];
1305}
amrex::Vector< std::string > BCNames
Definition REMORA.cpp:56
PlotfileType
plotfile format
#define Temp_comp
#define NC3D
#define Tracer_comp
#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.
static PlotfileType plotfile_type
Native or NetCDF plotfile output.
Definition REMORA.H:1386
void FillCoarsePatchPC(int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=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 an entire multifab by interpolating from the coarser level using the piecewise constant interpol...
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_EminusP
evaporation minus precipitation [kg/m^2/s], defined at rho-points
Definition REMORA.H:353
std::string riv_time_varname
Name of time field for river time.
Definition REMORA.H:1422
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:268
static amrex::Real fixed_dt
User specified fixed baroclinic time step.
Definition REMORA.H:1287
amrex::Real last_plot_file_time
Simulation time when we last output a plotfile.
Definition REMORA.H:1258
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.
void calculate_nodal_masks(int lev)
Calculate u-, v-, and psi-point masks based on rho-point masks after analytic initialization.
std::unique_ptr< NCTimeSeries > qair_data_from_file
Data container for specific humidity read from file.
Definition REMORA.H:1110
static amrex::Real previousCPUTimeUsed
Accumulator variable for CPU time used thusfar.
Definition REMORA.H:1482
amrex::Vector< std::unique_ptr< amrex::YAFluxRegister > > advflux_reg
array of flux registers for refluxing in multilevel
Definition REMORA.H:1224
std::unique_ptr< NCTimeSeries > sustr_data_from_file
Data container for u-component surface momentum flux read from file.
Definition REMORA.H:1100
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:405
void init_gls_vmix(int lev, SolverChoice solver_choice)
Initialize GLS variables.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:253
void set2DPlotVariables(const std::string &pp_plot_var_names_2d)
amrex::Vector< REMORAFillPatcher > FPr_v
Vector over levels of FillPatchers for v (3D)
Definition REMORA.H:1168
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:241
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
Definition REMORA.H:1097
void FillCoarsePatch(int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=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 an entire multifab by interpolating from the coarser level
void stretch_transform(int lev)
Calculate vertical stretched coordinates.
std::unique_ptr< NCTimeSeries > rain_data_from_file
Data container for precipitation rate read from file.
Definition REMORA.H:1118
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vwind
Wind in the v direction, defined at rho-points.
Definition REMORA.H:318
std::unique_ptr< ProblemBase > prob
Pointer to container of analytical functions for problem definition.
Definition REMORA.H:1197
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:389
void Construct_REMORAFillPatchers(int lev)
Construct FillPatchers.
Definition REMORA.cpp:361
static int sum_interval
Diagnostic sum output interval in number of steps.
Definition REMORA.H:1375
int history_count
Counter for which time index we are writing to in the netcdf history file.
Definition REMORA.H:1323
amrex::Real stop_time
Time to stop.
Definition REMORA.H:1273
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rain
precipitation rate [kg/m^2/s]
Definition REMORA.H:347
int do_substep
Whether to substep fine levels in time.
Definition REMORA.H:1296
void Evolve()
Advance solution to final time.
Definition REMORA.cpp:150
std::string bdry_time_varname
Default name of time field for boundary data.
Definition REMORA.H:1405
amrex::Real plotfile_fill_value
fill value for masked arrays in amrex plotfiles
Definition REMORA.H:1342
void ReadCheckpointFile()
read checkpoint file from disk
bool chunk_history_file
Whether to chunk netcdf history file.
Definition REMORA.H:1316
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:311
amrex::Real get_t_old(int lev) const
Accessor method for t_old to expose to outside classes.
Definition REMORA.cpp:1302
std::unique_ptr< NCTimeSeries > longwave_down_data_from_file
Data container for downward longwave radiation flux read from file.
Definition REMORA.H:1116
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:266
amrex::Vector< std::string > datalogname
Definition REMORA.H:1514
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:247
void init_only(int lev, amrex::Real time)
Init (NOT restart or regrid)
Definition REMORA.cpp:799
void init_set_vmix(int lev)
Initialize vertical mixing coefficients from file or analytic.
Definition REMORA.cpp:616
std::unique_ptr< NCTimeSeries > v_clim_data_from_file
Data container for v-velocity climatology data read from file.
Definition REMORA.H:1131
std::string clim_u_time_varname
Name of time field for u climatology data.
Definition REMORA.H:1414
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:587
amrex::Vector< REMORAFillPatcher > FPr_u
Vector over levels of FillPatchers for u (3D)
Definition REMORA.H:1166
std::string clim_temp_time_varname
Name of time field for temperature climatology data.
Definition REMORA.H:1420
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1204
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Holds info for dynamically generated tagging criteria.
Definition REMORA.H:1431
std::string nc_frc_file
NetCDF forcing file.
Definition REMORA.H:1396
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
Definition REMORA.H:276
std::unique_ptr< NCTimeSeriesRiver > river_source_transportbar
Data container for vertically integrated momentum transport in rivers.
Definition REMORA.H:1142
std::string clim_ubar_time_varname
Name of time field for ubar climatology data.
Definition REMORA.H:1410
std::unique_ptr< NCTimeSeries > u_clim_data_from_file
Data container for u-velocity climatology data read from file.
Definition REMORA.H:1129
std::string check_file
Checkpoint file prefix.
Definition REMORA.H:1309
static amrex::Real startCPUTime
Variable for CPU timing.
Definition REMORA.H:1480
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:234
amrex::Real start_time
Time of the start of the simulation, in seconds.
Definition REMORA.H:1276
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:245
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_uwind
Wind in the u direction, defined at rho-points.
Definition REMORA.H:316
static amrex::Real fixed_fast_dt
User specified fixed barotropic time step.
Definition REMORA.H:1289
int regrid_int
how often each level regrids the higher levels of refinement (after a level advances that many time s...
Definition REMORA.H:1299
amrex::Real check_int_time
Checkpoint output interval in seconds.
Definition REMORA.H:1313
void Define_REMORAFillPatchers(int lev)
Define FillPatchers.
Definition REMORA.cpp:410
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:278
amrex::Real plot_int_time
Plotfile output interval in seconds.
Definition REMORA.H:1307
amrex::Vector< int > num_files_at_level
how many netcdf input files specified at each level
Definition REMORA.H:1202
amrex::Vector< REMORAFillPatcher > FPr_vbar
Vector over levels of FillPatchers for vbar (2D)
Definition REMORA.H:1176
void AverageDownTo(int crse_lev)
more flexible version of AverageDown() that lets you average down across multiple levels
Definition REMORA.cpp:1276
int steps_per_history_file
Number of time steps per netcdf history file.
Definition REMORA.H:1321
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Called after every level 0 timestep.
Definition REMORA.cpp:224
int max_step
maximum number of steps
Definition REMORA.H:1271
amrex::Vector< amrex::MultiFab * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
Definition REMORA.H:238
std::unique_ptr< NCTimeSeries > svstr_data_from_file
Data container for v-component surface momentum flux read from file.
Definition REMORA.H:1102
amrex::Vector< int > num_boxes_at_level
how many boxes specified at each level by tagging criteria
Definition REMORA.H:1200
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:243
void refinement_criteria_setup()
Set refinement criteria.
int last_check_file_step
Step when we last output a checkpoint file.
Definition REMORA.H:1261
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:1412
amrex::Vector< int > nsubsteps
How many substeps on each level?
Definition REMORA.H:1209
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
Definition REMORA.H:1221
void ComputeDt()
a wrapper for estTimeStep()
void fill_3d_masks(int lev)
Copy maskr to all z levels.
std::unique_ptr< NCTimeSeries > EminusP_data_from_file
Data container for evaporation minus precipitation read from file.
Definition REMORA.H:1122
int plot_int
Plotfile output interval in iterations.
Definition REMORA.H:1305
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.
std::unique_ptr< NCTimeSeries > cloud_data_from_file
Data container for cloud cover fraction read from file.
Definition REMORA.H:1120
void InitData()
Initialize multilevel data.
Definition REMORA.cpp:256
void set3DPlotVariables(const std::string &pp_plot_var_names_3d)
amrex::Vector< int > istep
which step?
Definition REMORA.H:1207
void WriteCheckpointFile()
write checkpoint file to disk
std::string nc_clim_coeff_file
NetCDF climatology coefficient file.
Definition REMORA.H:1402
void setRecordDataInfo(int i, const std::string &filename)
Definition REMORA.H:1500
void set_analytic_vmix(int lev)
Set vertical mixing coefficients from analytic.
Definition REMORA.cpp:633
void init_flat_bathymetry(int lev)
Initialize flat bathymetry to value from problo.
Definition REMORA.cpp:673
std::unique_ptr< NCTimeSeries > temp_clim_data_from_file
Data container for temperature climatology data read from file.
Definition REMORA.H:1133
amrex::Vector< std::string > bdry_time_name_byvar
Name of time fields for boundary data.
Definition REMORA.H:1407
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
Definition REMORA.H:1380
void init_riv_pos_from_netcdf(int lev)
static amrex::Vector< std::string > nc_bdry_file
NetCDF boundary data.
Definition REMORA.H:51
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_r
Harmonic viscosity defined on the rho points (centers)
Definition REMORA.H:280
std::string nc_riv_file
Definition REMORA.H:1397
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:313
std::unique_ptr< NCTimeSeries > srflx_data_from_file
Data container for shortwave radiation flux read from file.
Definition REMORA.H:1114
REMORA()
Definition REMORA.cpp:64
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:1285
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:229
void set_bathymetry(int lev)
Initialize bathymetry from file or analytic.
Definition REMORA.cpp:498
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:236
std::unique_ptr< NCTimeSeries > ubar_clim_data_from_file
Data container for ubar climatology data read from file.
Definition REMORA.H:1125
amrex::Vector< REMORAFillPatcher > FPr_c
Vector over levels of FillPatchers for scalars.
Definition REMORA.H:1164
std::unique_ptr< NCTimeSeries > Tair_data_from_file
Data container for air temperature read from file.
Definition REMORA.H:1108
std::string clim_v_time_varname
Name of time field for v climatology data.
Definition REMORA.H:1416
amrex::Vector< REMORAFillPatcher > FPr_w
Vector over levels of FillPatchers for w.
Definition REMORA.H:1170
std::string clim_salt_time_varname
Name of time field for salinity climatology data.
Definition REMORA.H:1418
std::unique_ptr< NCTimeSeries > Uwind_data_from_file
Data container for u-direction wind read from file.
Definition REMORA.H:1104
std::unique_ptr< NCTimeSeries > Pair_data_from_file
Data container for air pressure read from file.
Definition REMORA.H:1112
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1211
void init_stretch_coeffs()
initialize and calculate stretch coefficients
void init_bdry_from_netcdf(int lev)
Boundary data initialization from NetCDF file.
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1336
void set_masks(int lev)
Initialize land-sea masks from file or analytic.
Definition REMORA.cpp:772
amrex::Vector< amrex::Vector< std::unique_ptr< NCTimeSeriesBoundary > > > boundary_series
Vector over BdyVars of boundary series data containers.
Definition REMORA.H:1145
int cf_set_width
Width for fixing values at coarse-fine interface.
Definition REMORA.H:1161
void ReadParameters()
read in some parameters from inputs file
Definition REMORA.cpp:1021
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:262
void sum_integrated_quantities(amrex::Real time)
Integrate conserved quantities for diagnostics.
static int total_nc_plot_file_step
Definition REMORA.H:1008
static int plot_staggered_vels
Whether to write the staggered velocities (not averaged to cell centers)
Definition REMORA.H:1383
static amrex::Vector< amrex::Vector< std::string > > nc_grid_file
NetCDF grid file.
Definition REMORA.H:53
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_longwave_down
Downward longwave radiation.
Definition REMORA.H:331
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
Definition REMORA.H:386
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:384
void set_zeta_to_Ztavg(int lev)
Set zeta components to be equal to time-averaged Zt_avg1.
bool expand_plotvars_to_unif_rr
whether plotfile variables should be expanded to a uniform refinement ratio
Definition REMORA.H:1339
int plot_file_on_restart
Whether to output a plotfile on restart from checkpoint.
Definition REMORA.H:1265
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:382
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:232
std::string frc_time_varname
Name of time field for forcing data.
Definition REMORA.H:1424
void set_wind(int lev)
Initialize or calculate wind speed from file or analytic.
Definition REMORA.cpp:702
amrex::Vector< REMORAFillPatcher > FPr_ubar
Vector over levels of FillPatchers for ubar (2D)
Definition REMORA.H:1174
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.
std::unique_ptr< NCTimeSeries > Vwind_data_from_file
Data container for v-direction wind read from file.
Definition REMORA.H:1106
void set_hmixcoef(int lev)
Initialize horizontal mixing coefficients.
Definition REMORA.cpp:647
amrex::Vector< std::unique_ptr< NCTimeSeriesRiver > > river_source_cons
Vector of data containers for scalar data in rivers.
Definition REMORA.H:1138
std::string nc_clim_his_file
NetCDF climatology history file.
Definition REMORA.H:1400
void timeStep(int lev, amrex::Real time, int iteration)
advance a level by dt, includes a recursive call for finer levels
std::unique_ptr< NCTimeSeriesRiver > river_source_transport
Data container for momentum transport in rivers.
Definition REMORA.H:1140
void AverageDown()
set covered coarse cells to be the average of overlying fine cells
Definition REMORA.cpp:1263
static int fixed_ndtfast_ratio
User specified, number of barotropic steps per baroclinic step.
Definition REMORA.H:1291
amrex::Real netcdf_fill_value
fill value for masked arrays in netcdf output
Definition REMORA.H:1344
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:1513
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
Definition REMORA.H:274
static amrex::Real cfl
CFL condition.
Definition REMORA.H:1283
void append3DPlotVariables(const std::string &pp_plot_var_names_3d)
static int verbose
Verbosity level of output.
Definition REMORA.H:1372
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_cloud
cloud cover fraction [0-1], defined at rho-points
Definition REMORA.H:351
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1303
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:308
std::unique_ptr< NCTimeSeries > salt_clim_data_from_file
Data container for salinity climatology data read from file.
Definition REMORA.H:1135
int check_int
Checkpoint output interval in iterations.
Definition REMORA.H:1311
void set_smflux(int lev)
Initialize or calculate surface momentum flux from file or analytic.
Definition REMORA.cpp:683
void WritePlotFile(int istep)
main driver for writing AMReX plotfiles
std::string restart_chkfile
If set, restart from this checkpoint file.
Definition REMORA.H:1279
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:264
int cf_width
Nudging width at coarse-fine interface.
Definition REMORA.H:1159
static amrex::Vector< amrex::Vector< std::string > > nc_init_file
NetCDF initialization file.
Definition REMORA.H:52
int last_plot_file_step
Step when we last output a plotfile.
Definition REMORA.H:1256
amrex::Vector< amrex::Real > t_old
old time at each level
Definition REMORA.H:1213
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_srflx
Shortwave radiation flux [W/m²], defined at rho-points.
Definition REMORA.H:327
amrex::Real last_check_file_time
Simulation time when we last output a checkpoint file.
Definition REMORA.H:1263
void append2DPlotVariables(const std::string &pp_plot_var_names_2d)
amrex::Vector< amrex::Real > dt
time step at each level
Definition REMORA.H:1215
static amrex::Real sum_per
Diagnostic sum output interval in time.
Definition REMORA.H:1377
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Pair
Air pressure [mb], defined at rho-points.
Definition REMORA.H:324
virtual ~REMORA()
Definition REMORA.cpp:145
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_qair
Specific humidity [kg/kg], defined at rho-points.
Definition REMORA.H:322
void restart()
Definition REMORA.cpp:457
std::unique_ptr< NCTimeSeries > vbar_clim_data_from_file
Data container for vbar climatology data read from file.
Definition REMORA.H:1127
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Tair
Air temperature [°C], defined at rho-points.
Definition REMORA.H:320
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
Definition REMORA.H:282
const char * buildInfoGetGitHash(int i)
HorizMixingType horiz_mixing_type
amrex::Vector< amrex::Real > tnu2
std::string longwave_netcdf_varname
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