REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_make_new_level.cpp
Go to the documentation of this file.
1/**
2 * \file REMORA_make_new_level.cpp
3 */
4
5#include <REMORA.H>
7
8#include <AMReX_buildInfo.H>
9
10using namespace amrex;
11
12/**
13 * Make a new level using provided BoxArray and DistributionMapping and
14 * fill with interpolated coarse level data (overrides the pure virtual function in AmrCore)
15 * regrid --> RemakeLevel (if level already existed)
16 * regrid --> MakeNewLevelFromCoarse (if adding new level)
17 *
18 * @param[in ] lev level to make
19 * @param[in ] time current time
20 * @param[in ] ba BoxArray for the level
21 * @param[in ] dm DistributionMapping for the level
22 */
23void
24REMORA::MakeNewLevelFromCoarse (int lev, Real time, const BoxArray& ba,
25 const DistributionMapping& dm)
26{
27 BoxList bl2d = ba.boxList();
28 for (auto& b : bl2d) {
29 b.setRange(2,0);
30 }
31 BoxArray ba2d(std::move(bl2d));
32
33 cons_new[lev] = new MultiFab(ba, dm, NCONS, cons_new[lev-1]->nGrowVect());
34 cons_old[lev] = new MultiFab(ba, dm, NCONS, cons_new[lev-1]->nGrowVect());
35
36 xvel_new[lev] = new MultiFab(convert(ba, IntVect(1,0,0)), dm, 1, xvel_new[lev-1]->nGrowVect());
37 xvel_old[lev] = new MultiFab(convert(ba, IntVect(1,0,0)), dm, 1, xvel_new[lev-1]->nGrowVect());
38
39 yvel_new[lev] = new MultiFab(convert(ba, IntVect(0,1,0)), dm, 1, yvel_new[lev-1]->nGrowVect());
40 yvel_old[lev] = new MultiFab(convert(ba, IntVect(0,1,0)), dm, 1, yvel_new[lev-1]->nGrowVect());
41
42 zvel_new[lev] = new MultiFab(convert(ba, IntVect(0,0,1)), dm, 1, zvel_new[lev-1]->nGrowVect());
43 zvel_old[lev] = new MultiFab(convert(ba, IntVect(0,0,1)), dm, 1, zvel_new[lev-1]->nGrowVect());
44
45 resize_stuff(lev);
46
47 vec_Zt_avg1[lev].reset(new MultiFab(ba2d ,dm,1,IntVect(NGROW+1,NGROW+1,0))); //2d, average of the free surface (zeta)
48 vec_hOfTheConfusingName[lev].reset(new MultiFab(ba2d ,dm,2,IntVect(NGROW+1,NGROW+1,0))); //2d, average of the free surface (zeta)
49 vec_ubar[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,3,IntVect(NGROW,NGROW,0)));
50 vec_vbar[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,3,IntVect(NGROW,NGROW,0)));
51
52 vec_ru[lev].reset(new MultiFab(convert(ba,IntVect(1,0,0)),dm,2,IntVect(NGROW,NGROW,0))); // RHS u (incl horizontal and vertical advection)
53 vec_rv[lev].reset(new MultiFab(convert(ba,IntVect(0,1,0)),dm,2,IntVect(NGROW,NGROW,0))); // RHS v
54
55 vec_ru2d[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,2,IntVect(NGROW,NGROW,0))); // RHS u for 2d
56 vec_rv2d[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,2,IntVect(NGROW,NGROW,0))); // RHS v for 2d
57
58 t_new[lev] = time;
59 t_old[lev] = time - 1.e200_rt;
60
61 init_masks(lev, ba, dm);
62
63 FillCoarsePatch(lev, time, cons_new[lev], cons_new[lev-1]);
64 FillCoarsePatch(lev, time, xvel_new[lev], xvel_new[lev-1]);
65 FillCoarsePatch(lev, time, yvel_new[lev], yvel_new[lev-1]);
66 FillCoarsePatch(lev, time, zvel_new[lev], zvel_new[lev-1]);
67
68 init_stuff(lev, ba, dm);
69
70 FillCoarsePatch(lev, time, vec_hOfTheConfusingName[lev].get(), vec_hOfTheConfusingName[lev-1].get());
71 FillCoarsePatch(lev, time, vec_Zt_avg1[lev].get(), vec_Zt_avg1[lev-1].get());
72 for (int icomp=0; icomp<3; icomp++) {
73 FillCoarsePatch(lev, time, vec_ubar[lev].get(), vec_ubar[lev-1].get(),icomp,false);
74 FillCoarsePatch(lev, time, vec_vbar[lev].get(), vec_vbar[lev-1].get(),icomp,false);
75 }
76 for (int icomp=0; icomp<2; icomp++) {
77 FillCoarsePatch(lev, time, vec_ru[lev].get(), vec_ru[lev-1].get(),icomp,false);
78 FillCoarsePatch(lev, time, vec_rv[lev].get(), vec_rv[lev-1].get(),icomp,false);
79 FillCoarsePatch(lev, time, vec_ru2d[lev].get(), vec_ru2d[lev-1].get(),icomp,false);
80 FillCoarsePatch(lev, time, vec_rv2d[lev].get(), vec_rv2d[lev-1].get(),icomp,false);
81 }
82
83
84 set_grid_scale(lev);
86
87 init_set_vmix(lev);
88 set_hmixcoef(lev);
89 set_coriolis(lev);
91 // Previously set smflux
92
93 // ********************************************************************************************
94 // If we are making a new level then the FillPatcher for this level hasn't been allocated yet
95 // ********************************************************************************************
96 if (cf_width >= 0) {
99 }
100
101#ifdef REMORA_USE_PARTICLES
102 // particleData.Redistribute();
103#endif
104}
105
106/**
107 * Remake an existing level using provided BoxArray and DistributionMapping and
108 * fill with existing fine and coarse data.
109 * overrides the pure virtual function in AmrCore
110 * @param[in ] lev level to make
111 * @param[in ] time current time
112 * @param[in ] ba BoxArray for the level
113 * @param[in ] dm DistributionMapping for the level
114 */
115void
116REMORA::RemakeLevel (int lev, Real time, const BoxArray& ba, const DistributionMapping& dm)
117{
118 BoxArray ba_old(cons_new[lev]->boxArray());
119 DistributionMapping dm_old(cons_new[lev]->DistributionMap());
120
121 BoxList bl2d = ba.boxList();
122 for (auto& b : bl2d) {
123 b.setRange(2,0);
124 }
125 BoxArray ba2d(std::move(bl2d));
126
127#if (NGROW==2)
128 int ngrow_state = ComputeGhostCells(solverChoice.spatial_order)+1;
129 int ngrow_vels = ComputeGhostCells(solverChoice.spatial_order)+1;
130 int ngrow_zeta = ComputeGhostCells(solverChoice.spatial_order)+1;
132 int ngrow_velbar = ComputeGhostCells(solverChoice.spatial_order);
133#else
134 int ngrow_state = ComputeGhostCells(solverChoice.spatial_order)+2;
135 int ngrow_vels = ComputeGhostCells(solverChoice.spatial_order)+2;
136 int ngrow_zeta = ComputeGhostCells(solverChoice.spatial_order)+2;
138 int ngrow_velbar = ComputeGhostCells(solverChoice.spatial_order)+1;
139#endif
140
141 MultiFab tmp_cons_new(ba, dm, NCONS, ngrow_state);
142 MultiFab tmp_cons_old(ba, dm, NCONS, ngrow_state);
143
144 MultiFab tmp_xvel_new(convert(ba, IntVect(1,0,0)), dm, 1, ngrow_vels);
145 MultiFab tmp_xvel_old(convert(ba, IntVect(1,0,0)), dm, 1, ngrow_vels);
146
147 MultiFab tmp_yvel_new(convert(ba, IntVect(0,1,0)), dm, 1, ngrow_vels);
148 MultiFab tmp_yvel_old(convert(ba, IntVect(0,1,0)), dm, 1, ngrow_vels);
149
150 MultiFab tmp_zvel_new(convert(ba, IntVect(0,0,1)), dm, 1, IntVect(ngrow_vels,ngrow_vels,0));
151 MultiFab tmp_zvel_old(convert(ba, IntVect(0,0,1)), dm, 1, IntVect(ngrow_vels,ngrow_vels,0));
152
153 MultiFab tmp_Zt_avg1_new(ba2d, dm, 1, IntVect(ngrow_zeta,ngrow_zeta,0));
154 MultiFab tmp_Zt_avg1_old(ba2d, dm, 1, IntVect(ngrow_zeta,ngrow_zeta,0));
155 MultiFab tmp_h(ba2d, dm, 2, IntVect(ngrow_h,ngrow_h,0));
156
157 MultiFab tmp_ubar_new(convert(ba2d, IntVect(1,0,0)), dm, 3, IntVect(ngrow_velbar,ngrow_velbar,0));
158 MultiFab tmp_ubar_old(convert(ba2d, IntVect(1,0,0)), dm, 3, IntVect(ngrow_velbar,ngrow_velbar,0));
159
160 MultiFab tmp_vbar_new(convert(ba2d, IntVect(0,1,0)), dm, 3, IntVect(ngrow_velbar,ngrow_velbar,0));
161 MultiFab tmp_vbar_old(convert(ba2d, IntVect(0,1,0)), dm, 3, IntVect(ngrow_velbar,ngrow_velbar,0));
162
163 MultiFab tmp_ru_new(convert(ba, IntVect(1,0,0)),dm,2,IntVect(NGROW,NGROW,0));
164 MultiFab tmp_rv_new(convert(ba, IntVect(0,1,0)),dm,2,IntVect(NGROW,NGROW,0));
165
166 MultiFab tmp_ru2d_new(convert(ba2d, IntVect(1,0,0)),dm,2,IntVect(NGROW,NGROW,0));
167 MultiFab tmp_rv2d_new(convert(ba2d, IntVect(0,1,0)),dm,2,IntVect(NGROW,NGROW,0));
168
169 init_masks(lev, ba, dm);
170
171 // This will fill the temporary MultiFabs with data from previous fine data as well as coarse where needed
172 FillPatch(lev, time, tmp_cons_new, cons_new, BCVars::cons_bc, BdyVars::t,0,true,false);
173 FillPatch(lev, time, tmp_xvel_new, xvel_new, BCVars::xvel_bc, BdyVars::u,0,true,false,0,0,0.0,tmp_xvel_new);
174 FillPatch(lev, time, tmp_yvel_new, yvel_new, BCVars::yvel_bc, BdyVars::v,0,true,false,0,0,0.0,tmp_yvel_new);
175 FillPatch(lev, time, tmp_zvel_new, zvel_new, BCVars::zvel_bc, BdyVars::null,0,true,false);
176
177 FillPatch(lev, time, tmp_h, GetVecOfPtrs(vec_hOfTheConfusingName), BCVars::cons_bc, BdyVars::null,0,false,false);
178 FillPatch(lev, time, tmp_h, GetVecOfPtrs(vec_hOfTheConfusingName), BCVars::cons_bc, BdyVars::null,1,false,false);
179 FillPatch(lev, time, tmp_Zt_avg1_new, GetVecOfPtrs(vec_Zt_avg1), BCVars::zeta_bc, BdyVars::null,0,true,false);
180 for (int icomp=0; icomp<3; icomp++) {
181 FillPatch(lev, time, tmp_ubar_new, GetVecOfPtrs(vec_ubar), BCVars::ubar_bc, BdyVars::ubar, icomp,false,false);
182 FillPatch(lev, time, tmp_vbar_new, GetVecOfPtrs(vec_vbar), BCVars::vbar_bc, BdyVars::vbar, icomp,false,false);
183 }
184 for (int icomp=0; icomp<2; icomp++) {
185 FillPatch(lev, time, tmp_ru_new, GetVecOfPtrs(vec_ru),BCVars::xvel_bc, BdyVars::null, icomp,false,false);
186 FillPatch(lev, time, tmp_rv_new, GetVecOfPtrs(vec_rv),BCVars::yvel_bc, BdyVars::null, icomp,false,false);
187 // These might want to have BCVars::ubar_bc and vbar_bc
188 FillPatch(lev, time, tmp_ru2d_new, GetVecOfPtrs(vec_ru2d),BCVars::xvel_bc, BdyVars::null, icomp,false,false);
189 FillPatch(lev, time, tmp_rv2d_new, GetVecOfPtrs(vec_rv2d),BCVars::yvel_bc, BdyVars::null, icomp,false,false);
190 }
191
192 MultiFab::Copy(tmp_cons_old,tmp_cons_new,0,0,NCONS,tmp_cons_new.nGrowVect());
193 MultiFab::Copy(tmp_xvel_old,tmp_xvel_new,0,0, 1,tmp_xvel_new.nGrowVect());
194 MultiFab::Copy(tmp_yvel_old,tmp_yvel_new,0,0, 1,tmp_yvel_new.nGrowVect());
195 MultiFab::Copy(tmp_zvel_old,tmp_zvel_new,0,0, 1,tmp_zvel_new.nGrowVect());
196
197 std::swap(tmp_cons_new, *cons_new[lev]);
198 std::swap(tmp_cons_old, *cons_old[lev]);
199 std::swap(tmp_xvel_new, *xvel_new[lev]);
200 std::swap(tmp_xvel_old, *xvel_old[lev]);
201 std::swap(tmp_yvel_new, *yvel_new[lev]);
202 std::swap(tmp_yvel_old, *yvel_old[lev]);
203 std::swap(tmp_zvel_new, *zvel_new[lev]);
204 std::swap(tmp_zvel_old, *zvel_old[lev]);
205 std::swap(tmp_Zt_avg1_new, *vec_Zt_avg1[lev]);
206 std::swap(tmp_h, *vec_hOfTheConfusingName[lev]);
207 std::swap(tmp_ubar_new, *vec_ubar[lev]);
208 std::swap(tmp_vbar_new, *vec_vbar[lev]);
209 std::swap(tmp_ru_new, *vec_ru[lev]);
210 std::swap(tmp_rv_new, *vec_rv[lev]);
211 std::swap(tmp_ru2d_new, *vec_ru2d[lev]);
212 std::swap(tmp_rv2d_new, *vec_rv2d[lev]);
213
214 t_new[lev] = time;
215 t_old[lev] = time - 1.e200_rt;
216
217 init_stuff(lev, ba, dm);
218
219 set_grid_scale(lev);
221
222 init_set_vmix(lev);
223 set_hmixcoef(lev);
224 set_coriolis(lev);
226 // Previously set smflux here
227
228 // We need to re-define the FillPatcher if the grids have changed
229 if (lev > 0 && cf_width >= 0) {
230 bool ba_changed = (ba != ba_old);
231 bool dm_changed = (dm != dm_old);
232 if (ba_changed || dm_changed) {
234 }
235 }
236
237#ifdef REMORA_USE_PARTICLES
238 particleData.Redistribute();
239#endif
240}
241
242/**
243 * Make a new level from scratch using provided BoxArray and DistributionMapping.
244 * This is called both for initialization and for restart
245 * (overrides the pure virtual function in AmrCore)
246 * main.cpp --> REMORA::InitData --> InitFromScratch --> MakeNewGrids --> MakeNewLevelFromScratch
247 * restart --> MakeNewGrids --> MakeNewLevelFromScratch
248 *
249 * @param[in ] lev level to make
250 * @param[in ] time current time
251 * @param[in ] ba BoxArray for the level
252 * @param[in ] dm DistributionMapping for the level
253 */
254void REMORA::MakeNewLevelFromScratch (int lev, Real time, const BoxArray& ba,
255 const DistributionMapping& dm)
256{
257 // Set BoxArray grids and DistributionMapping dmap in AMReX_AmrMesh.H class
258 SetBoxArray(lev, ba);
259 SetDistributionMap(lev, dm);
260
261 BoxList bl2d = ba.boxList();
262 for (auto& b : bl2d) {
263 b.setRange(2,0);
264 }
265 BoxArray ba2d(std::move(bl2d));
266
267 amrex::Print() << "GRIDS AT LEVEL " << lev << " ARE " << ba << std::endl;
268
269 // The number of ghost cells for density must be 1 greater than that for velocity
270 // so that we can go back in forth between velocity and momentum on all faces
271#if NGROW==2
272 int ngrow_state = ComputeGhostCells(solverChoice.spatial_order)+1;
273 int ngrow_vels = ComputeGhostCells(solverChoice.spatial_order)+1;
274#else
275 int ngrow_state = ComputeGhostCells(solverChoice.spatial_order)+2;
276 int ngrow_vels = ComputeGhostCells(solverChoice.spatial_order)+2;
277#endif
278
279 cons_old[lev] = new MultiFab(ba, dm, NCONS, ngrow_state);
280 cons_new[lev] = new MultiFab(ba, dm, NCONS, ngrow_state);
281
282 xvel_new[lev] = new MultiFab(convert(ba, IntVect(1,0,0)), dm, 1, ngrow_vels);
283 xvel_old[lev] = new MultiFab(convert(ba, IntVect(1,0,0)), dm, 1, ngrow_vels);
284
285 yvel_new[lev] = new MultiFab(convert(ba, IntVect(0,1,0)), dm, 1, ngrow_vels);
286 yvel_old[lev] = new MultiFab(convert(ba, IntVect(0,1,0)), dm, 1, ngrow_vels);
287
288 zvel_new[lev] = new MultiFab(convert(ba, IntVect(0,0,1)), dm, 1, IntVect(ngrow_vels,ngrow_vels,0));
289 zvel_old[lev] = new MultiFab(convert(ba, IntVect(0,0,1)), dm, 1, IntVect(ngrow_vels,ngrow_vels,0));
290
291 resize_stuff(lev);
292
293 vec_Zt_avg1[lev].reset(new MultiFab(ba2d ,dm,1,IntVect(NGROW+1,NGROW+1,0))); //2d, average of the free surface (zeta)
294 vec_hOfTheConfusingName[lev].reset(new MultiFab(ba2d ,dm,2,IntVect(NGROW+1,NGROW+1,0))); //2d, bathymetry
295 vec_ubar[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,3,IntVect(NGROW,NGROW,0)));
296 vec_vbar[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,3,IntVect(NGROW,NGROW,0)));
297
298 vec_ru[lev].reset(new MultiFab(convert(ba,IntVect(1,0,0)),dm,2,IntVect(NGROW,NGROW,0))); // RHS u (incl horizontal and vertical advection)
299 vec_rv[lev].reset(new MultiFab(convert(ba,IntVect(0,1,0)),dm,2,IntVect(NGROW,NGROW,0))); // RHS v
300
301 vec_ru2d[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,2,IntVect(NGROW,NGROW,0))); // RHS u (incl horizontal and vertical advection)
302 vec_rv2d[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,2,IntVect(NGROW,NGROW,0))); // RHS v
303
304 init_masks(lev, ba, dm);
305 init_stuff(lev, ba, dm);
306
307 init_only(lev, time);
308
309#ifdef REMORA_USE_PARTICLES
310 if (restart_chkfile.empty()) {
311 if (lev == 0) {
312 initializeTracers((ParGDBBase*)GetParGDB(),vec_z_phys_nd);
313 } else {
314 particleData.Redistribute();
315 }
316 }
317#endif
318}
319
320/**
321 * @param[in ] lev level do operate on
322 */
324{
325 vec_z_phys_nd.resize(lev+1);
326
327 vec_hOfTheConfusingName.resize(lev+1);
328 vec_Zt_avg1.resize(lev+1);
329 vec_s_r.resize(lev+1);
330 vec_s_w.resize(lev+1);
331 vec_Cs_r.resize(lev+1);
332 vec_Cs_w.resize(lev+1);
333 vec_z_w.resize(lev+1);
334 vec_z_r.resize(lev+1);
335 vec_Hz.resize(lev+1);
336 vec_Huon.resize(lev+1);
337 vec_Hvom.resize(lev+1);
338 vec_Akv.resize(lev+1);
339 vec_Akt.resize(lev+1);
340 vec_visc2_p.resize(lev+1);
341 vec_visc2_r.resize(lev+1);
342 vec_diff2.resize(lev+1);
343 vec_ru.resize(lev+1);
344 vec_rv.resize(lev+1);
345 vec_ru2d.resize(lev+1);
346 vec_rv2d.resize(lev+1);
347 vec_rufrc.resize(lev+1);
348 vec_rvfrc.resize(lev+1);
349 vec_sustr.resize(lev+1);
350 vec_svstr.resize(lev+1);
351 vec_btflx.resize(lev+1);
352 vec_stflx.resize(lev+1);
353 vec_btflux.resize(lev+1);
354 vec_stflux.resize(lev+1);
355 vec_lrflx.resize(lev+1);
356 vec_lhflx.resize(lev+1);
357 vec_shflx.resize(lev+1);
358 vec_rain.resize(lev+1);
359 vec_evap.resize(lev+1);
360 vec_rdrag.resize(lev+1);
361 vec_rdrag2.resize(lev+1);
362 vec_ZoBot.resize(lev+1);
363 vec_bustr.resize(lev+1);
364 vec_bvstr.resize(lev+1);
365 vec_uwind.resize(lev+1);
366 vec_vwind.resize(lev+1);
367 vec_alpha.resize(lev+1);
368 vec_beta.resize(lev+1);
369
370 vec_DU_avg1.resize(lev+1);
371 vec_DU_avg2.resize(lev+1);
372 vec_DV_avg1.resize(lev+1);
373 vec_DV_avg2.resize(lev+1);
374 vec_rubar.resize(lev+1);
375 vec_rvbar.resize(lev+1);
376 vec_rzeta.resize(lev+1);
377 vec_ubar.resize(lev+1);
378 vec_vbar.resize(lev+1);
379 vec_zeta.resize(lev+1);
380 vec_mskr.resize(lev+1);
381 vec_msku.resize(lev+1);
382 vec_mskv.resize(lev+1);
383 vec_mskp.resize(lev+1);
384 vec_sstore.resize(lev+1);
385
386 vec_pm.resize(lev+1);
387 vec_pn.resize(lev+1);
388 vec_fcor.resize(lev+1);
389
390 vec_xr.resize(lev+1);
391 vec_yr.resize(lev+1);
392 vec_xu.resize(lev+1);
393 vec_yu.resize(lev+1);
394 vec_xv.resize(lev+1);
395 vec_yv.resize(lev+1);
396 vec_xp.resize(lev+1);
397 vec_yp.resize(lev+1);
398
399 vec_rhoS.resize(lev+1);
400 vec_rhoA.resize(lev+1);
401 vec_bvf.resize(lev+1);
402
403 vec_tke.resize(lev+1);
404 vec_gls.resize(lev+1);
405 vec_Lscale.resize(lev+1);
406 vec_Akk.resize(lev+1);
407 vec_Akp.resize(lev+1);
408
409 vec_river_position.resize(lev+1);
410
411 if (lev==0) vec_nudg_coeff.resize(BdyVars::NumTypes);
412
413 vec_nudg_coeff[BdyVars::u].resize(lev+1);
414 vec_nudg_coeff[BdyVars::v].resize(lev+1);
415 vec_nudg_coeff[BdyVars::t].resize(lev+1);
416 vec_nudg_coeff[BdyVars::s].resize(lev+1);
417 vec_nudg_coeff[BdyVars::ubar].resize(lev+1);
418 vec_nudg_coeff[BdyVars::vbar].resize(lev+1);
419 vec_nudg_coeff[BdyVars::zeta].resize(lev+1);
420}
421
422/**
423 * @param[in ] lev level to operate on
424 * @param[in ] ba BoxArray for the level
425 * @param[in ] dm DistributionMapping for the level
426 */
427void REMORA::init_masks (int lev, const BoxArray& ba, const DistributionMapping& dm)
428{
429 BoxList bl2d = ba.boxList();
430 for (auto& b : bl2d) {
431 b.setRange(2,0);
432 }
433
434 BoxArray ba2d(std::move(bl2d));
435 vec_mskr[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW+1,NGROW+1,0)));
436 vec_msku[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,1,IntVect(NGROW+1,NGROW+1,0)));
437 vec_mskv[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,1,IntVect(NGROW+1,NGROW+1,0)));
438 vec_mskp[lev].reset(new MultiFab(convert(ba2d,IntVect(1,1,0)),dm,1,IntVect(NGROW+1,NGROW+1,0)));
439
440 vec_mskr[lev]->setVal(1.0_rt);
441 vec_msku[lev]->setVal(1.0_rt);
442 vec_mskv[lev]->setVal(1.0_rt);
443 vec_mskp[lev]->setVal(1.0_rt);
444}
445
446/**
447 * @param[in ] lev level to operate on
448 * @param[in ] ba BoxArray for the level
449 * @param[in ] dm DistributionMapping for the level
450 */
451void REMORA::init_stuff (int lev, const BoxArray& ba, const DistributionMapping& dm)
452{
453 // ********************************************************************************************
454 // Initialize the boundary conditions
455 // ********************************************************************************************
456 physbcs[lev] = std::make_unique<REMORAPhysBCFunct> (lev, geom[lev], domain_bcs_type, domain_bcs_type_d,
458
459 BoxList bl2d = ba.boxList();
460 for (auto& b : bl2d) {
461 b.setRange(2,0);
462 }
463 BoxArray ba2d(std::move(bl2d));
464
465 BoxList bl1d = ba.boxList();
466 for (auto& b : bl1d) {
467 b.setRange(0,0);
468 b.setRange(1,0);
469 }
470 BoxArray ba1d(std::move(bl1d));
471
472 BoxArray ba_nd(ba);
473 ba_nd.surroundingNodes();
474 BoxArray ba_w(ba);
475 ba_w.surroundingNodes(2);
476
477 vec_z_phys_nd[lev].reset (new MultiFab(ba_nd,dm,1,IntVect(NGROW,NGROW,1))); // z at psi points (nodes) MIGHT NEED NGROW+1
478
479 vec_s_r[lev].reset (new MultiFab(ba1d,dm,1,IntVect( 0, 0,0))); // scaled vertical coordinate [0,1] , transforms to z
480
481 vec_s_w[lev].reset (new MultiFab(convert(ba1d,IntVect(0,0,1)),dm,1,IntVect( 0, 0,0))); // scaled vertical coordinate at w-points [0,1] , transforms to z
482 //
483 vec_Cs_r[lev].reset (new MultiFab(ba1d,dm,1,IntVect( 0, 0,0)));
484 vec_Cs_w[lev].reset (new MultiFab(convert(ba1d,IntVect(0,0,1)),dm,1,IntVect( 0, 0,0)));
485
486 vec_z_w[lev].reset (new MultiFab(convert(ba,IntVect(0,0,1)),dm,1,IntVect(NGROW+1,NGROW+1,0))); // z at w points (cell faces)
487 vec_z_r[lev].reset (new MultiFab(ba,dm,1,IntVect(NGROW+1,NGROW+1,0))); // z at r points (cell center)
488 vec_Hz[lev].reset (new MultiFab(ba,dm,1,IntVect(NGROW+1,NGROW+1,NGROW+1))); // like in ROMS, thickness of cell in z
489
490 vec_Huon[lev].reset (new MultiFab(convert(ba,IntVect(1,0,0)),dm,1,IntVect(NGROW,NGROW,0))); // mass flux for u component
491 vec_Hvom[lev].reset (new MultiFab(convert(ba,IntVect(0,1,0)),dm,1,IntVect(NGROW,NGROW,0))); // mass flux for v component
492
493 vec_Akv[lev].reset (new MultiFab(convert(ba,IntVect(0,0,1)),dm,1,IntVect(NGROW,NGROW,0))); // vertical mixing coefficient (.in)
494 vec_Akt[lev].reset (new MultiFab(convert(ba,IntVect(0,0,1)),dm,NCONS,IntVect(NGROW,NGROW,0))); // vertical mixing coefficient (.in)
495
496 // check dimensionality
497 vec_visc2_p[lev].reset(new MultiFab(ba,dm,1,IntVect(NGROW,NGROW,0))); // harmonic viscosity at psi points -- difference to 3d?
498 vec_visc2_r[lev].reset(new MultiFab(ba,dm,1,IntVect(NGROW,NGROW,0))); // harmonic viscosity at rho points
499 vec_diff2[lev].reset(new MultiFab(ba,dm,NCONS,IntVect(NGROW,NGROW,0))); // harmonic diffusivity temperature/salt
500
501 //2d, (incl advection terms and surface/bottom stresses, integral over the whole column, k=0)
502 vec_rufrc[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,2,IntVect(NGROW,NGROW,0)));
503 vec_rvfrc[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,2,IntVect(NGROW,NGROW,0))); //2d, same as above but v
504
505 vec_sustr[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,1,IntVect(NGROW,NGROW,0))); //2d, surface stress
506 vec_svstr[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,1,IntVect(NGROW,NGROW,0))); //2d
507
509 //2d, linear drag coefficient [m/s], defined at rho, somehow related to rdrg
510 vec_rdrag[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
512 vec_rdrag2[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
513 }
514
517 vec_ZoBot[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
518 }
519
520 vec_bustr[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,1,IntVect(NGROW,NGROW,0))); //2d, bottom stress
521 vec_bvstr[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,1,IntVect(NGROW,NGROW,0)));
522
523 //all 2d -- all associated with the 2D advance
524 //2d DU: sum(height[incl free surface?] * u)
525 vec_DU_avg1[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,1,IntVect(NGROW,NGROW,0)));
526
527 //2d like above, but correct(or)?
528 vec_DU_avg2[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,1,IntVect(NGROW,NGROW,0)));
529
530 vec_DV_avg1[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,1,IntVect(NGROW,NGROW,0)));
531 vec_DV_avg2[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,1,IntVect(NGROW,NGROW,0)));
532
533 vec_rubar[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,4,IntVect(NGROW,NGROW,0))); // 2d RHS ubar
534 vec_rvbar[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,4,IntVect(NGROW,NGROW,0)));
535 vec_rzeta[lev].reset(new MultiFab(ba2d,dm,4,IntVect(NGROW,NGROW,0))); // 2d RHS zeta
536
537 // starts off kind of like a depth-averaged u, but exists at more points and more timesteps (b/c fast 2D update) than full u
538 vec_zeta[lev].reset(new MultiFab(ba2d,dm,3,IntVect(NGROW+1,NGROW+1,0))); // 2d free surface
539
540 vec_pm[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW+1,NGROW+2,0)));
541 vec_pn[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW+2,NGROW+1,0)));
542 vec_fcor[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW+1,NGROW+1,0)));
543
544 vec_xr[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW+1,NGROW+1,0)));
545 vec_yr[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW+1,NGROW+1,0)));
546
547 vec_xu[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,1,IntVect(NGROW,NGROW,0)));
548 vec_yu[lev].reset(new MultiFab(convert(ba2d,IntVect(1,0,0)),dm,1,IntVect(NGROW,NGROW,0)));
549 vec_xv[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,1,IntVect(NGROW,NGROW,0)));
550 vec_yv[lev].reset(new MultiFab(convert(ba2d,IntVect(0,1,0)),dm,1,IntVect(NGROW,NGROW,0)));
551 vec_xp[lev].reset(new MultiFab(convert(ba2d,IntVect(1,1,0)),dm,1,IntVect(NGROW,NGROW,0)));
552 vec_yp[lev].reset(new MultiFab(convert(ba2d,IntVect(1,1,0)),dm,1,IntVect(NGROW,NGROW,0)));
553
554
555 // tempstore, saltstore, etc
556 vec_sstore[lev].reset(new MultiFab(ba,dm,NCONS,IntVect(NGROW,NGROW,0)));
557
558 vec_rhoS[lev].reset(new MultiFab(ba,dm,1,IntVect(NGROW,NGROW,0)));
559 vec_rhoA[lev].reset(new MultiFab(ba,dm,1,IntVect(NGROW,NGROW,0)));
560 vec_bvf[lev].reset(new MultiFab(convert(ba,IntVect(0,0,1)),dm,1,IntVect(NGROW,NGROW,0)));
561
562 vec_tke[lev].reset(new MultiFab(convert(ba,IntVect(0,0,1)),dm,3,IntVect(NGROW,NGROW,0)));
563 vec_gls[lev].reset(new MultiFab(convert(ba,IntVect(0,0,1)),dm,3,IntVect(NGROW,NGROW,0)));
564 vec_Lscale[lev].reset(new MultiFab(convert(ba,IntVect(0,0,1)),dm,1,IntVect(NGROW,NGROW,0)));
565 vec_Akk[lev].reset(new MultiFab(convert(ba,IntVect(0,0,1)),dm,1,IntVect(NGROW,NGROW,0)));
566 vec_Akp[lev].reset(new MultiFab(convert(ba,IntVect(0,0,1)),dm,1,IntVect(NGROW,NGROW,0)));
567
568 // surface/bottom tracer fluxes for update
569 vec_stflx[lev].reset(new MultiFab(ba2d,dm,NCONS,IntVect(NGROW,NGROW,0)));
570 vec_btflx[lev].reset(new MultiFab(ba2d,dm,NCONS,IntVect(NGROW,NGROW,0)));
571 // surface/bottom tracer fluxes to be filled by inputs
572 vec_stflux[lev].reset(new MultiFab(ba2d,dm,NCONS,IntVect(NGROW,NGROW,0)));
573 vec_btflux[lev].reset(new MultiFab(ba2d,dm,NCONS,IntVect(NGROW,NGROW,0)));
574
576 vec_uwind[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0))); //2d, surface wind u
577 vec_vwind[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0))); //2d, surface wind v
578 vec_alpha[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
579 vec_beta[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
580 vec_lrflx[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
581 vec_lhflx[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
582 vec_shflx[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
583 vec_rain[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
584 vec_evap[lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
585 vec_lhflx[lev]->setVal(0.0_rt);
586 vec_shflx[lev]->setVal(0.0_rt);
587 vec_rain[lev]->setVal(solverChoice.rain);
588 }
589
591 vec_river_position[lev].reset(new iMultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
592 vec_river_position[lev]->setVal(-1);
593 }
594
595 vec_nudg_coeff[BdyVars::u][lev].reset(new MultiFab(ba,dm,1,IntVect(NGROW,NGROW,0)));
596 vec_nudg_coeff[BdyVars::v][lev].reset(new MultiFab(ba,dm,1,IntVect(NGROW,NGROW,0)));
597 vec_nudg_coeff[BdyVars::t][lev].reset(new MultiFab(ba,dm,1,IntVect(NGROW,NGROW,0)));
598 vec_nudg_coeff[BdyVars::s][lev].reset(new MultiFab(ba,dm,1,IntVect(NGROW,NGROW,0)));
599 vec_nudg_coeff[BdyVars::ubar][lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
600 vec_nudg_coeff[BdyVars::vbar][lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
601 vec_nudg_coeff[BdyVars::zeta][lev].reset(new MultiFab(ba2d,dm,1,IntVect(NGROW,NGROW,0)));
602
603 set_weights(lev);
604
605 vec_DU_avg1[lev]->setVal(0.0_rt);
606 vec_DU_avg2[lev]->setVal(0.0_rt);
607 vec_DV_avg1[lev]->setVal(0.0_rt);
608 vec_DV_avg2[lev]->setVal(0.0_rt);
609 vec_rubar[lev]->setVal(0.0_rt);
610 vec_rvbar[lev]->setVal(0.0_rt);
611 vec_rzeta[lev]->setVal(0.0_rt);
612
613 // Initialize these vars even if we aren't using GLS to
614 // avoid issues on e.g. checkpoint
615 vec_tke[lev]->setVal(solverChoice.gls_Kmin);
616 vec_gls[lev]->setVal(solverChoice.gls_Pmin);
617 vec_Lscale[lev]->setVal(0.0_rt);
618 vec_Akk[lev]->setVal(solverChoice.Akk_bak);
619 vec_Akp[lev]->setVal(solverChoice.Akp_bak);
620
621 vec_stflx[lev]->setVal(0.0_rt);
622 vec_btflx[lev]->setVal(0.0_rt);
623 vec_stflux[lev]->setVal(0.0_rt);
624 vec_btflux[lev]->setVal(0.0_rt);
625
626 // NOTE: Used to set vec_pm and vec_pn to 1e34 here to make foextrap work
627 // when init_type = real. However, this does not appear to be necessary so removing
628
629 // Set initial linear drag coefficient
631 vec_rdrag[lev]->setVal(solverChoice.rdrag);
633 vec_rdrag2[lev]->setVal(solverChoice.rdrag2);
634 }
635
638 vec_ZoBot[lev]->setVal(solverChoice.Zob);
639 }
640
641
642 // ********************************************************************************************
643 // Create the REMORAFillPatcher object
644 // ********************************************************************************************
645 if (lev > 0 && cf_width >= 0) {
648 }
649}
650
651/**
652 * Delete level data. Overrides the pure virtual function in AmrCore
653 *
654 * @param[in ] lev level to operate on
655 */
656void
658{
659 delete cons_new[lev]; delete xvel_new[lev]; delete yvel_new[lev]; delete zvel_new[lev];
660 delete cons_old[lev]; delete xvel_old[lev]; delete yvel_old[lev]; delete zvel_old[lev];
661}
662
663/**
664 * @param[in ] lev level to operate on
665 */
666void
668{
671 const auto dxi = Geom(lev).InvCellSize();
672 vec_pm[lev]->setVal(dxi[0]); vec_pm[lev]->FillBoundary(geom[lev].periodicity());
673 vec_pn[lev]->setVal(dxi[1]); vec_pn[lev]->FillBoundary(geom[lev].periodicity());
675 prob->init_analytic_grid_scale(lev, Geom(lev), solverChoice, *this, *vec_pm[lev].get(), *vec_pn[lev].get());
676 vec_pm[lev]->FillBoundary(geom[lev].periodicity());
677 vec_pn[lev]->FillBoundary(geom[lev].periodicity());
678 }
679
680 for ( MFIter mfi(*vec_xr[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
681 {
682 Array4<const Real> const& pm = vec_pm[lev]->const_array(mfi);
683 Array4<const Real> const& pn = vec_pn[lev]->const_array(mfi);
684 Array4<Real> const& xr = vec_xr[lev]->array(mfi);
685 Array4<Real> const& yr = vec_yr[lev]->array(mfi);
686 Array4<Real> const& xu = vec_xu[lev]->array(mfi);
687 Array4<Real> const& yu = vec_yu[lev]->array(mfi);
688 Array4<Real> const& xv = vec_xv[lev]->array(mfi);
689 Array4<Real> const& yv = vec_yv[lev]->array(mfi);
690 Array4<Real> const& xp = vec_xp[lev]->array(mfi);
691 Array4<Real> const& yp = vec_yp[lev]->array(mfi);
692
693 Box bx = mfi.growntilebox(IntVect(NGROW,NGROW,0));
694 ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int)
695 {
696 xr(i,j,0) = (i + 0.5_rt) / pm(i,j,0);
697 yr(i,j,0) = (j + 0.5_rt) / pn(i,j,0);
698 });
699
700 ParallelFor(grow(convert(bx,IntVect(1,0,0)),IntVect(-1,0,0)), [=] AMREX_GPU_DEVICE (int i, int j, int)
701 {
702 xu(i,j,0) = i / pm(i,j,0);
703 yu(i,j,0) = (j + 0.5_rt) / pn(i,j,0);
704 });
705
706 ParallelFor(grow(convert(bx,IntVect(0,1,0)),IntVect(0,-1,0)), [=] AMREX_GPU_DEVICE (int i, int j, int)
707 {
708 xv(i,j,0) = (i + 0.5_rt) / pm(i,j,0);
709 yv(i,j,0) = j / pn(i,j,0);
710 });
711
712 ParallelFor(grow(convert(bx,IntVect(1,1,0)),IntVect(-1,-1,0)), [=] AMREX_GPU_DEVICE (int i, int j, int)
713 {
714 xp(i,j,0) = i / pm(i,j,0);
715 yp(i,j,0) = j / pn(i,j,0);
716 });
717 }
718}
719
720/**
721 * @param[in ] lev level to operate on
722 */
723void
725{
726 std::unique_ptr<MultiFab>& mf_zeta = vec_zeta[lev];
727 std::unique_ptr<MultiFab>& mf_Zt_avg1 = vec_Zt_avg1[lev];
729 for ( MFIter mfi(*vec_zeta[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
730 {
731 Array4<Real> const& Zt_avg1 = (mf_Zt_avg1)->array(mfi);
732 Array4<const Real> const& evap = vec_evap[lev]->const_array(mfi);
733 Array4<const Real> const& rain = vec_rain[lev]->const_array(mfi);
734
735 Box bx2 = mfi.growntilebox(IntVect(NGROW,NGROW,0));// bx2.grow(IntVect(NGROW,NGROW,0));
736
737 Real cff = dt[lev] / rhow;
738
739 ParallelFor(bx2, [=] AMREX_GPU_DEVICE (int i, int j, int )
740 {
741 Zt_avg1(i,j,0) = Zt_avg1(i,j,0) - (evap(i,j,0) - rain(i,j,0)) * cff;
742 });
743 }
744 }
745 Gpu::streamSynchronize();
746
747 vec_Zt_avg1[lev]->FillBoundary(geom[lev].periodicity());
748
749 for ( MFIter mfi(*vec_zeta[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
750 {
751 Box bx3 = mfi.tilebox(); bx3.grow(IntVect(NGROW+1,NGROW+1,0));
752 Array4<Real> const& zeta = mf_zeta->array(mfi);
753 Array4<Real> const& Zt_avg1 = (mf_Zt_avg1)->array(mfi);
754
755 ParallelFor(bx3, 3, [=] AMREX_GPU_DEVICE (int i, int j, int , int n)
756 {
757 zeta(i,j,0,n) = Zt_avg1(i,j,0);
758 });
759 }
760}
761
762/**
763 * @param[in ] lev level to operate on
764 */
765void
767{
768 for ( MFIter mfi(*vec_mskr[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
769 {
770 Array4<const Real> const& mskr = vec_mskr[lev]->const_array(mfi);
771 Array4< Real> const& mskp = vec_mskp[lev]->array(mfi);
772
773 Box bx = mfi.tilebox(); bx.grow(IntVect(1,1,0)); bx.makeSlab(2,0);
774
775 Real cff1 = 1.0_rt;
776 Real cff2 = 2.0_rt;
777
778 ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int)
779 {
780 if ((mskr(i-1,j,0) > 0.5) and (mskr(i,j,0) > 0.5) and (mskr(i-1,j-1,0) > 0.5) and (mskr(i,j-1,0) > 0.5)) {
781 mskp(i,j,0) = 1.0_rt;
782 } else if ((mskr(i-1,j,0) < 0.5) and (mskr(i,j,0) > 0.5) and (mskr(i-1,j-1,0) > 0.5) and (mskr(i,j-1,0) > 0.5)) {
783 mskp(i,j,0) = cff1;
784 } else if ((mskr(i-1,j,0) > 0.5) and (mskr(i,j,0) < 0.5) and (mskr(i-1,j-1,0) > 0.5) and (mskr(i,j-1,0) > 0.5)) {
785 mskp(i,j,0) = cff1;
786 } else if ((mskr(i-1,j,0) > 0.5) and (mskr(i,j,0) > 0.5) and (mskr(i-1,j-1,0) < 0.5) and (mskr(i,j-1,0) > 0.5)) {
787 mskp(i,j,0) = cff1;
788 } else if ((mskr(i-1,j,0) > 0.5) and (mskr(i,j,0) > 0.5) and (mskr(i-1,j-1,0) > 0.5) and (mskr(i,j-1,0) < 0.5)) {
789 mskp(i,j,0) = cff1;
790 } else if ((mskr(i-1,j,0) > 0.5) and (mskr(i,j,0) < 0.5) and (mskr(i-1,j-1,0) > 0.5) and (mskr(i,j-1,0) < 0.5)) {
791 mskp(i,j,0) = cff2;
792 } else if ((mskr(i-1,j,0) < 0.5) and (mskr(i,j,0) > 0.5) and (mskr(i-1,j-1,0) < 0.5) and (mskr(i,j-1,0) > 0.5)) {
793 mskp(i,j,0) = cff2;
794 } else if ((mskr(i-1,j,0) > 0.5) and (mskr(i,j,0) > 0.5) and (mskr(i-1,j-1,0) < 0.5) and (mskr(i,j-1,0) < 0.5)) {
795 mskp(i,j,0) = cff2;
796 } else if ((mskr(i-1,j,0) < 0.5) and (mskr(i,j,0) < 0.5) and (mskr(i-1,j-1,0) > 0.5) and (mskr(i,j-1,0) > 0.5)) {
797 mskp(i,j,0) = cff2;
798 } else {
799 mskp(i,j,0) = 0.0_rt;
800 }
801
802 });
803 }
804}
constexpr amrex::Real rhow
#define NGROW
#define NCONS
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_evap
evaporation rate [kg/m^2/s]
Definition REMORA.H:310
amrex::Vector< amrex::BCRec > domain_bcs_type
vector (over BCVars) of BCRecs
Definition REMORA.H:1120
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Cs_r
Stretching coefficients at rho points.
Definition REMORA.H:271
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:360
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_btflux
Bottom tracer flux; input arrays.
Definition REMORA.H:305
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_ZoBot
Bottom roughness length [m], defined at rho points.
Definition REMORA.H:317
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_lrflx
longwave radiation
Definition REMORA.H:292
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yv
y_grid on v-points (2D)
Definition REMORA.H:375
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:217
virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Make a new level using provided BoxArray and DistributionMapping and fill with interpolated coarse le...
void stretch_transform(int lev)
Calculate vertical stretched coordinates.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vwind
Wind in the v direction, defined at rho-points.
Definition REMORA.H:289
std::unique_ptr< ProblemBase > prob
Pointer to container of analytical functions for problem definition.
Definition REMORA.H:1084
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:346
void Construct_REMORAFillPatchers(int lev)
Construct FillPatchers.
Definition REMORA.cpp:364
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rain
precipitation rate [kg/m^2/s]
Definition REMORA.H:308
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_tke
Turbulent kinetic energy.
Definition REMORA.H:403
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflx
Surface tracer flux; working arrays.
Definition REMORA.H:299
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_gls
Turbulent generic length scale.
Definition REMORA.H:405
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:282
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yp
y_grid on psi-points (2D)
Definition REMORA.H:380
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xr
x_grid on rho points (2D)
Definition REMORA.H:363
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:242
virtual void ClearLevel(int lev) override
Delete level data Overrides the pure virtual function in AmrCore.
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_FORCE_INLINE int ComputeGhostCells(const int &spatial_order)
Helper function to determine number of ghost cells.
Definition REMORA.H:1309
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sstore
additional scratch space for calculations on temp, salt, etc
Definition REMORA.H:383
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xv
x_grid on v-points (2D)
Definition REMORA.H:373
void init_only(int lev, amrex::Real time)
Init (NOT restart or regrid)
Definition REMORA.cpp:674
void init_set_vmix(int lev)
Initialize vertical mixing coefficients from file or analytic.
Definition REMORA.cpp:583
void set_grid_scale(int lev)
Set pm and pn arrays on level lev. Only works if using analytic initialization.
void set_coriolis(int lev)
Initialize Coriolis factor from file or analytic.
Definition REMORA.cpp:560
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Lscale
Vertical mixing turbulent length scale.
Definition REMORA.H:407
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hz
Width of cells in the vertical (z-) direction (3D, Hz in ROMS)
Definition REMORA.H:232
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
Definition REMORA.H:252
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< 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
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_uwind
Wind in the u direction, defined at rho-points.
Definition REMORA.H:287
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rufrc
u velocity RHS, integrated, including advection and bottom/surface stresses (2D)
Definition REMORA.H:246
void Define_REMORAFillPatchers(int lev)
Define FillPatchers.
Definition REMORA.cpp:412
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_shflx
sensible heat flux
Definition REMORA.H:296
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
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< std::unique_ptr< amrex::MultiFab > > vec_z_r
z coordinates at rho points (cell centers)
Definition REMORA.H:261
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:219
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lhflx
latent heat flux
Definition REMORA.H:294
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_s_w
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at w-points (cell faces)
Definition REMORA.H:268
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvf
Brunt-Vaisala frequency (3D)
Definition REMORA.H:390
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Remake an existing level using provided BoxArray and DistributionMapping and fill with existing fine ...
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskv
land/sea mask at y-faces (2D)
Definition REMORA.H:350
void init_masks(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Allocate MultiFabs for masks.
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
Definition REMORA.H:1108
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< 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
void update_mskp(int lev)
Set psi-point mask to be consistent with rho-point mask.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:284
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Huon
u-volume flux (3D)
Definition REMORA.H:234
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:212
void init_stuff(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Allocate MultiFabs for state and evolution variables.
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_Cs_w
Stretching coefficients at w points.
Definition REMORA.H:273
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
void resize_stuff(int lev)
Resize variable containers to accommodate data on levels 0 to max_lev.
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > vec_river_position
iMultiFab for river positions; contents are indices of rivers
Definition REMORA.H:1036
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akk
Turbulent kinetic energy vertical diffusion coefficient.
Definition REMORA.H:409
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag2
Quadratic drag coefficient [unitless], defined at rho points.
Definition REMORA.H:315
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:238
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xp
x_grid on psi-points (2D)
Definition REMORA.H:378
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 set_zeta_to_Ztavg(int lev)
Set zeta components to be equal to time-averaged Zt_avg1.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_alpha
Thermal expansion coefficient (3D)
Definition REMORA.H:392
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NCONS+8 > m_bc_extdir_vals
Array holding the Dirichlet values at walls which need them.
Definition REMORA.H:1128
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< std::unique_ptr< amrex::MultiFab > > vec_bustr
Bottom stress in the u direction.
Definition REMORA.H:320
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
void set_hmixcoef(int lev)
Initialize horizontal mixing coefficients.
Definition REMORA.cpp:612
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yu
y_grid on u-points (2D)
Definition REMORA.H:370
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvstr
Bottom stress in the v direction.
Definition REMORA.H:322
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_z_phys_nd
z coordinates at psi points (cell nodes)
Definition REMORA.H:276
void FillCoarsePatch(int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse, const int icomp=0, const bool fill_all=true)
fill an entire multifab by interpolating from the coarser level
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xu
x_grid on u-points (2D)
Definition REMORA.H:368
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > vec_nudg_coeff
Climatology nudging coefficients.
Definition REMORA.H:414
void set_weights(int lev)
Set weights for averaging 3D variables to 2D.
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_Akv
Vertical viscosity coefficient (3D)
Definition REMORA.H:250
virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Make a new level from scratch using provided BoxArray and DistributionMapping. Only used during initi...
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflux
Surface tracer flux; input arrays.
Definition REMORA.H:301
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag
Linear drag coefficient [m/s], defined at rho points.
Definition REMORA.H:313
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:279
std::string restart_chkfile
If set, restart from this checkpoint file.
Definition REMORA.H:1166
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:240
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_btflx
Bottom tracer flux; working arrays.
Definition REMORA.H:303
int cf_width
Nudging width at coarse-fine interface.
Definition REMORA.H:1046
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_beta
Saline contraction coefficient (3D)
Definition REMORA.H:394
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
amrex::Gpu::DeviceVector< amrex::BCRec > domain_bcs_type_d
GPU vector (over BCVars) of BCRecs.
Definition REMORA.H:1122
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_w
z coordinates at w points (faces between z-cells)
Definition REMORA.H:264
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yr
y_grid on rho points (2D)
Definition REMORA.H:365
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_s_r
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at rho points (cell centers)
Definition REMORA.H:266
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hvom
v-volume flux (3D)
Definition REMORA.H:236
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akp
Turbulent length scale vertical diffusion coefficient.
Definition REMORA.H:411
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
Definition REMORA.H:258
amrex::Real rdrag2
IC_BC_Type ic_bc_type
amrex::Real Akk_bak
BottomStressType bottom_stress_type
amrex::Real gls_Kmin
amrex::Real rdrag
VertMixingType vert_mixing_type
GridScaleType grid_scale_type
amrex::Real gls_Pmin
amrex::Real Akp_bak