REMORA
Energy Research and Forecasting: An Atmospheric Modeling Code
|
#include <REMORA.H>
Public Member Functions | |
REMORA () | |
virtual | ~REMORA () |
void | Evolve () |
virtual void | ErrorEst (int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override |
void | InitData () |
void | init_only (int lev, amrex::Real time) |
void | restart () |
void | post_timestep (int nstep, amrex::Real time, amrex::Real dt_lev) |
void | sum_integrated_quantities (amrex::Real time) |
amrex::Real | volWgtSumMF (int lev, const amrex::MultiFab &mf, int comp, bool local, bool finemask) |
bool | is_it_time_for_action (int nstep, amrex::Real time, amrex::Real dt, int action_interval, amrex::Real action_per) |
virtual void | MakeNewLevelFromCoarse (int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override |
virtual void | RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override |
virtual void | ClearLevel (int lev) override |
virtual void | MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override |
amrex::Real | estTimeStep (int lev) const |
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) |
amrex::MultiFab & | build_fine_mask (int lev) |
void | WritePlotFile (int which, amrex::Vector< std::string > plot_var_names) |
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< std::string > &varnames, amrex::Real time, const amrex::Vector< int > &level_steps, 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 |
void | WriteGenericPlotfileHeaderWithBathymetry (std::ostream &HeaderFile, int nlevels, const amrex::Vector< amrex::BoxArray > &bArray, const amrex::Vector< std::string > &varnames, amrex::Real time, const amrex::Vector< int > &level_steps, const std::string &versionName, const std::string &levelPrefix, const std::string &mfPrefix) const |
void | Advance (int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle) |
void | setup_step (int lev, amrex::Real time, amrex::Real dt_lev) |
void | advance_3d_ml (int lev, amrex::Real dt_lev) |
void | advance_2d_onestep (int lev, amrex::Real dt_lev, amrex::Real dtfast_lev, int my_iif, int nfast_counter) |
void | advance_2d (int lev, amrex::MultiFab const *mf_rhoS, amrex::MultiFab const *mf_rhoA, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, 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, const amrex::Real dtfast_lev, const bool predictor_2d_step, const bool first_2d_step, int my_iif, int &next_indx1) |
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, const int N, const amrex::Real dt_lev) |
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 int iic, const int nfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real dt_lev) |
void | prestep_t_advection (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, int iic, int ntfirst, int nrhs, int N, const amrex::Real dt_lev) |
void | rhs_t_3d (const amrex::Box &bx, const amrex::Box &gbx, 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, int nrhs, int nnew, int N, const amrex::Real dt_lev) |
void | rhs_uv_3d (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) |
void | rhs_uv_2d (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) |
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 const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &h, const int N) |
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 int nrhs, const int N) |
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 > &DC, 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) |
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 > &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) |
void | update_massflux_3d (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 int nnew) |
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 > &dxlen, const int nnew, const int N) |
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, int nrhs, int nnew, const amrex::Real dt_lev) |
void | t3dmix (const amrex::Box &bx, const amrex::Array4< amrex::Real > &state, 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::Real dt_lev, const int ncomp) |
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) |
void | set_2darrays (int lev) |
void | set_bathymetry (int lev) |
void | set_coriolis (int lev) |
void | stretch_transform (int lev) |
void | set_vmix (int lev) |
void | set_smflux (int lev, amrex::Real time) |
void | set_hmixcoef (int lev) |
void | set_drag (int lev) |
void | set_weights (int lev) |
void | FillPatch (int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bdy_var_type=BdyVars::null) |
void | fill_from_bdyfiles (amrex::MultiFab &mf_to_fill, const amrex::Real time, const int bdy_var_type) |
void | init_beta_plane_coriolis (int lev) |
void | writeJobInfo (const std::string &dir) const |
Static Public Member Functions | |
static void | writeBuildInfo (std::ostream &os) |
Public Attributes | |
std::string | pp_prefix {"remora"} |
amrex::Vector< amrex::MultiFab * > | cons_old |
amrex::Vector< amrex::MultiFab * > | xvel_old |
amrex::Vector< amrex::MultiFab * > | yvel_old |
amrex::Vector< amrex::MultiFab * > | zvel_old |
amrex::Vector< amrex::MultiFab * > | cons_new |
amrex::Vector< amrex::MultiFab * > | xvel_new |
amrex::Vector< amrex::MultiFab * > | yvel_new |
amrex::Vector< amrex::MultiFab * > | zvel_new |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_hOfTheConfusingName |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_Hz |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_Huon |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_Hvom |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_ru |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rv |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rufrc |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rvfrc |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_Akv |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_Akt |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_visc2_p |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_visc2_r |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_diff2 |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_x_r |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_y_r |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_z_r |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_z_w |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_s_r |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_z_phys_nd |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_Zt_avg1 |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_sustr |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_svstr |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rdrag |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_bustr |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_bvstr |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_DU_avg1 |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_DU_avg2 |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_DV_avg1 |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_DV_avg2 |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rubar |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rvbar |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rzeta |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_ubar |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_vbar |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_zeta |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_pm |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_pn |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_fcor |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_sstore |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_visc3d_r |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rhoS |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | vec_rhoA |
amrex::Vector< amrex::Real > | vec_weight1 |
amrex::Vector< amrex::Real > | vec_weight2 |
Private Member Functions | |
void | ReadParameters () |
void | AverageDown () |
void | init1DArrays () |
void | init_bcs () |
void | init_custom (int lev) |
void | init_stuff (int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) |
void | resize_stuff (int lev) |
void | AverageDownTo (int crse_lev) |
void | FillCoarsePatch (int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse) |
TimeInterpolatedData | GetDataAtTime (int lev, amrex::Real time) |
void | timeStep (int lev, amrex::Real time, int iteration) |
void | timeStepML (amrex::Real time, int iteration) |
void | initHSE () |
Initialize HSE. More... | |
void | initRayleigh () |
Initialize Rayleigh damping profiles. More... | |
void | ComputeDt () |
std::string | PlotFileName (int lev) const |
amrex::Vector< std::string > | PlotFileVarNames (amrex::Vector< std::string > plot_var_names) const |
void | setPlotVariables (const std::string &pp_plot_var_names, amrex::Vector< std::string > &plot_var_names) |
void | WriteCheckpointFile () const |
void | ReadCheckpointFile () |
void | InitializeFromFile () |
void | InitializeLevelFromData (int lev, const amrex::MultiFab &initial_data) |
void | post_update (amrex::MultiFab &state_mf, const amrex::Real time, const amrex::Geometry &geom) |
void | fill_rhs (amrex::MultiFab &rhs_mf, const amrex::MultiFab &state_mf, const amrex::Real time, const amrex::Geometry &geom) |
void | refinement_criteria_setup () |
AMREX_FORCE_INLINE int | ComputeGhostCells (const int &spatial_order) |
AMREX_FORCE_INLINE amrex::YAFluxRegister * | getAdvFluxReg (int lev) |
AMREX_FORCE_INLINE std::ostream & | DataLog (int i) |
AMREX_FORCE_INLINE int | NumDataLogs () noexcept |
amrex::Real | getCPUTime () const |
void | setRecordDataInfo (int i, const std::string &filename) |
const std::string | DataLogName (int i) const noexcept |
The filename of the ith datalog file. More... | |
Static Private Member Functions | |
static void | GotoNextLine (std::istream &is) |
Private Attributes | |
amrex::Vector< int > | num_boxes_at_level |
amrex::Vector< int > | num_files_at_level |
amrex::Vector< amrex::Vector< amrex::Box > > | boxes_at_level |
amrex::Vector< int > | istep |
amrex::Vector< int > | nsubsteps |
amrex::Vector< amrex::Real > | t_new |
amrex::Vector< amrex::Real > | t_old |
amrex::Vector< amrex::Real > | dt |
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > | physbcs |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | mapfac_m |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | mapfac_u |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | mapfac_v |
amrex::Vector< std::unique_ptr< amrex::MultiFab > > | sst |
amrex::Vector< amrex::YAFluxRegister * > | advflux_reg |
amrex::Vector< amrex::BCRec > | domain_bcs_type |
amrex::Gpu::DeviceVector< amrex::BCRec > | domain_bcs_type_d |
amrex::Array< std::string, 2 *AMREX_SPACEDIM > | domain_bc_type |
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NCONS > | m_bc_extdir_vals |
amrex::GpuArray< REMORA_BC, AMREX_SPACEDIM *2 > | phys_bc_type |
int | last_plot_file_step_1 |
int | last_plot_file_step_2 |
int | last_check_file_step |
int | plot_file_on_restart = 1 |
int | max_step = std::numeric_limits<int>::max() |
amrex::Real | stop_time = std::numeric_limits<amrex::Real>::max() |
std::string | restart_chkfile = "" |
int | nfast |
int | do_substep |
int | regrid_int = 2 |
std::string | plot_file_1 {"plt_1_"} |
std::string | plot_file_2 {"plt_2_"} |
int | plot_int_1 = -1 |
int | plot_int_2 = -1 |
std::string | check_file {"chk"} |
int | check_int = -1 |
amrex::Vector< std::string > | plot_var_names_1 |
amrex::Vector< std::string > | plot_var_names_2 |
const amrex::Vector< std::string > | cons_names {"temp", "salt", "scalar"} |
const amrex::Vector< std::string > | derived_names {"dpdx", "dpdy"} |
amrex::Vector< amrex::Real > | h_havg_density |
amrex::Vector< amrex::Real > | h_havg_temperature |
amrex::Vector< amrex::Real > | h_havg_pressure |
amrex::Gpu::DeviceVector< amrex::Real > | d_havg_density |
amrex::Gpu::DeviceVector< amrex::Real > | d_havg_temperature |
amrex::Gpu::DeviceVector< amrex::Real > | d_havg_pressure |
amrex::MultiFab | fine_mask |
amrex::Vector< std::unique_ptr< std::fstream > > | datalog |
amrex::Vector< std::string > | datalogname |
Static Private Attributes | |
static amrex::Real | cfl = 0.8 |
static amrex::Real | init_shrink = 1.0 |
static amrex::Real | change_max = 1.1 |
static amrex::Real | fixed_dt = -1.0 |
static amrex::Real | fixed_fast_dt = -1.0 |
static int | fixed_ndtfast_ratio = 0 |
static SolverChoice | solverChoice |
static int | verbose = 0 |
static bool | use_tracer_particles |
static int | sum_interval = -1 |
static amrex::Real | sum_per = -1.0 |
static PlotfileType | plotfile_type = PlotfileType::amrex |
static amrex::Vector< amrex::Vector< std::string > > | nc_init_file |
static amrex::Vector< amrex::Vector< std::string > > | nc_grid_file |
static std::string | nc_bdry_file |
static std::string | input_sounding_file |
static int | output_1d_column |
static int | column_interval |
static amrex::Real | column_per |
static amrex::Real | column_loc_x |
static amrex::Real | column_loc_y |
static std::string | column_file_name |
static int | output_bndry_planes |
static int | bndry_output_planes_interval |
static amrex::Real | bndry_output_planes_per |
static amrex::Real | bndry_output_planes_start_time |
static int | input_bndry_planes |
static amrex::Vector< amrex::AMRErrorTag > | ref_tags |
static amrex::Real | startCPUTime = 0.0 |
static amrex::Real | previousCPUTimeUsed = 0.0 |
REMORA::REMORA | ( | ) |
constructor - reads in parameters from inputs file
void REMORA::Advance | ( | int | lev, |
amrex::Real | time, | ||
amrex::Real | dt_lev, | ||
int | iteration, | ||
int | ncycle | ||
) |
advance a single level for a single time step
void REMORA::advance_2d | ( | int | lev, |
amrex::MultiFab const * | mf_rhoS, | ||
amrex::MultiFab const * | mf_rhoA, | ||
amrex::MultiFab * | mf_ru, | ||
amrex::MultiFab * | mf_rv, | ||
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, | ||
const amrex::Real | dtfast_lev, | ||
const bool | predictor_2d_step, | ||
const 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
Function that coordinates the evolution across levels – this calls Advance to do the actual advance at this level, then recursively calls itself at finer levels
[in] | lev | level of refinement (coarsest level is 0) |
[in] | mf_rhoS | |
[in] | mf_rhoA | |
[in] | mf_ru | |
[in] | mf_rv | |
[in,out] | mf_rufrc | |
[in,out] | mf_rvfrc | |
[in,out] | mf_Zt_avg1 | |
[in,out] | mf_DU_avg1 | |
[in,out] | mf_DU_avg2 | |
[in,out] | mf_DV_avg1 | |
[in,out] | mf_DV_avg2 | |
[in,out] | mf_rubar | |
[in,out] | mf_rvbar | |
[in,out] | mf_zeta | |
[in] | mf_h | |
[in,out] | mf_visc2_p | |
[in,out] | mf_visc2_f | |
[in] | dtfast_lev | |
[in] | predictor_2d_step | |
[in] | first_2d_step | |
[in] | my_iif | |
[in] | next_indx1 |
void REMORA::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
void REMORA::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, | ||
const int | N, | ||
const amrex::Real | dt_lev | ||
) |
Advance the 3D variables
void REMORA::advance_3d_ml | ( | int | lev, |
amrex::Real | dt_lev | ||
) |
3D advance on a single level
|
private |
set covered coarse cells to be the average of overlying fine cells
|
private |
more flexible version of AverageDown() that lets you average down across multiple levels
MultiFab & REMORA::build_fine_mask | ( | int | lev | ) |
Make mask to zero out covered cells
|
overridevirtual |
|
private |
a wrapper for estTimeStep()
|
inlineprivate |
void REMORA::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 | ||
) |
|
inlineprivate |
|
inlineprivatenoexcept |
The filename of the ith datalog file.
|
overridevirtual |
Tag cells for refinement
Function to tag cells for refinement – this overrides the pure virtual function in AmrCore
[in] | level | level of refinement (0 is coarsest leve) |
[out] | tags | array of tagged cells |
[in] | time | current time |
Real REMORA::estTimeStep | ( | int | lev | ) | const |
compute dt from CFL considerations
void REMORA::Evolve | ( | ) |
Advance solution to final time
Referenced by main().
void REMORA::fill_from_bdyfiles | ( | amrex::MultiFab & | mf_to_fill, |
const amrex::Real | time, | ||
const int | bdy_var_type | ||
) |
|
private |
|
private |
fill an entire multifab by interpolating from the coarser level this comes into play when a new level of refinement appears
void REMORA::FillPatch | ( | int | lev, |
amrex::Real | time, | ||
amrex::MultiFab & | mf_to_be_filled, | ||
amrex::Vector< amrex::MultiFab * > const & | mfs, | ||
const int | bdy_var_type = BdyVars::null |
||
) |
|
inlineprivate |
|
inlineprivate |
|
private |
utility to copy in data from old and/or new state into another multifab
|
staticprivate |
|
private |
|
private |
void REMORA::init_beta_plane_coriolis | ( | int | lev | ) |
|
private |
void REMORA::init_only | ( | int | lev, |
amrex::Real | time | ||
) |
Init (NOT restart or regrid)
|
private |
void REMORA::InitData | ( | ) |
Initialize multilevel data
Referenced by main().
|
private |
Initialize HSE.
|
private |
|
private |
|
private |
Initialize Rayleigh damping profiles.
bool REMORA::is_it_time_for_action | ( | int | nstep, |
amrex::Real | time, | ||
amrex::Real | dt, | ||
int | action_interval, | ||
amrex::Real | action_per | ||
) |
|
overridevirtual |
Make a new level using provided BoxArray and DistributionMapping and fill with interpolated coarse level data. Overrides the pure virtual function in AmrCore
|
overridevirtual |
Make a new level from scratch using provided BoxArray and DistributionMapping. Only used during initialization. Overrides the pure virtual function in AmrCore
|
inlineprivatenoexcept |
|
private |
get plotfile name
|
private |
void REMORA::post_timestep | ( | int | nstep, |
amrex::Real | time, | ||
amrex::Real | dt_lev | ||
) |
Called after every level 0 timestep
|
private |
void REMORA::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 int | iic, | ||
const int | nfirst, | ||
const int | nnew, | ||
int | nstp, | ||
int | nrhs, | ||
int | N, | ||
const amrex::Real | dt_lev | ||
) |
Wrapper function for prestep
prestep
[in] | mf_vold | |
[in,out] | mf_u | (looks like reset in update_vel <- prestep_uv_3d, so maybe just out |
[in,out] | mf_v | maybe just out |
[in,out] | mf_ru | |
[in,out] | mf_rv | |
[in] | S_old | |
[in,out] | S_new | |
[out] | mf_W | |
[none | ] mf_DC (temp) | |
[in] | mf_z_r | |
[in] | mf_z_w | |
[in] | mf_h | |
[in] | mf_sustr | |
[in] | mf_svstr | |
[in] | mf_bustr | |
[in] | mf_bvstr | |
[in] | iic | |
[in] | ntfirst | |
[in] | nnew | |
[in] | nstp | |
[in] | nrhs | |
[in] | N | |
[in] | dt_lev |
void REMORA::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 > & | DC, | ||
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 REMORA::prestep_t_advection | ( | 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, | ||
int | iic, | ||
int | ntfirst, | ||
int | nrhs, | ||
int | N, | ||
const amrex::Real | dt_lev | ||
) |
Prestep calculations for the tracers
void REMORA::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 int | nrhs, | ||
const int | N | ||
) |
|
private |
|
private |
read in some parameters from inputs file
|
private |
|
overridevirtual |
Remake an existing level using provided BoxArray and DistributionMapping and fill with existing fine and coarse data. Overrides the pure virtual function in AmrCore
void REMORA::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
|
private |
void REMORA::restart | ( | ) |
Restart
void REMORA::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 const > & | Hz, | ||
const amrex::Array4< amrex::Real const > & | z_w, | ||
const amrex::Array4< amrex::Real const > & | h, | ||
const int | N | ||
) |
rho_eos
[in] | bx | box for calculation |
[in] | state | state holds temp, salt |
[out] | rho | density |
[out] | rhoA | vertically-averaged density |
[out] | rhoS | density perturbation |
[in] | Hz | |
[in] | z_w | |
[in] | h | |
[in] | N |
void REMORA::rhs_t_3d | ( | const amrex::Box & | bx, |
const amrex::Box & | gbx, | ||
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, | ||
int | nrhs, | ||
int | nnew, | ||
int | N, | ||
const amrex::Real | dt_lev | ||
) |
rhs_t_3d
[in] | gbx | |
[in,out] | t | |
[in] | sstore | |
[in] | Huon | |
[in] | Hvom | |
[in] | Hz | |
[in] | pn | |
[in] | pm | |
[in] | W | |
[in,out] | FC | |
[in] | nrhs | |
[in] | nnew | |
[in] | N | |
[in] | dt_lev |
void REMORA::rhs_uv_2d | ( | 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_uv_2d
[in] | xbx | Box for operations on x-velocity |
[in] | ybx | Box for operations on y-velocity |
[in] | ubar | |
[in] | vbar | |
[out] | rhs_ubar | |
[out] | rhs_vbar | |
[in] | DUon | |
[in] | DVom | |
[in] | krhs |
void REMORA::rhs_uv_3d | ( | 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 | ||
) |
Calculation of the RHS
rhs_uv_2d
[in] | xbx | Box for operations on x-velocity |
[in] | ybx | Box for operations on y-velocity |
[in] | uold | |
[in] | vold | |
[out] | ru | |
[out] | rv | |
[out] | rufrc | |
[out] | rvfrc | |
[in] | sustr | |
[in] | svstr | |
[in] | bustr | |
[in] | bvstr | |
[in] | Huon | |
[in] | Hvom | |
[in] | pm | |
[in] | pn | |
[in] | W | |
[in,out] | FC | |
[in] | nrhs | |
[in] | N |
void REMORA::set_2darrays | ( | int | lev | ) |
void REMORA::set_bathymetry | ( | int | lev | ) |
void REMORA::set_coriolis | ( | int | lev | ) |
void REMORA::set_drag | ( | int | lev | ) |
void REMORA::set_hmixcoef | ( | int | lev | ) |
void REMORA::set_smflux | ( | int | lev, |
amrex::Real | time | ||
) |
void REMORA::set_vmix | ( | int | lev | ) |
void REMORA::set_weights | ( | int | lev | ) |
|
private |
set which variables and derived quantities go into plotfiles
|
inlineprivate |
void REMORA::setup_step | ( | int | lev, |
amrex::Real | time, | ||
amrex::Real | dt_lev | ||
) |
Set everything up for a step on a level
void REMORA::stretch_transform | ( | int | lev | ) |
void REMORA::sum_integrated_quantities | ( | amrex::Real | time | ) |
void REMORA::t3dmix | ( | const amrex::Box & | bx, |
const amrex::Array4< amrex::Real > & | state, | ||
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::Real | dt_lev, | ||
const int | ncomp | ||
) |
Harmonic diffusivity for tracers
|
private |
advance a level by dt, includes a recursive call for finer levels
|
private |
advance all levels by dt, loops over finer levels
void REMORA::update_massflux_3d | ( | 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 int | nnew | ||
) |
update_massflux_3d
[in] | bx | box on which to update |
[in] | ioff | offset in x-direction |
[in] | joff | offset in y-direction |
[in] | phi | u or v |
[out] | Hphi | H-weighted u or v |
[in] | Hz | weighting |
[in] | om_v_or_on_u | |
[in] | Dphi_avg1 | |
[in] | Dphi_avg2 | |
[in,out] | DC | |
[in,out] | FC | |
[in] | nnew | component of velocity |
void REMORA::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, | ||
int | nrhs, | ||
int | nnew, | ||
const amrex::Real | dt_lev | ||
) |
Harmonic viscosity
void REMORA::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 > & | dxlen, | ||
const int | nnew, | ||
const int | N | ||
) |
void REMORA::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 > & | 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 | ||
) |
Real REMORA::volWgtSumMF | ( | int | lev, |
const amrex::MultiFab & | mf, | ||
int | comp, | ||
bool | local, | ||
bool | finemask | ||
) |
|
static |
|
private |
void REMORA::WriteGenericPlotfileHeaderWithBathymetry | ( | std::ostream & | HeaderFile, |
int | nlevels, | ||
const amrex::Vector< amrex::BoxArray > & | bArray, | ||
const amrex::Vector< std::string > & | varnames, | ||
amrex::Real | time, | ||
const amrex::Vector< int > & | level_steps, | ||
const std::string & | versionName, | ||
const std::string & | levelPrefix, | ||
const std::string & | mfPrefix | ||
) | const |
void REMORA::writeJobInfo | ( | const std::string & | dir | ) | const |
void REMORA::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< std::string > & | varnames, | ||
amrex::Real | time, | ||
const amrex::Vector< int > & | level_steps, | ||
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 |
void REMORA::WritePlotFile | ( | int | which, |
amrex::Vector< std::string > | plot_var_names | ||
) |
write plotfile to disk
|
private |
Referenced by getAdvFluxReg().
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
private |
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
private |
amrex::Vector<amrex::MultiFab*> REMORA::cons_new |
amrex::Vector<amrex::MultiFab*> REMORA::cons_old |
|
private |
|
private |
|
private |
|
private |
Referenced by DataLog(), NumDataLogs(), and setRecordDataInfo().
|
private |
Referenced by DataLogName().
|
private |
|
private |
|
private |
Referenced by writeJobInfo().
|
private |
|
private |
|
private |
|
private |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
private |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
private |
|
private |
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
staticprivate |
std::string REMORA::pp_prefix {"remora"} |
|
staticprivate |
Referenced by getCPUTime().
|
staticprivate |
|
private |
|
private |
|
staticprivate |
|
private |
|
staticprivate |
Referenced by getCPUTime().
|
private |
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
staticprivate |
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_Akt |
Vertical diffusion coefficient (3D), set in .in file
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_Akv |
Vertical viscosity coefficient (3D), set in .in file
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_bustr |
Bottom stress in the u direction
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_bvstr |
Bottom stress in the v direction
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_diff2 |
Harmonic diffusivity for temperature / salinity
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_DU_avg1 |
time average of barotropic x velocity flux (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_DU_avg2 |
correct time average of barotropic x velocity flux for coupling (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_DV_avg1 |
time average of barotropic y velocity flux
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_DV_avg2 |
correct time average of barotropic y velocity flux for coupling (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_fcor |
coriolis factor (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_hOfTheConfusingName |
Bathymetry data (2D, positive valued, h in ROMS)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_Huon |
u-volume flux (3D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_Hvom |
v-volume flux (3D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_Hz |
Width of cells in the vertical (z-) direction (3D, Hz in ROMS)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_pm |
map factor (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_pn |
map factor (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rdrag |
Linear drag coefficient [m/s], defined at rho points
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rhoA |
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rhoS |
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_ru |
u velocity RHS (3D, includes horizontal and vertical advection)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rubar |
barotropic x velocity for the RHS (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rufrc |
u velocity RHS, integrated, including advection and bottom/surface stresses (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rv |
v velocity RHS (3D, includes horizontal and vertical advection)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rvbar |
barotropic y velocity for the RHS (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rvfrc |
v velocity RHS, integrated, including advection and bottom/surface stresses (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_rzeta |
free surface height for the RHS (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_s_r |
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at rho points (cell centers)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_sstore |
saltstore, tempstore, etc
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_sustr |
Surface stress in the u direction
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_svstr |
Surface stress in the v direction
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_ubar |
barotropic x velocity (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_vbar |
barotropic y velocity (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_visc2_p |
Harmonic viscosity devined on the psi points (corners of horizontal grid cells)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_visc2_r |
Harmonic viscosity devined on the rho points (centers)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_visc3d_r |
amrex::Vector<amrex::Real> REMORA::vec_weight1 |
amrex::Vector<amrex::Real> REMORA::vec_weight2 |
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_x_r |
x coordinates at rho points (cell centers)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_y_r |
y coordinates at rho points (cell centers)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_z_phys_nd |
z coordinates at psi points (cell nodes)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_z_r |
z coordinates at rho points (cell centers)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_z_w |
z coordinates at w points (faces between z-cells)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_zeta |
free surface height (2D)
amrex::Vector<std::unique_ptr<amrex::MultiFab> > REMORA::vec_Zt_avg1 |
Average of the free surface, zeta (2D)
|
staticprivate |
amrex::Vector<amrex::MultiFab*> REMORA::xvel_new |
amrex::Vector<amrex::MultiFab*> REMORA::xvel_old |
amrex::Vector<amrex::MultiFab*> REMORA::yvel_new |
amrex::Vector<amrex::MultiFab*> REMORA::yvel_old |
amrex::Vector<amrex::MultiFab*> REMORA::zvel_new |
amrex::Vector<amrex::MultiFab*> REMORA::zvel_old |