REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_Advance.cpp
Go to the documentation of this file.
1#include <REMORA.H>
2
3using namespace amrex;
4
5#ifdef REMORA_USE_FUNWAVE_FORT
7#endif
8
9/**
10 * @param[in] lev level of refinement
11 * @param[in] time simulation time at start of step
12 * @param[in] dt_lev baroclinic time step at level
13 * @param[in] iteration iteration in subcycling, if using
14 * @param[in] ncycle total number of subcycles, if using
15 */
16 void
17REMORA::Advance (int lev, Real time, Real dt_lev, int /*iteration*/, int /*ncycle*/)
18{
19 BL_PROFILE("REMORA::Advance()");
20
21 setup_step(lev, time, dt_lev);
22
24 {
25 int nfast_counter=nfast + 1;
26
27 //***************************************************
28 //Compute fast timestep from dt_lev and ratio
29 //***************************************************
30 Real dtfast_lev=dt_lev/Real(fixed_ndtfast_ratio);
31
32 //***************************************************
33 //Advance nfast_counter steps of the 2d integrator
34 //***************************************************
35 for (int my_iif = 0; my_iif < nfast_counter; my_iif++) {
36 advance_2d_onestep(lev, dt_lev, dtfast_lev, my_iif, nfast_counter);
37 }
38 }
39
40#ifdef REMORA_USE_FUNWAVE_FORT
41 MultiFab* mf_rhoS = vec_rhoS[lev].get();
42 for ( MFIter mfi(*mf_rhoS, TilingIfNotGPU()); mfi.isValid(); ++mfi )
43 {
44 Box bx = mfi.validbox();
45 int ims = bx.smallEnd(0);
46 int jms = bx.smallEnd(1);
47 int kms = bx.smallEnd(2);
48 int ime = bx.bigEnd(0);
49 int jme = bx.bigEnd(1);
50 int kme = bx.bigEnd(2);
51
52 Array4<Real> const& rho_salt = mf_rhoS->array(mfi);
53
54 funwave_advance_c(rho_salt.dataPtr(), ims, ime, jms, jme, kms, kme);
55 }
56#endif
57
58 //***************************************************
59 //Advance one step of the 3d integrator
60 //***************************************************
61 advance_3d_ml(lev, dt_lev);
62
63}
void funwave_advance_c(double *salt, int ims, int ime, int jms, int jme, int kms, int kme)
int nfast
Number of fast steps to take.
Definition REMORA.H:1180
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
void advance_3d_ml(int lev, amrex::Real dt_lev)
3D advance on a single level
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoS
density perturbation
Definition REMORA.H:386
void Advance(int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle)
advance a single level for a single time step
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1221
static int fixed_ndtfast_ratio
User specified, number of barotropic steps per baroclinic step.
Definition REMORA.H:1178
void setup_step(int lev, amrex::Real time, amrex::Real dt_lev)
Set everything up for a step on a level.