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
56/**
57 * constructor:
58 * - reads in parameters from inputs file
59 * - sizes multilevel arrays and data structures
60 * - initializes BCRec boundary condition object
61 */
63{
64 BL_PROFILE("REMORA::REMORA()");
65
66 if (ParallelDescriptor::IOProcessor()) {
67 const char* remora_hash = amrex::buildInfoGetGitHash(1);
68 const char* amrex_hash = amrex::buildInfoGetGitHash(2);
69 const char* buildgithash = amrex::buildInfoGetBuildGitHash();
70 const char* buildgitname = amrex::buildInfoGetBuildGitName();
71
72 if (strlen(remora_hash) > 0) {
73 amrex::Print() << "\n"
74 << "REMORA git hash: " << remora_hash << "\n";
75 }
76 if (strlen(amrex_hash) > 0) {
77 amrex::Print() << "AMReX git hash: " << amrex_hash << "\n";
78 }
79 if (strlen(buildgithash) > 0) {
80 amrex::Print() << buildgitname << " git hash: " << buildgithash << "\n";
81 }
82
83 amrex::Print() << "\n";
84 }
85
87
88 const std::string& pv3d = "plot_vars_3d"; set3DPlotVariables(pv3d);
89 const std::string& pv2d = "plot_vars_2d"; set2DPlotVariables(pv2d);
90
91 prob = amrex_probinit(geom[0].ProbLo(),geom[0].ProbHi());
92
93 // Geometry on all levels has been defined already.
94
95 // No valid BoxArray and DistributionMapping have been defined.
96 // But the arrays for them have been resized.
97
98 int nlevs_max = max_level + 1;
99
100 istep.resize(nlevs_max, 0);
101 nsubsteps.resize(nlevs_max, 1);
102 for (int lev = 1; lev <= max_level; ++lev) {
103 nsubsteps[lev] = do_substep ? MaxRefRatio(lev-1) : 1;
104 }
105
106 physbcs.resize(nlevs_max);
107
108 t_new.resize(nlevs_max, 0.0_rt);
109 t_old.resize(nlevs_max, -1.e100_rt);
110 dt.resize(nlevs_max, 1.e100_rt);
111
112 cons_new.resize(nlevs_max);
113 cons_old.resize(nlevs_max);
114 xvel_new.resize(nlevs_max);
115 xvel_old.resize(nlevs_max);
116 yvel_new.resize(nlevs_max);
117 yvel_old.resize(nlevs_max);
118 zvel_new.resize(nlevs_max);
119 zvel_old.resize(nlevs_max);
120
121 advflux_reg.resize(nlevs_max);
122
123 // Initialize tagging criteria for mesh refinement
125
126 IntVect cum_ref_ratio = IntVect(1,1,0);
127 cum_ref_ratios.push_back(cum_ref_ratio);
128 // We have already read in the ref_Ratio (via amr.ref_ratio =) but we need to enforce
129 // that there is no refinement in the vertical so we test on that here.
130 for (int lev = 0; lev < max_level; ++lev)
131 {
132 amrex::Print() << "Refinement ratio at level " << lev << " set to be " <<
133 ref_ratio[lev][0] << " " << ref_ratio[lev][1] << " " << ref_ratio[lev][2] << std::endl;
134
135 if (ref_ratio[lev][2] != 1)
136 {
137 amrex::Print() << "********************************************************************************" << std::endl;
138 amrex::Print() << "We don't allow refinement in the vertical -- make sure to set ref_ratio = 1 in z" << std::endl;
139 amrex::Print() << "It's possible you set amr.ref_ratio when you meant to set amr.ref_ratio_vect " << std::endl;
140 amrex::Print() << "********************************************************************************" << std::endl;
141 amrex::Abort();
142 }
143
144 cum_ref_ratio[0] *= ref_ratio[lev][0];
145 cum_ref_ratio[1] *= ref_ratio[lev][1];
146 cum_ref_ratios.push_back(cum_ref_ratio);
147 }
148}
149
150REMORA::REMORA (const amrex::RealBox& rb, int max_level_in, const amrex::Vector<int>& n_cell_in, int coord, const amrex::Vector<amrex::IntVect>& ref_ratio_in, const amrex::Array<int,AMREX_SPACEDIM>& is_per, std::string prefix)
151 : amrex::AmrCore (rb, max_level_in, n_cell_in, coord, ref_ratio_in, is_per)
152{
153 BL_PROFILE("REMORA::REMORA(explicit)");
154 pp_prefix = prefix;
155
156 if (ParallelDescriptor::IOProcessor()) {
157 const char* remora_hash = amrex::buildInfoGetGitHash(1);
158 const char* amrex_hash = amrex::buildInfoGetGitHash(2);
159 const char* buildgithash = amrex::buildInfoGetBuildGitHash();
160 const char* buildgitname = amrex::buildInfoGetBuildGitName();
161
162 if (strlen(remora_hash) > 0) {
163 amrex::Print() << "\n"
164 << "REMORA git hash: " << remora_hash << "\n";
165 }
166 if (strlen(amrex_hash) > 0) {
167 amrex::Print() << "AMReX git hash: " << amrex_hash << "\n";
168 }
169 if (strlen(buildgithash) > 0) {
170 amrex::Print() << buildgitname << " git hash: " << buildgithash << "\n";
171 }
172
173 amrex::Print() << "\n";
174 }
175
177
178 const std::string& pv3d = "plot_vars_3d"; set3DPlotVariables(pv3d);
179 const std::string& pv2d = "plot_vars_2d"; set2DPlotVariables(pv2d);
180
181 prob = amrex_probinit(geom[0].ProbLo(),geom[0].ProbHi());
182
183 int nlevs_max = max_level + 1;
184
185 istep.resize(nlevs_max, 0);
186 nsubsteps.resize(nlevs_max, 1);
187 for (int lev = 1; lev <= max_level; ++lev) {
188 nsubsteps[lev] = do_substep ? MaxRefRatio(lev-1) : 1;
189 }
190
191 physbcs.resize(nlevs_max);
192
193 t_new.resize(nlevs_max, 0.0_rt);
194 t_old.resize(nlevs_max, -1.e100_rt);
195 dt.resize(nlevs_max, 1.e100_rt);
196
197 cons_new.resize(nlevs_max);
198 cons_old.resize(nlevs_max);
199 xvel_new.resize(nlevs_max);
200 xvel_old.resize(nlevs_max);
201 yvel_new.resize(nlevs_max);
202 yvel_old.resize(nlevs_max);
203 zvel_new.resize(nlevs_max);
204 zvel_old.resize(nlevs_max);
205
206 advflux_reg.resize(nlevs_max);
207
209
210 for (int lev = 0; lev < max_level; ++lev)
211 {
212 amrex::Print() << "Refinement ratio at level " << lev << " set to be " <<
213 ref_ratio[lev][0] << " " << ref_ratio[lev][1] << " " << ref_ratio[lev][2] << std::endl;
214
215 if (ref_ratio[lev][2] != 1)
216 {
217 amrex::Print() << "********************************************************************************" << std::endl;
218 amrex::Print() << "We don't allow refinement in the vertical -- make sure to set ref_ratio = 1 in z" << std::endl;
219 amrex::Print() << "It's possible you set amr.ref_ratio when you meant to set amr.ref_ratio_vect " << std::endl;
220 amrex::Print() << "********************************************************************************" << std::endl;
221 amrex::Abort();
222 }
223 }
224}
225
227{
228}
229
230void
232{
233 cons_names.clear();
234 cons_names.reserve(ncons);
235 cons_names.emplace_back("temp");
236 cons_names.emplace_back("salt");
237 cons_names.emplace_back("tracer");
238 for (int i = 1; i < nscalar; ++i) {
239 cons_names.emplace_back("tracer_" + std::to_string(i));
240 }
241}
242
243void
245{
246 BL_PROFILE("REMORA::Evolve()");
247 Real cur_time = t_new[0];
248
249 // Take one coarse timestep by calling timeStep -- which recursively calls timeStep
250 // for finer levels (with or without subcycling)
251 for (int step = istep[0]; step < max_step && cur_time < stop_time; ++step)
252 {
253 amrex::Print() << "\nCoarse STEP " << step+1 << " starts ..." << std::endl;
254
255 ComputeDt();
256
257 int lev = 0;
258 int iteration = 1;
259 auto dEvolveTime0 = amrex::second();
260
261 if (max_level == 0) {
262 timeStep(lev, cur_time, iteration);
263 }
264 else {
265 timeStepML(cur_time, iteration);
266 }
267
268 cur_time += dt[0];
269
270 amrex::Print() << "Coarse STEP " << step+1 << " ends." << " TIME = " << cur_time
271 << " DT = " << dt[0] << std::endl;
272
273 if (verbose > 0)
274 {
275 auto dEvolveTime = amrex::second() - dEvolveTime0;
276 ParallelDescriptor::ReduceRealMax(dEvolveTime,ParallelDescriptor::IOProcessorNumber());
277 amrex::Print() << "Timestep time = " << dEvolveTime << " seconds." << '\n';
278 }
279
280 if ( (plot_int > 0 && (step+1 - last_plot_file_step) == plot_int ) ||
281 (plot_int_time > 0 && (cur_time >= (last_plot_file_time + plot_int_time))) )
282 {
283 last_plot_file_step = step+1;
284 last_plot_file_time = cur_time;
285 WritePlotFile(step+1);
287 }
288
289 if ((check_int > 0 && (step+1 - last_check_file_step) == check_int)
290 || (check_int_time > 0 && cur_time >= (last_check_file_time + check_int_time))) {
291 last_check_file_step = step+1;
292 last_check_file_time = cur_time;
294 }
295
296 post_timestep(step, cur_time, dt[0]);
297
298#ifdef AMREX_MEM_PROFILING
299 {
300 std::ostringstream ss;
301 ss << "[STEP " << step+1 << "]";
302 MemProfiler::report(ss.str());
303 }
304#endif
305
306 if (cur_time >= stop_time - 1.e-6*dt[0]) break;
307 }
308
309 if ( (plot_int > 0 || plot_int_time > 0.0) && istep[0] > last_plot_file_step)
310 {
313 }
314
315 if ((check_int > 0 || check_int_time > 0.0) && istep[0] > last_check_file_step) {
317 }
318}
319
320/**
321 * @param[in ] nstep which step we're on
322 * @param[in ] time current time
323 * @param[in ] dt_lev0 time step on level 0
324 */
325void
326REMORA::post_timestep (int nstep, Real time, Real dt_lev0)
327{
328 BL_PROFILE("REMORA::post_timestep()");
329
330#ifdef REMORA_USE_PARTICLES
331 particleData.Redistribute();
332#endif
333
335 {
336 for (int lev = finest_level-1; lev >= 0; lev--)
337 {
338 // This call refluxes from the lev/lev+1 interface onto lev
339 //getAdvFluxReg(lev+1)->Reflux(*cons_new[lev], 0, 0, NCONS);
340
341 // We need to do this before anything else because refluxing changes the
342 // values of coarse cells underneath fine grids with the assumption they'll
343 // be over-written by averaging down
344 //
345 AverageDownTo(lev);
346 }
347 }
348
349 if (is_it_time_for_action(nstep, time, dt_lev0, sum_interval, sum_per)) {
351 }
352}
353
354/**
355 * This is called from main.cpp and handles all initialization, whether from start or restart
356 */
357void
359{
360 BL_PROFILE("REMORA::InitData()");
361 // Initialize the start time for our CPU-time tracker
362 startCPUTime = Real(ParallelDescriptor::second());
363
364 // Map the words in the inputs file to BC types, then translate
365 // those types into what they mean for each variable
366 init_bcs();
367
368 // Init vertical stretching coeffs
370
373 last_plot_file_time = -1.0_rt;
374 last_check_file_time = -1.0_rt;
375
376 if (restart_chkfile == "") {
377 // start simulation from the beginning
378
379 InitFromScratch(start_time);
380
382 AverageDown();
383 }
384
385 } else { // Restart from a checkpoint
386
387 restart();
388
389 }
390#ifdef REMORA_USE_MOAB
391 InitMOABMesh();
392#endif
393 // Initialize flux registers (whether we start from scratch or restart)
395 advflux_reg[0] = nullptr;
396 for (int lev = 1; lev <= finest_level; lev++)
397 {
398 advflux_reg[lev].reset( new YAFluxRegister(grids[lev], grids[lev-1],
399 dmap[lev], dmap[lev-1],
400 geom[lev], geom[lev-1],
401 ref_ratio[lev-1], lev, ncons));
402 }
403 }
404
405 // Fill ghost cells/faces
406 for (int lev = 0; lev <= finest_level; ++lev)
407 {
408 if (lev > 0 && cf_width >= 0) {
410 }
411
412 if (restart_chkfile == "") {
413 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]);
414 FillPatch(lev, t_new[lev], *xvel_new[lev], xvel_new, xvel_bc(), BdyVars::u, 0, true, false,0,0,0.0,*xvel_new[lev]);
415 FillPatch(lev, t_new[lev], *yvel_new[lev], yvel_new, yvel_bc(), BdyVars::v, 0, true, false,0,0,0.0,*yvel_new[lev]);
416 FillPatch(lev, t_new[lev], *zvel_new[lev], zvel_new, zvel_bc(), BdyVars::null, 0, true, false);
417
418 // Copy from new into old just in case when initializing from scratch
419 int ngs = cons_new[lev]->nGrow();
420 int ngvel = xvel_new[lev]->nGrow();
421 MultiFab::Copy(*cons_old[lev],*cons_new[lev],0,0,ncons,ngs);
422 MultiFab::Copy(*xvel_old[lev],*xvel_new[lev],0,0,1,ngvel);
423 MultiFab::Copy(*yvel_old[lev],*yvel_new[lev],0,0,1,ngvel);
424 MultiFab::Copy(*zvel_old[lev],*zvel_new[lev],0,0,1,IntVect(ngvel,ngvel,0));
425 }
426 } // lev
427
428 // Check for additional plotting variables that are available after
429 // particle containers are setup.
430 const std::string& pv3d = "plot_vars_3d"; append3DPlotVariables(pv3d);
431 const std::string& pv2d = "plot_vars_2d"; append2DPlotVariables(pv2d);
432
433 if (restart_chkfile == "" && (check_int > 0 || check_int_time > 0.0_rt))
434 {
437 }
438
439 if ( (restart_chkfile == "") ||
441 {
442 if (plot_int > 0 || plot_int_time > 0.0)
443 {
444 int step0 = 0;
445 WritePlotFile(step0);
448 }
449 }
450
453 }
454
455 ComputeDt();
456
457}
458
459/**
460 * @param[in ] lev level to operate on
461 */
462void
464{
465 BL_PROFILE("REMORA::Construct_REMORAFillPatchers()");
466 amrex::Print() << ":::Construct_REMORAFillPatchers " << lev << std::endl;
467
468 auto& ba_fine = cons_new[lev ]->boxArray();
469 auto& ba_crse = cons_new[lev-1]->boxArray();
470 auto& dm_fine = cons_new[lev ]->DistributionMap();
471 auto& dm_crse = cons_new[lev-1]->DistributionMap();
472
473 BoxList bl2d_fine = ba_fine.boxList();
474 for (auto& b : bl2d_fine) {
475 b.setRange(2,0);
476 }
477 BoxArray ba2d_fine(std::move(bl2d_fine));
478
479 BoxList bl2d_crse = ba_crse.boxList();
480 for (auto& b : bl2d_crse) {
481 b.setRange(2,0);
482 }
483 BoxArray ba2d_crse(std::move(bl2d_crse));
484
485 int ncomp = cons_new[lev]->nComp();
486
487 FPr_c.emplace_back(ba_fine, dm_fine, geom[lev] ,
488 ba_crse, dm_crse, geom[lev-1],
489 -cf_width, -cf_set_width, ncomp, &cell_cons_interp);
490 FPr_u.emplace_back(convert(ba_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
491 convert(ba_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
492 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
493 FPr_v.emplace_back(convert(ba_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
494 convert(ba_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
495 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
496 FPr_w.emplace_back(convert(ba_fine, IntVect(0,0,1)), dm_fine, geom[lev] ,
497 convert(ba_crse, IntVect(0,0,1)), dm_crse, geom[lev-1],
498 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
499
500 FPr_ubar.emplace_back(convert(ba2d_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
501 convert(ba2d_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
502 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
503 FPr_vbar.emplace_back(convert(ba2d_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
504 convert(ba2d_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
505 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
506}
507
508/**
509 * @param[in ] lev level to operate on
510 */
511void
513{
514 BL_PROFILE("REMORA::Define_REMORAFillPatchers()");
515 amrex::Print() << ":::Define_REMORAFillPatchers " << lev << std::endl;
516
517 auto& ba_fine = cons_new[lev ]->boxArray();
518 auto& ba_crse = cons_new[lev-1]->boxArray();
519 auto& dm_fine = cons_new[lev ]->DistributionMap();
520 auto& dm_crse = cons_new[lev-1]->DistributionMap();
521
522 BoxList bl2d_fine = ba_fine.boxList();
523 for (auto& b : bl2d_fine) {
524 b.setRange(2,0);
525 }
526 BoxArray ba2d_fine(std::move(bl2d_fine));
527
528 BoxList bl2d_crse = ba_crse.boxList();
529 for (auto& b : bl2d_crse) {
530 b.setRange(2,0);
531 }
532 BoxArray ba2d_crse(std::move(bl2d_crse));
533
534
535 int ncomp = cons_new[lev]->nComp();
536
537 FPr_c[lev-1].Define(ba_fine, dm_fine, geom[lev] ,
538 ba_crse, dm_crse, geom[lev-1],
539 -cf_width, -cf_set_width, ncomp, &cell_cons_interp);
540 FPr_u[lev-1].Define(convert(ba_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
541 convert(ba_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
542 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
543 FPr_v[lev-1].Define(convert(ba_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
544 convert(ba_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
545 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
546 FPr_w[lev-1].Define(convert(ba_fine, IntVect(0,0,1)), dm_fine, geom[lev] ,
547 convert(ba_crse, IntVect(0,0,1)), dm_crse, geom[lev-1],
548 -cf_width, -cf_set_width, 1, &face_cons_linear_interp);
549
550 FPr_ubar[lev-1].Define(convert(ba2d_fine, IntVect(1,0,0)), dm_fine, geom[lev] ,
551 convert(ba2d_crse, IntVect(1,0,0)), dm_crse, geom[lev-1],
552 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
553 FPr_vbar[lev-1].Define(convert(ba2d_fine, IntVect(0,1,0)), dm_fine, geom[lev] ,
554 convert(ba2d_crse, IntVect(0,1,0)), dm_crse, geom[lev-1],
555 -cf_width, -cf_set_width, 3, &face_cons_linear_interp);
556}
557
558void
560{
561 BL_PROFILE("REMORA::restart()");
563
564 // We set this here so that we don't over-write the checkpoint file we just started from
566}
567
568/**
569 * @param[in ] lev level to operate on
570 */
571void
573{
574 BL_PROFILE("REMORA::set_zeta()");
575 if (lev==0) {
576 if (hires_init_level < 0) {
578 prob->init_analytic_zeta(lev, geom[lev], solverChoice, *this, *vec_zeta[lev]);
579 } else if (solverChoice.ic_type == IC_Type::netcdf) {
580#ifdef REMORA_USE_NETCDF
581 amrex::Print() << "Calling init_zeta_from_netcdf on level " << lev << std::endl;
583 amrex::Print() << "Sea surface height loaded from netcdf file \n " << std::endl;
584#endif
585 } else {
586 amrex::Abort("Unknown IC_Type");
587 }
588 } else {
590 }
591 vec_zeta[lev]->FillBoundary(geom[lev].periodicity());
592 } else {
593 // If our level is higher than the high resolution grid or initialization
594 // is analytic, interpolate from level below. Otherwise, copy over the bathymetry
595 // data that has been averaged down
596 if (lev > hires_init_level) {
597 Real dummy_time = 0.0_rt;
598 FillCoarsePatch(lev,dummy_time,vec_zeta[lev].get(), vec_zeta[lev-1].get(),BCVars::cons_bc);
599 } else {
601 vec_zeta[lev]->FillBoundary(geom[lev].periodicity());
602 }
603 }
604 set_zeta_average(lev);
605}
606
607/**
608 * @param[in ] lev level to operate on
609 */
610void
612{
613 BL_PROFILE("REMORA::bathymetry()");
614 // Only set bathymetry on level 0, and interpolate for finer levels
615 if (lev==0) {
618 // If grid data is not defined on a level > 0 (negative level) then
619 // initialize from low-resolution grid normally. Otherwise use high-resolution
620 // grid data averaged down to level 0
621 } else if (hires_grid_level < 0) {
623 prob->init_analytic_bathymetry(lev, geom[lev], solverChoice, *this, *vec_h[lev]);
624 } else if (solverChoice.ic_type == IC_Type::netcdf) {
625#ifdef REMORA_USE_NETCDF
626 amrex::Print() << "Calling init_bathymetry_from_netcdf " << std::endl;
628 amrex::Print() << "Bathymetry loaded from netcdf file \n " << std::endl;
629 amrex::Print() << "Calling init_grid_vars_from_netcdf " << std::endl;
631 amrex::Print() << "Grid variables loaded from netcdf file \n " << std::endl;
632#endif
633 } else {
634 amrex::Abort("Unknown IC_Type");
635 }
636 } else {
639 }
640 // Need FillBoundary to fill at grid-grid boundaries, and EnforcePeriodicity
641 // to make sure ghost cells in the domain corners are consistent.
642 vec_h[lev]->FillBoundary(geom[lev].periodicity());
643 vec_h[lev]->EnforcePeriodicity(geom[lev].periodicity());
644 } else {
645 // If our level is higher than the high resolution grid or initialization
646 // is analytic, interpolate from level below. Otherwise, copy over the bathymetry
647 // data that has been averaged down
648 if (lev > hires_grid_level) {
649 Real dummy_time = 0.0_rt;
650 FillCoarsePatch(lev,dummy_time,vec_h[lev].get(), vec_h[lev-1].get(),BCVars::cons_bc);
651 } else {
653 vec_h[lev]->FillBoundary(geom[lev].periodicity());
654 vec_h[lev]->EnforcePeriodicity(geom[lev].periodicity());
655 }
656 }
657 set_grid_scale(lev);
658}
659
660/**
661 * @param[in ] lev level to operate on
662 */
663void
665 Real dummy_time = 0.0_rt;
666 // Note: don't understand why the grow vector args aren't vec_h and then vec_h_full_domain
667 ParallelCopy(*vec_h[lev].get(), *vec_h_full_domain[lev].get(), 0, 0, 1,vec_h_full_domain[lev]->nGrowVect(),vec_h[lev]->nGrowVect());
668 ParallelCopy(*vec_h[lev].get(), *vec_h_full_domain[lev].get(), 0, 1, 1,vec_h_full_domain[lev]->nGrowVect(),vec_h[lev]->nGrowVect());
669 FillPatch(lev,dummy_time,*vec_h[lev],GetVecOfPtrs(vec_h),
671 BdyVars::null,0,false,true,1);
672 FillPatch(lev,dummy_time,*vec_h[lev],GetVecOfPtrs(vec_h),
674 BdyVars::null,1,false,true,1);
675}
676
677/**
678 * @param[in ] lev level to operate on
679 */
680void
682 Real dummy_time = 0.0_rt;
683 ParallelCopy(*vec_pm[lev].get(), *vec_pm_full_domain[lev].get(), 0, 0, 1,
684 vec_pm_full_domain[lev]->nGrowVect(),vec_pm[lev]->nGrowVect());
685 ParallelCopy(*vec_pn[lev].get(), *vec_pn_full_domain[lev].get(), 0, 0, 1,
686 vec_pn_full_domain[lev]->nGrowVect(),vec_pn[lev]->nGrowVect());
687 FillPatch(lev,dummy_time,*vec_pm[lev],GetVecOfPtrs(vec_pm),
689 BdyVars::null,0,false,true);
690 FillPatch(lev,dummy_time,*vec_pn[lev],GetVecOfPtrs(vec_pn),
692 BdyVars::null,0,false,true);
693}
694
695/**
696 * @param[in ] lev level to operate on
697 */
698void
700 ParallelCopy(*vec_zeta[lev].get(), *vec_zeta_full_domain[lev].get(), 0, 0, 1,
701 vec_zeta_full_domain[lev]->nGrowVect(),vec_zeta[lev]->nGrowVect());
702 FillPatch(lev, t_new[lev], *vec_zeta[lev], GetVecOfPtrs(vec_zeta), zeta_bc(), BdyVars::zeta,
703 0, false,false,0,0,0.0,*vec_zeta[lev]);
704}
705
706/**
707 * @param[in ] lev level to operate on
708 */
709void
711 ParallelCopy(*cons_new[lev], *vec_cons_full_domain[lev], 0, 0, ncons,
712 vec_cons_full_domain[lev]->nGrowVect(),cons_new[lev]->nGrowVect());
713 ParallelCopy(*xvel_new[lev], *vec_xvel_full_domain[lev], 0, 0, 1,
714 vec_xvel_full_domain[lev]->nGrowVect(),xvel_new[lev]->nGrowVect());
715 ParallelCopy(*yvel_new[lev], *vec_yvel_full_domain[lev], 0, 0, 1,
716 vec_yvel_full_domain[lev]->nGrowVect(),yvel_new[lev]->nGrowVect());
717
718 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]);
719 FillPatch(lev, t_new[lev], *xvel_new[lev], xvel_new, xvel_bc(), BdyVars::u, 0, true, false,0,0,0.0,*xvel_new[lev]);
720 FillPatch(lev, t_new[lev], *yvel_new[lev], yvel_new, yvel_bc(), BdyVars::v, 0, true, false,0,0,0.0,*yvel_new[lev]);
721}
722
723/**
724 * @param[in ] lev level to operate on
725 */
726void
728 BL_PROFILE("REMORA::set_coriolis()");
731 prob->init_analytic_coriolis(lev, geom[lev], solverChoice, *this, *vec_fcor[lev]);
734#ifdef REMORA_USE_NETCDF
736 if (lev == 0) {
737 amrex::Print() << "Calling init_coriolis_from_netcdf " << std::endl;
739 amrex::Print() << "Coriolis loaded from netcdf file \n" << std::endl;
740 } else {
741 Real dummy_time = 0.0_rt;
742 FillCoarsePatch(lev,dummy_time,vec_fcor[lev].get(), vec_fcor[lev-1].get(),BCVars::cons_bc);
743 }
744#endif
745 } else {
746 Abort("Don't know this coriolis_type!");
747 }
748
749 Real time = 0.0_rt;
750 FillPatch(lev, time, *vec_fcor[lev], GetVecOfPtrs(vec_fcor), foextrap_bc());
751 vec_fcor[lev]->EnforcePeriodicity(geom[lev].periodicity());
752 }
753}
754
755void
757 BL_PROFILE("REMORA::init_set_vmix()");
762 // The GLS initialization just sets the multifab to a value, so there's
763 // no need to call FillPatch here
764 } else {
765 Abort("Don't know this vertical mixing type");
766 }
767}
768
769/**
770 * @param[in ] lev level to operate on
771 */
772void
774 BL_PROFILE("REMORA::set_analytic_vmix()");
775 Real time = 0.0_rt;
776 vec_Akv[lev]->setVal(solverChoice.Akv_bak);
777 vec_Akt[lev]->setVal(solverChoice.Akt_bak);
778 prob->init_analytic_vmix(lev, geom[lev], solverChoice, *this,*vec_Akv[lev], *vec_Akt[lev]);
779 FillPatch(lev, time, *vec_Akv[lev], GetVecOfPtrs(vec_Akv), zvel_bc(), BdyVars::null,0,true,false);
780 for (int n = 0; n < ncons; n++) {
781 FillPatch(lev, time, *vec_Akt[lev], GetVecOfPtrs(vec_Akt), zvel_bc(), BdyVars::null,n,false,false);
782 }
783}
784
785/**
786 * @param[in ] lev level to operate on
787 */
788void
790{
792 prob->init_analytic_masks(lev,geom[lev], solverChoice, *this, *vec_mskr[lev]);
795#ifdef REMORA_USE_NETCDF
796 if (lev == 0) {
797 amrex::Print() << "Calling init_masks_from_netcdf level " << lev << std::endl;
799 amrex::Print() << "Masks loaded from netcdf file \n " << std::endl;
800 } else {
801 Real dummy_time = 0.0_rt;
802 FillCoarsePatchPC(lev, dummy_time, vec_mskr[lev].get(), vec_mskr[lev-1].get(),
803 foextrap_bc());
805 }
806#endif
807 }
808 fill_3d_masks(lev);
809}
810
811/**
812 * @param[in ] lev level to operate on
813 */
814void
816{
817 BL_PROFILE("REMORA::set_hmixcoef()");
818
819 // Optional AMR scaling: decrease coefficients on refined levels linearly
820 // with grid size (i.e., proportional to sqrt(cell area)). For a horizontal
821 // refinement ratio rx x ry, the effective scale factor is 1/sqrt(rx*ry).
822 Real lev_scale = 1.0_rt;
824 Real rf = 1.0_rt;
825 for (int l = 0; l < lev; ++l) {
826 rf *= std::sqrt(static_cast<Real>(ref_ratio[l][0]) * static_cast<Real>(ref_ratio[l][1]));
827 }
828 lev_scale = 1.0_rt / rf;
829 }
830
832 prob->init_analytic_hmix(lev, geom[lev], solverChoice,
833 *this, *vec_visc2_p[lev], *vec_visc2_r[lev], *vec_diff2[lev]);
834
836 vec_visc2_p[lev]->setVal(solverChoice.visc2 * lev_scale);
837 vec_visc2_r[lev]->setVal(solverChoice.visc2 * lev_scale);
838 for (int n = 0; n < ncons; n++) {
839 vec_diff2[lev]->setVal(solverChoice.tnu2[n] * lev_scale, n, 1);
840 }
841
842 // Scale harmonic viscosity and diffusivity by the grid size as ROMS
843 // does in Utility/ini_hmixcoef.F. Intended for curvilinear grids.
844 //
845 // Define the ROMS grid factor (grdscl):
846 // G(i,j) = sqrt( 1 / (pm(i,j) * pn(i,j)) )
847 // = sqrt(cell area)
848 // Gmax = max over grid of G(i,j)
849 //
850 // Then horizontal harmonic mixing coefficients are scaled as:
851 // nu(i,j) = nu0 * G(i,j) / Gmax
852 // kappa_n(i,j) = kappa0 * G(i,j) / Gmax
853 //
854 // where:
855 // nu0 = solverChoice.visc2
856 // kappa0 = solverChoice.tnu2[n]
857 //
858 // This makes mixing strongest where grid spacing is largest.
859 //
860 // NOTE: The normalization (Gmax) is computed over the entire grid (ignoring masks).
861 // Therefore, if the largest cell area occurs over land, the maximum over *wet* cells
862 // (or in masked output files) may be smaller than the user-specified value.
863
865
866 // ------------------------------------------------------------
867 // Step 1: Compute grdmax over entire grid
868 // ------------------------------------------------------------
869 vec_visc2_r[lev]->setVal(solverChoice.visc2);
870 vec_visc2_p[lev]->setVal(solverChoice.visc2);
871 for (int n = 0; n < ncons; n++) {
872 vec_diff2[lev]->setVal(solverChoice.tnu2[n], n, 1);
873 }
874
875 // NOTE: This must be GPU-safe. Do not dereference MultiFab data on host.
876 // Force the reduction to run in the GPU launch region if GPUs are enabled.
877 // (If the launch region is disabled at runtime, ReduceMax may fall back to
878 // a host path that can try to read device-only data.)
879 amrex::Gpu::LaunchSafeGuard lsg(true);
880 Real denom_min = amrex::ReduceMin(*vec_pm[lev], *vec_pn[lev], 0,
881 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
882 Array4<Real const> const& pm,
883 Array4<Real const> const& pn) -> Real
884 {
885 Real local_min = 1.0e200_rt;
886 amrex::Loop(bx, [=,&local_min] (int i, int j, int) noexcept
887 {
888 local_min = amrex::min(local_min, pm(i,j,0) * pn(i,j,0));
889 });
890 return local_min;
891 });
892
893 ParallelDescriptor::ReduceRealMin(denom_min);
894 if (denom_min <= 0.0_rt) {
895 Abort("scaled_to_grid: found non-positive pm*pn (grid metrics must be > 0)");
896 }
897
898 Real grdmax = amrex::ReduceMax(*vec_pm[lev], *vec_pn[lev], 0,
899 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
900 Array4<Real const> const& pm,
901 Array4<Real const> const& pn) -> Real
902 {
903 Real local_max = 0.0_rt;
904 amrex::Loop(bx, [=,&local_max] (int i, int j, int) noexcept
905 {
906 Real denom = pm(i,j,0) * pn(i,j,0);
907 if (denom > 0.0_rt) {
908 Real G = std::sqrt(1.0_rt / denom);
909 local_max = amrex::max(local_max, G);
910 }
911 });
912 return local_max;
913 });
914
915 ParallelDescriptor::ReduceRealMax(grdmax);
916 if (grdmax <= 0.0_rt) {
917 Abort("scaled_to_grid: grdmax <= 0");
918 }
919
920 // Optional AMR scaling: decrease coefficients on refined levels linearly
921 // with grid size (i.e., proportional to sqrt(cell area)). For a horizontal
922 // refinement ratio rx x ry, the effective scale factor is 1/sqrt(rx*ry).
923 lev_scale = 1.0_rt;
925 Real rf = 1.0_rt;
926 for (int l = 0; l < lev; ++l) {
927 rf *= std::sqrt(static_cast<Real>(ref_ratio[l][0]) * static_cast<Real>(ref_ratio[l][1]));
928 }
929 lev_scale = 1.0_rt / rf;
930 }
931
932 Real visc0 = solverChoice.visc2 * lev_scale;
933 Real cff = visc0 / grdmax;
934
935 // ------------------------------------------------------------
936 // Step 2: Set rho coefficients everywhere
937 // ------------------------------------------------------------
938 amrex::Gpu::DeviceVector<Real> diff0_d(ncons);
939 amrex::Gpu::copy(amrex::Gpu::hostToDevice,
940 solverChoice.tnu2.begin(), solverChoice.tnu2.begin() + ncons,
941 diff0_d.begin());
942 Real const* diff0_ptr = diff0_d.data();
943
944 for (MFIter mfi(*vec_visc2_r[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
945 {
946 const Box& bx = mfi.validbox();
947 auto pm = vec_pm[lev]->const_array(mfi);
948 auto pn = vec_pn[lev]->const_array(mfi);
949 auto visc2_r = vec_visc2_r[lev]->array(mfi);
950 auto diff2 = vec_diff2[lev]->array(mfi);
951
952 int ncons_local = ncons;
953 ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int) noexcept
954 {
955 Real denom = pm(i,j,0) * pn(i,j,0);
956 Real grdscl = (denom > 0.0_rt) ? std::sqrt(1.0_rt / denom) : 0.0_rt;
957 visc2_r(i,j,0) = cff * grdscl;
958
959 for (int n = 0; n < ncons_local; n++) {
960 diff2(i,j,0,n) = ((diff0_ptr[n] * lev_scale) / grdmax) * grdscl;
961 }
962 });
963 }
964
965 // Fill ghost cells for rho coefficients BEFORE psi averaging
966 Real time = 0.0_rt;
967 FillPatch(lev, time, *vec_visc2_r[lev], GetVecOfPtrs(vec_visc2_r), foextrap_periodic_bc());
968
969 // ------------------------------------------------------------
970 // Step 3: Psi coefficients = average of 4 surrounding rho
971 // ------------------------------------------------------------
972 for (MFIter mfi(*vec_visc2_p[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
973 {
974 const Box& bx = mfi.validbox();
975 auto visc2_p = vec_visc2_p[lev]->array(mfi);
976 auto visc2_r = vec_visc2_r[lev]->const_array(mfi);
977
978 ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int) noexcept
979 {
980 visc2_p(i,j,0) = 0.25_rt * (
981 visc2_r(i-1,j-1,0) +
982 visc2_r(i ,j-1,0) +
983 visc2_r(i-1,j ,0) +
984 visc2_r(i ,j ,0)
985 );
986 });
987 }
988
989 FillPatch(lev, time, *vec_visc2_p[lev], GetVecOfPtrs(vec_visc2_p), foextrap_periodic_bc());
990
991 // Diagnostics
992 // NOTE: coefficients are computed everywhere (including land). Output routines may later
993 // mask land points (e.g., to FillValue in NetCDF/plotfiles), and analysis tools may
994 // additionally apply mask_rho (setting land to 0). Report both conventions.
995 //
996 // Global (MPI-reduced) extrema over all valid cells (no ghost).
997 Real visc_min_all = vec_visc2_r[lev]->min(0,0,false);
998 Real visc_max_all = vec_visc2_r[lev]->max(0,0,false);
999
1000 // Global extrema over *wet* rho points only, k=0.
1001 amrex::Gpu::LaunchSafeGuard lsg_diag(true);
1002 Real visc_min_wet = amrex::ReduceMin(*vec_visc2_r[lev], *vec_mskr[lev], 0,
1003 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
1004 Array4<Real const> const& visc2,
1005 Array4<Real const> const& mskr) -> Real
1006 {
1007 Real local_min = 1.0e200_rt;
1008 amrex::Loop(bx, [=,&local_min] (int i, int j, int) noexcept
1009 {
1010 if (mskr(i,j,0) > 0.0_rt) {
1011 local_min = amrex::min(local_min, visc2(i,j,0));
1012 }
1013 });
1014 return local_min;
1015 });
1016 ParallelDescriptor::ReduceRealMin(visc_min_wet);
1017
1018 Real visc_max_wet = amrex::ReduceMax(*vec_visc2_r[lev], *vec_mskr[lev], 0,
1019 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
1020 Array4<Real const> const& visc2,
1021 Array4<Real const> const& mskr) -> Real
1022 {
1023 Real local_max = -1.0e200_rt;
1024 amrex::Loop(bx, [=,&local_max] (int i, int j, int) noexcept
1025 {
1026 if (mskr(i,j,0) > 0.0_rt) {
1027 local_max = amrex::max(local_max, visc2(i,j,0));
1028 }
1029 });
1030 return local_max;
1031 });
1032 ParallelDescriptor::ReduceRealMax(visc_max_wet);
1033
1034 // Mimic "apply mask_rho" convention (dry -> 0).
1035 Real visc_min_mask0 = amrex::ReduceMin(*vec_visc2_r[lev], *vec_mskr[lev], 0,
1036 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
1037 Array4<Real const> const& visc2,
1038 Array4<Real const> const& mskr) -> Real
1039 {
1040 Real local_min = 1.0e200_rt;
1041 amrex::Loop(bx, [=,&local_min] (int i, int j, int) noexcept
1042 {
1043 const Real v = (mskr(i,j,0) > 0.0_rt) ? visc2(i,j,0) : 0.0_rt;
1044 local_min = amrex::min(local_min, v);
1045 });
1046 return local_min;
1047 });
1048 ParallelDescriptor::ReduceRealMin(visc_min_mask0);
1049
1050 Real visc_max_mask0 = amrex::ReduceMax(*vec_visc2_r[lev], *vec_mskr[lev], 0,
1051 [=] AMREX_GPU_HOST_DEVICE (Box const& bx,
1052 Array4<Real const> const& visc2,
1053 Array4<Real const> const& mskr) -> Real
1054 {
1055 Real local_max = -1.0e200_rt;
1056 amrex::Loop(bx, [=,&local_max] (int i, int j, int) noexcept
1057 {
1058 const Real v = (mskr(i,j,0) > 0.0_rt) ? visc2(i,j,0) : 0.0_rt;
1059 local_max = amrex::max(local_max, v);
1060 });
1061 return local_max;
1062 });
1063 ParallelDescriptor::ReduceRealMax(visc_max_mask0);
1064 if (ParallelDescriptor::IOProcessor() && lev == 0)
1065 {
1066 Print() << "\nHorizontal mixing scaled by grid metric\n";
1067 Print() << "grdmax = " << grdmax << "\n";
1069 Print() << "AMR scaling (linear) lev_scale = " << lev_scale << "\n";
1070 }
1071 Print() << "visc2(all) min/max = "
1072 << visc_min_all << " / "
1073 << visc_max_all << "\n";
1074 Print() << "visc2(wet,k=0) min/max = "
1075 << visc_min_wet << " / "
1076 << visc_max_wet << "\n";
1077 Print() << "visc2(mask->0) min/max = "
1078 << visc_min_mask0 << " / "
1079 << visc_max_mask0 << "\n";
1080 }
1081
1082 } else {
1083 Abort("Don't know this horizontal mixing type");
1084 }
1085
1086 // Final FillPatch for all fields
1087 Real time = 0.0_rt;
1088 FillPatch(lev, time, *vec_visc2_p[lev], GetVecOfPtrs(vec_visc2_p), foextrap_periodic_bc());
1089 FillPatch(lev, time, *vec_visc2_r[lev], GetVecOfPtrs(vec_visc2_r), foextrap_periodic_bc());
1090 for (int n = 0; n < ncons; n++) {
1091 FillPatch(lev, time, *vec_diff2[lev], GetVecOfPtrs(vec_diff2),
1092 foextrap_periodic_bc(), BdyVars::null, n, false);
1093 }
1094}
1095
1096/**
1097 * @param[in ] lev level to operate on
1098 */
1099void
1101{
1102 BL_PROFILE("REMORA::init_flat_bathymetry()");
1103 vec_h[lev]->setVal(-geom[0].ProbLo()[2]);
1104}
1105
1106/**
1107 * @param[in ] lev level to operate on
1108 */
1109void
1111{
1112 BL_PROFILE("REMORA::set_smflux()");
1114 prob->init_analytic_smflux(lev, geom[lev], solverChoice, *this,*vec_sustr[lev], *vec_svstr[lev]);
1116#ifdef REMORA_USE_NETCDF
1117 sustr_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_sustr[lev].get(), geom, ref_ratio);
1118 svstr_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_svstr[lev].get(), geom, ref_ratio);
1119 FillPatch(lev, t_old[lev], *vec_sustr[lev], GetVecOfPtrs(vec_sustr), foextrap_periodic_bc(), BdyVars::null,0,false);
1120 FillPatch(lev, t_old[lev], *vec_svstr[lev], GetVecOfPtrs(vec_svstr), foextrap_periodic_bc(), BdyVars::null,0,false);
1121#endif
1122 }
1123}
1124
1125/**
1126 * @param[in ] lev level to operate on
1127 */
1128void
1130{
1131 BL_PROFILE("REMORA::set_wind()");
1132 const bool driver_has_uwind = driver_atmos_state_from_driver[0];
1133 const bool driver_has_vwind = driver_atmos_state_from_driver[1];
1134
1136 // The analytic wind initializer writes both components together, so only
1137 // invoke it when the driver has not already provided the wind pair.
1138 if (!(driver_has_uwind && driver_has_vwind)) {
1139 prob->init_analytic_wind(lev,geom[lev], solverChoice, *this, *vec_uwind[lev], *vec_vwind[lev]);
1140 }
1141 if (vec_uwind[lev] != nullptr) { vec_uwind[lev]->FillBoundary(geom[lev].periodicity()); }
1142 if (vec_vwind[lev] != nullptr) { vec_vwind[lev]->FillBoundary(geom[lev].periodicity()); }
1143 } else if (solverChoice.wind_type == WindType::netcdf) {
1144#ifdef REMORA_USE_NETCDF
1145 if (!driver_has_uwind) {
1146 Uwind_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_uwind[lev].get(), geom, ref_ratio);
1147 FillPatch(lev, t_old[lev], *vec_uwind[lev], GetVecOfPtrs(vec_uwind),
1149 } else {
1150 if (vec_uwind[lev] != nullptr) { vec_uwind[lev]->FillBoundary(geom[lev].periodicity()); }
1151 }
1152 if (!driver_has_vwind) {
1153 Vwind_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_vwind[lev].get(), geom, ref_ratio);
1154 FillPatch(lev, t_old[lev], *vec_vwind[lev], GetVecOfPtrs(vec_vwind),
1156 } else {
1157 if (vec_vwind[lev] != nullptr) { vec_vwind[lev]->FillBoundary(geom[lev].periodicity()); }
1158 }
1159
1160 // Conditionally update atmospheric fields if loaded from NetCDF
1162 Tair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_Tair[lev].get(), geom, ref_ratio);
1163 FillPatch(lev, t_old[lev], *vec_Tair[lev], GetVecOfPtrs(vec_Tair),
1165 } else if (vec_Tair[lev]) {
1166 vec_Tair[lev]->FillBoundary(geom[lev].periodicity());
1167 }
1169 qair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_qair[lev].get(), geom, ref_ratio);
1170 FillPatch(lev, t_old[lev], *vec_qair[lev], GetVecOfPtrs(vec_qair),
1172
1173 // Convert qair from percentage (0-100) to specific humidity (0-1) if needed
1175 vec_qair[lev]->mult(0.01);
1176
1177 // Update ghost cells after modification
1178 vec_qair[lev]->FillBoundary(geom[lev].periodicity());
1179 }
1180 } else if (vec_qair[lev]) {
1181 vec_qair[lev]->FillBoundary(geom[lev].periodicity());
1182 }
1184 Pair_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_Pair[lev].get(), geom, ref_ratio);
1185 FillPatch(lev, t_old[lev], *vec_Pair[lev], GetVecOfPtrs(vec_Pair),
1187 } else if (vec_Pair[lev]) {
1188 vec_Pair[lev]->FillBoundary(geom[lev].periodicity());
1189 }
1191 srflx_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_srflx[lev].get(), geom, ref_ratio);
1192 FillPatch(lev, t_old[lev], *vec_srflx[lev], GetVecOfPtrs(vec_srflx),
1194 } else if (vec_srflx[lev]) {
1195 vec_srflx[lev]->FillBoundary(geom[lev].periodicity());
1196 }
1198 longwave_down_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_longwave_down[lev].get(), geom, ref_ratio);
1199 FillPatch(lev, t_old[lev], *vec_longwave_down[lev], GetVecOfPtrs(vec_longwave_down),
1201 } else if (vec_longwave_down[lev]) {
1202 vec_longwave_down[lev]->FillBoundary(geom[lev].periodicity());
1203 }
1205 rain_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_rain[lev].get(), geom, ref_ratio);
1206 FillPatch(lev, t_old[lev], *vec_rain[lev], GetVecOfPtrs(vec_rain),
1208 } else if (vec_rain[lev]) {
1209 vec_rain[lev]->FillBoundary(geom[lev].periodicity());
1210 }
1212 cloud_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_cloud[lev].get(), geom, ref_ratio);
1213 FillPatch(lev, t_old[lev], *vec_cloud[lev], GetVecOfPtrs(vec_cloud),
1215 } else if (vec_cloud[lev]) {
1216 vec_cloud[lev]->FillBoundary(geom[lev].periodicity());
1217 }
1219 EminusP_data_from_file->update_interpolated_to_time(t_old[lev], lev, vec_EminusP[lev].get(), geom, ref_ratio);
1220 FillPatch(lev, t_old[lev], *vec_EminusP[lev], GetVecOfPtrs(vec_EminusP),
1222 } else if (vec_EminusP[lev]) {
1223 vec_EminusP[lev]->FillBoundary(geom[lev].periodicity());
1224 }
1225#endif
1226 } else {
1227 amrex::Abort("Unknown wind_type in REMORA::set_wind()");
1228 }
1229}
1230
1231/**
1232 * @param[in ] lev level to operate on
1233 * @param[in ] time current time for initialization
1234 */
1235void
1236REMORA::init_only (int lev, Real time)
1237{
1238 BL_PROFILE("REMORA::init_only()");
1239 t_new[lev] = time;
1240 t_old[lev] = time - 1.e200_rt;
1241
1242 cons_new[lev]->setVal(0.0_rt);
1243 xvel_new[lev]->setVal(0.0_rt);
1244 yvel_new[lev]->setVal(0.0_rt);
1245 zvel_new[lev]->setVal(0.0_rt);
1246
1247 xvel_old[lev]->setVal(0.0_rt);
1248 yvel_old[lev]->setVal(0.0_rt);
1249 zvel_old[lev]->setVal(0.0_rt);
1250
1251 vec_ru[lev]->setVal(0.0_rt);
1252 vec_rv[lev]->setVal(0.0_rt);
1253
1254 vec_ru2d[lev]->setVal(0.0_rt);
1255 vec_rv2d[lev]->setVal(0.0_rt);
1256
1258 set_grid_scale(lev);
1259 }
1260 set_masks(lev);
1261
1262#ifdef REMORA_USE_NETCDF
1265
1266 if (solverChoice.do_any_clim_nudg && lev == 0) {
1267 if (nc_clim_his_file.empty() || nc_clim_his_file[0].empty()) {
1268 amrex::Error("NetCDF climatology file name must be provided via input");
1269 }
1272 clim_ubar_time_varname, geom[lev].Domain(),vec_ubar[lev].get(),true,true));
1274 clim_ubar_time_varname, geom[lev].Domain(),vec_vbar[lev].get(),true,true));
1275 ubar_clim_data_from_file->Initialize();
1276 vbar_clim_data_from_file->Initialize();
1277 }
1279 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));
1280 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));
1281 u_clim_data_from_file->Initialize();
1282 v_clim_data_from_file->Initialize();
1283 }
1284 // Since the NCTimeSeries object isn't filling the cons_new MultiFab directly, we don't have to specify a component.
1285 // It just needs to know the shape of the MultiFab
1287 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));
1288 temp_clim_data_from_file->Initialize();
1289 }
1291 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));
1292 salt_clim_data_from_file->Initialize();
1293 }
1294 }
1295 }
1296
1298 amrex::Print() << "Calling init_bdry_from_netcdf at level " << lev << std::endl;
1300 amrex::Print() << "Boundary data loaded from netcdf file \n " << std::endl;
1301 }
1302
1303 // This will be a non-op if forcings specified analytically
1305 if (lev==0) {
1306 if (nc_frc_file.empty() || nc_frc_file[0].empty()) {
1307 amrex::Error("NetCDF forcing file name must be provided via input for winds");
1308 }
1309 Uwind_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Uwind", frc_time_varname, geom[lev].Domain(),vec_uwind[lev].get(), true, false));
1310 Vwind_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Vwind", frc_time_varname, geom[lev].Domain(),vec_vwind[lev].get(), true, false));
1311 Uwind_data_from_file->Initialize();
1312 Vwind_data_from_file->Initialize();
1313 } else {
1314 FillCoarsePatch(lev, time, vec_uwind[lev].get(), vec_uwind[lev-1].get(), foextrap_bc());
1315 FillCoarsePatch(lev, time, vec_vwind[lev].get(), vec_vwind[lev-1].get(), foextrap_bc());
1316 }
1318 if (lev==0) {
1319 if (nc_frc_file.empty() || nc_frc_file[0].empty()) {
1320 amrex::Error("NetCDF forcing file name must be provided via input for surface momentum fluxes");
1321 }
1322 sustr_data_from_file.reset(new NCTimeSeries(nc_frc_file, "sustr", frc_time_varname, geom[lev].Domain(),vec_sustr[lev].get(), true, false));
1323 svstr_data_from_file.reset(new NCTimeSeries(nc_frc_file, "svstr", frc_time_varname, geom[lev].Domain(),vec_svstr[lev].get(), true, false));
1324 sustr_data_from_file->Initialize();
1325 svstr_data_from_file->Initialize();
1326 } else {
1327 FillCoarsePatch(lev, time, vec_sustr[lev].get(), vec_sustr[lev-1].get(), foextrap_bc());
1328 FillCoarsePatch(lev, time, vec_svstr[lev].get(), vec_svstr[lev-1].get(), foextrap_bc());
1329 }
1330 }
1331
1332 // Conditionally load atmospheric forcing fields from NetCDF based on user flags
1333 if (lev==0) {
1335 Tair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Tair", frc_time_varname, geom[lev].Domain(),vec_Tair[lev].get(), true, false));
1336 Tair_data_from_file->Initialize();
1337 }
1339 qair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "qair", frc_time_varname, geom[lev].Domain(),vec_qair[lev].get(), true, false));
1340 qair_data_from_file->Initialize();
1341 }
1343 Pair_data_from_file.reset(new NCTimeSeries(nc_frc_file, "Pair", frc_time_varname, geom[lev].Domain(),vec_Pair[lev].get(), true, false));
1344 Pair_data_from_file->Initialize();
1345 }
1347 srflx_data_from_file.reset(new NCTimeSeries(nc_frc_file, "swrad", frc_time_varname, geom[lev].Domain(),vec_srflx[lev].get(), true, false));
1348 srflx_data_from_file->Initialize();
1349 }
1351 rain_data_from_file.reset(new NCTimeSeries(nc_frc_file, "rain", frc_time_varname, geom[lev].Domain(),vec_rain[lev].get(), true, false));
1352 rain_data_from_file->Initialize();
1353 }
1355 cloud_data_from_file.reset(new NCTimeSeries(nc_frc_file, "cloud", frc_time_varname, geom[lev].Domain(),vec_cloud[lev].get(), true, false));
1356 cloud_data_from_file->Initialize();
1357 }
1359 EminusP_data_from_file.reset(new NCTimeSeries(nc_frc_file, "EminusP", frc_time_varname, geom[lev].Domain(),vec_EminusP[lev].get(), true, false));
1360 EminusP_data_from_file->Initialize();
1361 }
1362 } else {
1364 FillCoarsePatch(lev, time, vec_Tair[lev].get(), vec_Tair[lev-1].get(), foextrap_bc());
1365 }
1367 FillCoarsePatch(lev, time, vec_qair[lev].get(), vec_qair[lev-1].get(), foextrap_bc());
1368 }
1370 FillCoarsePatch(lev, time, vec_Pair[lev].get(), vec_Pair[lev-1].get(), foextrap_bc());
1371 }
1373 FillCoarsePatch(lev, time, vec_srflx[lev].get(), vec_srflx[lev-1].get(), foextrap_bc());
1374 }
1376 FillCoarsePatch(lev, time, vec_rain[lev].get(), vec_rain[lev-1].get(), foextrap_bc());
1377 }
1379 FillCoarsePatch(lev, time, vec_cloud[lev].get(), vec_cloud[lev-1].get(), foextrap_bc());
1380 }
1382 FillCoarsePatch(lev, time, vec_EminusP[lev].get(), vec_EminusP[lev-1].get(), foextrap_bc());
1383 }
1384 }
1386 if (lev==0) {
1387 if (nc_frc_file.empty() || nc_frc_file[0].empty()) {
1388 amrex::Error("NetCDF forcing file name must be provided via input for longwave radiation");
1389 }
1391 geom[lev].Domain(), vec_longwave_down[lev].get(), true, false));
1392 longwave_down_data_from_file->Initialize();
1393 } else {
1394 FillCoarsePatch(lev, time, vec_longwave_down[lev].get(), vec_longwave_down[lev-1].get(), foextrap_bc());
1395 }
1396 }
1397
1398 // Only need to read in rivers on level 0
1399 // Will need to be on higher levels eventually
1400 if (solverChoice.do_rivers) {
1401 if (nc_riv_file.empty() || nc_riv_file[0].empty()) {
1402 amrex::Error("NetCDF river file name must be provided via input for rivers");
1403 }
1404 auto dom = geom[0].Domain();
1405 int nz = dom.length(2);
1406 river_source_cons.resize(ncons);
1409 river_source_cons[Salt_comp]->Initialize();
1410 }
1413 river_source_cons[Temp_comp]->Initialize();
1414 }
1417 river_source_cons[Tracer_comp]->Initialize();
1418 }
1419 river_source_transport.reset(new NCTimeSeriesRiver(nc_riv_file, "river_transport", riv_time_varname, nz));
1420 river_source_transport->Initialize();
1421 river_source_transportbar.reset(new NCTimeSeriesRiver(nc_riv_file, "river_transport", riv_time_varname, nz, 1));
1422 river_source_transportbar->Initialize();
1424 }
1425
1426 if (lev==0 and hires_grid_level > 0 and solverChoice.ic_type == IC_Type::netcdf) {
1427 amrex::Print() << "Reading high resolution bathymetry and grid data" << std::endl;
1431 amrex::Print() << "Done reading in high resolution bathymetry and grid data" << std::endl;
1432 }
1433 if (lev==0 and hires_init_level > 0 and solverChoice.ic_type == IC_Type::netcdf) {
1434 amrex::Print() << "Reading high resolution initial data" << std::endl;
1438 amrex::Print() << "Done reading in high resolution initial data" << std::endl;
1439 }
1440#else
1442 Abort("Not compiled with NetCDF, but selected boundary conditions require NetCDF");
1443 }
1444 if (solverChoice.do_rivers) {
1445 Abort("Not compiled with NetCDF, but using river sources requires NetCDF");
1446 }
1447#endif
1448
1449 if (lev==0 and hires_grid_level > 0 and solverChoice.ic_type == IC_Type::analytic) {
1452 }
1453
1454 if (lev==0 and hires_init_level > 0 and solverChoice.ic_type == IC_Type::analytic) {
1457 }
1458
1459 set_bathymetry(lev);
1460 set_zeta(lev);
1461 stretch_transform(lev);
1462
1463 if (lev==0 and hires_init_level > 0 and solverChoice.ic_type == IC_Type::analytic) {
1465 }
1466
1467 if (lev==0) {
1468 if (hires_init_level < 0) {
1470 init_analytic(lev);
1471 } else if (solverChoice.ic_type == IC_Type::netcdf) {
1472#ifdef REMORA_USE_NETCDF
1473 amrex::Print() << "Calling init_data_from_netcdf " << std::endl;
1475 set_zeta_to_Ztavg(lev);
1476 amrex::Print() << "Initial data loaded from netcdf file \n " << std::endl;
1477#endif
1478 } else {
1479 amrex::Abort("Unknown IC_Type");
1480 }
1481 } else {
1483 set_zeta_to_Ztavg(lev); // MAYBE???
1484 // Since set_grid_scale is usually called from init_analytic for analytic problems
1486 set_grid_scale(lev);
1487 }
1488 }
1489 } else {
1490 if (lev > hires_init_level) {
1492 FillCoarsePatch(lev, time, xvel_new[lev], xvel_new[lev-1], xvel_bc(), BdyVars::u);
1493 FillCoarsePatch(lev, time, yvel_new[lev], yvel_new[lev-1], yvel_bc(), BdyVars::v);
1494 FillCoarsePatch(lev, time, zvel_new[lev], zvel_new[lev-1], zvel_bc(), BdyVars::null);
1495 } else {
1497 set_zeta_to_Ztavg(lev); // MAYBE???
1499 // Since set_grid_scale is usually called from init_analytic for analytic problems
1500 set_grid_scale(lev);
1501 }
1502 }
1503 }
1504
1505 // Ensure that the face-based data are the same on both sides of a periodic domain.
1506 // The data associated with the lower grid ID is considered the correct value.
1507 xvel_new[lev]->OverrideSync(geom[lev].periodicity());
1508 yvel_new[lev]->OverrideSync(geom[lev].periodicity());
1509 zvel_new[lev]->OverrideSync(geom[lev].periodicity());
1510
1511 set_2darrays(lev);
1512
1513 init_set_vmix(lev);
1514 set_hmixcoef(lev);
1515 set_coriolis(lev);
1516
1517 // Previously set smflux here with OverrideSync:
1518// set_smflux(lev);
1519// prob->init_analytic_smflux(lev, geom[lev], solverChoice, *this, *vec_sustr[lev], *vec_svstr[lev]);
1520// vec_sustr[lev]->OverrideSync(geom[lev].periodicity());
1521// vec_svstr[lev]->OverrideSync(geom[lev].periodicity());
1522
1523}
1524
1525void
1527{
1528 BL_PROFILE("REMORA::ReadParameters()");
1529 {
1530 ParmParse pp; // Traditionally, max_step and stop_time do not have prefix, so allow it for now.
1531 bool noprefix_max_step = pp.queryAdd("max_step", max_step);
1532 bool noprefix_stop_time = pp.queryAdd("stop_time", stop_time);
1533 bool remora_max_step = pp.queryAdd("remora.max_step", max_step);
1534 bool remora_stop_time = pp.queryAdd("remora.stop_time", stop_time);
1535 if (remora_max_step and noprefix_max_step) {
1536 Abort("remora.max_step and max_step are both specified. Please use only one!");
1537 }
1538 if (remora_stop_time and noprefix_stop_time) {
1539 Abort("remora.stop_time and stop_time are both specified. Please use only one!");
1540 }
1541 }
1542
1543 ParmParse pp(pp_prefix);
1544
1545 // Common physics and simulation parameters
1546 pp.queryAdd("nscalar", nscalar);
1547 if (nscalar < 1) {
1548 amrex::Abort("remora.nscalar must be at least 1");
1549 }
1552
1553 pp.queryAdd("check_file", check_file);
1554 pp.queryAdd("check_int", check_int);
1555 pp.queryAdd("check_int_time", check_int_time);
1556 pp.queryAdd("expand_plotvars_to_unif_rr", expand_plotvars_to_unif_rr);
1557 pp.query("plotfile_fill_value", plotfile_fill_value);
1558 pp.query("netcdf_fill_value", netcdf_fill_value);
1559 pp.queryAdd("restart", restart_chkfile);
1560 pp.queryAdd("start_time", start_time);
1561
1562 num_boxes_at_level.resize(max_level + 1, 0);
1563 boxes_at_level.resize(max_level + 1);
1564 num_boxes_at_level[0] = 1;
1565 boxes_at_level[0].resize(1);
1566 boxes_at_level[0][0] = geom[0].Domain();
1567
1568 if (pp.contains("data_log")) {
1569 int num_datalogs = pp.countval("data_log");
1570 datalog.resize(num_datalogs);
1571 datalogname.resize(num_datalogs);
1572 pp.queryarr("data_log", datalogname, 0, num_datalogs);
1573 for (int i = 0; i < num_datalogs; i++)
1575 }
1576
1577 pp.queryAdd("v", verbose);
1578 pp.queryAdd("sum_interval", sum_interval);
1579 pp.queryAdd("sum_period", sum_per);
1580 pp.queryAdd("file_min_digits", file_min_digits);
1581
1582 if (file_min_digits < 0) {
1583 amrex::Abort("remora.file_min_digits must be non-negative");
1584 }
1585
1586 pp.queryAdd("cfl", cfl);
1587 pp.queryAdd("change_max", change_max);
1588 pp.queryAdd("fixed_dt", fixed_dt);
1589 pp.queryAdd("fixed_fast_dt", fixed_fast_dt);
1590 pp.queryAdd("fixed_ndtfast_ratio", fixed_ndtfast_ratio);
1591
1592 if (fixed_dt > 0. && fixed_fast_dt > 0. && fixed_ndtfast_ratio > 0) {
1594 amrex::Abort("Dt is over-specfied");
1595 }
1596 } else if (fixed_dt > 0. && fixed_fast_dt > 0. && fixed_ndtfast_ratio <= 0) {
1597 fixed_ndtfast_ratio = static_cast<int>(fixed_dt / fixed_fast_dt);
1598 }
1599 AMREX_ASSERT(cfl > 0. || fixed_dt > 0.);
1600
1601 num_files_at_level.resize(max_level + 1, 0);
1602 num_boxes_at_level.resize(max_level + 1, 0);
1603 boxes_at_level.resize(max_level + 1);
1604 num_boxes_at_level[0] = 1;
1605 boxes_at_level[0].resize(1);
1606 boxes_at_level[0][0] = geom[0].Domain();
1607
1608 pp.queryAdd("plot_file", plot_file_name);
1609 pp.queryAdd("plot_int", plot_int);
1610 pp.queryAdd("plot_int_time", plot_int_time);
1611 pp.queryAdd("plot_staggered_vels", plot_staggered_vels);
1612
1613 std::string plotfile_type_str = "amrex";
1614 pp.queryAdd("plotfile_type", plotfile_type_str);
1615 if (plotfile_type_str == "amrex") {
1617 } else if (plotfile_type_str == "netcdf" || plotfile_type_str == "NetCDF") {
1619#ifdef REMORA_USE_NETCDF
1620 pp.queryAdd("write_history_file",write_history_file);
1621 pp.queryAdd("chunk_history_file",chunk_history_file);
1622 pp.queryAdd("steps_per_history_file",steps_per_history_file);
1623 // Estimate size of domain for one timestep of netcdf
1624 auto dom = geom[0].Domain();
1625 int nx = dom.length(0) + 2;
1626 int ny = dom.length(1) + 2;
1627 int nz = dom.length(2);
1629 // Estimate number of steps that will fit into a 2GB file.
1630 steps_per_history_file = int((1.6e10 - NCH2D * nx * ny * 64.0_rt)
1631 / (nx * ny * 64.0_rt * (NC3D*nz + NC2D)));
1632 // If we calculate that a single step will exceed 2GB and the user has
1633 // requested automatic history file sizing, warn about a possible impending
1634 // error, and set steps_per_history_file = 1 to attempt output anyway.
1635 if (steps_per_history_file == 0) {
1636 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.");
1638 }
1639 } else if (write_history_file and !chunk_history_file) {
1640 // Estimate number of output steps we'll need
1641 int nt_out = int((max_step) / plot_int) + 1;
1642 Real est_hist_file_size = NCH2D * nx * ny * 64.0_rt + nt_out * nx * ny * 64.0_rt * (NC3D*nz + NC2D);
1643 if (est_hist_file_size > 1.6e10) {
1644 amrex::Warning("WARNING: NetCDF history file may be larger than 2GB limit. Consider setting remora.chunk_history_file=true");
1645 }
1646 }
1648 Print() << "NetCDF history files will have " << steps_per_history_file << " steps per file." << std::endl;
1649 }
1650#endif
1651 } else {
1652 amrex::Print() << "User selected plotfile_type = " << plotfile_type_str << std::endl;
1653 amrex::Abort("Dont know this plotfile_type");
1654 }
1655#ifndef REMORA_USE_NETCDF
1657 {
1658 amrex::Abort("Please compile with NetCDF in order to enable NetCDF plotfiles");
1659 }
1660
1661#endif
1662#ifdef REMORA_USE_NETCDF
1663 nc_init_file.resize(max_level+1);
1664 nc_grid_file.resize(max_level+1);
1665 num_files_at_level.resize(max_level + 1, 0);
1666
1667 boundary_series.resize(max_level+1);
1668
1669
1670 // NetCDF initialization files -- possibly multiple files at each of multiple levels
1671 // but we always have exactly one file at level 0
1672 for (int lev = 0; lev <= max_level; lev++)
1673 {
1674 const std::string nc_file_names = amrex::Concatenate("nc_init_file_",lev,1);
1675 const std::string nc_bathy_file_names = amrex::Concatenate("nc_grid_file_",lev,1);
1676
1677 if (pp.contains(nc_file_names.c_str()))
1678 {
1679 int num_files = pp.countval(nc_file_names.c_str());
1680 int num_bathy_files = pp.countval(nc_bathy_file_names.c_str());
1681 if (num_files != num_bathy_files) {
1682 amrex::Error("Must have same number of netcdf files for grid info as for solution");
1683 }
1684
1685 num_files_at_level[lev] = num_files;
1686 nc_init_file[lev].resize(num_files);
1687 nc_grid_file[lev].resize(num_files);
1688
1689 pp.queryarr(nc_file_names.c_str() , nc_init_file[lev] ,0,num_files);
1690 pp.queryarr(nc_bathy_file_names.c_str(), nc_grid_file[lev],0,num_files);
1691 }
1692 }
1693
1694 pp.queryAdd("nc_grid_file_hires", nc_grid_file_hires);
1695 pp.queryAdd("nc_init_file_hires", nc_init_file_hires);
1696
1697 // We only read boundary data at level 0
1698 pp.queryarr("nc_bdry_file", nc_bdry_file);
1699
1700 // Also only read forcings at level 0 (for now)
1701 if (pp.contains("nc_frc_file")) {
1702 int num_files = pp.countval("nc_frc_file");
1703 nc_frc_file.resize(num_files);
1704 pp.queryarr("nc_frc_file", nc_frc_file, 0, num_files);
1705 }
1706
1707 // Get river file
1708 if (pp.contains("nc_river_file")) {
1709 int num_files = pp.countval("nc_river_file");
1710 nc_riv_file.resize(num_files);
1711 pp.queryarr("nc_river_file", nc_riv_file, 0, num_files);
1712 }
1713
1714 // Read in file names for climatology history and nudging weights
1715 if (pp.contains("nc_clim_his_file")) {
1716 int num_files = pp.countval("nc_clim_his_file");
1717 nc_clim_his_file.resize(num_files);
1718 pp.queryarr("nc_clim_his_file", nc_clim_his_file, 0, num_files);
1719 }
1720 pp.queryAdd("nc_clim_coeff_file", nc_clim_coeff_file);
1721
1722 for (int i=0; i<BdyVars::NumTypes; i++) {
1723 bdry_time_name_byvar.push_back("");
1724 }
1725 pp.queryAdd("bdy_time_varname",bdry_time_varname);
1726 pp.queryAdd("bdy_temp_time_varname",bdry_time_name_byvar[BdyVars::t]);
1727 pp.queryAdd("bdy_salt_time_varname",bdry_time_name_byvar[BdyVars::s]);
1728 pp.queryAdd("bdy_u_time_varname",bdry_time_name_byvar[BdyVars::u]);
1729 pp.queryAdd("bdy_v_time_varname",bdry_time_name_byvar[BdyVars::v]);
1730 pp.queryAdd("bdy_ubar_time_varname",bdry_time_name_byvar[BdyVars::ubar]);
1731 pp.queryAdd("bdy_vbar_time_varname",bdry_time_name_byvar[BdyVars::vbar]);
1732 pp.queryAdd("bdy_zeta_time_varname",bdry_time_name_byvar[BdyVars::zeta]);
1733
1734 // If not specified per variable, populate with the default
1735 for (int i=0; i<BdyVars::NumTypes; i++) {
1736 if (bdry_time_name_byvar[i] == "") {
1738 }
1739 }
1740
1741 pp.queryAdd("frc_time_varname",frc_time_varname);
1742
1743 pp.queryAdd("riv_time_varname",riv_time_varname);
1744
1745 pp.queryAdd("clim_ubar_time_varname",clim_ubar_time_varname);
1746 pp.queryAdd("clim_vbar_time_varname",clim_vbar_time_varname);
1747 pp.queryAdd("clim_u_time_varname",clim_u_time_varname);
1748 pp.queryAdd("clim_v_time_varname",clim_v_time_varname);
1749 pp.queryAdd("clim_salt_time_varname",clim_salt_time_varname);
1750 pp.queryAdd("clim_temp_time_varname",clim_temp_time_varname);
1751
1752#endif
1753 pp.queryAdd("hires_grid_level", hires_grid_level);
1754 if (hires_grid_level > max_level) {
1755 amrex::Abort("hires_grid_level must be less than or equal to amr.max_level");
1756 }
1757 pp.queryAdd("hires_init_level", hires_init_level);
1758 if (hires_init_level > max_level) {
1759 amrex::Abort("hires_init_level must be less than or equal to amr.max_level");
1760 }
1761#ifdef REMORA_USE_PARTICLES
1762 readTracersParams();
1763#endif
1764
1765 {
1766 ParmParse pp_amr("amr");
1767 pp_amr.queryAdd("regrid_int", regrid_int);
1768 pp_amr.queryAdd("do_substep", do_substep);
1769 if (do_substep) {
1770 amrex::Abort("Time substepping is not yet implemented. amr.do_substep must be 0");
1771 }
1772
1773 }
1775
1776 // NOTE: This feature is not yet implemented because it will require passing x,y,z to prob functions.
1777 // Currently these are accessed by passing a pointer to the REMORA class. However, this requires the
1778 // coordinates at hires_init_level to already exist (and specifically for the hires_init_level level
1779 // to already be initialized), which is generally not the case. A solution is to create a separate
1780 // coordinates object that is passed to the prob functions instead of the REMORA object. Then x,y,z
1781 // coordinates can be calculated at any level without the corresponding level having been created.
1783 amrex::Abort("Cannot do high-resolution initialization for analytic initial conditions. Not yet implemented");
1784 }
1785
1786}
1787
1788
1789void
1791{
1792 BL_PROFILE("REMORA::AverageDown()");
1793 for (int lev = finest_level-1; lev >= 0; --lev)
1794 {
1795 AverageDownTo(lev);
1796 }
1797}
1798
1799/**
1800 * @param[in ] crse_lev level to average down to
1801 */
1802void
1804{
1805 BL_PROFILE("REMORA::AverageDownTo()");
1806 average_down(*cons_new[crse_lev+1], *cons_new[crse_lev],
1807 0, cons_new[crse_lev]->nComp(), refRatio(crse_lev));
1808 average_down(*vec_Zt_avg1[crse_lev+1].get(), *vec_Zt_avg1[crse_lev].get(),
1809 0, vec_Zt_avg1[crse_lev]->nComp(), refRatio(crse_lev));
1810
1811 Array<MultiFab*,AMREX_SPACEDIM> faces_crse;
1812 Array<MultiFab*,AMREX_SPACEDIM> faces_fine;
1813 faces_crse[0] = xvel_new[crse_lev];
1814 faces_crse[1] = yvel_new[crse_lev];
1815 faces_crse[2] = zvel_new[crse_lev];
1816
1817 faces_fine[0] = xvel_new[crse_lev+1];
1818 faces_fine[1] = yvel_new[crse_lev+1];
1819 faces_fine[2] = zvel_new[crse_lev+1];
1820
1821 average_down_faces(GetArrOfConstPtrs(faces_fine), faces_crse,
1822 refRatio(crse_lev),geom[crse_lev]);
1823 stretch_transform(crse_lev);
1824}
1825
1826/**
1827 * @param[in ] crse_lev level to average data down to
1828 * @param[inout] vec_mf vector over levels of multifabs containing data to average
1829 */
1830void
1831REMORA::average_down_with_grow_cells (int crse_lev, Vector<std::unique_ptr<MultiFab>>& vec_mf)
1832{
1833 auto const& crsema = vec_mf[crse_lev]->arrays();
1834 auto const& finema = vec_mf[crse_lev+1]->const_arrays();
1835 auto ref_ratio_crse = refRatio(crse_lev);
1836 auto index_type = (vec_mf[crse_lev]->boxArray().ixType()).toIntVect();
1837 auto nghost_crse = cum_ref_ratios[crse_lev] - index_type;
1838 if (index_type[0]==0 and index_type[1]==0) {
1839 ParallelFor(*vec_mf[crse_lev], nghost_crse, vec_mf[crse_lev]->nComp(),
1840 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1841 {
1842 amrex_avgdown(i,j,k,n,crsema[box_no],finema[box_no],0,0,ref_ratio_crse);
1843 });
1844 } else if (index_type[0]==1 and index_type[1]==0) {
1845 ParallelFor(*vec_mf[crse_lev], nghost_crse, 1,
1846 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1847 {
1848 amrex_avgdown_faces(i,j,k,n,crsema[box_no],finema[box_no],0,0,ref_ratio_crse,0);
1849 });
1850 } else if (index_type[0]==0 and index_type[1]==1) {
1851 ParallelFor(*vec_mf[crse_lev], nghost_crse, 1,
1852 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1853 {
1854 amrex_avgdown_faces(i,j,k,n,crsema[box_no],finema[box_no],0,0,ref_ratio_crse,1);
1855 });
1856 } else {
1857 amrex::Abort("Unexpected nodality in average_down_with_grow_cells");
1858 }
1859 Gpu::streamSynchronize();
1860}
1861
1862/**
1863 * @param[in ] lev level at which to get time
1864 */
1865amrex::Real REMORA::get_t_old(int lev) const
1866{
1867 return t_old[lev];
1868}
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:1544
std::string nc_grid_file_hires
Grid file for high resolution bathymetry.
Definition REMORA.H:1556
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:417
amrex::Vector< std::string > nc_riv_file
NetCDF river file(s)
Definition REMORA.H:1573
void set_grid_vars_averaged_down(int lev)
Set pm/pn by averaging down from higher-resolution grid.
Definition REMORA.cpp:681
std::string riv_time_varname
Name of time field for river time.
Definition REMORA.H:1598
int foextrap_periodic_bc() const noexcept
Definition REMORA.H:1145
amrex::Vector< std::string > nc_clim_his_file
NetCDF climatology history file(s)
Definition REMORA.H:1576
int ncons
Number of conserved scalars in the state (temperature + salt + passive scalars)
Definition REMORA.H:1437
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta_full_domain
high resolution initial free surface height (2D)
Definition REMORA.H:456
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:332
std::string nc_init_file_hires
Init file for high resolution.
Definition REMORA.H:1563
static amrex::Real fixed_dt
User specified fixed baroclinic time step.
Definition REMORA.H:1445
amrex::Real last_plot_file_time
Simulation time when we last output a plotfile.
Definition REMORA.H:1411
int zvel_bc() const noexcept
Definition REMORA.H:1140
void init_full_domain_zeta_from_analytic()
Initialize high resolution initial sea surface height from analytic functions.
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:1138
void set_zeta_averaged_down(int lev)
Copy over zeta data that has been averaged down from high res.
Definition REMORA.cpp:699
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:1263
static amrex::Real previousCPUTimeUsed
Accumulator variable for CPU time used thusfar.
Definition REMORA.H:1657
amrex::Vector< std::string > cons_names
Names of scalars for plotfile output.
Definition REMORA.H:1488
amrex::Vector< std::unique_ptr< amrex::YAFluxRegister > > advflux_reg
array of flux registers for refluxing in multilevel
Definition REMORA.H:1377
std::unique_ptr< NCTimeSeries > sustr_data_from_file
Data container for u-component surface momentum flux read from file.
Definition REMORA.H:1253
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:480
void allocate_init_full_domain()
Allocate multifabs for storing full-domain high resolution initial data.
void init_gls_vmix(int lev, SolverChoice solver_choice)
Initialize GLS variables.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xvel_full_domain
multilevel data container for high res initial x velocities (u in ROMS)
Definition REMORA.H:305
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:314
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:470
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:1321
void init_zeta_full_domain_from_netcdf()
Full-domain high res sea-surface height data initialization from NetCDF file.
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:294
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
Definition REMORA.H:1250
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:1271
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vwind
Wind in the v direction, defined at rho-points.
Definition REMORA.H:382
std::unique_ptr< ProblemBase > prob
Pointer to container of analytical functions for problem definition.
Definition REMORA.H:1350
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:459
void Construct_REMORAFillPatchers(int lev)
Construct FillPatchers.
Definition REMORA.cpp:463
void init_grid_vars_from_netcdf(int lev)
Grid variable initialization from NetCDF file.
static int sum_interval
Diagnostic sum output interval in number of steps.
Definition REMORA.H:1533
int history_count
Counter for which time index we are writing to in the netcdf history file.
Definition REMORA.H:1481
amrex::Real stop_time
Time to stop.
Definition REMORA.H:1426
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rain
precipitation rate [kg/m^2/s]
Definition REMORA.H:411
int do_substep
Whether to substep fine levels in time.
Definition REMORA.H:1454
void Evolve()
Advance solution to final time.
Definition REMORA.cpp:244
std::string bdry_time_varname
Default name of time field for boundary data.
Definition REMORA.H:1581
amrex::Real plotfile_fill_value
fill value for masked arrays in amrex plotfiles
Definition REMORA.H:1500
void ReadCheckpointFile()
read checkpoint file from disk
bool chunk_history_file
Whether to chunk netcdf history file.
Definition REMORA.H:1474
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:375
amrex::Real get_t_old(int lev) const
Accessor method for t_old to expose to outside classes.
Definition REMORA.cpp:1865
int yvel_bc() const noexcept
Definition REMORA.H:1139
std::unique_ptr< NCTimeSeries > longwave_down_data_from_file
Data container for downward longwave radiation flux read from file.
Definition REMORA.H:1269
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:330
amrex::Vector< std::string > datalogname
Definition REMORA.H:1689
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:300
void init_only(int lev, amrex::Real time)
Init (NOT restart or regrid)
Definition REMORA.cpp:1236
void init_set_vmix(int lev)
Initialize vertical mixing coefficients from file or analytic.
Definition REMORA.cpp:756
std::unique_ptr< NCTimeSeries > v_clim_data_from_file
Data container for v-velocity climatology data read from file.
Definition REMORA.H:1284
std::string clim_u_time_varname
Name of time field for u climatology data.
Definition REMORA.H:1590
void set_grid_scale(int lev)
Set pm and pn arrays and x/y coords on level lev.
void set_coriolis(int lev)
Initialize Coriolis factor from file or analytic.
Definition REMORA.cpp:727
int foextrap_bc() const noexcept
Definition REMORA.H:1146
amrex::Vector< REMORAFillPatcher > FPr_u
Vector over levels of FillPatchers for u (3D)
Definition REMORA.H:1319
std::string clim_temp_time_varname
Name of time field for temperature climatology data.
Definition REMORA.H:1596
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1357
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Holds info for dynamically generated tagging criteria.
Definition REMORA.H:1606
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
Definition REMORA.H:340
std::unique_ptr< NCTimeSeriesRiver > river_source_transportbar
Data container for vertically integrated momentum transport in rivers.
Definition REMORA.H:1295
std::array< bool, AtmosState::NumTypes > driver_atmos_state_from_driver
provenance flags for driver-supplied atmospheric forcing lanes
Definition REMORA.H:420
std::string clim_ubar_time_varname
Name of time field for ubar climatology data.
Definition REMORA.H:1586
std::unique_ptr< NCTimeSeries > u_clim_data_from_file
Data container for u-velocity climatology data read from file.
Definition REMORA.H:1282
std::string check_file
Checkpoint file prefix.
Definition REMORA.H:1467
static amrex::Real startCPUTime
Variable for CPU timing.
Definition REMORA.H:1655
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm_full_domain
horizontal scaling factor: 1 / dx (2D) on whole domain
Definition REMORA.H:474
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:287
amrex::Real start_time
Time of the start of the simulation, in seconds.
Definition REMORA.H:1429
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:298
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_uwind
Wind in the u direction, defined at rho-points.
Definition REMORA.H:380
static amrex::Real fixed_fast_dt
User specified fixed barotropic time step.
Definition REMORA.H:1447
int regrid_int
how often each level regrids the higher levels of refinement (after a level advances that many time s...
Definition REMORA.H:1457
amrex::Real check_int_time
Checkpoint output interval in seconds.
Definition REMORA.H:1471
void init_scalar_metadata()
Build runtime scalar names after nscalar is known.
Definition REMORA.cpp:231
int zeta_bc() const noexcept
Definition REMORA.H:1143
void Define_REMORAFillPatchers(int lev)
Define FillPatchers.
Definition REMORA.cpp:512
amrex::Vector< amrex::IntVect > cum_ref_ratios
Cumulative refinement ratio between level 0 and level i.
Definition REMORA.H:1568
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:342
amrex::Real plot_int_time
Plotfile output interval in seconds.
Definition REMORA.H:1465
amrex::Vector< int > num_files_at_level
how many netcdf input files specified at each level
Definition REMORA.H:1355
amrex::Vector< REMORAFillPatcher > FPr_vbar
Vector over levels of FillPatchers for vbar (2D)
Definition REMORA.H:1329
void AverageDownTo(int crse_lev)
more flexible version of AverageDown() that lets you average down across multiple levels
Definition REMORA.cpp:1803
int steps_per_history_file
Number of time steps per netcdf history file.
Definition REMORA.H:1479
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Called after every level 0 timestep.
Definition REMORA.cpp:326
int max_step
maximum number of steps
Definition REMORA.H:1424
amrex::Vector< amrex::MultiFab * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
Definition REMORA.H:291
std::unique_ptr< NCTimeSeries > svstr_data_from_file
Data container for v-component surface momentum flux read from file.
Definition REMORA.H:1255
amrex::Vector< std::string > nc_frc_file
NetCDF forcing file(s)
Definition REMORA.H:1571
amrex::Vector< int > num_boxes_at_level
how many boxes specified at each level by tagging criteria
Definition REMORA.H:1353
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:296
void refinement_criteria_setup()
Set refinement criteria.
int last_check_file_step
Step when we last output a checkpoint file.
Definition REMORA.H:1414
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:1588
amrex::Vector< int > nsubsteps
How many substeps on each level?
Definition REMORA.H:1362
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
Definition REMORA.H:1374
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:1275
int plot_int
Plotfile output interval in iterations.
Definition REMORA.H:1463
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:1273
void InitData()
Initialize multilevel data.
Definition REMORA.cpp:358
void set3DPlotVariables(const std::string &pp_plot_var_names_3d)
amrex::Vector< int > istep
which step?
Definition REMORA.H:1360
void WriteCheckpointFile()
write checkpoint file to disk
std::string nc_clim_coeff_file
NetCDF climatology coefficient file.
Definition REMORA.H:1578
void setRecordDataInfo(int i, const std::string &filename)
Definition REMORA.H:1675
void set_analytic_vmix(int lev)
Set vertical mixing coefficients from analytic.
Definition REMORA.cpp:773
void init_flat_bathymetry(int lev)
Initialize flat bathymetry to value from problo.
Definition REMORA.cpp:1100
std::unique_ptr< NCTimeSeries > temp_clim_data_from_file
Data container for temperature climatology data read from file.
Definition REMORA.H:1286
amrex::Vector< std::string > bdry_time_name_byvar
Name of time fields for boundary data.
Definition REMORA.H:1583
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
Definition REMORA.H:1538
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:344
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yvel_full_domain
multilevel data container for high res initial y velocities (v in ROMS)
Definition REMORA.H:307
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:377
std::unique_ptr< NCTimeSeries > srflx_data_from_file
Data container for shortwave radiation flux read from file.
Definition REMORA.H:1267
REMORA()
Definition REMORA.cpp:62
void set_zeta(int lev)
Initialize zeta from file or analytic.
Definition REMORA.cpp:572
static amrex::Real change_max
Fraction maximum change in subsequent time steps.
Definition REMORA.H:1443
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:282
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h_full_domain
Bathymetry data on the whole domain at each potential level.
Definition REMORA.H:317
void set_bathymetry(int lev)
Initialize bathymetry from file or analytic.
Definition REMORA.cpp:611
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:289
std::unique_ptr< NCTimeSeries > ubar_clim_data_from_file
Data container for ubar climatology data read from file.
Definition REMORA.H:1278
void init_full_domain_from_analytic()
Initialize high resolution initial problem data from analytic functions.
void init_data_full_domain_from_netcdf()
High resolution roblem initialization from NetCDF file.
amrex::Vector< REMORAFillPatcher > FPr_c
Vector over levels of FillPatchers for scalars.
Definition REMORA.H:1317
int hires_init_level
Which level the high resolution initialization data is at.
Definition REMORA.H:1561
std::unique_ptr< NCTimeSeries > Tair_data_from_file
Data container for air temperature read from file.
Definition REMORA.H:1261
int nscalar
Number of passive scalars carried in the state.
Definition REMORA.H:1435
std::string clim_v_time_varname
Name of time field for v climatology data.
Definition REMORA.H:1592
amrex::Vector< REMORAFillPatcher > FPr_w
Vector over levels of FillPatchers for w.
Definition REMORA.H:1323
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:1831
std::string clim_salt_time_varname
Name of time field for salinity climatology data.
Definition REMORA.H:1594
std::unique_ptr< NCTimeSeries > Uwind_data_from_file
Data container for u-direction wind read from file.
Definition REMORA.H:1257
std::unique_ptr< NCTimeSeries > Pair_data_from_file
Data container for air pressure read from file.
Definition REMORA.H:1265
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1364
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:1494
void set_masks(int lev)
Initialize land-sea masks from file or analytic.
Definition REMORA.cpp:789
amrex::Vector< amrex::Vector< std::unique_ptr< NCTimeSeriesBoundary > > > boundary_series
Vector over BdyVars of boundary series data containers.
Definition REMORA.H:1298
int cf_set_width
Width for fixing values at coarse-fine interface.
Definition REMORA.H:1314
void ReadParameters()
read in some parameters from inputs file
Definition REMORA.cpp:1526
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:326
void sum_integrated_quantities(amrex::Real time)
Integrate conserved quantities for diagnostics.
static int total_nc_plot_file_step
Definition REMORA.H:1153
static int plot_staggered_vels
Whether to write the staggered velocities (not averaged to cell centers)
Definition REMORA.H:1541
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:395
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
Definition REMORA.H:453
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:451
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:1497
int plot_file_on_restart
Whether to output a plotfile on restart from checkpoint.
Definition REMORA.H:1418
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:449
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:285
std::string frc_time_varname
Name of time field for forcing data.
Definition REMORA.H:1600
void set_wind(int lev)
Initialize or calculate wind speed from file or analytic.
Definition REMORA.cpp:1129
amrex::Vector< REMORAFillPatcher > FPr_ubar
Vector over levels of FillPatchers for ubar (2D)
Definition REMORA.H:1327
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.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_cons_full_domain
multilevel data container for high res initial data: temperature, salinity, passive tracer
Definition REMORA.H:303
std::unique_ptr< NCTimeSeries > Vwind_data_from_file
Data container for v-direction wind read from file.
Definition REMORA.H:1259
void init_bathymetry_full_domain_from_netcdf()
Full domain high-res bathymetry data initialization from NetCDF file.
void set_hmixcoef(int lev)
Initialize horizontal mixing coefficients.
Definition REMORA.cpp:815
amrex::Vector< std::unique_ptr< NCTimeSeriesRiver > > river_source_cons
Vector of data containers for scalar data in rivers.
Definition REMORA.H:1291
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:1293
void init_grid_vars_full_domain_from_netcdf()
Full domain high-res grid variable initialization from NetCDF file.
void AverageDown()
set covered coarse cells to be the average of overlying fine cells
Definition REMORA.cpp:1790
static int fixed_ndtfast_ratio
User specified, number of barotropic steps per baroclinic step.
Definition REMORA.H:1449
amrex::Real netcdf_fill_value
fill value for masked arrays in netcdf output
Definition REMORA.H:1502
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn_full_domain
horizontal scaling factor: 1 / dy (2D) on whole domain
Definition REMORA.H:476
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:472
amrex::Vector< std::unique_ptr< std::fstream > > datalog
Definition REMORA.H:1688
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
Definition REMORA.H:338
static amrex::Real cfl
CFL condition.
Definition REMORA.H:1441
void append3DPlotVariables(const std::string &pp_plot_var_names_3d)
void allocate_bathymetry_grid_vars_full_domain()
Allocate multifabs for storing full-domain bathymetry and grid vars data.
void set_init_data_averaged_down(int lev)
Problem initialization from averaged-down high resolution data.
Definition REMORA.cpp:710
static int verbose
Verbosity level of output.
Definition REMORA.H:1530
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_cloud
cloud cover fraction [0-1], defined at rho-points
Definition REMORA.H:415
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1461
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:372
std::unique_ptr< NCTimeSeries > salt_clim_data_from_file
Data container for salinity climatology data read from file.
Definition REMORA.H:1288
int check_int
Checkpoint output interval in iterations.
Definition REMORA.H:1469
void set_smflux(int lev)
Initialize or calculate surface momentum flux from file or analytic.
Definition REMORA.cpp:1110
void WritePlotFile(int istep)
main driver for writing AMReX plotfiles
std::string restart_chkfile
If set, restart from this checkpoint file.
Definition REMORA.H:1432
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:328
int cf_width
Nudging width at coarse-fine interface.
Definition REMORA.H:1312
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:1409
amrex::Vector< amrex::Real > t_old
old time at each level
Definition REMORA.H:1366
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_srflx
Shortwave radiation flux [W/m²], defined at rho-points.
Definition REMORA.H:391
amrex::Real last_check_file_time
Simulation time when we last output a checkpoint file.
Definition REMORA.H:1416
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:664
amrex::Vector< amrex::Real > dt
time step at each level
Definition REMORA.H:1368
static amrex::Real sum_per
Diagnostic sum output interval in time.
Definition REMORA.H:1535
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Pair
Air pressure [mb], defined at rho-points.
Definition REMORA.H:388
virtual ~REMORA()
Definition REMORA.cpp:226
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_qair
Specific humidity [kg/kg], defined at rho-points.
Definition REMORA.H:386
int hires_grid_level
Which level the high resolution bathymetry is at.
Definition REMORA.H:1554
void restart()
Definition REMORA.cpp:559
std::unique_ptr< NCTimeSeries > vbar_clim_data_from_file
Data container for vbar climatology data read from file.
Definition REMORA.H:1280
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Tair
Air temperature [°C], defined at rho-points.
Definition REMORA.H:384
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
Definition REMORA.H:346
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