REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_TimeStep.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] time simulation time at start of step
8 * @param[in] iteration iteration in subcycling, if using
9 */
10void
11REMORA::timeStep (int lev, Real time, int iteration)
12{
13 if (regrid_int > 0) // We may need to regrid
14 {
15 // help keep track of whether a level was already regridded
16 // from a coarser level call to regrid
17 static Vector<int> last_regrid_step(max_level+1, 0);
18
19 // regrid changes level "lev+1" so we don't regrid on max_level
20 // also make sure we don't regrid fine levels again if
21 // it was taken care of during a coarser regrid
22 if (lev < max_level && istep[lev] > last_regrid_step[lev])
23 {
24 if (istep[lev] % regrid_int == 0)
25 {
26 // regrid could add newly refine levels (if finest_level < max_level)
27 // so we save the previous finest level index
28 int old_finest = finest_level;
29 regrid(lev, time);
30
31#ifdef REMORA_USE_PARTICLES
32 if (finest_level != old_finest) {
33 particleData.Redistribute();
34 }
35#endif
36
37 // Mark that we have regridded this level already
38 for (int k = lev; k <= finest_level; ++k) {
39 last_regrid_step[k] = istep[k];
40 }
41
42 // If there are newly created levels, set the time step
43 for (int k = old_finest+1; k <= finest_level; ++k) {
44 dt[k] = dt[k-1] / MaxRefRatio(k-1);
45 }
46 }
47 }
48 }
49 // Update what we call "old" and "new" time
50 t_old[lev] = t_new[lev];
51 t_new[lev] += dt[lev];
52
53 if (Verbose()) {
54 amrex::Print() << "[Level " << lev << " step " << istep[lev]+1 << "] ";
55 amrex::Print() << "ADVANCE from time = " << t_old[lev] << " to " << t_new[lev]
56 << " with dt = " << dt[lev] << std::endl;
57 }
58
59 // We must swap the pointers so the previous step's "new" is now this step's "old"
60 std::swap(cons_old[lev], cons_new[lev]);
61 std::swap(xvel_old[lev], xvel_new[lev]);
62 std::swap(yvel_old[lev], yvel_new[lev]);
63 std::swap(zvel_old[lev], zvel_new[lev]);
64
65 // Advance a single level for a single time step
66 Advance(lev, time, dt[lev], iteration, nsubsteps[lev]);
67
68 ++istep[lev];
69
70 if (Verbose())
71 {
72 amrex::Print() << "[Level " << lev << " step " << istep[lev] << "] ";
73 amrex::Print() << "Advanced " << CountCells(lev) << " cells" << std::endl;
74 }
75
76 if (lev < finest_level)
77 {
78 // recursive call for next-finer level
79 for (int i = 1; i <= nsubsteps[lev+1]; ++i)
80 {
81 timeStep(lev+1, time+(i-1)*dt[lev+1], i);
82 }
83
85 AverageDownTo(lev); // average lev+1 down to lev
86 }
87 }
88}
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:217
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:223
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:210
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:221
int regrid_int
how often each level regrids the higher levels of refinement (after a level advances that many time s...
Definition REMORA.H:1186
void AverageDownTo(int crse_lev)
more flexible version of AverageDown() that lets you average down across multiple levels
Definition REMORA.cpp:1069
amrex::Vector< amrex::MultiFab * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
Definition REMORA.H:214
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:219
amrex::Vector< int > nsubsteps
How many substeps on each level?
Definition REMORA.H:1096
amrex::Vector< int > istep
which step?
Definition REMORA.H:1094
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:212
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< amrex::Real > t_new
new time at each level
Definition REMORA.H:1098
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1221
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:208
void timeStep(int lev, amrex::Real time, int iteration)
advance a level by dt, includes a recursive call for finer levels
amrex::Vector< amrex::Real > t_old
old time at each level
Definition REMORA.H:1100
amrex::Vector< amrex::Real > dt
time step at each level
Definition REMORA.H:1102
CouplingType coupling_type