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 cons_names.clear();
153 cons_names.reserve(ncons);
154 cons_names.emplace_back("temp");
155 cons_names.emplace_back("salt");
156 cons_names.emplace_back("tracer");
157 for (int i = 1; i < nscalar; ++i) {
158 cons_names.emplace_back("tracer_" + std::to_string(i));
159 }
160}
161
162void
164{
165 BL_PROFILE("REMORA::Evolve()");
166 Real cur_time = t_new[0];
167
168 // Take one coarse timestep by calling timeStep -- which recursively calls timeStep
169 // for finer levels (with or without subcycling)
170 for (int step = istep[0]; step < max_step && cur_time < stop_time; ++step)
171 {
172 amrex::Print() << "\nCoarse STEP " << step+1 << " starts ..." << std::endl;
173
174 ComputeDt();
175
176 int lev = 0;
177 int iteration = 1;
178
179 if (max_level == 0) {
180 timeStep(lev, cur_time, iteration);
181 }
182 else {
183 timeStepML(cur_time, iteration);
184 }
185
186 cur_time += dt[0];
187
188 amrex::Print() << "Coarse STEP " << step+1 << " ends." << " TIME = " << cur_time
189 << " DT = " << dt[0] << std::endl;
190
191 if ( (plot_int > 0 && (step+1 - last_plot_file_step) == plot_int ) ||
192 (plot_int_time > 0 && (cur_time >= (last_plot_file_time + plot_int_time))) )
193 {
194 last_plot_file_step = step+1;
195 last_plot_file_time = cur_time;
196 WritePlotFile(step+1);
198 }
199
200 if ((check_int > 0 && (step+1 - last_check_file_step) == check_int)
201 || (check_int_time > 0 && cur_time >= (last_check_file_time + check_int_time))) {
202 last_check_file_step = step+1;
203 last_check_file_time = cur_time;
205 }
206
207 post_timestep(step, cur_time, dt[0]);
208
209#ifdef AMREX_MEM_PROFILING
210 {
211 std::ostringstream ss;
212 ss << "[STEP " << step+1 << "]";
213 MemProfiler::report(ss.str());
214 }
215#endif
216
217 if (cur_time >= stop_time - 1.e-6*dt[0]) break;
218 }
219
220 if ( (plot_int > 0 || plot_int_time > 0.0) && istep[0] > last_plot_file_step)
221 {
224 }
225
226 if ((check_int > 0 || check_int_time > 0.0) && istep[0] > last_check_file_step) {
228 }
229}
230
231/**
232 * @param[in ] nstep which step we're on
233 * @param[in ] time current time
234 * @param[in ] dt_lev0 time step on level 0
235 */
236void
237REMORA::post_timestep (int nstep, Real time, Real dt_lev0)
238{
239 BL_PROFILE("REMORA::post_timestep()");
240
241#ifdef REMORA_USE_PARTICLES
242 particleData.Redistribute();
243#endif
244
246 {
247 for (int lev = finest_level-1; lev >= 0; lev--)
248 {
249 // This call refluxes from the lev/lev+1 interface onto lev
250 //getAdvFluxReg(lev+1)->Reflux(*cons_new[lev], 0, 0, NCONS);
251
252 // We need to do this before anything else because refluxing changes the
253 // values of coarse cells underneath fine grids with the assumption they'll
254 // be over-written by averaging down
255 //
256 AverageDownTo(lev);
257 }
258 }
259
260 if (is_it_time_for_action(nstep, time, dt_lev0, sum_interval, sum_per)) {
262 }
263}
264
265/**
266 * This is called from main.cpp and handles all initialization, whether from start or restart
267 */
268void
270{
271 BL_PROFILE("REMORA::InitData()");
272 // Initialize the start time for our CPU-time tracker
273 startCPUTime = Real(ParallelDescriptor::second());
274
275 // Map the words in the inputs file to BC types, then translate
276 // those types into what they mean for each variable
277 init_bcs();
278
279 // Init vertical stretching coeffs
281
284 last_plot_file_time = -1.0_rt;
285 last_check_file_time = -1.0_rt;
286
287 if (restart_chkfile == "") {
288 // start simulation from the beginning
289
290 InitFromScratch(start_time);
291
293 AverageDown();
294 }
295
296 } else { // Restart from a checkpoint
297
298 restart();
299
300 }
301#ifdef REMORA_USE_MOAB
302 InitMOABMesh();
303#endif
304 // Initialize flux registers (whether we start from scratch or restart)
306 advflux_reg[0] = nullptr;
307 for (int lev = 1; lev <= finest_level; lev++)
308 {
309 advflux_reg[lev].reset( new YAFluxRegister(grids[lev], grids[lev-1],
310 dmap[lev], dmap[lev-1],
311 geom[lev], geom[lev-1],
312 ref_ratio[lev-1], lev, ncons));
313 }
314 }
315
316 // Fill ghost cells/faces
317 for (int lev = 0; lev <= finest_level; ++lev)
318 {
319 if (lev > 0 && cf_width >= 0) {
321 }
322
323 if (restart_chkfile == "") {
324 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]);
325 FillPatch(lev, t_new[lev], *xvel_new[lev], xvel_new, xvel_bc(), BdyVars::u, 0, true, false,0,0,0.0,*xvel_new[lev]);
326 FillPatch(lev, t_new[lev], *yvel_new[lev], yvel_new, yvel_bc(), BdyVars::v, 0, true, false,0,0,0.0,*yvel_new[lev]);
327 FillPatch(lev, t_new[lev], *zvel_new[lev], zvel_new, zvel_bc(), BdyVars::null, 0, true, false);
328
329 // Copy from new into old just in case when initializing from scratch
330 int ngs = cons_new[lev]->nGrow();
331 int ngvel = xvel_new[lev]->nGrow();
332 MultiFab::Copy(*cons_old[lev],*cons_new[lev],0,0,ncons,ngs);
333 MultiFab::Copy(*xvel_old[lev],*xvel_new[lev],0,0,1,ngvel);
334 MultiFab::Copy(*yvel_old[lev],*yvel_new[lev],0,0,1,ngvel);
335 MultiFab::Copy(*zvel_old[lev],*zvel_new[lev],0,0,1,IntVect(ngvel,ngvel,0));
336 }
337 } // lev
338
339 // Check for additional plotting variables that are available after
340 // particle containers are setup.
341 const std::string& pv3d = "plot_vars_3d"; append3DPlotVariables(pv3d);
342 const std::string& pv2d = "plot_vars_2d"; append2DPlotVariables(pv2d);
343
344 if (restart_chkfile == "" && (check_int > 0 || check_int_time > 0.0_rt))
345 {
348 }
349
350 if ( (restart_chkfile == "") ||
352 {
353 if (plot_int > 0 || plot_int_time > 0.0)
354 {
355 int step0 = 0;
356 WritePlotFile(step0);
359 }
360 }
361
364 }
365
366 ComputeDt();
367
368}
369
370/**
371 * @param[in ] lev level to operate on
372 */
373void
375{
376 BL_PROFILE("REMORA::Construct_REMORAFillPatchers()");
377 amrex::Print() << ":::Construct_REMORAFillPatchers " << lev << std::endl;
378
379 auto& ba_fine = cons_new[lev ]->boxArray();
380 auto& ba_crse = cons_new[lev-1]->boxArray();
381 auto& dm_fine = cons_new[lev ]->DistributionMap();
382 auto& dm_crse = cons_new[lev-1]->DistributionMap();
383
384 BoxList bl2d_fine = ba_fine.boxList();
385 for (auto& b : bl2d_fine) {
386 b.setRange(2,0);
387 }
388 BoxArray ba2d_fine(std::move(bl2d_fine));
389
390 BoxList bl2d_crse = ba_crse.boxList();
391 for (auto& b : bl2d_crse) {
392 b.setRange(2,0);
393 }
394 BoxArray ba2d_crse(std::move(bl2d_crse));
395
396 int ncomp = cons_new[lev]->nComp();
397
398 FPr_c.emplace_back(ba_fine, dm_fine, geom[lev] ,
399 ba_crse, dm_crse, geom[lev-1],
400 -cf_width, -cf_set_width, ncomp, &cell_cons_interp);
401 FPr_u.emplace_back(convert(ba_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
402 convert(ba_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
403 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
404 FPr_v.emplace_back(convert(ba_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
405 convert(ba_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
406 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
407 FPr_w.emplace_back(convert(ba_fine, IntVect(0,0,1)), dm_fine, geom[lev] ,
408 convert(ba_crse, IntVect(0,0,1)), dm_crse, geom[lev-1],
409 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
410
411 FPr_ubar.emplace_back(convert(ba2d_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
412 convert(ba2d_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
413 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
414 FPr_vbar.emplace_back(convert(ba2d_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
415 convert(ba2d_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
416 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
417}
418
419/**
420 * @param[in ] lev level to operate on
421 */
422void
424{
425 BL_PROFILE("REMORA::Define_REMORAFillPatchers()");
426 amrex::Print() << ":::Define_REMORAFillPatchers " << lev << std::endl;
427
428 auto& ba_fine = cons_new[lev ]->boxArray();
429 auto& ba_crse = cons_new[lev-1]->boxArray();
430 auto& dm_fine = cons_new[lev ]->DistributionMap();
431 auto& dm_crse = cons_new[lev-1]->DistributionMap();
432
433 BoxList bl2d_fine = ba_fine.boxList();
434 for (auto& b : bl2d_fine) {
435 b.setRange(2,0);
436 }
437 BoxArray ba2d_fine(std::move(bl2d_fine));
438
439 BoxList bl2d_crse = ba_crse.boxList();
440 for (auto& b : bl2d_crse) {
441 b.setRange(2,0);
442 }
443 BoxArray ba2d_crse(std::move(bl2d_crse));
444
445
446 int ncomp = cons_new[lev]->nComp();
447
448 FPr_c[lev-1].Define(ba_fine, dm_fine, geom[lev] ,
449 ba_crse, dm_crse, geom[lev-1],
450 -cf_width, -cf_set_width, ncomp, &cell_cons_interp);
451 FPr_u[lev-1].Define(convert(ba_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
452 convert(ba_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
453 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
454 FPr_v[lev-1].Define(convert(ba_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
455 convert(ba_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
456 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
457 FPr_w[lev-1].Define(convert(ba_fine, IntVect(0,0,1)), dm_fine, geom[lev] ,
458 convert(ba_crse, IntVect(0,0,1)), dm_crse, geom[lev-1],
459 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
460
461 FPr_ubar[lev-1].Define(convert(ba2d_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
462 convert(ba2d_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
463 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
464 FPr_vbar[lev-1].Define(convert(ba2d_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
465 convert(ba2d_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
466 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
467}
468
469void
471{
472 BL_PROFILE("REMORA::restart()");
474
475 // We set this here so that we don't over-write the checkpoint file we just started from
477}
478
479/**
480 * @param[in ] lev level to operate on
481 */
482void
484{
485 BL_PROFILE("REMORA::set_zeta()");
487 prob->init_analytic_zeta(lev, geom[lev], solverChoice, *this, *vec_zeta[lev]);
488
489#ifdef REMORA_USE_NETCDF
490 } else if (solverChoice.ic_type == IC_Type::netcdf) {
491 if (lev == 0) {
492 amrex::Print() << "Calling init_zeta_from_netcdf on level " << lev << std::endl;
494 amrex::Print() << "Sea surface height loaded from netcdf file \n " << std::endl;
495 } else {
496 Real dummy_time = 0.0_rt;
497 FillCoarsePatch(lev,dummy_time,vec_zeta[lev].get(), vec_zeta[lev-1].get(),BCVars::cons_bc);
498 }
499#endif
500 } else {
501 Abort("Don't know this ic_type!");
502 }
503 vec_zeta[lev]->FillBoundary(geom[lev].periodicity());
504 set_zeta_average(lev);
505}
506
507/**
508 * @param[in ] lev level to operate on
509 */
510void
512{
513 BL_PROFILE("REMORA::bathymetry()");
514 // Only set bathymetry on level 0, and interpolate for finer levels
515 if (lev==0) {
518 // If grid data is not defined on a level > 0 (negative level) then
519 // initialize from low-resolution grid normally. Otherwise use high-resolution
520 // grid data averaged down to level 0
521 } else if (hires_grid_level < 0) {
523 prob->init_analytic_bathymetry(lev, geom[lev], solverChoice, *this, *vec_h[lev]);
524 } else if (solverChoice.ic_type == IC_Type::netcdf) {
525#ifdef REMORA_USE_NETCDF
526 amrex::Print() << "Calling init_bathymetry_from_netcdf " << std::endl;
528 amrex::Print() << "Bathymetry loaded from netcdf file \n " << std::endl;
529#endif
530 } else {
531 amrex::Abort("Unknown IC_Type");
532 }
533 } else {
535 }
536 // Need FillBoundary to fill at grid-grid boundaries, and EnforcePeriodicity
537 // to make sure ghost cells in the domain corners are consistent.
538 vec_h[lev]->FillBoundary(geom[lev].periodicity());
539 vec_h[lev]->EnforcePeriodicity(geom[lev].periodicity());
540 } else {
541 // If our level is higher than the high resolution grid or initialization
542 // is analytic, interpolate from level below. Otherwise, copy over the bathymetry
543 // data that has been averaged down
544 if (lev > hires_grid_level) {
545 Real dummy_time = 0.0_rt;
546 FillCoarsePatch(lev,dummy_time,vec_h[lev].get(), vec_h[lev-1].get(),BCVars::cons_bc);
547 } else {
549 vec_h[lev]->FillBoundary(geom[lev].periodicity());
550 vec_h[lev]->EnforcePeriodicity(geom[lev].periodicity());
551 }
553 set_grid_scale(lev);
554 }
555 }
556}
557
558/**
559 * @param[in ] lev level to operate on
560 */
561void
563 Real dummy_time = 0.0_rt;
564 ParallelCopy(*vec_h[lev].get(), *vec_h_full_domain[lev].get(), 0, 0, 1,IntVect(1,1,0),IntVect(1,1,0));
565 ParallelCopy(*vec_h[lev].get(), *vec_h_full_domain[lev].get(), 0, 1, 1,IntVect(1,1,0),IntVect(1,1,0));
566 FillPatch(lev,dummy_time,*vec_h[lev],GetVecOfPtrs(vec_h),
568 BdyVars::null,0,false,true,1);
569 FillPatch(lev,dummy_time,*vec_h[lev],GetVecOfPtrs(vec_h),
571 BdyVars::null,1,false,true,1);
572}
573
574/**
575 * @param[in ] lev level to operate on
576 */
577void
579 BL_PROFILE("REMORA::set_coriolis()");
582 prob->init_analytic_coriolis(lev, geom[lev], solverChoice, *this, *vec_fcor[lev]);
585#ifdef REMORA_USE_NETCDF
587 if (lev == 0) {
588 amrex::Print() << "Calling init_coriolis_from_netcdf " << std::endl;
590 amrex::Print() << "Coriolis loaded from netcdf file \n" << std::endl;
591 } else {
592 Real dummy_time = 0.0_rt;
593 FillCoarsePatch(lev,dummy_time,vec_fcor[lev].get(), vec_fcor[lev-1].get(),BCVars::cons_bc);
594 }
595#endif
596 } else {
597 Abort("Don't know this coriolis_type!");
598 }
599
600 Real time = 0.0_rt;
601 FillPatch(lev, time, *vec_fcor[lev], GetVecOfPtrs(vec_fcor), foextrap_bc());
602 vec_fcor[lev]->EnforcePeriodicity(geom[lev].periodicity());
603 }
604}
605
606void
608 BL_PROFILE("REMORA::init_set_vmix()");
613 // The GLS initialization just sets the multifab to a value, so there's
614 // no need to call FillPatch here
615 } else {
616 Abort("Don't know this vertical mixing type");
617 }
618}
619
620/**
621 * @param[in ] lev level to operate on
622 */
623void
625 BL_PROFILE("REMORA::set_analytic_vmix()");
626 Real time = 0.0_rt;
627 vec_Akv[lev]->setVal(solverChoice.Akv_bak);
628 vec_Akt[lev]->setVal(solverChoice.Akt_bak);
629 prob->init_analytic_vmix(lev, geom[lev], solverChoice, *this,*vec_Akv[lev], *vec_Akt[lev]);
630 FillPatch(lev, time, *vec_Akv[lev], GetVecOfPtrs(vec_Akv), zvel_bc(), BdyVars::null,0,true,false);
631 for (int n = 0; n < ncons; n++) {
632 FillPatch(lev, time, *vec_Akt[lev], GetVecOfPtrs(vec_Akt), zvel_bc(), BdyVars::null,n,false,false);
633 }
634}
635
636/**
637 * @param[in ] lev level to operate on
638 */
639void
641{
643 prob->init_analytic_masks(lev,geom[lev], solverChoice, *this, *vec_mskr[lev]);
646#ifdef REMORA_USE_NETCDF
647 if (lev == 0) {
648 amrex::Print() << "Calling init_masks_from_netcdf level " << lev << std::endl;
650 amrex::Print() << "Masks loaded from netcdf file \n " << std::endl;
651 } else {
652 Real dummy_time = 0.0_rt;
653 FillCoarsePatchPC(lev, dummy_time, vec_mskr[lev].get(), vec_mskr[lev-1].get(),
654 foextrap_bc());
656 }
657#endif
658 }
659 fill_3d_masks(lev);
660}
661
662/**
663 * @param[in ] lev level to operate on
664 */
665void
667{
668 BL_PROFILE("REMORA::set_hmixcoef()");
669
670 // Optional AMR scaling: decrease coefficients on refined levels linearly
671 // with grid size (i.e., proportional to sqrt(cell area)). For a horizontal
672 // refinement ratio rx x ry, the effective scale factor is 1/sqrt(rx*ry).
673 Real lev_scale = 1.0_rt;
675 Real rf = 1.0_rt;
676 for (int l = 0; l < lev; ++l) {
677 rf *= std::sqrt(static_cast<Real>(ref_ratio[l][0]) * static_cast<Real>(ref_ratio[l][1]));
678 }
679 lev_scale = 1.0_rt / rf;
680 }
681
683 prob->init_analytic_hmix(lev, geom[lev], solverChoice,
684 *this, *vec_visc2_p[lev], *vec_visc2_r[lev], *vec_diff2[lev]);
685
687 vec_visc2_p[lev]->setVal(solverChoice.visc2 * lev_scale);
688 vec_visc2_r[lev]->setVal(solverChoice.visc2 * lev_scale);
689 for (int n = 0; n < ncons; n++) {
690 vec_diff2[lev]->setVal(solverChoice.tnu2[n] * lev_scale, n, 1);
691 }
692
693 // Scale harmonic viscosity and diffusivity by the grid size as ROMS
694 // does in Utility/ini_hmixcoef.F. Intended for curvilinear grids.
695 //
696 // Define the ROMS grid factor (grdscl):
697 // G(i,j) = sqrt( 1 / (pm(i,j) * pn(i,j)) )
698 // = sqrt(cell area)
699 // Gmax = max over grid of G(i,j)
700 //
701 // Then horizontal harmonic mixing coefficients are scaled as:
702 // nu(i,j) = nu0 * G(i,j) / Gmax
703 // kappa_n(i,j) = kappa0 * G(i,j) / Gmax
704 //
705 // where:
706 // nu0 = solverChoice.visc2
707 // kappa0 = solverChoice.tnu2[n]
708 //
709 // This makes mixing strongest where grid spacing is largest.
710 //
711 // NOTE: The normalization (Gmax) is computed over the entire grid (ignoring masks).
712 // Therefore, if the largest cell area occurs over land, the maximum over *wet* cells
713 // (or in masked output files) may be smaller than the user-specified value.
714
716
717 // ------------------------------------------------------------
718 // Step 1: Compute grdmax over entire grid
719 // ------------------------------------------------------------
720 vec_visc2_r[lev]->setVal(solverChoice.visc2);
721 vec_visc2_p[lev]->setVal(solverChoice.visc2);
722 for (int n = 0; n < ncons; n++) {
723 vec_diff2[lev]->setVal(solverChoice.tnu2[n], n, 1);
724 }
725
726 // NOTE: This must be GPU-safe. Do not dereference MultiFab data on host.
727 // Force the reduction to run in the GPU launch region if GPUs are enabled.
728 // (If the launch region is disabled at runtime, ReduceMax may fall back to
729 // a host path that can try to read device-only data.)
730 amrex::Gpu::LaunchSafeGuard lsg(true);
731 Real denom_min = amrex::ReduceMin(*vec_pm[lev], *vec_pn[lev], 0,
732 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
733 Array4<Real const> const& pm,
734 Array4<Real const> const& pn) -> Real
735 {
736 Real local_min = 1.0e200_rt;
737 amrex::Loop(bx, [=,&local_min] (int i, int j, int) noexcept
738 {
739 local_min = amrex::min(local_min, pm(i,j,0) * pn(i,j,0));
740 });
741 return local_min;
742 });
743
744 ParallelDescriptor::ReduceRealMin(denom_min);
745 if (denom_min <= 0.0_rt) {
746 Abort("scaled_to_grid: found non-positive pm*pn (grid metrics must be > 0)");
747 }
748
749 Real grdmax = amrex::ReduceMax(*vec_pm[lev], *vec_pn[lev], 0,
750 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
751 Array4<Real const> const& pm,
752 Array4<Real const> const& pn) -> Real
753 {
754 Real local_max = 0.0_rt;
755 amrex::Loop(bx, [=,&local_max] (int i, int j, int) noexcept
756 {
757 Real denom = pm(i,j,0) * pn(i,j,0);
758 if (denom > 0.0_rt) {
759 Real G = std::sqrt(1.0_rt / denom);
760 local_max = amrex::max(local_max, G);
761 }
762 });
763 return local_max;
764 });
765
766 ParallelDescriptor::ReduceRealMax(grdmax);
767 if (grdmax <= 0.0_rt) {
768 Abort("scaled_to_grid: grdmax <= 0");
769 }
770
771 // Optional AMR scaling: decrease coefficients on refined levels linearly
772 // with grid size (i.e., proportional to sqrt(cell area)). For a horizontal
773 // refinement ratio rx x ry, the effective scale factor is 1/sqrt(rx*ry).
774 lev_scale = 1.0_rt;
776 Real rf = 1.0_rt;
777 for (int l = 0; l < lev; ++l) {
778 rf *= std::sqrt(static_cast<Real>(ref_ratio[l][0]) * static_cast<Real>(ref_ratio[l][1]));
779 }
780 lev_scale = 1.0_rt / rf;
781 }
782
783 Real visc0 = solverChoice.visc2 * lev_scale;
784 Real cff = visc0 / grdmax;
785
786 // ------------------------------------------------------------
787 // Step 2: Set rho coefficients everywhere
788 // ------------------------------------------------------------
789 amrex::Gpu::DeviceVector<Real> diff0_d(ncons);
790 amrex::Gpu::copy(amrex::Gpu::hostToDevice,
791 solverChoice.tnu2.begin(), solverChoice.tnu2.begin() + ncons,
792 diff0_d.begin());
793 Real const* diff0_ptr = diff0_d.data();
794
795 for (MFIter mfi(*vec_visc2_r[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
796 {
797 const Box& bx = mfi.validbox();
798 auto pm = vec_pm[lev]->const_array(mfi);
799 auto pn = vec_pn[lev]->const_array(mfi);
800 auto visc2_r = vec_visc2_r[lev]->array(mfi);
801 auto diff2 = vec_diff2[lev]->array(mfi);
802
803 int ncons_local = ncons;
804 ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int) noexcept
805 {
806 Real denom = pm(i,j,0) * pn(i,j,0);
807 Real grdscl = (denom > 0.0_rt) ? std::sqrt(1.0_rt / denom) : 0.0_rt;
808 visc2_r(i,j,0) = cff * grdscl;
809
810 for (int n = 0; n < ncons_local; n++) {
811 diff2(i,j,0,n) = ((diff0_ptr[n] * lev_scale) / grdmax) * grdscl;
812 }
813 });
814 }
815
816 // Fill ghost cells for rho coefficients BEFORE psi averaging
817 Real time = 0.0_rt;
818 FillPatch(lev, time, *vec_visc2_r[lev], GetVecOfPtrs(vec_visc2_r), foextrap_periodic_bc());
819
820 // ------------------------------------------------------------
821 // Step 3: Psi coefficients = average of 4 surrounding rho
822 // ------------------------------------------------------------
823 for (MFIter mfi(*vec_visc2_p[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
824 {
825 const Box& bx = mfi.validbox();
826 auto visc2_p = vec_visc2_p[lev]->array(mfi);
827 auto visc2_r = vec_visc2_r[lev]->const_array(mfi);
828
829 ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int) noexcept
830 {
831 visc2_p(i,j,0) = 0.25_rt * (
832 visc2_r(i-1,j-1,0) +
833 visc2_r(i ,j-1,0) +
834 visc2_r(i-1,j ,0) +
835 visc2_r(i ,j ,0)
836 );
837 });
838 }
839
840 FillPatch(lev, time, *vec_visc2_p[lev], GetVecOfPtrs(vec_visc2_p), foextrap_periodic_bc());
841
842 // Diagnostics
843 // NOTE: coefficients are computed everywhere (including land). Output routines may later
844 // mask land points (e.g., to FillValue in NetCDF/plotfiles), and analysis tools may
845 // additionally apply mask_rho (setting land to 0). Report both conventions.
846 //
847 // Global (MPI-reduced) extrema over all valid cells (no ghost).
848 Real visc_min_all = vec_visc2_r[lev]->min(0,0,false);
849 Real visc_max_all = vec_visc2_r[lev]->max(0,0,false);
850
851 // Global extrema over *wet* rho points only, k=0.
852 amrex::Gpu::LaunchSafeGuard lsg_diag(true);
853 Real visc_min_wet = amrex::ReduceMin(*vec_visc2_r[lev], *vec_mskr[lev], 0,
854 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
855 Array4<Real const> const& visc2,
856 Array4<Real const> const& mskr) -> Real
857 {
858 Real local_min = 1.0e200_rt;
859 amrex::Loop(bx, [=,&local_min] (int i, int j, int) noexcept
860 {
861 if (mskr(i,j,0) > 0.0_rt) {
862 local_min = amrex::min(local_min, visc2(i,j,0));
863 }
864 });
865 return local_min;
866 });
867 ParallelDescriptor::ReduceRealMin(visc_min_wet);
868
869 Real visc_max_wet = amrex::ReduceMax(*vec_visc2_r[lev], *vec_mskr[lev], 0,
870 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
871 Array4<Real const> const& visc2,
872 Array4<Real const> const& mskr) -> Real
873 {
874 Real local_max = -1.0e200_rt;
875 amrex::Loop(bx, [=,&local_max] (int i, int j, int) noexcept
876 {
877 if (mskr(i,j,0) > 0.0_rt) {
878 local_max = amrex::max(local_max, visc2(i,j,0));
879 }
880 });
881 return local_max;
882 });
883 ParallelDescriptor::ReduceRealMax(visc_max_wet);
884
885 // Mimic "apply mask_rho" convention (dry -> 0).
886 Real visc_min_mask0 = amrex::ReduceMin(*vec_visc2_r[lev], *vec_mskr[lev], 0,
887 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
888 Array4<Real const> const& visc2,
889 Array4<Real const> const& mskr) -> Real
890 {
891 Real local_min = 1.0e200_rt;
892 amrex::Loop(bx, [=,&local_min] (int i, int j, int) noexcept
893 {
894 const Real v = (mskr(i,j,0) > 0.0_rt) ? visc2(i,j,0) : 0.0_rt;
895 local_min = amrex::min(local_min, v);
896 });
897 return local_min;
898 });
899 ParallelDescriptor::ReduceRealMin(visc_min_mask0);
900
901 Real visc_max_mask0 = amrex::ReduceMax(*vec_visc2_r[lev], *vec_mskr[lev], 0,
902 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
903 Array4<Real const> const& visc2,
904 Array4<Real const> const& mskr) -> Real
905 {
906 Real local_max = -1.0e200_rt;
907 amrex::Loop(bx, [=,&local_max] (int i, int j, int) noexcept
908 {
909 const Real v = (mskr(i,j,0) > 0.0_rt) ? visc2(i,j,0) : 0.0_rt;
910 local_max = amrex::max(local_max, v);
911 });
912 return local_max;
913 });
914 ParallelDescriptor::ReduceRealMax(visc_max_mask0);
915 if (ParallelDescriptor::IOProcessor() && lev == 0)
916 {
917 Print() << "\nHorizontal mixing scaled by grid metric\n";
918 Print() << "grdmax = " << grdmax << "\n";
920 Print() << "AMR scaling (linear) lev_scale = " << lev_scale << "\n";
921 }
922 Print() << "visc2(all) min/max = "
923 << visc_min_all << " / "
924 << visc_max_all << "\n";
925 Print() << "visc2(wet,k=0) min/max = "
926 << visc_min_wet << " / "
927 << visc_max_wet << "\n";
928 Print() << "visc2(mask->0) min/max = "
929 << visc_min_mask0 << " / "
930 << visc_max_mask0 << "\n";
931 }
932
933 } else {
934 Abort("Don't know this horizontal mixing type");
935 }
936
937 // Final FillPatch for all fields
938 Real time = 0.0_rt;
939 FillPatch(lev, time, *vec_visc2_p[lev], GetVecOfPtrs(vec_visc2_p), foextrap_periodic_bc());
940 FillPatch(lev, time, *vec_visc2_r[lev], GetVecOfPtrs(vec_visc2_r), foextrap_periodic_bc());
941 for (int n = 0; n < ncons; n++) {
942 FillPatch(lev, time, *vec_diff2[lev], GetVecOfPtrs(vec_diff2),
944 }
945}
946
947/**
948 * @param[in ] lev level to operate on
949 */
950void
952{
953 BL_PROFILE("REMORA::init_flat_bathymetry()");
954 vec_h[lev]->setVal(-geom[0].ProbLo()[2]);
955}
956
957/**
958 * @param[in ] lev level to operate on
959 */
960void
962{
963 BL_PROFILE("REMORA::set_smflux()");
965 prob->init_analytic_smflux(lev, geom[lev], solverChoice, *this,*vec_sustr[lev], *vec_svstr[lev]);
967#ifdef REMORA_USE_NETCDF
968 sustr_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_sustr[lev].get(), geom, ref_ratio);
969 svstr_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_svstr[lev].get(), geom, ref_ratio);
970 FillPatch(lev, t_old[lev], *vec_sustr[lev], GetVecOfPtrs(vec_sustr), foextrap_periodic_bc(), BdyVars::null,0,false);
971 FillPatch(lev, t_old[lev], *vec_svstr[lev], GetVecOfPtrs(vec_svstr), foextrap_periodic_bc(), BdyVars::null,0,false);
972#endif
973 }
974}
975
976/**
977 * @param[in ] lev level to operate on
978 */
979void
981{
982 BL_PROFILE("REMORA::set_wind()");
984 prob->init_analytic_wind(lev,geom[lev], solverChoice, *this, *vec_uwind[lev], *vec_vwind[lev]);
986#ifdef REMORA_USE_NETCDF
987 Uwind_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_uwind[lev].get(), geom, ref_ratio);
988 Vwind_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_vwind[lev].get(), geom, ref_ratio);
989 FillPatch(lev, t_old[lev], *vec_uwind[lev], GetVecOfPtrs(vec_uwind),
991 FillPatch(lev, t_old[lev], *vec_vwind[lev], GetVecOfPtrs(vec_vwind),
993
994 // Conditionally update atmospheric fields if loaded from NetCDF
996 Tair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_Tair[lev].get(), geom, ref_ratio);
997 FillPatch(lev, t_old[lev], *vec_Tair[lev], GetVecOfPtrs(vec_Tair),
999 }
1001 qair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_qair[lev].get(), geom, ref_ratio);
1002 FillPatch(lev, t_old[lev], *vec_qair[lev], GetVecOfPtrs(vec_qair),
1004
1005 // Convert qair from percentage (0-100) to specific humidity (0-1) if needed
1007 vec_qair[lev]->mult(0.01);
1008
1009 // Update ghost cells after modification
1010 vec_qair[lev]->FillBoundary(geom[lev].periodicity());
1011 }
1012 }
1014 Pair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_Pair[lev].get(), geom, ref_ratio);
1015 FillPatch(lev, t_old[lev], *vec_Pair[lev], GetVecOfPtrs(vec_Pair),
1017 }
1019 srflx_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_srflx[lev].get(), geom, ref_ratio);
1020 FillPatch(lev, t_old[lev], *vec_srflx[lev], GetVecOfPtrs(vec_srflx),
1022 }
1024 longwave_down_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_longwave_down[lev].get(), geom, ref_ratio);
1025 FillPatch(lev, t_old[lev], *vec_longwave_down[lev], GetVecOfPtrs(vec_longwave_down),
1027 }
1029 rain_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_rain[lev].get(), geom, ref_ratio);
1030 FillPatch(lev, t_old[lev], *vec_rain[lev], GetVecOfPtrs(vec_rain),
1032 }
1034 cloud_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_cloud[lev].get(), geom, ref_ratio);
1035 FillPatch(lev, t_old[lev], *vec_cloud[lev], GetVecOfPtrs(vec_cloud),
1037 }
1039 EminusP_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_EminusP[lev].get(), geom, ref_ratio);
1040 FillPatch(lev, t_old[lev], *vec_EminusP[lev], GetVecOfPtrs(vec_EminusP),
1042 }
1043#endif
1044 }
1045}
1046
1047/**
1048 * @param[in ] lev level to operate on
1049 * @param[in ] time current time for initialization
1050 */
1051void
1052REMORA::init_only (int lev, Real time)
1053{
1054 BL_PROFILE("REMORA::init_only()");
1055 t_new[lev] = time;
1056 t_old[lev] = time - 1.e200_rt;
1057
1058 cons_new[lev]->setVal(0.0_rt);
1059 xvel_new[lev]->setVal(0.0_rt);
1060 yvel_new[lev]->setVal(0.0_rt);
1061 zvel_new[lev]->setVal(0.0_rt);
1062
1063 xvel_old[lev]->setVal(0.0_rt);
1064 yvel_old[lev]->setVal(0.0_rt);
1065 zvel_old[lev]->setVal(0.0_rt);
1066
1067 vec_ru[lev]->setVal(0.0_rt);
1068 vec_rv[lev]->setVal(0.0_rt);
1069
1070 vec_ru2d[lev]->setVal(0.0_rt);
1071 vec_rv2d[lev]->setVal(0.0_rt);
1072
1074 set_grid_scale(lev);
1075 }
1076 set_masks(lev);
1077
1078#ifdef REMORA_USE_NETCDF
1081
1082 if (solverChoice.do_any_clim_nudg && lev == 0) {
1083 if (nc_clim_his_file.empty() || nc_clim_his_file[0].empty()) {
1084 amrex::Error("NetCDF climatology file name must be provided via input");
1085 }
1088 clim_ubar_time_varname, geom[lev].Domain(),vec_ubar[lev].get(),true,true));
1090 clim_ubar_time_varname, geom[lev].Domain(),vec_vbar[lev].get(),true,true));
1091 ubar_clim_data_from_file->Initialize();
1092 vbar_clim_data_from_file->Initialize();
1093 }
1095 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));
1096 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));
1097 u_clim_data_from_file->Initialize();
1098 v_clim_data_from_file->Initialize();
1099 }
1100 // Since the NCTimeSeries object isn't filling the cons_new MultiFab directly, we don't have to specify a component.
1101 // It just needs to know the shape of the MultiFab
1103 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));
1104 temp_clim_data_from_file->Initialize();
1105 }
1107 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));
1108 salt_clim_data_from_file->Initialize();
1109 }
1110 }
1111 }
1112
1114 amrex::Print() << "Calling init_bdry_from_netcdf at level " << lev << std::endl;
1116 amrex::Print() << "Boundary data loaded from netcdf file \n " << std::endl;
1117 }
1118
1119 // This will be a non-op if forcings specified analytically
1121 if (lev==0) {
1122 if (nc_frc_file.empty() || nc_frc_file[0].empty()) {
1123 amrex::Error("NetCDF forcing file name must be provided via input for winds");
1124 }
1125 Uwind_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Uwind", frc_time_varname, geom[lev].Domain(),vec_uwind[lev].get(), true, false));
1126 Vwind_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Vwind", frc_time_varname, geom[lev].Domain(),vec_vwind[lev].get(), true, false));
1127 Uwind_data_from_file->Initialize();
1128 Vwind_data_from_file->Initialize();
1129 } else {
1130 FillCoarsePatch(lev, time, vec_uwind[lev].get(), vec_uwind[lev-1].get(), foextrap_bc());
1131 FillCoarsePatch(lev, time, vec_vwind[lev].get(), vec_vwind[lev-1].get(), foextrap_bc());
1132 }
1134 if (lev==0) {
1135 if (nc_frc_file.empty() || nc_frc_file[0].empty()) {
1136 amrex::Error("NetCDF forcing file name must be provided via input for surface momentum fluxes");
1137 }
1138 sustr_data_from_file.reset(new NCTimeSeries(nc_frc_file, "sustr", frc_time_varname, geom[lev].Domain(),vec_sustr[lev].get(), true, false));
1139 svstr_data_from_file.reset(new NCTimeSeries(nc_frc_file, "svstr", frc_time_varname, geom[lev].Domain(),vec_svstr[lev].get(), true, false));
1140 sustr_data_from_file->Initialize();
1141 svstr_data_from_file->Initialize();
1142 } else {
1143 FillCoarsePatch(lev, time, vec_sustr[lev].get(), vec_sustr[lev-1].get(), foextrap_bc());
1144 FillCoarsePatch(lev, time, vec_svstr[lev].get(), vec_svstr[lev-1].get(), foextrap_bc());
1145 }
1146 }
1147
1148 // Conditionally load atmospheric forcing fields from NetCDF based on user flags
1149 if (lev==0) {
1151 Tair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Tair", frc_time_varname, geom[lev].Domain(),vec_Tair[lev].get(), true, false));
1152 Tair_data_from_file->Initialize();
1153 }
1155 qair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "qair", frc_time_varname, geom[lev].Domain(),vec_qair[lev].get(), true, false));
1156 qair_data_from_file->Initialize();
1157 }
1159 Pair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Pair", frc_time_varname, geom[lev].Domain(),vec_Pair[lev].get(), true, false));
1160 Pair_data_from_file->Initialize();
1161 }
1163 srflx_data_from_file.reset(new NCTimeSeries(nc_frc_file, "swrad", frc_time_varname, geom[lev].Domain(),vec_srflx[lev].get(), true, false));
1164 srflx_data_from_file->Initialize();
1165 }
1167 rain_data_from_file.reset(new NCTimeSeries(nc_frc_file, "rain", frc_time_varname, geom[lev].Domain(),vec_rain[lev].get(), true, false));
1168 rain_data_from_file->Initialize();
1169 }
1171 cloud_data_from_file.reset(new NCTimeSeries(nc_frc_file, "cloud", frc_time_varname, geom[lev].Domain(),vec_cloud[lev].get(), true, false));
1172 cloud_data_from_file->Initialize();
1173 }
1175 EminusP_data_from_file.reset(new NCTimeSeries(nc_frc_file, "EminusP", frc_time_varname, geom[lev].Domain(),vec_EminusP[lev].get(), true, false));
1176 EminusP_data_from_file->Initialize();
1177 }
1178 } else {
1180 FillCoarsePatch(lev, time, vec_Tair[lev].get(), vec_Tair[lev-1].get(), foextrap_bc());
1181 }
1183 FillCoarsePatch(lev, time, vec_qair[lev].get(), vec_qair[lev-1].get(), foextrap_bc());
1184 }
1186 FillCoarsePatch(lev, time, vec_Pair[lev].get(), vec_Pair[lev-1].get(), foextrap_bc());
1187 }
1189 FillCoarsePatch(lev, time, vec_srflx[lev].get(), vec_srflx[lev-1].get(), foextrap_bc());
1190 }
1192 FillCoarsePatch(lev, time, vec_rain[lev].get(), vec_rain[lev-1].get(), foextrap_bc());
1193 }
1195 FillCoarsePatch(lev, time, vec_cloud[lev].get(), vec_cloud[lev-1].get(), foextrap_bc());
1196 }
1198 FillCoarsePatch(lev, time, vec_EminusP[lev].get(), vec_EminusP[lev-1].get(), foextrap_bc());
1199 }
1200 }
1202 if (lev==0) {
1203 if (nc_frc_file.empty() || nc_frc_file[0].empty()) {
1204 amrex::Error("NetCDF forcing file name must be provided via input for longwave radiation");
1205 }
1207 geom[lev].Domain(), vec_longwave_down[lev].get(), true, false));
1208 longwave_down_data_from_file->Initialize();
1209 } else {
1210 FillCoarsePatch(lev, time, vec_longwave_down[lev].get(), vec_longwave_down[lev-1].get(), foextrap_bc());
1211 }
1212 }
1213
1214 // Only need to read in rivers on level 0
1215 // Will need to be on higher levels eventually
1216 if (solverChoice.do_rivers) {
1217 if (nc_riv_file.empty() || nc_riv_file[0].empty()) {
1218 amrex::Error("NetCDF river file name must be provided via input for rivers");
1219 }
1220 auto dom = geom[0].Domain();
1221 int nz = dom.length(2);
1222 river_source_cons.resize(ncons);
1225 river_source_cons[Salt_comp]->Initialize();
1226 }
1229 river_source_cons[Temp_comp]->Initialize();
1230 }
1233 river_source_cons[Tracer_comp]->Initialize();
1234 }
1235 river_source_transport.reset(new NCTimeSeriesRiver(nc_riv_file, "river_transport", riv_time_varname, nz));
1236 river_source_transport->Initialize();
1237 river_source_transportbar.reset(new NCTimeSeriesRiver(nc_riv_file, "river_transport", riv_time_varname, nz, 1));
1238 river_source_transportbar->Initialize();
1240 }
1241
1242 if (lev==0 and hires_grid_level > 0 and solverChoice.ic_type == IC_Type::netcdf) {
1243 amrex::Print() << "Reading high resolution bathymetry data" << std::endl;
1246 amrex::Print() << "Done reading in high resolution bathymetry data" << std::endl;
1247 }
1248#else
1250 Abort("Not compiled with NetCDF, but selected boundary conditions require NetCDF");
1251 }
1252 if (solverChoice.do_rivers) {
1253 Abort("Not compiled with NetCDF, but using river sources requires NetCDF");
1254 }
1255#endif
1256
1257 if (lev==0 and hires_grid_level > 0 and solverChoice.ic_type == IC_Type::analytic) {
1260 }
1261
1262 set_bathymetry(lev);
1263 set_zeta(lev);
1264 stretch_transform(lev);
1265
1267 if (lev==0) {
1269 {
1270 init_analytic(lev);
1271#ifdef REMORA_USE_NETCDF
1272 } else if (solverChoice.ic_type == IC_Type::netcdf) {
1273 amrex::Print() << "Calling init_data_from_netcdf " << std::endl;
1275 set_zeta_to_Ztavg(lev);
1276 amrex::Print() << "Initial data loaded from netcdf file \n " << std::endl;
1277
1278#endif
1279 } else {
1280 Abort("Need to specify ic_type");
1281 }
1282 } else {
1284 FillCoarsePatch(lev, time, xvel_new[lev], xvel_new[lev-1], xvel_bc(), BdyVars::u);
1285 FillCoarsePatch(lev, time, yvel_new[lev], yvel_new[lev-1], yvel_bc(), BdyVars::v);
1286 FillCoarsePatch(lev, time, zvel_new[lev], zvel_new[lev-1], zvel_bc(), BdyVars::null);
1287 }
1290 {
1291 init_analytic(lev);
1292#ifdef REMORA_USE_NETCDF
1293 } else if (solverChoice.ic_type == IC_Type::netcdf) {
1294 amrex::Print() << "Calling init_data_from_netcdf on level " << lev << std::endl;
1296 set_zeta_to_Ztavg(lev);
1297 amrex::Print() << "Initial data loaded from netcdf file on level " << lev << std::endl;
1298
1299#endif
1300 } else {
1301 Abort("Need to specify ic_type");
1302 }
1303 } else {
1304 amrex::Abort("Need to specify T init procedure");
1305 }
1306
1307 // Ensure that the face-based data are the same on both sides of a periodic domain.
1308 // The data associated with the lower grid ID is considered the correct value.
1309 xvel_new[lev]->OverrideSync(geom[lev].periodicity());
1310 yvel_new[lev]->OverrideSync(geom[lev].periodicity());
1311 zvel_new[lev]->OverrideSync(geom[lev].periodicity());
1312
1313 set_2darrays(lev);
1314
1315 init_set_vmix(lev);
1316 set_hmixcoef(lev);
1317 set_coriolis(lev);
1318
1319 // Previously set smflux here with OverrideSync:
1320// set_smflux(lev);
1321// prob->init_analytic_smflux(lev, geom[lev], solverChoice, *this, *vec_sustr[lev], *vec_svstr[lev]);
1322// vec_sustr[lev]->OverrideSync(geom[lev].periodicity());
1323// vec_svstr[lev]->OverrideSync(geom[lev].periodicity());
1324
1325}
1326
1327void
1329{
1330 BL_PROFILE("REMORA::ReadParameters()");
1331 {
1332 ParmParse pp; // Traditionally, max_step and stop_time do not have prefix, so allow it for now.
1333 bool noprefix_max_step = pp.queryAdd("max_step", max_step);
1334 bool noprefix_stop_time = pp.queryAdd("stop_time", stop_time);
1335 bool remora_max_step = pp.queryAdd("remora.max_step", max_step);
1336 bool remora_stop_time = pp.queryAdd("remora.stop_time", stop_time);
1337 if (remora_max_step and noprefix_max_step) {
1338 Abort("remora.max_step and max_step are both specified. Please use only one!");
1339 }
1340 if (remora_stop_time and noprefix_stop_time) {
1341 Abort("remora.stop_time and stop_time are both specified. Please use only one!");
1342 }
1343 }
1344
1345 ParmParse pp(pp_prefix);
1346 ParmParse pp_amr("amr");
1347 {
1348 pp.queryAdd("nscalar", nscalar);
1349 if (nscalar < 1) {
1350 amrex::Abort("remora.nscalar must be at least 1");
1351 }
1354
1355 pp_amr.queryAdd("regrid_int", regrid_int);
1356 pp.queryAdd("check_file", check_file);
1357 pp.queryAdd("check_int", check_int);
1358 pp_amr.queryAdd("check_int", check_int);
1359 pp.queryAdd("check_int_time", check_int_time);
1360 pp_amr.queryAdd("check_int_time", check_int_time);
1361
1362 pp.queryAdd("expand_plotvars_to_unif_rr", expand_plotvars_to_unif_rr);
1363
1364 pp.query("plotfile_fill_value", plotfile_fill_value);
1365 pp.query("netcdf_fill_value", netcdf_fill_value);
1366
1367 pp.queryAdd("restart", restart_chkfile);
1368 pp_amr.queryAdd("restart", restart_chkfile);
1369 pp.queryAdd("start_time",start_time);
1370
1371 if (pp.contains("data_log"))
1372 {
1373 int num_datalogs = pp.countval("data_log");
1374 datalog.resize(num_datalogs);
1375 datalogname.resize(num_datalogs);
1376 pp.queryarr("data_log",datalogname,0,num_datalogs);
1377 for (int i = 0; i < num_datalogs; i++)
1379 }
1380
1381 // Verbosity
1382 pp.queryAdd("v", verbose);
1383
1384 // Frequency of diagnostic output
1385 pp.queryAdd("sum_interval", sum_interval);
1386 pp.queryAdd("sum_period" , sum_per);
1387 pp.queryAdd("file_min_digits", file_min_digits);
1388
1389 if (file_min_digits < 0) {
1390 amrex::Abort("remora.file_min_digits must be non-negative");
1391 }
1392
1393 // Time step controls
1394 pp.queryAdd("cfl", cfl);
1395 pp.queryAdd("change_max", change_max);
1396
1397 pp.queryAdd("fixed_dt", fixed_dt);
1398 pp.queryAdd("fixed_fast_dt", fixed_fast_dt);
1399
1400 pp.queryAdd("fixed_ndtfast_ratio", fixed_ndtfast_ratio);
1401
1402 // If all three are specified, they must be consistent
1403 if (fixed_dt > 0. && fixed_fast_dt > 0. && fixed_ndtfast_ratio > 0)
1404 {
1406 {
1407 amrex::Abort("Dt is over-specfied");
1408 }
1409 }
1410 // If two are specified, initialize fixed_ndtfast_ratio
1411 else if (fixed_dt > 0. && fixed_fast_dt > 0. && fixed_ndtfast_ratio <= 0)
1412 {
1413 fixed_ndtfast_ratio = static_cast<int>(fixed_dt / fixed_fast_dt);
1414 }
1415
1416 AMREX_ASSERT(cfl > 0. || fixed_dt > 0.);
1417
1418 pp_amr.queryAdd("do_substep", do_substep);
1419 if (do_substep) {
1420 amrex::Abort("Time substepping is not yet implemented. amr.do_substep must be 0");
1421 }
1422
1423 // We use this to keep track of how many boxes we read in from WRF initialization
1424 num_files_at_level.resize(max_level+1,0);
1425
1426 // We use this to keep track of how many boxes are specified thru the refinement indicators
1427 num_boxes_at_level.resize(max_level+1,0);
1428 boxes_at_level.resize(max_level+1);
1429
1430 // We always have exactly one file at level 0
1431 num_boxes_at_level[0] = 1;
1432 boxes_at_level[0].resize(1);
1433 boxes_at_level[0][0] = geom[0].Domain();
1434
1435 // Plotfile name and frequency
1436 pp.queryAdd("plot_file", plot_file_name);
1437 pp.queryAdd("plot_int", plot_int);
1438 pp.queryAdd("plot_int_time", plot_int_time);
1439
1440 // Should we plot the staggered face velocities (without averaging to cell centers)
1441 pp.queryAdd("plot_staggered_vels", plot_staggered_vels);
1442
1443 // Output format
1444 std::string plotfile_type_str = "amrex";
1445 pp.queryAdd("plotfile_type", plotfile_type_str);
1446 if (plotfile_type_str == "amrex") {
1448 } else if (plotfile_type_str == "netcdf" || plotfile_type_str == "NetCDF") {
1450#ifdef REMORA_USE_NETCDF
1451 pp.queryAdd("write_history_file",write_history_file);
1452 pp.queryAdd("chunk_history_file",chunk_history_file);
1453 pp.queryAdd("steps_per_history_file",steps_per_history_file);
1454 // Estimate size of domain for one timestep of netcdf
1455 auto dom = geom[0].Domain();
1456 int nx = dom.length(0) + 2;
1457 int ny = dom.length(1) + 2;
1458 int nz = dom.length(2);
1460 // Estimate number of steps that will fit into a 2GB file.
1461 steps_per_history_file = int((1.6e10 - NCH2D * nx * ny * 64.0_rt)
1462 / (nx * ny * 64.0_rt * (NC3D*nz + NC2D)));
1463 // If we calculate that a single step will exceed 2GB and the user has
1464 // requested automatic history file sizing, warn about a possible impending
1465 // error, and set steps_per_history_file = 1 to attempt output anyway.
1466 if (steps_per_history_file == 0) {
1467 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.");
1469 }
1470 } else if (write_history_file and !chunk_history_file) {
1471 // Estimate number of output steps we'll need
1472 int nt_out = int((max_step) / plot_int) + 1;
1473 Real est_hist_file_size = NCH2D * nx * ny * 64.0_rt + nt_out * nx * ny * 64.0_rt * (NC3D*nz + NC2D);
1474 if (est_hist_file_size > 1.6e10) {
1475 amrex::Warning("WARNING: NetCDF history file may be larger than 2GB limit. Consider setting remora.chunk_history_file=true");
1476 }
1477 }
1479 Print() << "NetCDF history files will have " << steps_per_history_file << " steps per file." << std::endl;
1480 }
1481#endif
1482 } else {
1483 amrex::Print() << "User selected plotfile_type = " << plotfile_type_str << std::endl;
1484 amrex::Abort("Dont know this plotfile_type");
1485 }
1486#ifndef REMORA_USE_NETCDF
1488 {
1489 amrex::Abort("Please compile with NetCDF in order to enable NetCDF plotfiles");
1490 }
1491
1492#endif
1493#ifdef REMORA_USE_NETCDF
1494 nc_init_file.resize(max_level+1);
1495 nc_grid_file.resize(max_level+1);
1496
1497 boundary_series.resize(max_level+1);
1498
1499 // NetCDF initialization files -- possibly multiple files at each of multiple levels
1500 // but we always have exactly one file at level 0
1501 for (int lev = 0; lev <= max_level; lev++)
1502 {
1503 const std::string nc_file_names = amrex::Concatenate("nc_init_file_",lev,1);
1504 const std::string nc_bathy_file_names = amrex::Concatenate("nc_grid_file_",lev,1);
1505
1506 if (pp.contains(nc_file_names.c_str()))
1507 {
1508 int num_files = pp.countval(nc_file_names.c_str());
1509 int num_bathy_files = pp.countval(nc_bathy_file_names.c_str());
1510 if (num_files != num_bathy_files) {
1511 amrex::Error("Must have same number of netcdf files for grid info as for solution");
1512 }
1513
1514 num_files_at_level[lev] = num_files;
1515 nc_init_file[lev].resize(num_files);
1516 nc_grid_file[lev].resize(num_files);
1517
1518 pp.queryarr(nc_file_names.c_str() , nc_init_file[lev] ,0,num_files);
1519 pp.queryarr(nc_bathy_file_names.c_str(), nc_grid_file[lev],0,num_files);
1520 }
1521 }
1522
1523 pp.queryAdd("nc_grid_file_hires", nc_grid_file_hires);
1524
1525 // We only read boundary data at level 0
1526 pp.queryarr("nc_bdry_file", nc_bdry_file);
1527
1528 // Also only read forcings at level 0 (for now)
1529 if (pp.contains("nc_frc_file")) {
1530 int num_files = pp.countval("nc_frc_file");
1531 nc_frc_file.resize(num_files);
1532 pp.queryarr("nc_frc_file", nc_frc_file, 0, num_files);
1533 }
1534
1535 // Get river file
1536 if (pp.contains("nc_river_file")) {
1537 int num_files = pp.countval("nc_river_file");
1538 nc_riv_file.resize(num_files);
1539 pp.queryarr("nc_river_file", nc_riv_file, 0, num_files);
1540 }
1541
1542 // Read in file names for climatology history and nudging weights
1543 if (pp.contains("nc_clim_his_file")) {
1544 int num_files = pp.countval("nc_clim_his_file");
1545 nc_clim_his_file.resize(num_files);
1546 pp.queryarr("nc_clim_his_file", nc_clim_his_file, 0, num_files);
1547 }
1548 pp.queryAdd("nc_clim_coeff_file", nc_clim_coeff_file);
1549
1550 for (int i=0; i<BdyVars::NumTypes; i++) {
1551 bdry_time_name_byvar.push_back("");
1552 }
1553 pp.queryAdd("bdy_time_varname",bdry_time_varname);
1554 pp.queryAdd("bdy_temp_time_varname",bdry_time_name_byvar[BdyVars::t]);
1555 pp.queryAdd("bdy_salt_time_varname",bdry_time_name_byvar[BdyVars::s]);
1556 pp.queryAdd("bdy_u_time_varname",bdry_time_name_byvar[BdyVars::u]);
1557 pp.queryAdd("bdy_v_time_varname",bdry_time_name_byvar[BdyVars::v]);
1558 pp.queryAdd("bdy_ubar_time_varname",bdry_time_name_byvar[BdyVars::ubar]);
1559 pp.queryAdd("bdy_vbar_time_varname",bdry_time_name_byvar[BdyVars::vbar]);
1560 pp.queryAdd("bdy_zeta_time_varname",bdry_time_name_byvar[BdyVars::zeta]);
1561
1562 // If not specified per variable, populate with the default
1563 for (int i=0; i<BdyVars::NumTypes; i++) {
1564 if (bdry_time_name_byvar[i] == "") {
1566 }
1567 }
1568
1569 pp.queryAdd("frc_time_varname",frc_time_varname);
1570
1571 pp.queryAdd("riv_time_varname",riv_time_varname);
1572
1573 pp.queryAdd("clim_ubar_time_varname",clim_ubar_time_varname);
1574 pp.queryAdd("clim_vbar_time_varname",clim_vbar_time_varname);
1575 pp.queryAdd("clim_u_time_varname",clim_u_time_varname);
1576 pp.queryAdd("clim_v_time_varname",clim_v_time_varname);
1577 pp.queryAdd("clim_salt_time_varname",clim_salt_time_varname);
1578 pp.queryAdd("clim_temp_time_varname",clim_temp_time_varname);
1579
1580#endif
1581 pp.queryAdd("hires_grid_level", hires_grid_level);
1582 if (hires_grid_level > max_level) {
1583 amrex::Abort("hires_grid_level must be less than or equal to amr.max_level");
1584 }
1585
1586#ifdef REMORA_USE_PARTICLES
1587 readTracersParams();
1588#endif
1589 }
1590
1592}
1593
1594void
1596{
1597 BL_PROFILE("REMORA::AverageDown()");
1598 for (int lev = finest_level-1; lev >= 0; --lev)
1599 {
1600 AverageDownTo(lev);
1601 }
1602}
1603
1604/**
1605 * @param[in ] crse_lev level to average down to
1606 */
1607void
1609{
1610 BL_PROFILE("REMORA::AverageDownTo()");
1611 average_down(*cons_new[crse_lev+1], *cons_new[crse_lev],
1612 0, cons_new[crse_lev]->nComp(), refRatio(crse_lev));
1613 average_down(*vec_Zt_avg1[crse_lev+1].get(), *vec_Zt_avg1[crse_lev].get(),
1614 0, vec_Zt_avg1[crse_lev]->nComp(), refRatio(crse_lev));
1615
1616 Array<MultiFab*,AMREX_SPACEDIM> faces_crse;
1617 Array<MultiFab*,AMREX_SPACEDIM> faces_fine;
1618 faces_crse[0] = xvel_new[crse_lev];
1619 faces_crse[1] = yvel_new[crse_lev];
1620 faces_crse[2] = zvel_new[crse_lev];
1621
1622 faces_fine[0] = xvel_new[crse_lev+1];
1623 faces_fine[1] = yvel_new[crse_lev+1];
1624 faces_fine[2] = zvel_new[crse_lev+1];
1625
1626 average_down_faces(GetArrOfConstPtrs(faces_fine), faces_crse,
1627 refRatio(crse_lev),geom[crse_lev]);
1628 stretch_transform(crse_lev);
1629}
1630
1631/**
1632 * @param[in ] crse_lev level to average data down to
1633 * @param[inout] vec_mf vector over levels of multifabs containing data to average
1634 */
1635void
1636REMORA::average_down_with_grow_cells(int crse_lev, Vector<std::unique_ptr<MultiFab>>& vec_mf)
1637{
1638 auto const& crsema = vec_mf[crse_lev]->arrays();
1639 auto const& finema = vec_mf[crse_lev+1]->const_arrays();
1640 auto nghost_crse = cum_ref_ratios[crse_lev];
1641 auto ref_ratio_crse = refRatio(crse_lev);
1642 ParallelFor(*vec_mf[crse_lev], nghost_crse, 1,
1643 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1644 {
1645 amrex_avgdown(i,j,k,n,crsema[box_no],finema[box_no],0,0,ref_ratio_crse);
1646 });
1647 Gpu::streamSynchronize();
1648}
1649
1650/**
1651 * @param[in ] lev level at which to get time
1652 */
1653amrex::Real REMORA::get_t_old(int lev) const
1654{
1655 return t_old[lev];
1656}
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
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:1453
std::string nc_grid_file_hires
Grid file for high resolution bathymetry.
Definition REMORA.H:1465
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:357
amrex::Vector< std::string > nc_riv_file
NetCDF river file(s)
Definition REMORA.H:1474
void allocate_bathymetry_full_domain()
Allocate multifabs for storing full-domain bathymetry data.
std::string riv_time_varname
Name of time field for river time.
Definition REMORA.H:1499
int foextrap_periodic_bc() const noexcept
Definition REMORA.H:1059
amrex::Vector< std::string > nc_clim_his_file
NetCDF climatology history file(s)
Definition REMORA.H:1477
int ncons
Number of conserved scalars in the state (temperature + salt + passive scalars)
Definition REMORA.H:1346
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:272
static amrex::Real fixed_dt
User specified fixed baroclinic time step.
Definition REMORA.H:1354
amrex::Real last_plot_file_time
Simulation time when we last output a plotfile.
Definition REMORA.H:1320
int zvel_bc() const noexcept
Definition REMORA.H:1054
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.
int xvel_bc() const noexcept
Definition REMORA.H:1052
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:1172
static amrex::Real previousCPUTimeUsed
Accumulator variable for CPU time used thusfar.
Definition REMORA.H:1558
amrex::Vector< std::string > cons_names
Names of scalars for plotfile output.
Definition REMORA.H:1397
amrex::Vector< std::unique_ptr< amrex::YAFluxRegister > > advflux_reg
array of flux registers for refluxing in multilevel
Definition REMORA.H:1286
std::unique_ptr< NCTimeSeries > sustr_data_from_file
Data container for u-component surface momentum flux read from file.
Definition REMORA.H:1162
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:409
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:254
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:404
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:1230
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:242
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
Definition REMORA.H:1159
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:1180
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vwind
Wind in the v direction, defined at rho-points.
Definition REMORA.H:322
std::unique_ptr< ProblemBase > prob
Pointer to container of analytical functions for problem definition.
Definition REMORA.H:1259
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:393
void Construct_REMORAFillPatchers(int lev)
Construct FillPatchers.
Definition REMORA.cpp:374
static int sum_interval
Diagnostic sum output interval in number of steps.
Definition REMORA.H:1442
int history_count
Counter for which time index we are writing to in the netcdf history file.
Definition REMORA.H:1390
amrex::Real stop_time
Time to stop.
Definition REMORA.H:1335
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rain
precipitation rate [kg/m^2/s]
Definition REMORA.H:351
int do_substep
Whether to substep fine levels in time.
Definition REMORA.H:1363
void Evolve()
Advance solution to final time.
Definition REMORA.cpp:163
std::string bdry_time_varname
Default name of time field for boundary data.
Definition REMORA.H:1482
amrex::Real plotfile_fill_value
fill value for masked arrays in amrex plotfiles
Definition REMORA.H:1409
void ReadCheckpointFile()
read checkpoint file from disk
bool chunk_history_file
Whether to chunk netcdf history file.
Definition REMORA.H:1383
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:315
amrex::Real get_t_old(int lev) const
Accessor method for t_old to expose to outside classes.
Definition REMORA.cpp:1653
int yvel_bc() const noexcept
Definition REMORA.H:1053
std::unique_ptr< NCTimeSeries > longwave_down_data_from_file
Data container for downward longwave radiation flux read from file.
Definition REMORA.H:1178
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:270
amrex::Vector< std::string > datalogname
Definition REMORA.H:1590
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:248
void init_only(int lev, amrex::Real time)
Init (NOT restart or regrid)
Definition REMORA.cpp:1052
void init_set_vmix(int lev)
Initialize vertical mixing coefficients from file or analytic.
Definition REMORA.cpp:607
std::unique_ptr< NCTimeSeries > v_clim_data_from_file
Data container for v-velocity climatology data read from file.
Definition REMORA.H:1193
std::string clim_u_time_varname
Name of time field for u climatology data.
Definition REMORA.H:1491
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:578
int foextrap_bc() const noexcept
Definition REMORA.H:1060
amrex::Vector< REMORAFillPatcher > FPr_u
Vector over levels of FillPatchers for u (3D)
Definition REMORA.H:1228
std::string clim_temp_time_varname
Name of time field for temperature climatology data.
Definition REMORA.H:1497
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1266
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Holds info for dynamically generated tagging criteria.
Definition REMORA.H:1507
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
Definition REMORA.H:280
std::unique_ptr< NCTimeSeriesRiver > river_source_transportbar
Data container for vertically integrated momentum transport in rivers.
Definition REMORA.H:1204
std::string clim_ubar_time_varname
Name of time field for ubar climatology data.
Definition REMORA.H:1487
std::unique_ptr< NCTimeSeries > u_clim_data_from_file
Data container for u-velocity climatology data read from file.
Definition REMORA.H:1191
std::string check_file
Checkpoint file prefix.
Definition REMORA.H:1376
static amrex::Real startCPUTime
Variable for CPU timing.
Definition REMORA.H:1556
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:235
amrex::Real start_time
Time of the start of the simulation, in seconds.
Definition REMORA.H:1338
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:246
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_uwind
Wind in the u direction, defined at rho-points.
Definition REMORA.H:320
static amrex::Real fixed_fast_dt
User specified fixed barotropic time step.
Definition REMORA.H:1356
int regrid_int
how often each level regrids the higher levels of refinement (after a level advances that many time s...
Definition REMORA.H:1366
amrex::Real check_int_time
Checkpoint output interval in seconds.
Definition REMORA.H:1380
void init_scalar_metadata()
Build runtime scalar names after nscalar is known.
Definition REMORA.cpp:150
void Define_REMORAFillPatchers(int lev)
Define FillPatchers.
Definition REMORA.cpp:423
amrex::Vector< amrex::IntVect > cum_ref_ratios
Cumulative refinement ratio between level 0 and level i.
Definition REMORA.H:1469
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:282
amrex::Real plot_int_time
Plotfile output interval in seconds.
Definition REMORA.H:1374
amrex::Vector< int > num_files_at_level
how many netcdf input files specified at each level
Definition REMORA.H:1264
amrex::Vector< REMORAFillPatcher > FPr_vbar
Vector over levels of FillPatchers for vbar (2D)
Definition REMORA.H:1238
void AverageDownTo(int crse_lev)
more flexible version of AverageDown() that lets you average down across multiple levels
Definition REMORA.cpp:1608
int steps_per_history_file
Number of time steps per netcdf history file.
Definition REMORA.H:1388
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Called after every level 0 timestep.
Definition REMORA.cpp:237
int max_step
maximum number of steps
Definition REMORA.H:1333
amrex::Vector< amrex::MultiFab * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
Definition REMORA.H:239
std::unique_ptr< NCTimeSeries > svstr_data_from_file
Data container for v-component surface momentum flux read from file.
Definition REMORA.H:1164
amrex::Vector< std::string > nc_frc_file
NetCDF forcing file(s)
Definition REMORA.H:1472
amrex::Vector< int > num_boxes_at_level
how many boxes specified at each level by tagging criteria
Definition REMORA.H:1262
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:244
void refinement_criteria_setup()
Set refinement criteria.
int last_check_file_step
Step when we last output a checkpoint file.
Definition REMORA.H:1323
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:1489
amrex::Vector< int > nsubsteps
How many substeps on each level?
Definition REMORA.H:1271
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
Definition REMORA.H:1283
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:1184
int plot_int
Plotfile output interval in iterations.
Definition REMORA.H:1372
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:1182
void InitData()
Initialize multilevel data.
Definition REMORA.cpp:269
void set3DPlotVariables(const std::string &pp_plot_var_names_3d)
amrex::Vector< int > istep
which step?
Definition REMORA.H:1269
void WriteCheckpointFile()
write checkpoint file to disk
std::string nc_clim_coeff_file
NetCDF climatology coefficient file.
Definition REMORA.H:1479
void setRecordDataInfo(int i, const std::string &filename)
Definition REMORA.H:1576
void set_analytic_vmix(int lev)
Set vertical mixing coefficients from analytic.
Definition REMORA.cpp:624
void init_flat_bathymetry(int lev)
Initialize flat bathymetry to value from problo.
Definition REMORA.cpp:951
std::unique_ptr< NCTimeSeries > temp_clim_data_from_file
Data container for temperature climatology data read from file.
Definition REMORA.H:1195
amrex::Vector< std::string > bdry_time_name_byvar
Name of time fields for boundary data.
Definition REMORA.H:1484
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
Definition REMORA.H:1447
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:284
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:317
std::unique_ptr< NCTimeSeries > srflx_data_from_file
Data container for shortwave radiation flux read from file.
Definition REMORA.H:1176
REMORA()
Definition REMORA.cpp:64
void set_zeta(int lev)
Initialize zeta from file or analytic.
Definition REMORA.cpp:483
static amrex::Real change_max
Fraction maximum change in subsequent time steps.
Definition REMORA.H:1352
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:230
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h_full_domain
Bathymetry data on the whole domain at each potential level.
Definition REMORA.H:257
void set_bathymetry(int lev)
Initialize bathymetry from file or analytic.
Definition REMORA.cpp:511
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:237
std::unique_ptr< NCTimeSeries > ubar_clim_data_from_file
Data container for ubar climatology data read from file.
Definition REMORA.H:1187
amrex::Vector< REMORAFillPatcher > FPr_c
Vector over levels of FillPatchers for scalars.
Definition REMORA.H:1226
std::unique_ptr< NCTimeSeries > Tair_data_from_file
Data container for air temperature read from file.
Definition REMORA.H:1170
int nscalar
Number of passive scalars carried in the state.
Definition REMORA.H:1344
std::string clim_v_time_varname
Name of time field for v climatology data.
Definition REMORA.H:1493
amrex::Vector< REMORAFillPatcher > FPr_w
Vector over levels of FillPatchers for w.
Definition REMORA.H:1232
void average_down_with_grow_cells(int lev, amrex::Vector< std::unique_ptr< amrex::MultiFab > > &mf)
Average down from level lev+1 to lev in mf, including grow cells.
Definition REMORA.cpp:1636
std::string clim_salt_time_varname
Name of time field for salinity climatology data.
Definition REMORA.H:1495
std::unique_ptr< NCTimeSeries > Uwind_data_from_file
Data container for u-direction wind read from file.
Definition REMORA.H:1166
std::unique_ptr< NCTimeSeries > Pair_data_from_file
Data container for air pressure read from file.
Definition REMORA.H:1174
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1273
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:1403
void set_masks(int lev)
Initialize land-sea masks from file or analytic.
Definition REMORA.cpp:640
amrex::Vector< amrex::Vector< std::unique_ptr< NCTimeSeriesBoundary > > > boundary_series
Vector over BdyVars of boundary series data containers.
Definition REMORA.H:1207
int cf_set_width
Width for fixing values at coarse-fine interface.
Definition REMORA.H:1223
void ReadParameters()
read in some parameters from inputs file
Definition REMORA.cpp:1328
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:266
void sum_integrated_quantities(amrex::Real time)
Integrate conserved quantities for diagnostics.
static int total_nc_plot_file_step
Definition REMORA.H:1067
static int plot_staggered_vels
Whether to write the staggered velocities (not averaged to cell centers)
Definition REMORA.H:1450
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:335
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
Definition REMORA.H:390
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:388
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:1406
int plot_file_on_restart
Whether to output a plotfile on restart from checkpoint.
Definition REMORA.H:1327
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:386
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:233
std::string frc_time_varname
Name of time field for forcing data.
Definition REMORA.H:1501
void set_wind(int lev)
Initialize or calculate wind speed from file or analytic.
Definition REMORA.cpp:980
amrex::Vector< REMORAFillPatcher > FPr_ubar
Vector over levels of FillPatchers for ubar (2D)
Definition REMORA.H:1236
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:1168
void init_bathymetry_full_domain_from_netcdf()
Full domain bathymetry data initialization from NetCDF file.
void set_hmixcoef(int lev)
Initialize horizontal mixing coefficients.
Definition REMORA.cpp:666
amrex::Vector< std::unique_ptr< NCTimeSeriesRiver > > river_source_cons
Vector of data containers for scalar data in rivers.
Definition REMORA.H:1200
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:1202
void AverageDown()
set covered coarse cells to be the average of overlying fine cells
Definition REMORA.cpp:1595
static int fixed_ndtfast_ratio
User specified, number of barotropic steps per baroclinic step.
Definition REMORA.H:1358
amrex::Real netcdf_fill_value
fill value for masked arrays in netcdf output
Definition REMORA.H:1411
void timeStepML(amrex::Real time, int iteration)
advance all levels by dt, loops over finer levels
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:406
amrex::Vector< std::unique_ptr< std::fstream > > datalog
Definition REMORA.H:1589
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
Definition REMORA.H:278
static amrex::Real cfl
CFL condition.
Definition REMORA.H:1350
void append3DPlotVariables(const std::string &pp_plot_var_names_3d)
static int verbose
Verbosity level of output.
Definition REMORA.H:1439
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_cloud
cloud cover fraction [0-1], defined at rho-points
Definition REMORA.H:355
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1370
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:312
std::unique_ptr< NCTimeSeries > salt_clim_data_from_file
Data container for salinity climatology data read from file.
Definition REMORA.H:1197
int check_int
Checkpoint output interval in iterations.
Definition REMORA.H:1378
void set_smflux(int lev)
Initialize or calculate surface momentum flux from file or analytic.
Definition REMORA.cpp:961
void WritePlotFile(int istep)
main driver for writing AMReX plotfiles
std::string restart_chkfile
If set, restart from this checkpoint file.
Definition REMORA.H:1341
void init_clim_nudg_coeff(int lev)
Wrapper to initialize climatology nudging coefficient.
void init_bathymetry_full_domain_from_analytic()
Full domain bathymetry data initialization from analytic.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:268
int cf_width
Nudging width at coarse-fine interface.
Definition REMORA.H:1221
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:1318
amrex::Vector< amrex::Real > t_old
old time at each level
Definition REMORA.H:1275
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_srflx
Shortwave radiation flux [W/m²], defined at rho-points.
Definition REMORA.H:331
amrex::Real last_check_file_time
Simulation time when we last output a checkpoint file.
Definition REMORA.H:1325
void append2DPlotVariables(const std::string &pp_plot_var_names_2d)
void set_bathymetry_averaged_down(int lev)
Copy over bathymetry data that has been averaged down from high resolution input netcdf file.
Definition REMORA.cpp:562
amrex::Vector< amrex::Real > dt
time step at each level
Definition REMORA.H:1277
static amrex::Real sum_per
Diagnostic sum output interval in time.
Definition REMORA.H:1444
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Pair
Air pressure [mb], defined at rho-points.
Definition REMORA.H:328
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:326
int hires_grid_level
Which level the high resolution bathymetry is at.
Definition REMORA.H:1463
void restart()
Definition REMORA.cpp:470
std::unique_ptr< NCTimeSeries > vbar_clim_data_from_file
Data container for vbar climatology data read from file.
Definition REMORA.H:1189
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Tair
Air temperature [°C], defined at rho-points.
Definition REMORA.H:324
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
Definition REMORA.H:286
static constexpr int cons_bc
static constexpr int Temp_bc_comp
const char * buildInfoGetGitHash(int i)
HorizMixingType horiz_mixing_type
amrex::Real Akv_bak
amrex::Vector< amrex::Real > tnu2
std::string longwave_netcdf_varname
amrex::Vector< int > do_rivers_cons
ScaledToGridAMRScaling scaled_to_grid_amr_scaling
amrex::Real Akt_bak
amrex::Real visc2
void init_params(int ncons)
read in and initialize parameters
SMFluxType smflux_type
VertMixingType vert_mixing_type
CouplingType coupling_type