REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_setup_step.cpp
Go to the documentation of this file.
1#include <REMORA.H>
2
3using namespace amrex;
4
5/**
6 * @param[in ] lev level to operate on
7 * @param[in ] time time at start of step
8 * @param[in ] dt_lev time step at level
9 */
10void
11REMORA::setup_step (int lev, Real time, Real dt_lev)
12{
13 BL_PROFILE("REMORA::setup_step()");
14
15 MultiFab& S_old = *cons_old[lev];
16 MultiFab& S_new = *cons_new[lev];
17
18 MultiFab& U_old = *xvel_old[lev];
19 MultiFab& V_old = *yvel_old[lev];
20 MultiFab& W_old = *zvel_old[lev];
21
22 MultiFab& U_new = *xvel_new[lev];
23 MultiFab& V_new = *yvel_new[lev];
24 MultiFab& W_new = *zvel_new[lev];
25
26 [[maybe_unused]] int nvars = S_old.nComp();
27
28 // Fill ghost cells/faces at old time
29 FillPatchNoBC(lev, time, *cons_old[lev], cons_old, BdyVars::t);
30 FillPatchNoBC(lev, time, *xvel_old[lev], xvel_old, BdyVars::u);
31 FillPatchNoBC(lev, time, *yvel_old[lev], yvel_old, BdyVars::v);
33
34 ////////// //pre_step3d corrections to boundaries
35
36 const BoxArray& ba = S_old.boxArray();
37 const DistributionMapping& dm = S_old.DistributionMap();
38
39 const int nrhs = 0;
40 int nstp = 0;
41
42 //-----------------------------------------------------------------------
43 // Time step momentum equation
44 //-----------------------------------------------------------------------
45
46 //Only used locally, probably should be rearranged into FArrayBox declaration
47
48 MultiFab mf_DC(ba,dm,1,IntVect(NGROW,NGROW,NGROW-1)); //2d missing j coordinate
49 MultiFab mf_logdrg_tmp(ba,dm,1,IntVect(NGROW,NGROW,0));
50 MultiFab mf_rho(ba,dm,1,IntVect(NGROW,NGROW,0));
51
52 MultiFab* mf_z_r = vec_z_r[lev].get();
53 MultiFab* mf_z_w = vec_z_w[lev].get();
54 MultiFab* mf_h = vec_h[lev].get();
55 MultiFab* mf_pm = vec_pm[lev].get();
56 MultiFab* mf_pn = vec_pn[lev].get();
57 MultiFab* mf_fcor = vec_fcor[lev].get();
58
59 MultiFab* mf_gls = vec_gls[lev].get();
60 MultiFab* mf_tke = vec_tke[lev].get();
61
62 //Consider passing these into the advance function or renaming relevant things
63
64 std::unique_ptr<MultiFab>& mf_rhoS = vec_rhoS[lev];
65 std::unique_ptr<MultiFab>& mf_rhoA = vec_rhoA[lev];
66 std::unique_ptr<MultiFab>& mf_bvf = vec_bvf[lev];
67 std::unique_ptr<MultiFab>& mf_ru = vec_ru[lev];
68 std::unique_ptr<MultiFab>& mf_rv = vec_rv[lev];
69 std::unique_ptr<MultiFab>& mf_rufrc = vec_rufrc[lev];
70 std::unique_ptr<MultiFab>& mf_rvfrc = vec_rvfrc[lev];
71 std::unique_ptr<MultiFab>& mf_sustr = vec_sustr[lev];
72 std::unique_ptr<MultiFab>& mf_svstr = vec_svstr[lev];
73 std::unique_ptr<MultiFab>& mf_rdrag = vec_rdrag[lev];
74 std::unique_ptr<MultiFab>& mf_rdrag2 = vec_rdrag2[lev];
75 std::unique_ptr<MultiFab>& mf_ZoBot = vec_ZoBot[lev];
76 std::unique_ptr<MultiFab>& mf_bustr = vec_bustr[lev];
77 std::unique_ptr<MultiFab>& mf_bvstr = vec_bvstr[lev];
78
79 std::unique_ptr<MultiFab>& mf_mskr = vec_mskr[lev];
80 std::unique_ptr<MultiFab>& mf_msku = vec_msku[lev];
81 std::unique_ptr<MultiFab>& mf_mskv = vec_mskv[lev];
82 std::unique_ptr<MultiFab>& mf_mskp = vec_mskp[lev];
83
84 std::unique_ptr<MultiFab>& mf_visc2_p = vec_visc2_p[lev];
85 std::unique_ptr<MultiFab>& mf_visc2_r = vec_visc2_r[lev];
86
87 // We need to set these because otherwise in the first call to remora_advance we may
88 // read uninitialized data on ghost values in setting the bc's on the velocities
89 mf_rho.setVal(0.e34_rt,IntVect(AMREX_D_DECL(NGROW-1,NGROW-1,0)));
90 mf_rhoS->setVal(0.e34_rt,IntVect(AMREX_D_DECL(NGROW-1,NGROW-1,0)));
91 mf_rhoA->setVal(0.e34_rt,IntVect(AMREX_D_DECL(NGROW-1,NGROW-1,0)));
92 mf_DC.setVal(0);
93
94 FillPatchNoBC(lev, time, *cons_new[lev], cons_new, BdyVars::t);
95 FillPatchNoBC(lev, time, *xvel_new[lev], xvel_new, BdyVars::u);
96 FillPatchNoBC(lev, time, *yvel_new[lev], yvel_new, BdyVars::v);
97
98 mf_rufrc->setVal(0);
99 mf_rvfrc->setVal(0);
100
101 int iic = istep[lev];
102 int ntfirst = 0;
103 if(iic==ntfirst) {
104 MultiFab::Copy(S_new,S_old,0,0,S_new.nComp(),S_new.nGrowVect());
105 MultiFab::Copy(U_new,U_old,0,0,U_new.nComp(),U_new.nGrowVect());
106 MultiFab::Copy(V_new,V_old,0,0,V_new.nComp(),V_new.nGrowVect());
107 MultiFab::Copy(W_new,W_old,0,0,W_new.nComp(),W_new.nGrowVect());
108 }
109
110 // If we're not doing bulk fluxes, set surface momentum fluxes directly.
111 // Otherwise, calculate them from winds, so those need to be set
113 set_smflux(lev);
114 } else {
115 set_wind(lev);
116 }
117
118 auto N = Geom(lev).Domain().size()[2]-1; // Number of vertical "levs" aka, NZ
119
120 const auto dlo = amrex::lbound(Geom(lev).Domain());
121 const auto dhi = amrex::ubound(Geom(lev).Domain());
122
123 for ( MFIter mfi(S_new, TilingIfNotGPU()); mfi.isValid(); ++mfi )
124 {
125 Array4<Real const> const& h = vec_h[lev]->const_array(mfi);
126 Array4<Real const> const& Hz = vec_Hz[lev]->const_array(mfi);
127 Array4<Real > const& Huon = vec_Huon[lev]->array(mfi);
128 Array4<Real > const& Hvom = vec_Hvom[lev]->array(mfi);
129
130 Array4<Real const> const& z_w = mf_z_w->const_array(mfi);
131 Array4<Real const> const& z_r = mf_z_r->const_array(mfi);
132 Array4<Real const> const& uold = U_old.const_array(mfi);
133 Array4<Real const> const& vold = V_old.const_array(mfi);
134 Array4<Real > const& rho = mf_rho.array(mfi);
135 Array4<Real > const& rhoA = mf_rhoA->array(mfi);
136 Array4<Real > const& rhoS = mf_rhoS->array(mfi);
137 Array4<Real > const& bvf = mf_bvf->array(mfi);
138 Array4<Real > const& alpha = (solverChoice.bulk_fluxes) ? vec_alpha[lev]->array(mfi) : Array4<Real>();
139 Array4<Real > const& beta = (solverChoice.bulk_fluxes) ? vec_beta[lev]->array(mfi) : Array4<Real>();
140
141 Array4<Real const> const& pm = mf_pm->const_array(mfi);
142 Array4<Real const> const& pn = mf_pn->const_array(mfi);
143
144 Array4<Real const> const& mskr = mf_mskr->const_array(mfi);
145
146 Box bx = mfi.tilebox();
147 Box gbx1 = mfi.growntilebox(IntVect(NGROW-1,NGROW-1,0));
148 Box gbx2 = mfi.growntilebox(IntVect(NGROW,NGROW,0));
149 Box ugbx2 = mfi.grownnodaltilebox(0,IntVect(NGROW,NGROW,0));
150 Box vgbx2 = mfi.grownnodaltilebox(1,IntVect(NGROW,NGROW,0));
151
152 Box bxD = bx;
153 bxD.makeSlab(2,0);
154 Box gbx1D = gbx1;
155 gbx1D.makeSlab(2,0);
156 Box gbx2D = gbx2;
157 gbx2D.makeSlab(2,0);
158
159 //
160 //-----------------------------------------------------------------------
161 // Compute horizontal mass fluxes, Hz*u/n and Hz*v/m (set_massflux_3d)
162 //-----------------------------------------------------------------------
163 //
164 ParallelFor(ugbx2, [=] AMREX_GPU_DEVICE (int i, int j, int k)
165 {
166 Real on_u = 2.0_rt / (pn(i-1,j,0)+pn(i,j,0));
167 Huon(i,j,k)=0.5_rt*(Hz(i,j,k)+Hz(i-1,j,k))*uold(i,j,k)* on_u;
168 });
169
170 ParallelFor(vgbx2, [=] AMREX_GPU_DEVICE (int i, int j, int k)
171 {
172 Real om_v= 2.0_rt / (pm(i,j-1,0)+pm(i,j,0));
173 Hvom(i,j,k)=0.5_rt*(Hz(i,j,k)+Hz(i,j-1,k))*vold(i,j,k)* om_v;
174 });
175
176 Array4<Real const> const& state_old = S_old.const_array(mfi);
177 rho_eos(gbx2,state_old,rho,rhoA,rhoS,bvf,alpha,beta,Hz,z_w,z_r,h,mskr,N);
178 }
179
180 const Real Cdb_min = solverChoice.Cdb_min;
181 const Real Cdb_max = solverChoice.Cdb_max;
182
184 amrex::Abort("remora.longwave_netcdf_is_net=true requires remora.longwave_down_from_netcdf=true");
185 }
186
188 amrex::Abort("remora.longwave_down=true currently requires remora.longwave_down_from_netcdf=true");
189 }
190
191 MultiFab* lw_ptr = nullptr;
192
194 lw_ptr = vec_longwave_down[lev].get();
196 bulk_fluxes(lev, cons_old[lev],vec_uwind[lev].get(),vec_vwind[lev].get(),
197 vec_Tair[lev].get(),vec_qair[lev].get(),vec_Pair[lev].get(),
198 vec_srflx[lev].get(),
199 lw_ptr,
200 vec_evap[lev].get(),
201 vec_sustr[lev].get(),vec_svstr[lev].get(),vec_stflux[lev].get(),
202 vec_lrflx[lev].get(),vec_lhflx[lev].get(),vec_shflx[lev].get(),N);
203 vec_evap[lev]->FillBoundary(geom[lev].periodicity());
204 }
205
207 for ( MFIter mfi(S_new, TilingIfNotGPU()); mfi.isValid(); ++mfi )
208 {
209 Array4<Real > const& stflx = vec_stflx[lev]->array(mfi);
210 Array4<Real > const& btflx = vec_btflx[lev]->array(mfi);
211 Array4<Real const> const& stflux = vec_stflux[lev]->array(mfi);
212 Array4<Real const> const& btflux = vec_btflux[lev]->array(mfi);
213 Box gbx2 = mfi.growntilebox(IntVect(NGROW,NGROW,0));
214 Box gbx2D = gbx2;
215 gbx2D.makeSlab(2,0);
216 ParallelFor(gbx2D, [=] AMREX_GPU_DEVICE (int i, int j, int ) {
217 stflx(i,j,0,Temp_comp) = stflux(i,j,0,Temp_comp);
218 btflx(i,j,0,Temp_comp) = btflux(i,j,0,Temp_comp);
219 });
220 }
221 }
223 for ( MFIter mfi(S_new, TilingIfNotGPU()); mfi.isValid(); ++mfi )
224 {
225 Array4<Real > const& stflx = vec_stflx[lev]->array(mfi);
226 Array4<Real > const& btflx = vec_btflx[lev]->array(mfi);
227 Array4<Real const> const& stflux = vec_stflux[lev]->const_array(mfi);
228 Array4<Real const> const& salt_old = S_old.const_array(mfi,Salt_comp);
229 Box gbx2 = mfi.growntilebox(IntVect(NGROW,NGROW,0));
230 Box gbx2D = gbx2;
231 gbx2D.makeSlab(2,0);
232 ParallelFor(gbx2D, [=] AMREX_GPU_DEVICE (int i, int j, int ) {
233 stflx(i,j,0,Salt_comp) = stflux(i,j,0,Salt_comp) * salt_old(i,j,N);
234 // The fact that this is btflx on the RHS matches what's in ROMS even though
235 // it's weird -- if it's non-zero, does that mean that it will run away since
236 // it's always getting multiplied by salt?
237 btflx(i,j,0,Salt_comp) = btflx(i,j,0,Salt_comp) * salt_old(i,j,0);
238 });
239 }
240 }
241
242 for ( MFIter mfi(S_new, TilingIfNotGPU()); mfi.isValid(); ++mfi )
243 {
244 Array4<Real > const& bustr = mf_bustr->array(mfi);
245 Array4<Real > const& bvstr = mf_bvstr->array(mfi);
246 Array4<Real const> const& uold = U_old.const_array(mfi);
247 Array4<Real const> const& vold = V_old.const_array(mfi);
248 Array4<Real > const& logdrg_tmp = mf_logdrg_tmp.array(mfi);
249 Array4<Real const> const& z_r = mf_z_r->const_array(mfi);
250 Array4<Real const> const& z_w = mf_z_w->const_array(mfi);
251
252 Box gbx2 = mfi.growntilebox(IntVect(NGROW,NGROW,0));
253 Box gbx2D = gbx2;
254 gbx2D.makeSlab(2,0);
255 Box ubx1 = mfi.grownnodaltilebox(0,IntVect(NGROW-1,NGROW-1,0));
256 Box ubx1D = ubx1;
257 ubx1D.makeSlab(2,0);
258 Box vbx1 = mfi.grownnodaltilebox(1,IntVect(NGROW-1,NGROW-1,0));
259 Box vbx1D = vbx1;
260 vbx1D.makeSlab(2,0);
261 // Set bottom stress as defined in set_vbx.F
263 Array4<Real const> const& rdrag = mf_rdrag->const_array(mfi);
264 ParallelFor(ubx1D, [=] AMREX_GPU_DEVICE (int i, int j, int )
265 {
266 bustr(i,j,0) = 0.5_rt * (rdrag(i-1,j,0)+rdrag(i,j,0))*(uold(i,j,0));
267 });
268 ParallelFor(vbx1D, [=] AMREX_GPU_DEVICE (int i, int j, int )
269 {
270 bvstr(i,j,0) = 0.5_rt * (rdrag(i,j-1,0)+rdrag(i,j,0))*(vold(i,j,0));
271 });
273 Array4<Real const> const& rdrag2 = mf_rdrag2->const_array(mfi);
274 ParallelFor(ubx1D, [=] AMREX_GPU_DEVICE (int i, int j, int )
275 {
276 Real avg_v = 0.25_rt * (vold(i,j,0) + vold(i,j+1,0) + vold(i-1,j,0) + vold(i-1,j+1,0));
277 Real vel_mag = std::sqrt(uold(i,j,0)*uold(i,j,0) + avg_v * avg_v);
278 bustr(i,j,0) = 0.5_rt * (rdrag2(i-1,j,0) + rdrag2(i,j,0)) * uold(i,j,0) * vel_mag;
279 });
280 ParallelFor(vbx1D, [=] AMREX_GPU_DEVICE (int i, int j, int )
281 {
282 Real avg_u = 0.25_rt * (uold(i,j,0) + uold(i+1,j,0) + uold(i,j-1,0) + uold(i+1,j-1,0));
283 Real vel_mag = std::sqrt(avg_u * avg_u + vold(i,j,0) * vold(i,j,0));
284 bvstr(i,j,0) = 0.5_rt * (rdrag2(i,j-1,0) + rdrag2(i,j,0)) * vold(i,j,0) * vel_mag;
285 });
287 Array4<Real const> const& ZoBot = mf_ZoBot->const_array(mfi);
288 ParallelFor(gbx2D, [=] AMREX_GPU_DEVICE (int i, int j, int )
289 {
290 Real logz = 1.0_rt / (std::log((z_r(i,j,0) - z_w(i,j,0)) / ZoBot(i,j,0)));
291 Real cff = vonKar * vonKar * logz * logz;
292 logdrg_tmp(i,j,0) = std::min(Cdb_max,std::max(Cdb_min,cff));
293 });
294 ParallelFor(ubx1D, [=] AMREX_GPU_DEVICE (int i, int j, int )
295 {
296 Real avg_v = 0.25_rt * (vold(i,j,0) + vold(i,j+1,0) + vold(i-1,j,0) + vold(i-1,j+1,0));
297 Real vel_mag = std::sqrt(uold(i,j,0)*uold(i,j,0) + avg_v * avg_v);
298 bustr(i,j,0) = 0.5_rt * (logdrg_tmp(i-1,j,0)+logdrg_tmp(i,j,0)) * uold(i,j,0) * vel_mag;
299 });
300 ParallelFor(vbx1D, [=] AMREX_GPU_DEVICE (int i, int j, int )
301 {
302 Real avg_u = 0.25_rt * (uold(i,j,0) + uold(i+1,j,0) + uold(i,j-1,0) + uold(i+1,j-1,0));
303 Real vel_mag = std::sqrt(avg_u * avg_u + vold(i,j,0) * vold(i,j,0));
304 bvstr(i,j,0) = 0.5_rt * (logdrg_tmp(i,j-1,0) + logdrg_tmp(i,j,0)) * vold(i,j,0) * vel_mag;
305 });
306 }
307 }
308 FillPatch(lev, time, *vec_bustr[lev].get(), GetVecOfPtrs(vec_bustr), BCVars::u2d_simple_bc, BdyVars::null,0,true,false);
309 FillPatch(lev, time, *vec_bvstr[lev].get(), GetVecOfPtrs(vec_bvstr), BCVars::v2d_simple_bc, BdyVars::null,0,true,false);
310
312 // Update Akv if using analytic mixing
314 }
315
317
318 MultiFab mf_W(convert(ba,IntVect(0,0,1)),dm,1,IntVect(NGROW+1,NGROW+1,0));
319 mf_W.setVal(0.0_rt);
320
321
323 const int nnew = 0;
324 prestep(lev, U_old, V_old, U_new, V_new,
325 mf_ru.get(), mf_rv.get(),
326 S_old, S_new, mf_W,
327 mf_DC, mf_z_r, mf_z_w, mf_h, mf_pm, mf_pn,
328 mf_sustr.get(), mf_svstr.get(), mf_bustr.get(), mf_bvstr.get(),
329 mf_msku.get(), mf_mskv.get(),
330 iic, ntfirst, nnew, nstp, nrhs, N, dt_lev);
331 }
332
333 // We use FillBoundary not FillPatch here since mf_W is single-level scratch space
334 mf_W.FillBoundary(geom[lev].periodicity());
335 (*physbcs[lev])(mf_W,*mf_mskr.get(),0,1,mf_W.nGrowVect(),t_new[lev],BCVars::zvel_bc);
336
337#ifdef REMORA_USE_NETCDF
338 // Get u and v climatology if we're going to do nudging
340 u_clim_data_from_file->update_interpolated_to_time(t_new[lev], lev, xvel_new[lev], geom, ref_ratio);
341 v_clim_data_from_file->update_interpolated_to_time(t_new[lev], lev, yvel_new[lev], geom, ref_ratio);
342 }
343#endif
344
345 for ( MFIter mfi(S_old, TilingIfNotGPU()); mfi.isValid(); ++mfi )
346 {
347 Array4<Real const> const& Hz = vec_Hz[lev]->const_array(mfi);
348 Array4<Real> const& Huon = vec_Huon[lev]->array(mfi);
349 Array4<Real> const& Hvom = vec_Hvom[lev]->array(mfi);
350 Array4<Real> const& z_r = (mf_z_r)->array(mfi);
351 Array4<Real> const& z_w = (mf_z_w)->array(mfi);
352 Array4<Real const> const& uold = U_old.const_array(mfi);
353 Array4<Real const> const& vold = V_old.const_array(mfi);
354 Array4<Real> const& u = U_new.array(mfi);
355 Array4<Real> const& v = V_new.array(mfi);
356 Array4<Real> const& rho = (mf_rho).array(mfi);
357 Array4<Real> const& ru = (mf_ru)->array(mfi);
358 Array4<Real> const& rv = (mf_rv)->array(mfi);
359 Array4<Real> const& rufrc = (mf_rufrc)->array(mfi);
360 Array4<Real> const& rvfrc = (mf_rvfrc)->array(mfi);
361 Array4<Real> const& W = (mf_W).array(mfi);
362 Array4<Real> const& sustr = (mf_sustr)->array(mfi);
363 Array4<Real> const& svstr = (mf_svstr)->array(mfi);
364 Array4<Real> const& bustr = (mf_bustr)->array(mfi);
365 Array4<Real> const& bvstr = (mf_bvstr)->array(mfi);
366 Array4<Real> const& visc2_p = (mf_visc2_p)->array(mfi);
367 Array4<Real> const& visc2_r = (mf_visc2_r)->array(mfi);
368
369 Array4<Real const> const& pm = mf_pm->const_array(mfi);
370 Array4<Real const> const& pn = mf_pn->const_array(mfi);
371 Array4<Real const> const& fcor = mf_fcor->const_array(mfi);
372
373 Array4<Real const> const& msku = mf_msku->const_array(mfi);
374 Array4<Real const> const& mskv = mf_mskv->const_array(mfi);
375 Array4<Real const> const& mskp = mf_mskp->const_array(mfi);
376
377 Box bx = mfi.tilebox();
378
379 Box tbxp1 = bx;
380 Box tbxp2 = bx;
381 Box xbx = mfi.nodaltilebox(0);
382 Box ybx = mfi.nodaltilebox(1);
383 Box xbx_adj = mfi.nodaltilebox(0);
384 Box ybx_adj = mfi.nodaltilebox(1);
385
386 auto xbx_lo = lbound(xbx_adj);
387 auto xbx_hi = ubound(xbx_adj);
388 auto ybx_lo = lbound(ybx_adj);
389 auto ybx_hi = ubound(ybx_adj);
390
391 if (xbx_lo.x == dlo.x) {
392 xbx_adj.growLo(0,-1);
393 } else if (xbx_hi.x == dhi.x) {
394 xbx_adj.growHi(0,-1);
395 }
396
397 if (ybx_lo.y == dlo.y) {
398 ybx_adj.growLo(1,-1);
399 } else if (ybx_hi.y == dhi.y) {
400 ybx_adj.growHi(1,-1);
401 }
402
403 Box gbx1 = mfi.growntilebox(IntVect(NGROW-1,NGROW-1,0));
404 Box gbx2 = mfi.growntilebox(IntVect(NGROW,NGROW,0));
405
406 Box utbx = mfi.nodaltilebox(0);
407 Box vtbx = mfi.nodaltilebox(1);
408
409 tbxp1.grow(IntVect(NGROW-1,NGROW-1,0));
410 tbxp2.grow(IntVect(NGROW,NGROW,0));
411
412 Box bxD = bx;
413 bxD.makeSlab(2,0);
414 Box gbx1D = gbx1;
415 gbx1D.makeSlab(2,0);
416 Box gbx2D = gbx2;
417 gbx2D.makeSlab(2,0);
418
419 Box tbxp1D = tbxp1;
420 tbxp1D.makeSlab(2,0);
421 Box tbxp2D = tbxp2;
422 tbxp2D.makeSlab(2,0);
423
424 FArrayBox fab_FC(surroundingNodes(tbxp2,2),1,amrex::The_Async_Arena()); //3D
425 auto FC=fab_FC.array();
426
427 FArrayBox fab_fomn(tbxp2D,1,amrex::The_Async_Arena());
428 auto fomn=fab_fomn.array();
429
431 ParallelFor(tbxp2D, [=] AMREX_GPU_DEVICE (int i, int j, int )
432 {
433 fomn(i,j,0) = fcor(i,j,0)*(1.0_rt/(pm(i,j,0)*pn(i,j,0)));
434 });
435 }
436
437 ParallelFor(gbx2, [=] AMREX_GPU_DEVICE (int i, int j, int k)
438 {
439 FC(i,j,k)=0.0_rt;
440 });
441
442 prsgrd(tbxp1,gbx1,utbx,vtbx,ru,rv,pn,pm,rho,FC,Hz,z_r,z_w,msku,mskv,nrhs,N);
443
444 // Apply mixing to temperature and, if use_salt, salt
445 int ncomp = solverChoice.use_salt ? 2 : 1;
446 Array4<Real> const& s_arr = S_new.array(mfi);
447 Array4<Real> const& s_arr_rhs = S_old.array(mfi);
448 Array4<Real> const& diff2_arr = vec_diff2[lev]->array(mfi);
449
450 t3dmix(bx, s_arr, s_arr_rhs, diff2_arr, Hz, pm, pn, msku, mskv, dt_lev, ncomp);
451
452 Array4<Real> const& diff2_arr_scalar = vec_diff2[lev]->array(mfi,Tracer_comp);
453 t3dmix(bx, S_new.array(mfi,Tracer_comp), S_old.array(mfi,Tracer_comp), diff2_arr_scalar, Hz, pm, pn, msku, mskv, dt_lev, 1);
454
456 //-----------------------------------------------------------------------
457 // coriolis
458 //-----------------------------------------------------------------------
459 //
460 // ru, rv updated
461 // In ROMS, coriolis is the first (un-ifdefed) thing to happen in rhs3d_tile, which gets called after t3dmix
462 coriolis(xbx, ybx, uold, vold, ru, rv, Hz, fomn, nrhs, nrhs);
463 }
464
465#ifdef REMORA_USE_NETCDF
467 Array4<const Real> const& uclim = u_clim_data_from_file->get_interpolated_mf(lev)->const_array(mfi);
468 Array4<const Real> const& vclim = v_clim_data_from_file->get_interpolated_mf(lev)->const_array(mfi);
469 Array4<const Real> const& u_nudg_coeff = vec_nudg_coeff[BdyVars::u][lev]->const_array(mfi);
470 Array4<const Real> const& v_nudg_coeff = vec_nudg_coeff[BdyVars::v][lev]->const_array(mfi);
471 // These boxes are set to match ROMS
472 apply_clim_nudg(xbx_adj, 1, 0, ru, uold, uclim, u_nudg_coeff, Hz, pm, pn);
473 apply_clim_nudg(ybx_adj, 0, 1, rv, vold, vclim, v_nudg_coeff, Hz, pm, pn);
474 }
475#endif
476
477 ////rufrc from 3d is set to ru, then the wind stress (and bottom stress) is added, then the mixing is added
478 //rufrc=ru+sustr*om_u*on_u
479
480 rhs_uv_3d(lev, xbx, ybx, uold, vold, ru, rv, rufrc, rvfrc,
481 sustr, svstr, bustr, bvstr, Huon, Hvom,
482 pm, pn, W, FC, nrhs, N);
483
485 const int nnew = 0;
486 uv3dmix(xbx, ybx, u, v, uold, vold, rufrc, rvfrc, visc2_p, visc2_r, Hz, pm, pn, mskp, nrhs, nnew, dt_lev);
487 }
488 } // MFIter
489
490 int nnew = (iic +1)% 2;
491 nstp = iic % 2;
493 gls_prestep(lev, mf_gls, mf_tke, mf_W, mf_msku.get(), mf_mskv.get(),
494 nstp, nnew, iic, ntfirst, N, dt_lev);
495 }
496 nstp = 0;
497
498 // Commenting out for now, but not sure it's necessary
499 //FillPatch(lev, time, *cons_old[lev], cons_old, BCVars::cons_bc, BdyVars::t);
500 //FillPatch(lev, time, *cons_new[lev], cons_new, BCVars::cons_bc, BdyVars::t);
501 FillPatch(lev, time, *vec_sstore[lev], GetVecOfPtrs(vec_sstore), BCVars::cons_bc, BdyVars::t,0,true,true,0,0,dt_lev,*cons_old[lev]);
502
503 // Don't actually want to apply boundary conditions here
504 vec_Huon[lev]->FillBoundary(geom[lev].periodicity());
505 vec_Hvom[lev]->FillBoundary(geom[lev].periodicity());
506}
constexpr amrex::Real vonKar
#define NGROW
#define Temp_comp
#define Tracer_comp
#define Salt_comp
void prsgrd(const amrex::Box &bx, const amrex::Box &gbx, const amrex::Box &utbx, const amrex::Box &vtbx, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &rho, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const int nrhs, const int N)
Calculate pressure gradient.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_evap
evaporation rate [kg/m^2/s]
Definition REMORA.H:349
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:405
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_btflux
Bottom tracer flux; input arrays.
Definition REMORA.H:344
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:253
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:400
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ZoBot
Bottom roughness length [m], defined at rho points.
Definition REMORA.H:360
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lrflx
longwave radiation
Definition REMORA.H:329
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:241
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vwind
Wind in the v direction, defined at rho-points.
Definition REMORA.H:318
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:389
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_tke
Turbulent kinetic energy.
Definition REMORA.H:448
void rho_eos(const amrex::Box &bx, const amrex::Array4< amrex::Real const > &state, const amrex::Array4< amrex::Real > &rho, const amrex::Array4< amrex::Real > &rhoA, const amrex::Array4< amrex::Real > &rhoS, const amrex::Array4< amrex::Real > &bvf, const amrex::Array4< amrex::Real > &alpha, const amrex::Array4< amrex::Real > &beta, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &mskr, const int N)
Wrapper around equation of state calculation.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflx
Surface tracer flux; working arrays.
Definition REMORA.H:338
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_gls
Turbulent generic length scale.
Definition REMORA.H:450
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:311
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:247
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sstore
additional scratch space for calculations on temp, salt, etc
Definition REMORA.H:428
std::unique_ptr< NCTimeSeries > v_clim_data_from_file
Data container for v-velocity climatology data read from file.
Definition REMORA.H:1131
void rhs_uv_3d(int lev, const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real > &rufrc, const amrex::Array4< amrex::Real > &rvfrc, const amrex::Array4< amrex::Real const > &sustr, const amrex::Array4< amrex::Real const > &svstr, const amrex::Array4< amrex::Real const > &bustr, const amrex::Array4< amrex::Real const > &bvstr, const amrex::Array4< amrex::Real const > &Huon, const amrex::Array4< amrex::Real const > &Hvom, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &W, const amrex::Array4< amrex::Real > &FC, int nrhs, int N)
RHS terms for 3D momentum.
void t3dmix(const amrex::Box &bx, const amrex::Array4< amrex::Real > &state, const amrex::Array4< amrex::Real > &state_rhs, const amrex::Array4< amrex::Real const > &diff2, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Real dt_lev, const int ncomp)
Harmonic diffusivity for tracers.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hz
Width of cells in the vertical (z-) direction (3D, Hz in ROMS)
Definition REMORA.H:256
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_msku
land/sea mask at x-faces (2D)
Definition REMORA.H:391
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvfrc
v velocity RHS, integrated, including advection and bottom/surface stresses (2D)
Definition REMORA.H:272
std::unique_ptr< NCTimeSeries > u_clim_data_from_file
Data container for u-velocity climatology data read from file.
Definition REMORA.H:1129
void FillPatchNoBC(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true)
Fill a new MultiFab by copying in phi from valid region and filling ghost cells without applying boun...
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:234
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:245
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_uwind
Wind in the u direction, defined at rho-points.
Definition REMORA.H:316
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rufrc
u velocity RHS, integrated, including advection and bottom/surface stresses (2D)
Definition REMORA.H:270
void gls_prestep(int lev, amrex::MultiFab *mf_gls, amrex::MultiFab *mf_tke, amrex::MultiFab &mf_W, amrex::MultiFab *mf_msku, amrex::MultiFab *mf_mskv, const int nstp, const int nnew, const int iic, const int ntfirst, const int N, const amrex::Real dt_lev)
Prestep for GLS calculation.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_shflx
sensible heat flux
Definition REMORA.H:335
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:278
amrex::Vector< amrex::MultiFab * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
Definition REMORA.H:238
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_r
z coordinates at rho points (cell centers)
Definition REMORA.H:285
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:243
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lhflx
latent heat flux
Definition REMORA.H:333
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskp
land/sea mask at cell corners (2D)
Definition REMORA.H:395
void prestep(int lev, amrex::MultiFab &mf_uold, amrex::MultiFab &mf_vold, amrex::MultiFab &mf_u, amrex::MultiFab &mf_v, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, amrex::MultiFab &S_old, amrex::MultiFab &S_new, amrex::MultiFab &mf_W, amrex::MultiFab &mf_DC, const amrex::MultiFab *mf_z_r, const amrex::MultiFab *mf_z_w, const amrex::MultiFab *mf_h, const amrex::MultiFab *mf_pm, const amrex::MultiFab *mf_pn, const amrex::MultiFab *mf_sustr, const amrex::MultiFab *mf_svstr, const amrex::MultiFab *mf_bustr, const amrex::MultiFab *mf_bvstr, const amrex::MultiFab *mf_msku, const amrex::MultiFab *mf_mskv, const int iic, const int nfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real dt_lev)
Wrapper function for prestep.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvf
Brunt-Vaisala frequency (3D)
Definition REMORA.H:435
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskv
land/sea mask at y-faces (2D)
Definition REMORA.H:393
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
Definition REMORA.H:1221
void FillPatch(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true, const int n_not_fill=0, const int icomp_calc=0, const amrex::Real dt=amrex::Real(0.0), const amrex::MultiFab &mf_calc=amrex::MultiFab())
Fill a new MultiFab by copying in phi from valid region and filling ghost cells.
amrex::Vector< int > istep
which step?
Definition REMORA.H:1207
void apply_clim_nudg(const amrex::Box &bx, int ioff, int joff, const amrex::Array4< amrex::Real > &var, const amrex::Array4< amrex::Real const > &var_old, const amrex::Array4< amrex::Real const > &var_clim, const amrex::Array4< amrex::Real const > &clim_coeff, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Real dt_lev=amrex::Real(0.0))
Apply climatology nudging.
void uv3dmix(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real > &u, const amrex::Array4< amrex::Real > &v, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &rufrc, const amrex::Array4< amrex::Real > &rvfrc, const amrex::Array4< amrex::Real const > &visc2_p, const amrex::Array4< amrex::Real const > &visc2_r, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &mskp, int nrhs, int nnew, const amrex::Real dt_lev)
Harmonic viscosity.
void set_analytic_vmix(int lev)
Set vertical mixing coefficients from analytic.
Definition REMORA.cpp:633
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoS
density perturbation
Definition REMORA.H:431
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_r
Harmonic viscosity defined on the rho points (centers)
Definition REMORA.H:280
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:313
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Huon
u-volume flux (3D)
Definition REMORA.H:258
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:236
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoA
vertically-averaged density
Definition REMORA.H:433
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1211
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1336
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag2
Quadratic drag coefficient [unitless], defined at rho points.
Definition REMORA.H:358
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:262
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_longwave_down
Downward longwave radiation.
Definition REMORA.H:331
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:437
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:232
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bustr
Bottom stress in the u direction.
Definition REMORA.H:363
void set_wind(int lev)
Initialize or calculate wind speed from file or analytic.
Definition REMORA.cpp:702
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvstr
Bottom stress in the v direction.
Definition REMORA.H:365
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > vec_nudg_coeff
Climatology nudging coefficients.
Definition REMORA.H:459
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:402
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflux
Surface tracer flux; input arrays.
Definition REMORA.H:340
void setup_step(int lev, amrex::Real time, amrex::Real dt_lev)
Set everything up for a step on a level.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag
Linear drag coefficient [m/s], defined at rho points.
Definition REMORA.H:356
void set_smflux(int lev)
Initialize or calculate surface momentum flux from file or analytic.
Definition REMORA.cpp:683
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:264
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_btflx
Bottom tracer flux; working arrays.
Definition REMORA.H:342
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_beta
Saline contraction coefficient (3D)
Definition REMORA.H:439
void bulk_fluxes(int lev, amrex::MultiFab *mf_cons, amrex::MultiFab *mf_uwind, amrex::MultiFab *mf_vwind, amrex::MultiFab *mf_Tair, amrex::MultiFab *mf_qair, amrex::MultiFab *mf_Pair, amrex::MultiFab *mf_srflx, amrex::MultiFab *mf_longwave_down, amrex::MultiFab *mf_evap, amrex::MultiFab *mf_sustr, amrex::MultiFab *mf_svstr, amrex::MultiFab *mf_stflux, amrex::MultiFab *mf_lrflx, amrex::MultiFab *mf_lhflx, amrex::MultiFab *mf_shflx, const int N)
Calculate bulk temperature, salinity, wind fluxes.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_srflx
Shortwave radiation flux [W/m²], defined at rho-points.
Definition REMORA.H:327
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Pair
Air pressure [mb], defined at rho-points.
Definition REMORA.H:324
void coriolis(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &fomn, int nrhs, int nr)
Calculate Coriolis terms.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_qair
Specific humidity [kg/kg], defined at rho-points.
Definition REMORA.H:322
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_w
z coordinates at w points (faces between z-cells)
Definition REMORA.H:288
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Tair
Air temperature [°C], defined at rho-points.
Definition REMORA.H:320
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hvom
v-volume flux (3D)
Definition REMORA.H:260
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
Definition REMORA.H:282
amrex::Real Cdb_min
BottomStressType bottom_stress_type
VertMixingType vert_mixing_type
amrex::Real Cdb_max