REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_advance_2d_onestep.cpp
Go to the documentation of this file.
1#include <REMORA.H>
2
3using namespace amrex;
4
5/**
6 * @param[in] lev level of refinement
7 * @param[in] dt_lev baroclinic time step at level
8 * @param[in] dtfast_lev barotropic time step at level
9 * @param[in] my_iif current barotropic time step
10 * @param[in] nfast_counter how many bartropicc time steps to do per baroclinic step
11 */
12void REMORA::advance_2d_onestep (int lev, Real /*dt_lev*/, Real dtfast_lev, int my_iif, int nfast_counter)
13{
14 bool first_2d_step=(my_iif==0);
15
16 //Predictor
17 bool predictor_2d_step=true;
18 int next_indx1 = 0;
19 advance_2d(lev,
20 vec_rhoS[lev].get() , vec_rhoA[lev].get(),
21 vec_ru2d[lev].get() , vec_rv2d[lev].get(),
22 vec_rufrc[lev].get(), vec_rvfrc[lev].get(),
23 vec_Zt_avg1[lev].get(),
24 vec_DU_avg1[lev], vec_DU_avg2[lev],
25 vec_DV_avg1[lev], vec_DV_avg2[lev],
26 vec_rubar[lev], vec_rvbar[lev], vec_rzeta[lev],
27 vec_ubar[lev], vec_vbar[lev], vec_zeta[lev].get(),
28 vec_hOfTheConfusingName[lev].get(),
29 vec_pm[lev].get(), vec_pn[lev].get(),
30 vec_fcor[lev].get(),
31 vec_visc2_p[lev].get(), vec_visc2_r[lev].get(),
32 vec_mskr[lev].get(), vec_msku[lev].get(), vec_mskv[lev].get(),
33 vec_mskp[lev].get(),
34 dtfast_lev, predictor_2d_step, first_2d_step, my_iif, next_indx1);
35
36 //Corrector. Skip it on last fast step
37 predictor_2d_step=false;
38 if (my_iif < nfast_counter - 1) {
39
40 // These are needed to pass ctests
41
42 advance_2d(lev,
43 vec_rhoS[lev].get(), vec_rhoA[lev].get(),
44 vec_ru2d[lev].get(), vec_rv2d[lev].get(),
45 vec_rufrc[lev].get(), vec_rvfrc[lev].get(),
46 vec_Zt_avg1[lev].get(),
47 vec_DU_avg1[lev], vec_DU_avg2[lev],
48 vec_DV_avg1[lev], vec_DV_avg2[lev],
49 vec_rubar[lev], vec_rvbar[lev], vec_rzeta[lev],
50 vec_ubar[lev], vec_vbar[lev], vec_zeta[lev].get(),
51 vec_hOfTheConfusingName[lev].get(),
52 vec_pm[lev].get(), vec_pn[lev].get(),
53 vec_fcor[lev].get(),
54 vec_visc2_p[lev].get(), vec_visc2_r[lev].get(),
55 vec_mskr[lev].get(), vec_msku[lev].get(), vec_mskv[lev].get(),
56 vec_mskp[lev].get(),
57 dtfast_lev, predictor_2d_step, first_2d_step, my_iif, next_indx1);
58 }
59}
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:244
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:360
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_hOfTheConfusingName
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:229
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rubar
barotropic x velocity for the RHS (2D)
Definition REMORA.H:333
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:355
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg2
correct time average of barotropic x velocity flux for coupling (2D)
Definition REMORA.H:327
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:346
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:242
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_msku
land/sea mask at x-faces (2D)
Definition REMORA.H:348
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvfrc
v velocity RHS, integrated, including advection and bottom/surface stresses (2D)
Definition REMORA.H:248
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rufrc
u velocity RHS, integrated, including advection and bottom/surface stresses (2D)
Definition REMORA.H:246
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_p
Harmonic viscosity defined on the psi points (corners of horizontal grid cells)
Definition REMORA.H:254
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvbar
barotropic y velocity for the RHS (2D)
Definition REMORA.H:335
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< std::unique_ptr< amrex::MultiFab > > vec_mskp
land/sea mask at cell corners (2D)
Definition REMORA.H:352
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskv
land/sea mask at y-faces (2D)
Definition REMORA.H:350
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoS
density perturbation
Definition REMORA.H:386
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_r
Harmonic viscosity defined on the rho points (centers)
Definition REMORA.H:256
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoA
vertically-averaged density
Definition REMORA.H:388
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg1
time average of barotropic y velocity flux
Definition REMORA.H:329
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
Definition REMORA.H:343
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:341
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg1
time average of barotropic x velocity flux (2D)
Definition REMORA.H:325
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.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:339
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg2
correct time average of barotropic y velocity flux for coupling (2D)
Definition REMORA.H:331
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rzeta
free surface height for the RHS (2D)
Definition REMORA.H:337
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:357
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:279