REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_TimeStepML.cpp
Go to the documentation of this file.
1#include <REMORA.H>
2
3using namespace amrex;
4
5/**
6 * @param[in] time simulation time at start of step
7 * @param[in] iteration iteration in subcycling, if using
8 */
9void
10REMORA::timeStepML (Real time, int /*iteration*/)
11{
12 if (time == 0.0_rt && solverChoice.init_l1ad_T) {
13 average_down(*cons_new[1], *cons_new[0],
14 0, cons_new[0]->nComp(), refRatio(0));
16 }
17
18 // HACK HACK so lev is defined and compiler won't complain, but always say regrid_int=-1
19 for (int lev=0; lev <= finest_level;lev++) {
20 if (regrid_int > 0) // We may need to regrid
21 {
22 // help keep track of whether a level was already regridded
23 // from a coarser level call to regrid
24 static Vector<int> last_regrid_step(max_level+1, 0);
25
26 // regrid changes level "lev+1" so we don't regrid on max_level
27 // also make sure we don't regrid fine levels again if
28 // it was taken care of during a coarser regrid
29 if (lev < max_level && istep[lev] > last_regrid_step[lev])
30 {
31 if (istep[lev] % regrid_int == 0)
32 {
33 // regrid could add newly refine levels (if finest_level < max_level)
34 // so we save the previous finest level index
35 int old_finest = finest_level;
36 regrid(lev, time);
37
38 // Mark that we have regridded this level already
39 for (int k = lev; k <= finest_level; ++k) {
40 last_regrid_step[k] = istep[k];
41 }
42
43 // If there are newly created levels, set the time step
44 for (int k = old_finest+1; k <= finest_level; ++k) {
45 dt[k] = dt[k-1] / nsubsteps[lev];
46 }
47 }
48 }
49 }
50 }
51
53
54 for (int lev=0; lev <= finest_level;lev++)
55 {
56 // Update what we call "old" and "new" time
57 t_old[lev] = t_new[lev];
58 t_new[lev] += dt[lev];
59
60 if (Verbose()) {
61 amrex::Print() << "[Level " << lev << " step " << istep[lev]+1 << "] ";
62 amrex::Print() << "ADVANCE from time = " << t_old[lev] << " to " << t_new[lev]
63 << " with dt = " << dt[lev] << std::endl;
64 }
65
66 // We must swap the pointers so the previous step's "new" is now this step's "old"
67 std::swap(cons_old[lev], cons_new[lev]);
68 std::swap(xvel_old[lev], xvel_new[lev]);
69 std::swap(yvel_old[lev], yvel_new[lev]);
70 std::swap(zvel_old[lev], zvel_new[lev]);
71
72 setup_step(lev, time, dt[lev]);
73
74 // **************************************************************************************
75 // Register old and new coarse data if we are at a level less than the finest level
76 // **************************************************************************************
77 if (lev < finest_level)
78 {
79 if (cf_width > 0) {
80 // We must fill the ghost cells of these so that the parallel copy works correctly
81 cons_old[lev]->FillBoundary(geom[lev].periodicity());
82 cons_new[lev]->FillBoundary(geom[lev].periodicity());
83 FPr_c[lev].RegisterCoarseData({cons_old[lev], cons_new[lev]}, {time, time + dt[lev]});
84 }
85
86 if (cf_width >= 0) {
87 // We must fill the ghost cells of these so that the parallel copy works correctly
88 xvel_old[lev]->FillBoundary(geom[lev].periodicity());
89 xvel_new[lev]->FillBoundary(geom[lev].periodicity());
90 FPr_u[lev].RegisterCoarseData({xvel_old[lev], xvel_new[lev]}, {time, time + dt[lev]});
91
92 yvel_old[lev]->FillBoundary(geom[lev].periodicity());
93 yvel_new[lev]->FillBoundary(geom[lev].periodicity());
94 FPr_v[lev].RegisterCoarseData({yvel_old[lev], yvel_new[lev]}, {time, time + dt[lev]});
95
96 zvel_old[lev]->FillBoundary(geom[lev].periodicity());
97 zvel_new[lev]->FillBoundary(geom[lev].periodicity());
98 FPr_w[lev].RegisterCoarseData({zvel_old[lev], zvel_new[lev]}, {time, time + dt[lev]});
99 }
100 }
101 }
102
104 {
105 int nfast_counter=nfast + 1;
106
107 for (int my_iif = 0; my_iif < nfast_counter; my_iif++) {
108 //Compute fast timestep from dt[lev] and ratio
109 for (int lev=0; lev <= finest_level; lev++)
110 {
111 Real dtfast_lev=dt[lev]/Real(fixed_ndtfast_ratio);
112 advance_2d_onestep(lev, dt[lev], dtfast_lev, my_iif, nfast_counter);
113 // **************************************************************************************
114 // Register old and new coarse data if we are at a level less than the finest level
115 // **************************************************************************************
116 if (lev < finest_level)
117 {
118 if (cf_width >= 0) {
119 // We must fill the ghost cells of these so that the parallel copy works correctly
120 vec_ubar[lev]->FillBoundary(geom[lev].periodicity());
121 FPr_ubar[lev].RegisterCoarseData({vec_ubar[lev].get(), vec_ubar[lev].get()}, {time, time + dt[lev]});
122
123 vec_vbar[lev]->FillBoundary(geom[lev].periodicity());
124 FPr_vbar[lev].RegisterCoarseData({vec_vbar[lev].get(), vec_vbar[lev].get()}, {time, time + dt[lev]});
125 }
126 }
127 } // my_iif
128 } // lev
129 } // use_barotropic
130
131 for (int lev=0; lev <= finest_level; lev++) {
132 advance_3d_ml(lev, dt[lev]);
133 ++istep[lev];
134
135 if (Verbose())
136 {
137 amrex::Print() << "[Level " << lev << " step " << istep[lev] << "] ";
138 amrex::Print() << "Advanced " << CountCells(lev) << " cells" << std::endl;
139 }
140 // **************************************************************************************
141 // Register old and new coarse data if we are at a level less than the finest level
142 // **************************************************************************************
143 if (lev < finest_level)
144 {
145 if (cf_width > 0) {
146 // We must fill the ghost cells of these so that the parallel copy works correctly
147 cons_old[lev]->FillBoundary(geom[lev].periodicity());
148 cons_new[lev]->FillBoundary(geom[lev].periodicity());
149 FPr_c[lev].RegisterCoarseData({cons_old[lev], cons_new[lev]}, {time, time + dt[lev]});
150 }
151
152 if (cf_width >= 0) {
153 // We must fill the ghost cells of these so that the parallel copy works correctly
154 xvel_old[lev]->FillBoundary(geom[lev].periodicity());
155 xvel_new[lev]->FillBoundary(geom[lev].periodicity());
156 FPr_u[lev].RegisterCoarseData({xvel_old[lev], xvel_new[lev]}, {time, time + dt[lev]});
157
158 yvel_old[lev]->FillBoundary(geom[lev].periodicity());
159 yvel_new[lev]->FillBoundary(geom[lev].periodicity());
160 FPr_v[lev].RegisterCoarseData({yvel_old[lev], yvel_new[lev]}, {time, time + dt[lev]});
161
162 zvel_old[lev]->FillBoundary(geom[lev].periodicity());
163 zvel_new[lev]->FillBoundary(geom[lev].periodicity());
164 FPr_w[lev].RegisterCoarseData({zvel_old[lev], zvel_new[lev]}, {time, time + dt[lev]});
165 }
166 }
167 FillPatch(lev, t_new[lev], *xvel_new[lev], xvel_new, BCVars::xvel_bc, BdyVars::u,0,true,true);
168 FillPatch(lev, t_new[lev], *yvel_new[lev], yvel_new, BCVars::yvel_bc, BdyVars::v,0,true,true);
169 FillPatch(lev, t_new[lev], *zvel_new[lev], zvel_new, BCVars::zvel_bc, BdyVars::null,0,true,true);
170 }
171
173
175 for (int lev=0; lev <= finest_level-1; lev++) {
176 AverageDownTo(lev); // average lev+1 down to lev
177 }
178 }
179}
int nfast
Number of fast steps to take.
Definition REMORA.H:1180
void scale_rhs_vars()
Scale RHS momentum variables by 1/cell area, needed before FillPatch to different levels.
amrex::Vector< REMORAFillPatcher > FPr_v
Vector over levels of FillPatchers for v (3D)
Definition REMORA.H:1055
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< REMORAFillPatcher > FPr_u
Vector over levels of FillPatchers for u (3D)
Definition REMORA.H:1053
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 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)
Definition REMORA.H:1063
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
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::Vector< int > istep
which step?
Definition REMORA.H:1094
void advance_3d_ml(int lev, amrex::Real dt_lev)
3D advance on a single level
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:212
amrex::Vector< REMORAFillPatcher > FPr_c
Vector over levels of FillPatchers for scalars.
Definition REMORA.H:1051
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.
Definition REMORA.H:1057
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< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:341
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:339
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:208
amrex::Vector< REMORAFillPatcher > FPr_ubar
Vector over levels of FillPatchers for ubar (2D)
Definition REMORA.H:1061
static int fixed_ndtfast_ratio
User specified, number of barotropic steps per baroclinic step.
Definition REMORA.H:1178
void timeStepML(amrex::Real time, int iteration)
advance all levels by dt, loops over finer levels
void setup_step(int lev, amrex::Real time, amrex::Real dt_lev)
Set everything up for a step on a level.
void WritePlotFile()
main driver for writing AMReX plotfiles
int cf_width
Nudging width at coarse-fine interface.
Definition REMORA.H:1046
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