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