12#include <AMReX_AmrCore.H>
13#include <AMReX_BCRec.H>
14#include <AMReX_InterpFaceRegister.H>
16#include <AMReX_ParallelDescriptor.H>
17#include <AMReX_ParmParse.H>
18#include <AMReX_MultiFabUtil.H>
19#include <AMReX_FillPatchUtil.H>
20#include <AMReX_VisMF.H>
21#include <AMReX_PhysBCFunct.H>
22#include <AMReX_YAFluxRegister.H>
23#include <AMReX_ErrorList.H>
25#ifdef AMREX_MEM_PROFILING
26#include <AMReX_MemProfiler.H>
39#ifdef REMORA_USE_PARTICLES
43#ifdef REMORA_USE_NETCDF
57#include <AMReX_Lazy.H>
61using amrex::MPI_COMM_WORLD;
79#ifdef REMORA_USE_NETCDF
89 :
public amrex::AmrCore
100 virtual void ErrorEst (
int lev, amrex::TagBoxArray& tags, amrex::Real time,
int ngrow)
override;
106 void init_only (
int lev, amrex::Real time);
112 void post_timestep (
int nstep, amrex::Real time, amrex::Real dt_lev);
122 const amrex::MultiFab& mf,
int comp,
bool local,
bool finemask);
126 int action_interval, amrex::Real action_per);
133 const amrex::DistributionMapping& dm)
override;
139 virtual void RemakeLevel (
int lev, amrex::Real time,
const amrex::BoxArray& ba,
140 const amrex::DistributionMapping& dm)
override;
152 const amrex::DistributionMapping& dm)
override;
174 amrex::MultiFab& source,
175 const amrex::Geometry fine_geom,
176 const amrex::Real
dt,
const amrex::Real time);
187 const amrex::Vector<const amrex::MultiFab*> &mf,
188 const amrex::Vector<const amrex::MultiFab*> &mf_nd,
189 const amrex::Vector<const amrex::MultiFab*> &mf_u,
190 const amrex::Vector<const amrex::MultiFab*> &mf_v,
191 const amrex::Vector<const amrex::MultiFab*> &mf_w,
192 const amrex::Vector<std::string> &varnames,
193 const amrex::Vector<amrex::Geometry>& my_geom,
195 const amrex::Vector<int> &level_steps,
196 const amrex::Vector<amrex::IntVect>& rr,
197 const std::string &versionName =
"HyperCLaw-V1.1",
198 const std::string &levelPrefix =
"Level_",
199 const std::string &mfPrefix =
"Cell",
200 const amrex::Vector<std::string>& extra_dirs = amrex::Vector<std::string>())
const;
206 const amrex::Vector<amrex::BoxArray> &bArray,
207 const amrex::Vector<std::string> &varnames,
208 const amrex::Vector<amrex::Geometry>& my_geom,
210 const amrex::Vector<int> &level_steps,
211 const amrex::Vector<amrex::IntVect>& rr,
212 const std::string &versionName,
213 const std::string &levelPrefix,
214 const std::string &mfPrefix)
const;
241 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_h;
244 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Hz;
246 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Huon;
248 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Hvom;
250 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_ru;
252 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rv;
254 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_ru2d;
256 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rv2d;
258 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rufrc;
260 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rvfrc;
262 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Akv;
264 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Akt;
270 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_diff2;
273 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_z_r;
276 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_z_w;
278 amrex::Gpu::DeviceVector<amrex::Real>
s_r;
280 amrex::Gpu::DeviceVector<amrex::Real>
s_w;
288 amrex::Gpu::DeviceVector<amrex::Real>
Cs_r;
290 amrex::Gpu::DeviceVector<amrex::Real>
Cs_w;
299 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_sustr;
301 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_svstr;
304 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_uwind;
306 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_vwind;
309 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_lrflx;
311 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_lhflx;
313 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_shflx;
316 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_stflx;
320 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_btflx;
325 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rain;
327 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_evap;
330 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rdrag;
334 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_ZoBot;
337 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_bustr;
339 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_bvstr;
350 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rubar;
352 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rvbar;
354 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rzeta;
356 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_ubar;
358 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_vbar;
360 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_zeta;
363 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_mskr;
365 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_msku;
367 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_mskv;
369 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_mskp;
372 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_pm;
374 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_pn;
377 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_fcor;
380 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_xr;
382 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_yr;
385 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_xu;
387 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_yu;
390 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_xv;
392 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_yv;
395 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_xp;
397 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_yp;
403 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rhoS;
405 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rhoA;
407 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_bvf;
409 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_alpha;
411 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_beta;
420 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_tke;
422 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_gls;
426 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Akk;
428 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Akp;
434 void Advance (
int lev, amrex::Real time, amrex::Real dt_lev,
int iteration,
int ncycle);
437 void setup_step(
int lev, amrex::Real time, amrex::Real dt_lev);
443 void advance_2d_onestep (
int lev, amrex::Real dt_lev, amrex::Real dtfast_lev,
int my_iif,
int nfast_counter);
447 amrex::MultiFab
const* mf_rhoS,
448 amrex::MultiFab
const* mf_rhoA,
449 amrex::MultiFab * mf_ru2d,
450 amrex::MultiFab * mf_rv2d,
451 amrex::MultiFab * mf_rufrc,
452 amrex::MultiFab * mf_rvfrc,
453 amrex::MultiFab * mf_Zt_avg1,
454 std::unique_ptr<amrex::MultiFab>& mf_DU_avg1,
455 std::unique_ptr<amrex::MultiFab>& mf_DU_avg2,
456 std::unique_ptr<amrex::MultiFab>& mf_DV_avg1,
457 std::unique_ptr<amrex::MultiFab>& mf_DV_avg2,
458 std::unique_ptr<amrex::MultiFab>& mf_rubar,
459 std::unique_ptr<amrex::MultiFab>& mf_rvbar,
460 std::unique_ptr<amrex::MultiFab>& mf_rzeta,
461 std::unique_ptr<amrex::MultiFab>& mf_ubar,
462 std::unique_ptr<amrex::MultiFab>& mf_vbar,
463 amrex::MultiFab * mf_zeta,
464 amrex::MultiFab
const* mf_h,
465 amrex::MultiFab
const* mf_pm,
466 amrex::MultiFab
const* mf_pn,
467 amrex::MultiFab
const* mf_fcor,
468 amrex::MultiFab
const* mf_visc2_p,
469 amrex::MultiFab
const* mf_visc2_r,
470 amrex::MultiFab
const* mf_mskr,
471 amrex::MultiFab
const* mf_msku,
472 amrex::MultiFab
const* mf_mskv,
473 amrex::MultiFab
const* mf_mskp,
474 amrex::Real dtfast_lev,
475 bool predictor_2d_step,
476 bool first_2d_step,
int my_iif,
481 amrex::MultiFab& mf_cons,
482 amrex::MultiFab& mf_u , amrex::MultiFab& mf_v,
483 amrex::MultiFab* mf_sstore,
484 amrex::MultiFab* mf_ru , amrex::MultiFab* mf_rv,
485 std::unique_ptr<amrex::MultiFab>& mf_DU_avg1,
486 std::unique_ptr<amrex::MultiFab>& mf_DU_avg2,
487 std::unique_ptr<amrex::MultiFab>& mf_DV_avg1,
488 std::unique_ptr<amrex::MultiFab>& mf_DV_avg2,
489 std::unique_ptr<amrex::MultiFab>& mf_ubar,
490 std::unique_ptr<amrex::MultiFab>& mf_vbar,
491 std::unique_ptr<amrex::MultiFab>& mf_Akv,
492 std::unique_ptr<amrex::MultiFab>& mf_Akt,
493 std::unique_ptr<amrex::MultiFab>& mf_Hz,
494 std::unique_ptr<amrex::MultiFab>& mf_Huon,
495 std::unique_ptr<amrex::MultiFab>& mf_Hvom,
496 std::unique_ptr<amrex::MultiFab>& mf_z_w,
497 amrex::MultiFab
const* mf_h,
498 amrex::MultiFab
const* mf_pm,
499 amrex::MultiFab
const* mf_pn,
500 amrex::MultiFab
const* mf_mskr,
501 amrex::MultiFab
const* mf_msku,
502 amrex::MultiFab
const* mf_mskv,
504 const amrex::Real dt_lev);
508 amrex::MultiFab* mf_cons,
509 amrex::MultiFab* mf_uwind,
510 amrex::MultiFab* mf_vwind,
511 amrex::MultiFab* mf_evap,
512 amrex::MultiFab* mf_sustr,
513 amrex::MultiFab* mf_svstr,
514 amrex::MultiFab* mf_stflux,
515 amrex::MultiFab* mf_lrflx,
516 amrex::MultiFab* mf_lhflx,
517 amrex::MultiFab* mf_shflx,
522 amrex::MultiFab& mf_uold, amrex::MultiFab& mf_vold,
523 amrex::MultiFab& mf_u, amrex::MultiFab& mf_v,
524 amrex::MultiFab* mf_ru,
525 amrex::MultiFab* mf_rv,
526 amrex::MultiFab& S_old,
527 amrex::MultiFab& S_new,
528 amrex::MultiFab& mf_W, amrex::MultiFab& mf_DC,
530 const amrex::MultiFab* mf_z_r,
531 const amrex::MultiFab* mf_z_w,
532 const amrex::MultiFab* mf_h,
533 const amrex::MultiFab* mf_pm,
534 const amrex::MultiFab* mf_pn,
535 const amrex::MultiFab* mf_sustr,
536 const amrex::MultiFab* mf_svstr,
537 const amrex::MultiFab* mf_bustr,
538 const amrex::MultiFab* mf_bvstr,
539 const amrex::MultiFab* mf_msku,
540 const amrex::MultiFab* mf_mskv,
541 const int iic,
const int nfirst,
542 const int nnew,
int nstp,
int nrhs,
543 int N,
const amrex::Real dt_lev);
547 const amrex::Box& gbx,
548 const amrex::Array4<amrex::Real >& tempold,
549 const amrex::Array4<amrex::Real >& tempcache,
550 const amrex::Array4<amrex::Real >& Hz,
551 const amrex::Array4<amrex::Real >& Huon,
552 const amrex::Array4<amrex::Real >& Hvom,
553 const amrex::Array4<amrex::Real >& W,
554 const amrex::Array4<amrex::Real >& DC,
555 const amrex::Array4<amrex::Real >& FC,
556 const amrex::Array4<amrex::Real >& sstore,
557 const amrex::Array4<amrex::Real const>& z_w,
558 const amrex::Array4<amrex::Real const>& h,
559 const amrex::Array4<amrex::Real const>& pm,
560 const amrex::Array4<amrex::Real const>& pn,
561 const amrex::Array4<amrex::Real const>& msku,
562 const amrex::Array4<amrex::Real const>& mskv,
563 const amrex::Array4<int const>& river_pos,
564 const amrex::Array4<amrex::Real const>& river_source,
565 int iic,
int ntfirst,
int nrhs,
int N,
566 const amrex::Real dt_lev);
570 const amrex::Box& bx,
571 const amrex::Array4<amrex::Real >& t,
572 const amrex::Array4<amrex::Real const>& tempstore,
573 const amrex::Array4<amrex::Real const>& Huon,
574 const amrex::Array4<amrex::Real const>& Hvom,
575 const amrex::Array4<amrex::Real const>& Hz,
576 const amrex::Array4<amrex::Real const>& pn,
577 const amrex::Array4<amrex::Real const>& pm,
578 const amrex::Array4<amrex::Real const>& W,
579 const amrex::Array4<amrex::Real >& FC,
580 const amrex::Array4<amrex::Real const>& mskr,
581 const amrex::Array4<amrex::Real const>& msku,
582 const amrex::Array4<amrex::Real const>& mskv,
583 const amrex::Array4<int const>& river_pos,
584 const amrex::Array4<amrex::Real const>& river_source,
585 int nrhs,
int nnew,
int N,
const amrex::Real dt_lev);
589 const amrex::Box& xbx,
590 const amrex::Box& ybx,
591 const amrex::Array4<amrex::Real const>& uold,
592 const amrex::Array4<amrex::Real const>& vold,
593 const amrex::Array4<amrex::Real >& ru,
594 const amrex::Array4<amrex::Real >& rv,
595 const amrex::Array4<amrex::Real >& rufrc,
596 const amrex::Array4<amrex::Real >& rvfrc,
597 const amrex::Array4<amrex::Real const>& sustr,
598 const amrex::Array4<amrex::Real const>& svstr,
599 const amrex::Array4<amrex::Real const>& bustr,
600 const amrex::Array4<amrex::Real const>& bvstr,
601 const amrex::Array4<amrex::Real const>& Huon,
602 const amrex::Array4<amrex::Real const>& Hvom,
603 const amrex::Array4<amrex::Real const>& pm,
604 const amrex::Array4<amrex::Real const>& pn,
605 const amrex::Array4<amrex::Real const>& W,
606 const amrex::Array4<amrex::Real >& FC,
611 const amrex::Box& xbx,
612 const amrex::Box& ybx,
613 const amrex::Array4<amrex::Real const>& uold,
614 const amrex::Array4<amrex::Real const>& vold,
615 const amrex::Array4<amrex::Real >& ru,
616 const amrex::Array4<amrex::Real >& rv,
617 const amrex::Array4<amrex::Real const>& Duon,
618 const amrex::Array4<amrex::Real const>& Dvom,
622 void rho_eos (
const amrex::Box& bx,
623 const amrex::Array4<amrex::Real const>& state,
624 const amrex::Array4<amrex::Real >& rho,
625 const amrex::Array4<amrex::Real >& rhoA,
626 const amrex::Array4<amrex::Real >& rhoS,
627 const amrex::Array4<amrex::Real >& bvf,
628 const amrex::Array4<amrex::Real >& alpha,
629 const amrex::Array4<amrex::Real >& beta,
630 const amrex::Array4<amrex::Real const>& Hz,
631 const amrex::Array4<amrex::Real const>& z_w,
632 const amrex::Array4<amrex::Real const>& z_r,
633 const amrex::Array4<amrex::Real const>& h,
634 const amrex::Array4<amrex::Real const>& mskr,
638 void lin_eos (
const amrex::Box& bx,
639 const amrex::Array4<amrex::Real const>& state,
640 const amrex::Array4<amrex::Real >& rho,
641 const amrex::Array4<amrex::Real >& rhoA,
642 const amrex::Array4<amrex::Real >& rhoS,
643 const amrex::Array4<amrex::Real >& bvf,
644 const amrex::Array4<amrex::Real const>& Hz,
645 const amrex::Array4<amrex::Real const>& z_w,
646 const amrex::Array4<amrex::Real const>& z_r,
647 const amrex::Array4<amrex::Real const>& h,
648 const amrex::Array4<amrex::Real const>& mskr,
653 const amrex::Array4<amrex::Real const>& state,
654 const amrex::Array4<amrex::Real >& rho,
655 const amrex::Array4<amrex::Real >& rhoA,
656 const amrex::Array4<amrex::Real >& rhoS,
657 const amrex::Array4<amrex::Real >& bvf,
658 const amrex::Array4<amrex::Real >& alpha,
659 const amrex::Array4<amrex::Real >& beta,
660 const amrex::Array4<amrex::Real const>& Hz,
661 const amrex::Array4<amrex::Real const>& z_w,
662 const amrex::Array4<amrex::Real const>& z_r,
663 const amrex::Array4<amrex::Real const>& h,
664 const amrex::Array4<amrex::Real const>& mskr,
668 void prsgrd (
const amrex::Box& bx,
669 const amrex::Box& gbx,
670 const amrex::Box& utbx,
671 const amrex::Box& vtbx,
672 const amrex::Array4<amrex::Real >& ru,
673 const amrex::Array4<amrex::Real >& rv,
674 const amrex::Array4<amrex::Real const>& pn,
675 const amrex::Array4<amrex::Real const>& pm,
676 const amrex::Array4<amrex::Real const>& rho,
677 const amrex::Array4<amrex::Real >& FC,
678 const amrex::Array4<amrex::Real const>& Hz,
679 const amrex::Array4<amrex::Real const>& z_r,
680 const amrex::Array4<amrex::Real const>& z_w,
681 const amrex::Array4<amrex::Real const>& msku,
682 const amrex::Array4<amrex::Real const>& mskv,
683 const int nrhs,
const int N);
688 const amrex::Box& gbx,
689 const int ioff,
const int joff,
690 const amrex::Array4<amrex::Real >& vel,
691 const amrex::Array4<amrex::Real const>& vel_old,
692 const amrex::Array4<amrex::Real >& rvel,
693 const amrex::Array4<amrex::Real const>& Hz,
694 const amrex::Array4<amrex::Real const>& Akv,
695 const amrex::Array4<amrex::Real >& FC,
696 const amrex::Array4<amrex::Real const>& sstr,
697 const amrex::Array4<amrex::Real const>& bstr,
698 const amrex::Array4<amrex::Real const>& z_r,
699 const amrex::Array4<amrex::Real const>& pm,
700 const amrex::Array4<amrex::Real const>& pn,
701 const int iic,
const int ntfirst,
const int nnew,
int nstp,
int nrhs,
int N,
702 const amrex::Real lambda,
const amrex::Real dt_lev);
706 const int ioff,
const int joff,
707 const amrex::Array4<amrex::Real >& phi,
708 const amrex::Array4<amrex::Real const>& Hz,
709 const amrex::Array4<amrex::Real >& Hzk,
710 const amrex::Array4<amrex::Real >& AK,
711 const amrex::Array4<amrex::Real const>& Akv,
712 const amrex::Array4<amrex::Real >& BC,
713 const amrex::Array4<amrex::Real >& DC,
714 const amrex::Array4<amrex::Real >& FC,
715 const amrex::Array4<amrex::Real >& CF,
716 const int nnew,
const int N,
717 const amrex::Real dt_lev);
721 const int ioff,
const int joff,
722 const amrex::Array4<amrex::Real >& phi,
723 const amrex::Array4<amrex::Real >& phibar,
724 const amrex::Array4<amrex::Real >& Hphi,
725 const amrex::Array4<amrex::Real const>& Hz,
726 const amrex::Array4<amrex::Real const>& pm_or_pn,
727 const amrex::Array4<amrex::Real const>& Dphi1,
728 const amrex::Array4<amrex::Real const>& Dphi2,
729 const amrex::Array4<amrex::Real >& DC,
730 const amrex::Array4<amrex::Real >& FC,
731 const amrex::Array4<amrex::Real const>& msk,
736 const int ioff,
const int joff,
737 const amrex::Array4<amrex::Real >& phi,
738 const amrex::Array4<amrex::Real const>& Hz,
739 const amrex::Array4<amrex::Real const>& Dphi_avg1,
740 const amrex::Array4<amrex::Real >& DC,
741 const amrex::Array4<amrex::Real >& CF,
742 const amrex::Array4<amrex::Real const>& pm_or_pn,
743 const amrex::Array4<amrex::Real const>& msk,
744 const int nnew,
const int N);
747 void uv3dmix (
const amrex::Box& xbx,
748 const amrex::Box& ybx,
749 const amrex::Array4<amrex::Real >& u,
750 const amrex::Array4<amrex::Real >& v,
751 const amrex::Array4<amrex::Real const>& uold,
752 const amrex::Array4<amrex::Real const>& vold,
753 const amrex::Array4<amrex::Real >& rufrc,
754 const amrex::Array4<amrex::Real >& rvfrc,
755 const amrex::Array4<amrex::Real const>& visc2_p,
756 const amrex::Array4<amrex::Real const>& visc2_r,
757 const amrex::Array4<amrex::Real const>& Hz,
758 const amrex::Array4<amrex::Real const>& pm,
759 const amrex::Array4<amrex::Real const>& pn,
760 const amrex::Array4<amrex::Real const>& mskp,
762 const amrex::Real dt_lev);
765 void t3dmix (
const amrex::Box& bx,
766 const amrex::Array4<amrex::Real >& state,
767 const amrex::Array4<amrex::Real >& state_rhs,
768 const amrex::Array4<amrex::Real const>& diff2,
769 const amrex::Array4<amrex::Real const>& Hz,
770 const amrex::Array4<amrex::Real const>& pm,
771 const amrex::Array4<amrex::Real const>& pn,
772 const amrex::Array4<amrex::Real const>& msku,
773 const amrex::Array4<amrex::Real const>& mskv,
774 const amrex::Real dt_lev,
778 void coriolis (
const amrex::Box& xbx,
779 const amrex::Box& ybx,
780 const amrex::Array4<amrex::Real const>& uold,
781 const amrex::Array4<amrex::Real const>& vold,
782 const amrex::Array4<amrex::Real >& ru,
783 const amrex::Array4<amrex::Real >& rv,
784 const amrex::Array4<amrex::Real const>& Hz,
785 const amrex::Array4<amrex::Real const>& fomn,
791 const amrex::Array4<amrex::Real >& var,
792 const amrex::Array4<amrex::Real const>& var_old,
793 const amrex::Array4<amrex::Real const>& var_clim,
794 const amrex::Array4<amrex::Real const>& clim_coeff,
795 const amrex::Array4<amrex::Real const>& Hz,
796 const amrex::Array4<amrex::Real const>& pm,
797 const amrex::Array4<amrex::Real const>& pn,
798 const amrex::Real dt_lev = amrex::Real(0.0));
828 void gls_prestep (
int lev, amrex::MultiFab* mf_gls, amrex::MultiFab* mf_tke,
829 amrex::MultiFab& mf_W,
830 amrex::MultiFab* mf_msku, amrex::MultiFab* mf_mskv,
831 const int nstp,
const int nnew,
const int iic,
const int ntfirst,
832 const int N,
const amrex::Real dt_lev);
834 void gls_corrector (
int lev, amrex::MultiFab* mf_gls, amrex::MultiFab* mf_tke,
835 amrex::MultiFab& mf_W, amrex::MultiFab* mf_Akv, amrex::MultiFab* mf_Akt,
836 amrex::MultiFab* mf_Akk, amrex::MultiFab* mf_Akp,
837 amrex::MultiFab* mf_mskr,
838 amrex::MultiFab* mf_msku, amrex::MultiFab* mf_mskv,
839 const int nstp,
const int nnew,
840 const int N,
const amrex::Real dt_lev);
869 void FillPatch (
int lev, amrex::Real time,
870 amrex::MultiFab& mf_to_be_filled,
871 amrex::Vector<amrex::MultiFab*>
const& mfs,
875 const bool fill_all=
true,
876 const bool fill_set=
true,
877 const int n_not_fill=0,
878 const int icomp_calc=0,
879 const amrex::Real
dt = amrex::Real(0.0),
880 const amrex::MultiFab& mf_calc = amrex::MultiFab());
884 amrex::MultiFab& mf_to_be_filled,
885 amrex::Vector<amrex::MultiFab*>
const& mfs,
888 const bool fill_all=
true,
889 const bool fill_set=
true);
892 const amrex::MultiFab& mf_mask,
893 const amrex::Real time,
895 const int bdy_var_type,
896 const int icomp_to_fill,
897 const int icomp_calc = 0,
898 const amrex::MultiFab& mf_calc = amrex::MultiFab(),
899 const amrex::Real = amrex::Real(0.0));
904 amrex::MultiFab* mf_crse,
908 const bool fill_all =
true,
909 const int n_not_fill=0,
910 const int icomp_calc=0,
911 const amrex::Real
dt = amrex::Real(0.0),
912 const amrex::MultiFab& mf_calc = amrex::MultiFab());
915 void FillCoarsePatch (
int lev, amrex::Real time, amrex::MultiFab* mf_fine,
916 amrex::MultiFab* mf_crse,
920 const bool fill_all =
true,
921 const int n_not_fill=0,
922 const int icomp_calc=0,
923 const amrex::Real
dt = amrex::Real(0.0),
924 const amrex::MultiFab& mf_calc = amrex::MultiFab());
929 amrex::MultiFab* mf_crse,
933 const bool fill_all =
true,
934 const int n_not_fill=0,
935 const int icomp_calc=0,
936 const amrex::Real
dt = amrex::Real(0.0),
937 const amrex::MultiFab& mf_calc = amrex::MultiFab(),
938 amrex::Interpolater* mapper =
nullptr);
973#ifdef REMORA_USE_NETCDF
996 void init_stuff (
int lev,
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm);
999 void init_masks (
int lev,
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm);
1017 void timeStep (
int lev, amrex::Real time,
int iteration);
1020 void timeStepML (amrex::Real time,
int iteration);
1034#ifdef REMORA_USE_NETCDF
1044 const std::string& name,
1045 bool set_ghost =
false)
const;
1049 const std::string &name,
1050 int coordinatorProc = amrex::ParallelDescriptor::IOProcessorNumber(),
1051 int allow_empty_mf = 0);
1100#ifdef REMORA_USE_MOAB
1102 void InitMOABMesh();
1143 std::unique_ptr<ProblemBase>
prob =
nullptr;
1161 amrex::Vector<amrex::Real>
dt;
1167 amrex::Vector<std::unique_ptr<REMORAPhysBCFunct>>
physbcs;
1219 amrex::Real
stop_time = std::numeric_limits<amrex::Real>::max();
1274 const amrex::Vector<std::string>
cons_names {
"temp",
"salt",
"scalar"};
1285#ifdef REMORA_USE_PARTICLES
1287 ParticleData particleData;
1291 bool m_use_tracer_particles;
1293 bool m_use_hydro_particles;
1296 void readTracersParams();
1299 void initializeTracers ( amrex::ParGDBBase*,
1300 const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& );
1303 void evolveTracers(
int lev,
1305 amrex::Vector<amrex::MultiFab const*>& flow_vel,
1306 const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& );
1327 static amrex::Vector<amrex::Vector<std::string>>
nc_init_file;
1329 static amrex::Vector<amrex::Vector<std::string>>
nc_grid_file;
1378 int nGhostCells = 0;
1379 switch (spatial_order) {
1396 amrex::Error(
"Must specify spatial order to be 2,3,4,5 or 6");
1433 int numCores = amrex::ParallelDescriptor::NProcs();
1435 numCores = numCores * omp_get_max_threads();
1439 numCores * (amrex::Real(amrex::ParallelDescriptor::second()) -
startCPUTime) +
1447 if (amrex::ParallelDescriptor::IOProcessor())
1449 datalog[i] = std::make_unique<std::fstream>();
1450 datalog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1452 amrex::FileOpenFailed(filename);
1455 amrex::ParallelDescriptor::Barrier(
"REMORA::setRecordDataInfo");
1458 amrex::Vector<std::unique_ptr<std::fstream> >
datalog;
1471 static void print_usage(MPI_Comm , std::ostream& );
1472 static void print_error(MPI_Comm ,
const std::string& msg);
1480 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1483 using namespace amrex;
1488 Real
x = std::pow(1.0_rt-15.0_rt*ZoL,0.25_rt);
1489 Real psik = 2.0_rt*std::log(0.5_rt*(1.0_rt+
x))+
1490 std::log(0.5_rt*(1.0_rt+
x*
x))-
1491 2.0_rt*std::atan(
x)+0.5_rt*
PI;
1493 Real sqrt3 = std::sqrt(3.0_rt);
1494 Real
y = std::pow(1.0_rt-10.15_rt*ZoL,1.0_rt/3.0_rt);
1495 Real psic = 1.5_rt*std::log(1.0_rt/3.0_rt*(1.0_rt+
y+
y*
y))-
1496 sqrt3*std::atan((1.0_rt+2.0_rt*
y)/sqrt3)+
PI/sqrt3;
1498 Real Zol2 = ZoL*ZoL;
1499 Real Fw = Zol2/(1.0_rt+Zol2);
1501 psiu = (1.0_rt-Fw)*psik+Fw*psic;
1504 Real cff=std::min(50.0_rt,0.35_rt*ZoL);
1505 psiu = -((1.0_rt+ZoL)+0.6667_rt*(ZoL-14.28_rt)/
1506 std::exp(cff)+8.525_rt);
1512 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1515 using namespace amrex;
1520 Real
x = std::pow(1.0_rt-15.0_rt*ZoL,0.5_rt);
1521 Real psik = 2.0_rt*std::log(0.5_rt*(1.0_rt+
x));
1523 Real sqrt3=std::sqrt(3.0_rt);
1524 Real
y=std::pow(1.0_rt-34.15_rt*ZoL,1.0_rt/3.0_rt);
1525 Real psic=1.5_rt*std::log(1.0_rt/3.0_rt*(1.0_rt+
y+
y*
y))-
1526 sqrt3*std::atan((1.0_rt+2.0_rt*
y)/sqrt3)+
PI/sqrt3;
1530 Real Fw=ZoL2/(1.0_rt+ZoL2);
1531 psit = (1.0_rt-Fw)*psik+Fw*psic;
1534 Real cff=std::min(50.0_rt,0.35_rt*ZoL);
1535 psit=-(std::pow(1.0_rt+2.0_rt*ZoL,1.5_rt)+
1536 0.6667_rt*(ZoL-14.28_rt)/std::exp(cff)+8.525_rt);
PlotfileType
plotfile format
A class to hold and interpolate time series data read from a NetCDF file.
A class to hold and interpolate time series data read from a NetCDF file.
Class that stores all relevant simulation state data with methods for time stepping.
static PlotfileType plotfile_type
Native or NetCDF plotfile output.
NCTimeSeriesRiver * river_source_transportbar
Data container for vertically integrated momentum transport in rivers.
const amrex::Vector< std::string > cons_names
Names of scalars for plotfile output.
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...
std::string riv_time_varname
Name of time field for river time.
static void GotoNextLine(std::istream &is)
utility to skip to next line in Header
int nfast
Number of fast steps to take.
void WriteNCMultiFab(const amrex::FabArray< amrex::FArrayBox > &fab, const std::string &name, bool set_ghost=false) const
Write MultiFab in NetCDF format.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
void prsgrd(const amrex::Box &bx, const amrex::Box &gbx, const amrex::Box &utbx, const amrex::Box &vtbx, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &rho, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const int nrhs, const int N)
Calculate pressure gradient.
static amrex::Real fixed_dt
User specified fixed baroclinic time step.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_evap
evaporation rate [kg/m^2/s]
amrex::Real last_plot_file_time
Simulation time when we last output a plotfile.
void update_massflux_3d(int lev, const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real > &phibar, const amrex::Array4< amrex::Real > &Hphi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm_or_pn, const amrex::Array4< amrex::Real const > &Dphi1, const amrex::Array4< amrex::Real const > &Dphi2, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &msk, const int nnew)
Correct mass flux.
void scale_rhs_vars()
Scale RHS momentum variables by 1/cell area, needed before FillPatch to different levels.
amrex::Vector< NCTimeSeriesRiver * > river_source_cons
Vector of data containers for scalar data in rivers.
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.
amrex::Vector< amrex::BCRec > domain_bcs_type
vector (over BCVars) of BCRecs
bool set_bcs_by_var
whether to set boundary conditions by variable rather than just by side
void calculate_nodal_masks(int lev)
Calculate u-, v-, and psi-point masks based on rho-point masks after analytic initialization.
static amrex::Real previousCPUTimeUsed
Accumulator variable for CPU time used thusfar.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_btflux
Bottom tracer flux; input arrays.
void init_gls_vmix(int lev, SolverChoice solver_choice)
Initialize GLS variables.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rubar
barotropic x velocity for the RHS (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
Bathymetry data (2D, positive valued, h in ROMS)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
amrex::MultiFab fine_mask
Mask that zeroes out values on a coarse level underlying grids on the next finest level.
amrex::Vector< REMORAFillPatcher > FPr_v
Vector over levels of FillPatchers for v (3D)
amrex::Vector< std::string > plot_var_names
Names of variables to output to AMReX plotfile.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ZoBot
Bottom roughness length [m], defined at rho points.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg2
correct time average of barotropic x velocity flux for coupling (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lrflx
longwave radiation
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yv
y_grid on v-points (2D)
static void print_error(MPI_Comm, const std::string &msg)
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Make a new level using provided BoxArray and DistributionMapping and fill with interpolated coarse le...
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.
amrex::Gpu::DeviceVector< amrex::Real > s_w
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at w-points (cell faces)
amrex::GpuArray< amrex::GpuArray< REMORA_BC, AMREX_SPACEDIM *2 >, BCVars::NumTypes > phys_bc_type
Array holding the "physical" boundary condition types (e.g. "inflow")
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vwind
Wind in the v direction, defined at rho-points.
std::unique_ptr< ProblemBase > prob
Pointer to container of analytical functions for problem definition.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
void Construct_REMORAFillPatchers(int lev)
Construct FillPatchers.
void advance_3d(int lev, amrex::MultiFab &mf_cons, amrex::MultiFab &mf_u, amrex::MultiFab &mf_v, amrex::MultiFab *mf_sstore, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, std::unique_ptr< amrex::MultiFab > &mf_DU_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg2, std::unique_ptr< amrex::MultiFab > &mf_DV_avg1, std::unique_ptr< amrex::MultiFab > &mf_DV_avg2, std::unique_ptr< amrex::MultiFab > &mf_ubar, std::unique_ptr< amrex::MultiFab > &mf_vbar, std::unique_ptr< amrex::MultiFab > &mf_Akv, std::unique_ptr< amrex::MultiFab > &mf_Akt, std::unique_ptr< amrex::MultiFab > &mf_Hz, std::unique_ptr< amrex::MultiFab > &mf_Huon, std::unique_ptr< amrex::MultiFab > &mf_Hvom, std::unique_ptr< amrex::MultiFab > &mf_z_w, amrex::MultiFab const *mf_h, amrex::MultiFab const *mf_pm, amrex::MultiFab const *mf_pn, amrex::MultiFab const *mf_mskr, amrex::MultiFab const *mf_msku, amrex::MultiFab const *mf_mskv, const int N, const amrex::Real dt_lev)
Advance the 3D variables.
static int sum_interval
Diagnostic sum output interval in number of steps.
int history_count
Counter for which time index we are writing to in the netcdf history file.
amrex::Real stop_time
Time to stop.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rain
precipitation rate [kg/m^2/s]
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_tke
Turbulent kinetic energy.
void rho_eos(const amrex::Box &bx, const amrex::Array4< amrex::Real const > &state, const amrex::Array4< amrex::Real > &rho, const amrex::Array4< amrex::Real > &rhoA, const amrex::Array4< amrex::Real > &rhoS, const amrex::Array4< amrex::Real > &bvf, const amrex::Array4< amrex::Real > &alpha, const amrex::Array4< amrex::Real > &beta, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &mskr, const int N)
Wrapper around equation of state calculation.
static void print_summary(std::ostream &)
int do_substep
Whether to substep fine levels in time.
NCTimeSeries * vbar_clim_data_from_file
Data container for vbar climatology data read from file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflx
Surface tracer flux; working arrays.
void Evolve()
Advance solution to final time.
NCTimeSeries * Uwind_data_from_file
Data container for u-direction wind read from file.
std::string bdry_time_varname
Default name of time field for boundary data.
void ReadCheckpointFile()
read checkpoint file from disk
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_gls
Turbulent generic length scale.
bool chunk_history_file
Whether to chunk netcdf history file.
void writeJobInfo(const std::string &dir) const
Write job info to stdout.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
amrex::Real get_t_old(int lev) const
Accessor method for t_old to expose to outside classes.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yp
y_grid on psi-points (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xr
x_grid on rho points (2D)
NCTimeSeries * v_clim_data_from_file
Data container for v-velocity climatology data read from file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
virtual void ClearLevel(int lev) override
Delete level data Overrides the pure virtual function in AmrCore.
amrex::Vector< std::string > datalogname
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
AMREX_FORCE_INLINE int ComputeGhostCells(const int &spatial_order)
Helper function to determine number of ghost cells.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sstore
additional scratch space for calculations on temp, salt, etc
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xv
x_grid on v-points (2D)
void init_only(int lev, amrex::Real time)
Init (NOT restart or regrid)
void init_set_vmix(int lev)
Initialize vertical mixing coefficients from file or analytic.
std::string clim_u_time_varname
Name of time field for u climatology data.
void rhs_uv_3d(int lev, const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real > &rufrc, const amrex::Array4< amrex::Real > &rvfrc, const amrex::Array4< amrex::Real const > &sustr, const amrex::Array4< amrex::Real const > &svstr, const amrex::Array4< amrex::Real const > &bustr, const amrex::Array4< amrex::Real const > &bvstr, const amrex::Array4< amrex::Real const > &Huon, const amrex::Array4< amrex::Real const > &Hvom, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &W, const amrex::Array4< amrex::Real > &FC, int nrhs, int N)
RHS terms for 3D momentum.
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.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Lscale
Vertical mixing turbulent length scale.
amrex::Vector< REMORAFillPatcher > FPr_u
Vector over levels of FillPatchers for u (3D)
std::string clim_temp_time_varname
Name of time field for temperature climatology data.
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
void prestep_diffusion(const amrex::Box &bx, const amrex::Box &gbx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &vel, const amrex::Array4< amrex::Real const > &vel_old, const amrex::Array4< amrex::Real > &rvel, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &Akv, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &sstr, const amrex::Array4< amrex::Real const > &bstr, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const int iic, const int ntfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real lambda, const amrex::Real dt_lev)
Update velocities or tracers with diffusion/viscosity as the last part of the prestep.
void vert_mean_3d(const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &Dphi_avg1, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &CF, const amrex::Array4< amrex::Real const > &pm_or_pn, const amrex::Array4< amrex::Real const > &msk, const int nnew, const int N)
Adjust 3D momentum variables based on vertical mean momentum.
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Holds info for dynamically generated tagging criteria.
std::string nc_frc_file
NetCDF forcing file.
void t3dmix(const amrex::Box &bx, const amrex::Array4< amrex::Real > &state, const amrex::Array4< amrex::Real > &state_rhs, const amrex::Array4< amrex::Real const > &diff2, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Real dt_lev, const int ncomp)
Harmonic diffusivity for tracers.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hz
Width of cells in the vertical (z-) direction (3D, Hz in ROMS)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_msku
land/sea mask at x-faces (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvfrc
v velocity RHS, integrated, including advection and bottom/surface stresses (2D)
std::string clim_ubar_time_varname
Name of time field for ubar climatology data.
amrex::Vector< NCTimeSeriesBoundary * > boundary_series
Vector over BdyVars of boundary series data containers.
void WriteNCPlotFile(int istep)
Write plotfile using NetCDF (wrapper)
std::string check_file
Checkpoint file prefix.
void FillPatchNoBC(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true)
Fill a new MultiFab by copying in phi from valid region and filling ghost cells without applying boun...
static amrex::Real startCPUTime
Variable for CPU timing.
void setPlotVariables(const std::string &pp_plot_var_names)
NCTimeSeries * u_clim_data_from_file
Data container for u-velocity climatology data read from file.
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
amrex::Real start_time
Time of the start of the simulation, in seconds.
NCTimeSeries * ubar_clim_data_from_file
Data container for ubar climatology data read from file.
void init_data_from_netcdf(int lev)
Problem initialization from NetCDF file.
void WriteNCPlotFile_which(int lev, int which_subdomain, bool write_header, ncutils::NCFile &ncf, bool is_history)
Write a particular NetCDF plotfile.
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)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_uwind
Wind in the u direction, defined at rho-points.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rufrc
u velocity RHS, integrated, including advection and bottom/surface stresses (2D)
static amrex::Real fixed_fast_dt
User specified fixed barotropic time step.
int regrid_int
how often each level regrids the higher levels of refinement (after a level advances that many time s...
amrex::Real check_int_time
Checkpoint output interval in seconds.
amrex::Real bdy_time_interval
Interval between boundary data times.
void init_bdry_from_netcdf()
Boundary data initialization from NetCDF file.
void gls_prestep(int lev, amrex::MultiFab *mf_gls, amrex::MultiFab *mf_tke, amrex::MultiFab &mf_W, amrex::MultiFab *mf_msku, amrex::MultiFab *mf_mskv, const int nstp, const int nnew, const int iic, const int ntfirst, const int N, const amrex::Real dt_lev)
Prestep for GLS calculation.
amrex::Real start_bdy_time
Start time in the time series of boundary data.
NCTimeSeries * sustr_data_from_file
Data container for u-component surface momentum flux read from file.
amrex::Vector< amrex::Real > vec_weight2
Weights for calculating avg2 in 2D advance.
amrex::Vector< int > bdy_index
Container to connect boundary data being read in boundary condition containers.
amrex::Gpu::DeviceVector< amrex::Real > s_r
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at rho points (cell centers)
void Define_REMORAFillPatchers(int lev)
Define FillPatchers.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_shflx
sensible heat flux
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_p
Harmonic viscosity defined on the psi points (corners of horizontal grid cells)
void set_drag(int lev)
Initialize or calculate bottom drag.
amrex::MultiFab & build_fine_mask(int lev)
Make mask to zero out covered cells (for mesh refinement)
amrex::Real plot_int_time
Plotfile output interval in seconds.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvbar
barotropic y velocity for the RHS (2D)
amrex::Vector< int > num_files_at_level
how many netcdf input files specified at each level
void advance_2d_onestep(int lev, amrex::Real dt_lev, amrex::Real dtfast_lev, int my_iif, int nfast_counter)
2D advance, one predictor/corrector step
amrex::Vector< REMORAFillPatcher > FPr_vbar
Vector over levels of FillPatchers for vbar (2D)
void AverageDownTo(int crse_lev)
more flexible version of AverageDown() that lets you average down across multiple levels
int steps_per_history_file
Number of time steps per netcdf history file.
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Called after every level 0 timestep.
int max_step
maximum number of steps
amrex::Vector< amrex::MultiFab * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_r
z coordinates at rho points (cell centers)
amrex::Vector< int > num_boxes_at_level
how many boxes specified at each level by tagging criteria
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lhflx
latent heat flux
void refinement_criteria_setup()
Set refinement criteria.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskp
land/sea mask at cell corners (2D)
void remora_advance(int level, amrex::MultiFab &cons_old, amrex::MultiFab &cons_new, amrex::MultiFab &xvel_old, amrex::MultiFab &yvel_old, amrex::MultiFab &zvel_old, amrex::MultiFab &xvel_new, amrex::MultiFab &yvel_new, amrex::MultiFab &zvel_new, amrex::MultiFab &source, const amrex::Geometry fine_geom, const amrex::Real dt, const amrex::Real time)
Interface for advancing the data at one level by one "slow" timestep.
void WriteMultiLevelPlotfileWithBathymetry(const std::string &plotfilename, int nlevels, const amrex::Vector< const amrex::MultiFab * > &mf, const amrex::Vector< const amrex::MultiFab * > &mf_nd, const amrex::Vector< const amrex::MultiFab * > &mf_u, const amrex::Vector< const amrex::MultiFab * > &mf_v, const amrex::Vector< const amrex::MultiFab * > &mf_w, const amrex::Vector< std::string > &varnames, const amrex::Vector< amrex::Geometry > &my_geom, amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &rr, const std::string &versionName="HyperCLaw-V1.1", const std::string &levelPrefix="Level_", const std::string &mfPrefix="Cell", const amrex::Vector< std::string > &extra_dirs=amrex::Vector< std::string >()) const
write out particular data to an AMReX plotfile
int last_check_file_step
Step when we last output a checkpoint file.
void prestep(int lev, amrex::MultiFab &mf_uold, amrex::MultiFab &mf_vold, amrex::MultiFab &mf_u, amrex::MultiFab &mf_v, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, amrex::MultiFab &S_old, amrex::MultiFab &S_new, amrex::MultiFab &mf_W, amrex::MultiFab &mf_DC, const amrex::MultiFab *mf_z_r, const amrex::MultiFab *mf_z_w, const amrex::MultiFab *mf_h, const amrex::MultiFab *mf_pm, const amrex::MultiFab *mf_pn, const amrex::MultiFab *mf_sustr, const amrex::MultiFab *mf_svstr, const amrex::MultiFab *mf_bustr, const amrex::MultiFab *mf_bvstr, const amrex::MultiFab *mf_msku, const amrex::MultiFab *mf_mskv, const int iic, const int nfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real dt_lev)
Wrapper function for prestep.
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.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvf
Brunt-Vaisala frequency (3D)
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Remake an existing level using provided BoxArray and DistributionMapping and fill with existing fine ...
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskv
land/sea mask at y-faces (2D)
void init_masks(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Allocate MultiFabs for masks.
amrex::Vector< int > nsubsteps
How many substeps on each level?
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
void ComputeDt()
a wrapper for estTimeStep()
static void writeBuildInfo(std::ostream &os)
Write build info to os.
virtual void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Tag cells for refinement.
int plot_int
Plotfile output interval in iterations.
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.
amrex::Real getCPUTime() const
Get CPU time used.
void InitData()
Initialize multilevel data.
amrex::Vector< int > istep
which step?
void apply_clim_nudg(const amrex::Box &bx, int ioff, int joff, const amrex::Array4< amrex::Real > &var, const amrex::Array4< amrex::Real const > &var_old, const amrex::Array4< amrex::Real const > &var_clim, const amrex::Array4< amrex::Real const > &clim_coeff, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Real dt_lev=amrex::Real(0.0))
Apply climatology nudging.
void uv3dmix(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real > &u, const amrex::Array4< amrex::Real > &v, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &rufrc, const amrex::Array4< amrex::Real > &rvfrc, const amrex::Array4< amrex::Real const > &visc2_p, const amrex::Array4< amrex::Real const > &visc2_r, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &mskp, int nrhs, int nnew, const amrex::Real dt_lev)
Harmonic viscosity.
void WriteCheckpointFile()
write checkpoint file to disk
std::string nc_clim_coeff_file
NetCDF climatology coefficient file.
void lin_eos(const amrex::Box &bx, const amrex::Array4< amrex::Real const > &state, const amrex::Array4< amrex::Real > &rho, const amrex::Array4< amrex::Real > &rhoA, const amrex::Array4< amrex::Real > &rhoS, const amrex::Array4< amrex::Real > &bvf, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &mskr, const int N)
Calculate density and related quantities from linear equation of state.
void setRecordDataInfo(int i, const std::string &filename)
void advance_3d_ml(int lev, amrex::Real dt_lev)
3D advance on a single level
const std::string DataLogName(int i) const noexcept
The filename of the ith datalog file.
void set_analytic_vmix(int lev)
Set vertical mixing coefficients from analytic.
void mask_arrays_for_write(int lev, amrex::Real fill_value, amrex::Real fill_where)
Mask data arrays before writing output.
void init_flat_bathymetry(int lev)
Initialize flat bathymetry to value from problo.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoS
density perturbation
void WriteGenericPlotfileHeaderWithBathymetry(std::ostream &HeaderFile, int nlevels, const amrex::Vector< amrex::BoxArray > &bArray, const amrex::Vector< std::string > &varnames, const amrex::Vector< amrex::Geometry > &my_geom, amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &rr, const std::string &versionName, const std::string &levelPrefix, const std::string &mfPrefix) const
write out header data for an AMReX plotfile
void fill_from_bdyfiles(amrex::MultiFab &mf_to_fill, const amrex::MultiFab &mf_mask, const amrex::Real time, const int bccomp, const int bdy_var_type, const int icomp_to_fill, const int icomp_calc=0, const amrex::MultiFab &mf_calc=amrex::MultiFab(), const amrex::Real=amrex::Real(0.0))
Fill boundary data from netcdf file.
amrex::Vector< std::string > bdry_time_name_byvar
Name of time fields for boundary data.
void rhs_uv_2d(int lev, const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &Duon, const amrex::Array4< amrex::Real const > &Dvom, const int nrhs)
RHS terms for 2D momentum.
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
void init_riv_pos_from_netcdf(int lev)
static amrex::Vector< std::string > nc_bdry_file
NetCDF boundary data.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_r
Harmonic viscosity defined on the rho points (centers)
void update_mskp(int lev)
Set psi-point mask to be consistent with rho-point mask.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Huon
u-volume flux (3D)
void FillBdyCCVels(int lev, amrex::MultiFab &mf_cc_vel)
Fill the physical boundary conditions for cell-centered velocity (diagnostic only)
void set_zeta(int lev)
Initialize zeta from file or analytic.
static amrex::Real change_max
Fraction maximum change in subsequent time steps.
void calc_stretch_coeffs()
calculate vertical stretch coefficients
void bulk_fluxes(int lev, amrex::MultiFab *mf_cons, amrex::MultiFab *mf_uwind, amrex::MultiFab *mf_vwind, amrex::MultiFab *mf_evap, amrex::MultiFab *mf_sustr, amrex::MultiFab *mf_svstr, amrex::MultiFab *mf_stflux, amrex::MultiFab *mf_lrflx, amrex::MultiFab *mf_lhflx, amrex::MultiFab *mf_shflx, const int N)
Calculate bulk temperature, salinity, wind fluxes.
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.
static void print_usage(MPI_Comm, std::ostream &)
std::string pp_prefix
default prefix for input file parameters
void set_bathymetry(int lev)
Initialize bathymetry from file or analytic.
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
void init_stuff(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Allocate MultiFabs for state and evolution variables.
TimeInterpolatedData GetDataAtTime(int lev, amrex::Real time)
utility to copy in data from old and/or new state into another multifab
void Advance(int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle)
advance a single level for a single time step
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoA
vertically-averaged density
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg1
time average of barotropic y velocity flux
amrex::Vector< REMORAFillPatcher > FPr_c
Vector over levels of FillPatchers for scalars.
void rhs_t_3d(int lev, const amrex::Box &bx, const amrex::Array4< amrex::Real > &t, const amrex::Array4< amrex::Real const > &tempstore, const amrex::Array4< amrex::Real const > &Huon, const amrex::Array4< amrex::Real const > &Hvom, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &W, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &mskr, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Array4< int const > &river_pos, const amrex::Array4< amrex::Real const > &river_source, int nrhs, int nnew, int N, const amrex::Real dt_lev)
RHS terms for tracer.
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real bulk_psiu(amrex::Real ZoL)
Evaluate stability function psi for wind speed.
static void print_banner(MPI_Comm, std::ostream &)
std::string clim_v_time_varname
Name of time field for v climatology data.
void scale_rhs_vars_inv()
Scale RHS momentum variables by cell area, needed after FillPatch to different levels.
amrex::Vector< REMORAFillPatcher > FPr_w
Vector over levels of FillPatchers for w.
amrex::Vector< amrex::YAFluxRegister * > advflux_reg
array of flux registers for refluxing in multilevel
std::string clim_salt_time_varname
Name of time field for salinity climatology data.
amrex::Gpu::DeviceVector< amrex::Real > Cs_r
Stretching coefficients at rho points.
NCTimeSeries * Vwind_data_from_file
Data container for v-direction wind read from file.
amrex::Vector< amrex::Real > t_new
new time at each level
NCTimeSeriesRiver * river_source_transport
Data container for momentum transport in rivers.
void init_stretch_coeffs()
initialize and calculate stretch coefficients
NCTimeSeries * temp_clim_data_from_file
Data container for temperature climatology data read from file.
static SolverChoice solverChoice
Container for algorithmic choices.
void resize_stuff(int lev)
Resize variable containers to accommodate data on levels 0 to max_lev.
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > vec_river_position
iMultiFab for river positions; contents are indices of rivers
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akk
Turbulent kinetic energy vertical diffusion coefficient.
void set_masks(int lev)
Initialize land-sea masks from file or analytic.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag2
Quadratic drag coefficient [unitless], defined at rho points.
int cf_set_width
Width for fixing values at coarse-fine interface.
NCTimeSeries * salt_clim_data_from_file
Data container for salinity climatology data read from file.
void ReadParameters()
read in some parameters from inputs file
static void print_tpls(std::ostream &)
void vert_visc_3d(const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real > &Hzk, const amrex::Array4< amrex::Real > &AK, const amrex::Array4< amrex::Real const > &Akv, const amrex::Array4< amrex::Real > &BC, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real > &CF, const int nnew, const int N, const amrex::Real dt_lev)
Calculate effects of vertical viscosity or diffusivity.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
NCTimeSeries * svstr_data_from_file
Data container for v-component surface momentum flux read from file.
void sum_integrated_quantities(amrex::Real time)
Integrate conserved quantities for diagnostics.
static int total_nc_plot_file_step
amrex::Vector< amrex::Real > vec_weight1
Weights for calculating avg1 in 2D advance.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xp
x_grid on psi-points (2D)
static int plot_staggered_vels
Whether to write the staggered velocities (not averaged to cell centers)
static amrex::Vector< amrex::Vector< std::string > > nc_grid_file
NetCDF grid file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
amrex::Gpu::DeviceVector< int > river_direction
Vector over rivers of river direction: 0: u-face; 1: v-face; 2: w-face.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
void ReadNCMultiFab(amrex::FabArray< amrex::FArrayBox > &fab, const std::string &name, int coordinatorProc=amrex::ParallelDescriptor::IOProcessorNumber(), int allow_empty_mf=0)
Read MultiFab in NetCDF format.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg1
time average of barotropic x velocity flux (2D)
amrex::Gpu::DeviceVector< amrex::Real > Cs_w
Stretching coefficients at w points.
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
void advance_2d(int lev, amrex::MultiFab const *mf_rhoS, amrex::MultiFab const *mf_rhoA, amrex::MultiFab *mf_ru2d, amrex::MultiFab *mf_rv2d, amrex::MultiFab *mf_rufrc, amrex::MultiFab *mf_rvfrc, amrex::MultiFab *mf_Zt_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg2, std::unique_ptr< amrex::MultiFab > &mf_DV_avg1, std::unique_ptr< amrex::MultiFab > &mf_DV_avg2, std::unique_ptr< amrex::MultiFab > &mf_rubar, std::unique_ptr< amrex::MultiFab > &mf_rvbar, std::unique_ptr< amrex::MultiFab > &mf_rzeta, std::unique_ptr< amrex::MultiFab > &mf_ubar, std::unique_ptr< amrex::MultiFab > &mf_vbar, amrex::MultiFab *mf_zeta, amrex::MultiFab const *mf_h, amrex::MultiFab const *mf_pm, amrex::MultiFab const *mf_pn, amrex::MultiFab const *mf_fcor, amrex::MultiFab const *mf_visc2_p, amrex::MultiFab const *mf_visc2_r, amrex::MultiFab const *mf_mskr, amrex::MultiFab const *mf_msku, amrex::MultiFab const *mf_mskv, amrex::MultiFab const *mf_mskp, amrex::Real dtfast_lev, bool predictor_2d_step, bool first_2d_step, int my_iif, int &next_indx1)
Perform a 2D predictor (predictor_2d_step=True) or corrector (predictor_2d_step=False) step.
int plot_file_on_restart
Whether to output a plotfile on restart from checkpoint.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_alpha
Thermal expansion coefficient (3D)
void set_2darrays(int lev)
Set 2D momentum arrays from 3D momentum.
void FillCoarsePatchMap(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(), amrex::Interpolater *mapper=nullptr)
fill an entire multifab by interpolating from the coarser level, explicitly specifying interpolator t...
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NCONS+8 > m_bc_extdir_vals
Array holding the Dirichlet values at walls which need them.
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)
void InitializeLevelFromData(int lev, const amrex::MultiFab &initial_data)
Initialize the new-time data at a level from the initial_data MultiFab.
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive scalar
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bustr
Bottom stress in the u direction.
AMREX_FORCE_INLINE int NumDataLogs() noexcept
amrex::Real volWgtSumMF(int lev, const amrex::MultiFab &mf, int comp, bool local, bool finemask)
Perform the volume-weighted sum.
std::string frc_time_varname
Name of time field for forcing data.
void set_wind(int lev)
Initialize or calculate wind speed from file or analytic.
void init_clim_nudg_coeff_from_netcdf(int lev)
Climatology nudging coefficient initialization from NetCDF file.
amrex::Vector< REMORAFillPatcher > FPr_ubar
Vector over levels of FillPatchers for ubar (2D)
bool is_it_time_for_action(int nstep, amrex::Real time, amrex::Real dt, int action_interval, amrex::Real action_per)
Decide if it is time to take an action.
void appendPlotVariables(const std::string &pp_plot_var_names)
std::string PlotFileName(int lev) const
get plotfile name
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg2
correct time average of barotropic y velocity flux for coupling (2D)
void set_hmixcoef(int lev)
Initialize horizontal mixing coefficients.
amrex::Array< std::string, 2 *AMREX_SPACEDIM > domain_bc_type
Array of strings describing domain boundary conditions.
amrex::Real estTimeStep(int lev) const
compute dt from CFL considerations
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yu
y_grid on u-points (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvstr
Bottom stress in the v direction.
std::string nc_clim_his_file
NetCDF climatology history file.
void timeStep(int lev, amrex::Real time, int iteration)
advance a level by dt, includes a recursive call for finer levels
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rzeta
free surface height for the RHS (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_phys_nd
z coordinates at psi points (cell nodes)
AMREX_FORCE_INLINE amrex::YAFluxRegister * getAdvFluxReg(int lev)
Get flux register.
void AverageDown()
set covered coarse cells to be the average of overlying fine cells
static int fixed_ndtfast_ratio
User specified, number of barotropic steps per baroclinic step.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xu
x_grid on u-points (2D)
void timeStepML(amrex::Real time, int iteration)
advance all levels by dt, loops over finer levels
amrex::GpuArray< amrex::GpuArray< bool, AMREX_SPACEDIM *2 >, BdyVars::NumTypes+1 > phys_bc_need_data
These are flags that indicate whether we need to read in boundary data from file.
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > vec_nudg_coeff
Climatology nudging coefficients.
void set_weights(int lev)
Set weights for averaging 3D variables to 2D.
void prestep_t_advection(int lev, const amrex::Box &tbx, const amrex::Box &gbx, const amrex::Array4< amrex::Real > &tempold, const amrex::Array4< amrex::Real > &tempcache, const amrex::Array4< amrex::Real > &Hz, const amrex::Array4< amrex::Real > &Huon, const amrex::Array4< amrex::Real > &Hvom, const amrex::Array4< amrex::Real > &W, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real > &sstore, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Array4< int const > &river_pos, const amrex::Array4< amrex::Real const > &river_source, int iic, int ntfirst, int nrhs, int N, const amrex::Real dt_lev)
Prestep advection calculations for the tracers.
AMREX_FORCE_INLINE std::ostream & DataLog(int i)
Helper function for IO stream.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
amrex::Vector< std::unique_ptr< std::fstream > > datalog
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
static amrex::Real cfl
CFL condition.
virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Make a new level from scratch using provided BoxArray and DistributionMapping. Only used during initi...
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflux
Surface tracer flux; input arrays.
static int verbose
Verbosity level of output.
void setup_step(int lev, amrex::Real time, amrex::Real dt_lev)
Set everything up for a step on a level.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag
Linear drag coefficient [m/s], defined at rho points.
void nonlin_eos(const amrex::Box &bx, const amrex::Array4< amrex::Real const > &state, const amrex::Array4< amrex::Real > &rho, const amrex::Array4< amrex::Real > &rhoA, const amrex::Array4< amrex::Real > &rhoS, const amrex::Array4< amrex::Real > &bvf, const amrex::Array4< amrex::Real > &alpha, const amrex::Array4< amrex::Real > &beta, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &mskr, const int N)
Calculate density and related quantities from nonlinear equation of state.
std::string plot_file_name
Plotfile prefix.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
int check_int
Checkpoint output interval in iterations.
void WritePlotFile()
main driver for writing AMReX plotfiles
void set_smflux(int lev)
Initialize or calculate surface momentum flux from file or analytic.
std::string restart_chkfile
If set, restart from this checkpoint file.
void init_clim_nudg_coeff(int lev)
Wrapper to initialize climatology nudging coefficient.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_btflx
Bottom tracer flux; working arrays.
int cf_width
Nudging width at coarse-fine interface.
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real bulk_psit(amrex::Real ZoL)
Evaluate stability function psi for moisture and heat.
static amrex::Vector< amrex::Vector< std::string > > nc_init_file
NetCDF initialization file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_beta
Saline contraction coefficient (3D)
int last_plot_file_step
Step when we last output a plotfile.
const amrex::Vector< std::string > derived_names
Names of derived fields for plotfiles.
amrex::Vector< amrex::Real > t_old
old time at each level
amrex::Real last_check_file_time
Simulation time when we last output a checkpoint file.
amrex::Vector< amrex::Real > dt
time step at each level
static amrex::Real sum_per
Diagnostic sum output interval in time.
void InitializeFromFile()
Read the file passed to remora.restart and use it as an initial condition for the current simulation.
void coriolis(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &fomn, int nrhs, int nr)
Calculate Coriolis terms.
amrex::Gpu::DeviceVector< amrex::BCRec > domain_bcs_type_d
GPU vector (over BCVars) of BCRecs.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_w
z coordinates at w points (faces between z-cells)
void convert_inv_days_to_inv_s(amrex::MultiFab *)
Convert data in a multifab from inverse days to inverse seconds.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yr
y_grid on rho points (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hvom
v-volume flux (3D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akp
Turbulent length scale vertical diffusion coefficient.
void gls_corrector(int lev, amrex::MultiFab *mf_gls, amrex::MultiFab *mf_tke, amrex::MultiFab &mf_W, amrex::MultiFab *mf_Akv, amrex::MultiFab *mf_Akt, amrex::MultiFab *mf_Akk, amrex::MultiFab *mf_Akp, amrex::MultiFab *mf_mskr, amrex::MultiFab *mf_msku, amrex::MultiFab *mf_mskv, const int nstp, const int nnew, const int N, const amrex::Real dt_lev)
Corrector step for GLS calculation.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
@ CellConservativeQuartic
@ CellConservativeProtected