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