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