Function for computing the slow RHS for the evolution equations for the density, potential temperature and momentum.
69 BL_PROFILE_REGION(
"erf_make_mom_sources()");
71 Box domain(geom.Domain());
72 const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
106 if (solverChoice.
terrain_type == TerrainType::ImmersedForcing) {
108 amrex::Error(
" Currently forest canopy cannot be used with immersed forcing");
118 auto cosphi = solverChoice.
cosphi;
119 auto sinphi = solverChoice.
sinphi;
151 Real rhoUA_target{0};
152 Real rhoVA_target{0};
159 Table1D<Real> dptr_r_plane, dptr_u_plane, dptr_v_plane;
160 TableData<Real, 1> r_plane_tab, u_plane_tab, v_plane_tab;
163 enforce_massflux_x || enforce_massflux_y))
166 const int u_offset = 1;
167 const int v_offset = 1;
176 r_ave.compute_averages(
ZDir(), r_ave.field());
178 int ncell = r_ave.ncell_line();
179 Gpu::HostVector< Real> r_plane_h(ncell);
180 Gpu::DeviceVector< Real> r_plane_d(ncell);
182 r_ave.line_average(
Rho_comp, r_plane_h);
184 Gpu::copyAsync(Gpu::hostToDevice, r_plane_h.begin(), r_plane_h.end(), r_plane_d.begin());
186 Real* dptr_r = r_plane_d.data();
188 Box tdomain = domain; tdomain.grow(2,ng_c[2]);
189 r_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
191 dptr_r_plane = r_plane_tab.table();
192 ParallelFor(ncell, [=] AMREX_GPU_DEVICE (
int k) noexcept
194 dptr_r_plane(k-
offset) = dptr_r[k];
198 IntVect ng_u = S_data[
IntVars::xmom].nGrowVect(); ng_u[2] = u_offset;
201 IntVect ng_v = S_data[
IntVars::ymom].nGrowVect(); ng_v[2] = v_offset;
204 u_ave.compute_averages(
ZDir(), u_ave.field());
205 v_ave.compute_averages(
ZDir(), v_ave.field());
207 int u_ncell = u_ave.ncell_line();
208 int v_ncell = v_ave.ncell_line();
209 Gpu::HostVector< Real> u_plane_h(u_ncell), v_plane_h(v_ncell);
210 Gpu::DeviceVector< Real> u_plane_d(u_ncell), v_plane_d(v_ncell);
212 u_ave.line_average(0, u_plane_h);
213 v_ave.line_average(0, v_plane_h);
215 Gpu::copyAsync(Gpu::hostToDevice, u_plane_h.begin(), u_plane_h.end(), u_plane_d.begin());
216 Gpu::copyAsync(Gpu::hostToDevice, v_plane_h.begin(), v_plane_h.end(), v_plane_d.begin());
218 Real* dptr_u = u_plane_d.data();
219 Real* dptr_v = v_plane_d.data();
221 Box udomain = domain; udomain.grow(2,ng_u[2]);
222 Box vdomain = domain; vdomain.grow(2,ng_v[2]);
223 u_plane_tab.resize({udomain.smallEnd(2)}, {udomain.bigEnd(2)});
224 v_plane_tab.resize({vdomain.smallEnd(2)}, {vdomain.bigEnd(2)});
226 dptr_u_plane = u_plane_tab.table();
227 ParallelFor(u_ncell, [=] AMREX_GPU_DEVICE (
int k) noexcept
229 dptr_u_plane(k-u_offset) = dptr_u[k];
232 dptr_v_plane = v_plane_tab.table();
233 ParallelFor(v_ncell, [=] AMREX_GPU_DEVICE (
int k) noexcept
235 dptr_v_plane(k-v_offset) = dptr_v[k];
239 if (enforce_massflux_x || enforce_massflux_y) {
240 Real Lx = geom.ProbHi(0) - geom.ProbLo(0);
241 Real Ly = geom.ProbHi(1) - geom.ProbLo(1);
243 if (solverChoice.
mesh_type == MeshType::ConstantDz) {
245 rhoUA = std::accumulate(u_plane_h.begin() + u_offset + massflux_klo,
246 u_plane_h.begin() + u_offset + massflux_khi+1, 0.0);
247 rhoVA = std::accumulate(v_plane_h.begin() + v_offset + massflux_klo,
248 v_plane_h.begin() + v_offset + massflux_khi+1, 0.0);
249 rhoUA_target = std::accumulate(r_plane_h.begin() +
offset + massflux_klo,
250 r_plane_h.begin() +
offset + massflux_khi+1, 0.0);
251 rhoVA_target = rhoUA_target;
253 rhoUA *= geom.CellSize(2) * Ly;
254 rhoVA *= geom.CellSize(2) * Lx;
255 rhoUA_target *= geom.CellSize(2) * Ly;
256 rhoVA_target *= geom.CellSize(2) * Lx;
258 }
else if (solverChoice.
mesh_type == MeshType::StretchedDz) {
260 for (
int k=massflux_klo; k < massflux_khi; ++k) {
261 rhoUA += u_plane_h[k + u_offset] * stretched_dz_h[k];
262 rhoVA += v_plane_h[k + v_offset] * stretched_dz_h[k];
263 rhoUA_target += r_plane_h[k +
offset] * stretched_dz_h[k];
265 rhoVA_target = rhoUA_target;
274 rhoUA_target *= U_target;
275 rhoVA_target *= V_target;
277 Print() <<
"Integrated mass flux : " << rhoUA <<
" " << rhoVA
278 <<
" (target: " << rhoUA_target <<
" " << rhoVA_target <<
")"
286 for ( MFIter mfi(S_data[
IntVars::cons]); mfi.isValid(); ++mfi)
288 Box tbx = mfi.nodaltilebox(0);
289 Box tby = mfi.nodaltilebox(1);
290 Box tbz = mfi.nodaltilebox(2);
291 if (tbz.bigEnd(2) == domain.bigEnd(2)+1) tbz.growHi(2,-1);
293 const Array4<const Real>& cell_data = S_data[
IntVars::cons].array(mfi);
294 const Array4<const Real>& rho_u = S_data[
IntVars::xmom].array(mfi);
295 const Array4<const Real>& rho_v = S_data[
IntVars::ymom].array(mfi);
296 const Array4<const Real>& rho_w = S_data[
IntVars::zmom].array(mfi);
298 const Array4<const Real>& u =
xvel.array(mfi);
299 const Array4<const Real>& v =
yvel.array(mfi);
300 const Array4<const Real>& w = wvel.array(mfi);
302 const Array4< Real>& xmom_src_arr = xmom_src.array(mfi);
303 const Array4< Real>& ymom_src_arr = ymom_src.array(mfi);
304 const Array4< Real>& zmom_src_arr = zmom_src.array(mfi);
306 const Array4<const Real>& r0 = r_hse.const_array(mfi);
308 const Array4<const Real>& f_drag_arr = (forest_drag) ? forest_drag->const_array(mfi) :
309 Array4<const Real>{};
310 const Array4<const Real>& t_blank_arr = (terrain_blank) ? terrain_blank->const_array(mfi) :
311 Array4<const Real>{};
313 const Array4<const Real>& cphi_arr = (cosPhi_mf) ? cosPhi_mf->const_array(mfi) :
314 Array4<const Real>{};
315 const Array4<const Real>& sphi_arr = (sinPhi_mf) ? sinPhi_mf->const_array(mfi) :
316 Array4<const Real>{};
318 const Array4<const Real>& z_nd_arr = z_phys_nd->const_array(mfi);
319 const Array4<const Real>& z_cc_arr = z_phys_cc->const_array(mfi);
325 if (use_coriolis && is_slow_step) {
326 if(solverChoice.
init_type == InitType::HindCast) {
327 const Array4<const Real>& latlon_arr = (*forecast_state_at_lev)[4].array(mfi);
328 ParallelFor(tbx, tby, tbz,
329 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
331 Real rho_v_loc = 0.25 * (rho_v(i,j+1,k) + rho_v(i,j,k) + rho_v(i-1,j+1,k) + rho_v(i-1,j,k));
332 Real rho_w_loc = 0.25 * (rho_w(i,j,k+1) + rho_w(i,j,k) + rho_w(i,j-1,k+1) + rho_w(i,j-1,k));
333 Real latitude = latlon_arr(i,j,k,0);
334 Real sphi_loc = std::sin(latitude*
PI/180.0);
335 Real cphi_loc = std::cos(latitude*
PI/180.0);
336 xmom_src_arr(i, j, k) += coriolis_factor * (rho_v_loc * sphi_loc - rho_w_loc * cphi_loc);
338 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
339 Real rho_u_loc = 0.25 * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j-1,k) + rho_u(i,j-1,k));
340 Real latitude = latlon_arr(i,j,k,0);
341 Real sphi_loc = std::sin(latitude*
PI/180.0);
342 ymom_src_arr(i, j, k) += -coriolis_factor * rho_u_loc * sphi_loc;
344 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
345 Real rho_u_loc = 0.25 * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j,k-1) + rho_u(i,j,k-1));
346 Real latitude = latlon_arr(i,j,k,0);
347 Real cphi_loc = std::cos(latitude*
PI/180.0);
348 zmom_src_arr(i, j, k) += coriolis_factor * rho_u_loc * cphi_loc;
351 else if (var_coriolis && has_lat_lon) {
352 ParallelFor(tbx, tby, tbz,
353 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
355 Real rho_v_loc = 0.25 * (rho_v(i,j+1,k) + rho_v(i,j,k) + rho_v(i-1,j+1,k) + rho_v(i-1,j,k));
356 Real rho_w_loc = 0.25 * (rho_w(i,j,k+1) + rho_w(i,j,k) + rho_w(i-1,j,k+1) + rho_w(i-1,j,k));
357 Real sphi_loc = 0.5 * (sphi_arr(i,j,0) + sphi_arr(i-1,j,0));
358 Real cphi_loc = 0.5 * (cphi_arr(i,j,0) + cphi_arr(i-1,j,0));
359 xmom_src_arr(i, j, k) += coriolis_factor * (rho_v_loc * sphi_loc - rho_w_loc * cphi_loc);
361 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
362 Real rho_u_loc = 0.25 * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j-1,k) + rho_u(i,j-1,k));
363 Real sphi_loc = 0.5 * (sphi_arr(i,j,0) + sphi_arr(i,j-1,0));
364 ymom_src_arr(i, j, k) += -coriolis_factor * rho_u_loc * sphi_loc;
366 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
367 Real rho_u_loc = 0.25 * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j,k-1) + rho_u(i,j,k-1));
368 zmom_src_arr(i, j, k) += coriolis_factor * rho_u_loc * cphi_arr(i,j,0);
371 ParallelFor(tbx, tby, tbz,
372 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
374 Real rho_v_loc = 0.25 * (rho_v(i,j+1,k) + rho_v(i,j,k) + rho_v(i-1,j+1,k) + rho_v(i-1,j,k));
375 Real rho_w_loc = 0.25 * (rho_w(i,j,k+1) + rho_w(i,j,k) + rho_w(i-1,j,k+1) + rho_w(i-1,j,k));
376 xmom_src_arr(i, j, k) += coriolis_factor * (rho_v_loc * sinphi - rho_w_loc * cosphi);
378 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
379 Real rho_u_loc = 0.25 * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j-1,k) + rho_u(i,j-1,k));
380 ymom_src_arr(i, j, k) += -coriolis_factor * rho_u_loc * sinphi;
382 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) {
383 Real rho_u_loc = 0.25 * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j,k-1) + rho_u(i,j,k-1));
384 zmom_src_arr(i, j, k) += coriolis_factor * rho_u_loc * cosphi;
394 if ( (is_slow_step && !use_Rayleigh_fast_uv) || (!is_slow_step && use_Rayleigh_fast_uv)) {
395 if (rayleigh_damp_U) {
396 ParallelFor(tbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
399 Real uu = rho_u(i,j,k) / rho_on_u_face;
400 Real sinesq = d_sinesq_at_lev[k];
401 xmom_src_arr(i, j, k) -= dampcoef*sinesq * (uu -
ubar[k]) * rho_on_u_face;
405 if (rayleigh_damp_V) {
406 ParallelFor(tby, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
409 Real vv = rho_v(i,j,k) / rho_on_v_face;
410 Real sinesq = d_sinesq_at_lev[k];
411 ymom_src_arr(i, j, k) -= dampcoef*sinesq * (vv -
vbar[k]) * rho_on_v_face;
416 if ( (is_slow_step && !use_Rayleigh_fast_w) || (!is_slow_step && use_Rayleigh_fast_w)) {
417 if (rayleigh_damp_W) {
418 ParallelFor(tbz, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
421 Real ww = rho_w(i,j,k) / rho_on_w_face;
422 Real sinesq = d_sinesq_stag_at_lev[k];
423 zmom_src_arr(i, j, k) -= dampcoef*sinesq * (ww -
wbar[k]) * rho_on_w_face;
432 ParallelFor(tbx, tby, tbz,
433 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
436 xmom_src_arr(i, j, k) += rho_on_u_face * abl_geo_forcing[0];
438 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
441 ymom_src_arr(i, j, k) += rho_on_v_face * abl_geo_forcing[1];
443 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
446 zmom_src_arr(i, j, k) += rho_on_w_face * abl_geo_forcing[2];
453 if (geo_wind_profile && is_slow_step) {
454 ParallelFor(tbx, tby,
455 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
458 xmom_src_arr(i, j, k) -= coriolis_factor * rho_on_u_face * dptr_v_geos[k] * sinphi;
460 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
463 ymom_src_arr(i, j, k) += coriolis_factor * rho_on_v_face * dptr_u_geos[k] * sinphi;
473 ParallelFor(tbx, tby,
474 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
476 Real dzInv = 0.5*dxInv[2];
478 Real z_xf_lo = 0.25 * ( z_nd_arr(i,j,k ) + z_nd_arr(i,j+1,k )
479 + z_nd_arr(i,j,k-1) + z_nd_arr(i,j+1,k-1) );
480 Real z_xf_hi = 0.25 * ( z_nd_arr(i,j,k+1) + z_nd_arr(i,j+1,k+1)
481 + z_nd_arr(i,j,k+2) + z_nd_arr(i,j+1,k+2) );
482 dzInv = 1.0 / (z_xf_hi - z_xf_lo);
484 Real rho_on_u_face = 0.5 * ( cell_data(i,j,k,
nr) + cell_data(i-1,j,k,
nr) );
485 Real U_hi = dptr_u_plane(k+1) / dptr_r_plane(k+1);
486 Real U_lo = dptr_u_plane(k-1) / dptr_r_plane(k-1);
487 Real wbar_xf = 0.5 * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
488 xmom_src_arr(i, j, k) -= rho_on_u_face * wbar_xf * (U_hi - U_lo) * dzInv;
490 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
492 Real dzInv = 0.5*dxInv[2];
494 Real z_yf_lo = 0.25 * ( z_nd_arr(i,j,k ) + z_nd_arr(i+1,j,k )
495 + z_nd_arr(i,j,k-1) + z_nd_arr(i+1,j,k-1) );
496 Real z_yf_hi = 0.25 * ( z_nd_arr(i,j,k+1) + z_nd_arr(i+1,j,k+1)
497 + z_nd_arr(i,j,k+2) + z_nd_arr(i+1,j,k+2) );
498 dzInv = 1.0 / (z_yf_hi - z_yf_lo);
500 Real rho_on_v_face = 0.5 * ( cell_data(i,j,k,
nr) + cell_data(i,j-1,k,
nr) );
501 Real V_hi = dptr_v_plane(k+1) / dptr_r_plane(k+1);
502 Real V_lo = dptr_v_plane(k-1) / dptr_r_plane(k-1);
503 Real wbar_yf = 0.5 * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
504 ymom_src_arr(i, j, k) -= rho_on_v_face * wbar_yf * (V_hi - V_lo) * dzInv;
507 ParallelFor(tbx, tby,
508 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
510 Real dzInv = 0.5*dxInv[2];
512 Real z_xf_lo = 0.25 * ( z_nd_arr(i,j,k ) + z_nd_arr(i,j+1,k )
513 + z_nd_arr(i,j,k-1) + z_nd_arr(i,j+1,k-1) );
514 Real z_xf_hi = 0.25 * ( z_nd_arr(i,j,k+1) + z_nd_arr(i,j+1,k+1)
515 + z_nd_arr(i,j,k+2) + z_nd_arr(i,j+1,k+2) );
516 dzInv = 1.0 / (z_xf_hi - z_xf_lo);
518 Real U_hi = dptr_u_plane(k+1) / dptr_r_plane(k+1);
519 Real U_lo = dptr_u_plane(k-1) / dptr_r_plane(k-1);
520 Real wbar_xf = 0.5 * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
521 xmom_src_arr(i, j, k) -= wbar_xf * (U_hi - U_lo) * dzInv;
523 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
525 Real dzInv = 0.5*dxInv[2];
527 Real z_yf_lo = 0.25 * ( z_nd_arr(i,j,k ) + z_nd_arr(i+1,j,k )
528 + z_nd_arr(i,j,k-1) + z_nd_arr(i+1,j,k-1) );
529 Real z_yf_hi = 0.25 * ( z_nd_arr(i,j,k+1) + z_nd_arr(i+1,j,k+1)
530 + z_nd_arr(i,j,k+2) + z_nd_arr(i+1,j,k+2) );
531 dzInv = 1.0 / (z_yf_hi - z_yf_lo);
533 Real V_hi = dptr_v_plane(k+1) / dptr_r_plane(k+1);
534 Real V_lo = dptr_v_plane(k-1) / dptr_r_plane(k-1);
535 Real wbar_yf = 0.5 * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
536 ymom_src_arr(i, j, k) -= wbar_yf * (V_hi - V_lo) * dzInv;
555 for (
int nt = 1; nt < n_sounding_times; nt++) {
558 if (itime_n == n_sounding_times-1) {
561 itime_np1 = itime_n+1;
564 coeff_n =
Real(1.0) - coeff_np1;
569 const Real* u_inp_sound_n = input_sounding_data.
U_inp_sound_d[itime_n].dataPtr();
570 const Real* u_inp_sound_np1 = input_sounding_data.
U_inp_sound_d[itime_np1].dataPtr();
571 const Real* v_inp_sound_n = input_sounding_data.
V_inp_sound_d[itime_n].dataPtr();
572 const Real* v_inp_sound_np1 = input_sounding_data.
V_inp_sound_d[itime_np1].dataPtr();
573 ParallelFor(tbx, tby,
574 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
576 Real nudge_u = (coeff_n*u_inp_sound_n[k] + coeff_np1*u_inp_sound_np1[k]) - (dptr_u_plane(k)/dptr_r_plane(k));
578 xmom_src_arr(i, j, k) += cell_data(i, j, k,
nr) * nudge_u;
580 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
582 Real nudge_v = (coeff_n*v_inp_sound_n[k] + coeff_np1*v_inp_sound_np1[k]) - (dptr_v_plane(k)/dptr_r_plane(k));
584 ymom_src_arr(i, j, k) += cell_data(i, j, k,
nr) * nudge_v;
593 const Array4<const Real>& mf_ux = mapfac[MapFac::ux]->const_array(mfi);
594 const Array4<const Real>& mf_uy = mapfac[MapFac::uy]->const_array(mfi);
595 const Array4<const Real>& mf_vx = mapfac[MapFac::vx]->const_array(mfi);
596 const Array4<const Real>& mf_vy = mapfac[MapFac::vy]->const_array(mfi);
598 u, cell_data, xmom_src_arr, mf_ux, mf_uy);
600 v, cell_data, ymom_src_arr, mf_vx, mf_vy);
611 z_cc_arr, xmom_src_arr, ymom_src_arr,
612 rho_u, rho_v, d_sponge_ptrs_at_lev);
617 xmom_src_arr, ymom_src_arr, zmom_src_arr, rho_u, rho_v, rho_w,
618 r0, z_nd_arr, z_cc_arr);
623 const Array4<const Real>& rho_u_forecast_state = (*forecast_state_at_lev)[
IntVars::xmom].array(mfi);
624 const Array4<const Real>& rho_v_forecast_state = (*forecast_state_at_lev)[
IntVars::ymom].array(mfi);
625 const Array4<const Real>& rho_w_forecast_state = (*forecast_state_at_lev)[
IntVars::zmom].array(mfi);
626 const Array4<const Real>& cons_forecast_state = (*forecast_state_at_lev)[
IntVars::cons].array(mfi);
628 xmom_src_arr, ymom_src_arr, zmom_src_arr,
630 rho_u_forecast_state, rho_v_forecast_state, rho_w_forecast_state,
631 cons_forecast_state);
634 const Array4<const Real>& surface_state_arr = (*surface_state_at_lev).array(mfi);
636 xmom_src_arr, ymom_src_arr,
647 ((is_slow_step && !use_canopy_fast) || (!is_slow_step && use_canopy_fast))) {
648 ParallelFor(tbx, tby, tbz,
649 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
651 const Real ux = u(i, j, k);
652 const Real uy = 0.25 * ( v(i, j , k ) + v(i-1, j , k )
653 + v(i, j+1, k ) + v(i-1, j+1, k ) );
654 const Real uz = 0.25 * ( w(i, j , k ) + w(i-1, j , k )
655 + w(i, j , k+1) + w(i-1, j , k+1) );
656 const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
657 const Real f_drag = 0.5 * (f_drag_arr(i, j, k) + f_drag_arr(i-1, j, k));
658 xmom_src_arr(i, j, k) -= f_drag * ux * windspeed;
660 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
662 const Real ux = 0.25 * ( u(i , j , k ) + u(i , j-1, k )
663 + u(i+1, j , k ) + u(i+1, j-1, k ) );
664 const Real uy = v(i, j, k);
665 const Real uz = 0.25 * ( w(i , j , k ) + w(i , j-1, k )
666 + w(i , j , k+1) + w(i , j-1, k+1) );
667 const amrex::Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
668 const Real f_drag = 0.5 * (f_drag_arr(i, j, k) + f_drag_arr(i, j-1, k));
669 ymom_src_arr(i, j, k) -= f_drag * uy * windspeed;
671 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
673 const amrex::Real ux = 0.25 * ( u(i , j , k ) + u(i+1, j , k )
674 + u(i , j , k-1) + u(i+1, j , k-1) );
675 const amrex::Real uy = 0.25 * ( v(i , j , k ) + v(i , j+1, k )
676 + v(i , j , k-1) + v(i , j+1, k-1) );
678 const amrex::Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
679 const Real f_drag = 0.5 * (f_drag_arr(i, j, k) + f_drag_arr(i, j, k-1));
680 zmom_src_arr(i, j, k) -= f_drag * uz * windspeed;
686 if (solverChoice.
terrain_type == TerrainType::ImmersedForcing &&
687 ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast))) {
689 const Real* dx_arr = geom.CellSize();
690 const Real dx_x = dx_arr[0];
691 const Real dx_y = dx_arr[1];
692 const Real dx_z = dx_arr[2];
695 const Real drag_coefficient = alpha_m / std::pow(dx_x*dx_y*dx_z, 1./3.);
697 const Real U_s = 1.0;
708 ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
710 const Real ux = u(i, j, k);
711 const Real uy = 0.25 * ( v(i, j , k ) + v(i-1, j , k )
712 + v(i, j+1, k ) + v(i-1, j+1, k ) );
713 const Real uz = 0.25 * ( w(i, j , k ) + w(i-1, j , k )
714 + w(i, j , k+1) + w(i-1, j , k+1) );
715 const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
716 const Real t_blank = 0.5 * (t_blank_arr(i, j, k) + t_blank_arr(i-1, j, k));
717 const Real t_blank_above = 0.5 * (t_blank_arr(i, j, k+1) + t_blank_arr(i-1, j, k+1));
718 const Real CdM = std::min(drag_coefficient / (windspeed + tiny), 1000.0);
721 if ((t_blank > 0 && (t_blank_above == 0.0)) && l_use_most) {
723 const Real ux2r = u(i, j, k+1) ;
724 const Real uy2r = 0.25 * ( v(i, j , k+1) + v(i-1, j , k+1)
725 + v(i, j+1, k+1) + v(i-1, j+1, k+1) ) ;
726 const Real h_windspeed2r = std::sqrt(ux2r * ux2r + uy2r * uy2r);
730 const Real rho_xface_below = 0.5 * ( cell_data(i,j,k-1,
Rho_comp) + cell_data(i-1,j,k-1,
Rho_comp) );
732 const Real theta_surf = theta_xface_below;
736 Real ustar = h_windspeed2r * kappa / (std::log(1.5 * dx_z / z0) - psi_m);
737 Real tflux = (tflux_in != 1e-8) ? tflux_in : -(theta_xface - theta_surf) * ustar * kappa / (std::log(1.5 * dx_z / z0) - psi_h);
738 Real Olen = (Olen_in != 1e-8) ? Olen_in : -ustar * ustar * ustar * theta_xface / (kappa * ggg * tflux + tiny);
739 Real zeta = 1.5 * dx_z / Olen;
744 ustar = h_windspeed2r * kappa / (std::log(1.5 * dx_z / z0) - psi_m);
747 if (!(ustar > 0.0 && !std::isnan(ustar))) { ustar = 0.0; }
748 if (!(ustar < 2.0 && !std::isnan(ustar))) { ustar = 2.0; }
749 if (psi_m > std::log(0.5 * dx_z / z0)) { psi_m = std::log(0.5 * dx_z / z0); }
752 const Real uTarget = ustar / kappa * (std::log(0.5 * dx_z / z0) - psi_m);
753 Real uxTarget = uTarget * ux2r / (tiny + h_windspeed2r);
754 const Real bc_forcing_x = -(uxTarget - ux);
755 xmom_src_arr(i, j, k) -= (1-t_blank) * rho_xface * CdM * U_s * bc_forcing_x;
757 xmom_src_arr(i, j, k) -= t_blank * rho_xface * CdM * ux * windspeed;
760 ParallelFor(tby, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
762 const Real ux = 0.25 * ( u(i , j , k ) + u(i , j-1, k )
763 + u(i+1, j , k ) + u(i+1, j-1, k ) );
764 const Real uy = v(i, j, k);
765 const Real uz = 0.25 * ( w(i , j , k ) + w(i , j-1, k )
766 + w(i , j , k+1) + w(i , j-1, k+1) );
767 const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
768 const Real t_blank = 0.5 * (t_blank_arr(i, j, k) + t_blank_arr(i, j-1, k));
769 const Real t_blank_above = 0.5 * (t_blank_arr(i, j, k+1) + t_blank_arr(i, j-1, k+1));
770 const Real CdM = std::min(drag_coefficient / (windspeed + tiny), 1000.0);
773 if ((t_blank > 0 && (t_blank_above == 0.0)) && l_use_most) {
775 const Real ux2r = 0.25 * ( u(i , j , k+1) + u(i , j-1, k+1)
776 + u(i+1, j , k+1) + u(i+1, j-1, k+1) );
777 const Real uy2r = v(i, j, k+1) ;
778 const Real h_windspeed2r = std::sqrt(ux2r * ux2r + uy2r * uy2r);
782 const Real rho_yface_below = 0.5 * ( cell_data(i,j,k-1,
Rho_comp) + cell_data(i,j-1,k-1,
Rho_comp) );
784 const Real theta_surf = theta_yface_below;
788 Real ustar = h_windspeed2r * kappa / (std::log(1.5 * dx_z / z0) - psi_m);
789 Real tflux = (tflux_in != 1e-8) ? tflux_in : -(theta_yface - theta_surf) * ustar * kappa / (std::log(1.5 * dx_z / z0) - psi_h);
790 Real Olen = (Olen_in != 1e-8) ? Olen_in : -ustar * ustar * ustar * theta_yface / (kappa * ggg * tflux + tiny);
791 Real zeta = 1.5 * dx_z / Olen;
796 ustar = h_windspeed2r * kappa / (std::log(1.5 * dx_z / z0) - psi_m);
799 if (!(ustar > 0.0 && !std::isnan(ustar))) { ustar = 0.0; }
800 if (!(ustar < 2.0 && !std::isnan(ustar))) { ustar = 2.0; }
801 if (psi_m > std::log(0.5 * dx_z / z0)) { psi_m = std::log(0.5 * dx_z / z0); }
804 const Real uTarget = ustar / kappa * (std::log(0.5 * dx_z / z0) - psi_m);
805 Real uyTarget = uTarget * uy2r / (tiny + h_windspeed2r);
806 const Real bc_forcing_y = -(uyTarget - uy);
807 ymom_src_arr(i, j, k) -= (1 - t_blank) * rho_yface * CdM * U_s * bc_forcing_y;
809 ymom_src_arr(i, j, k) -= t_blank * rho_yface * CdM * uy * windspeed;
812 ParallelFor(tbz, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
814 const Real ux = 0.25 * ( u(i , j , k ) + u(i+1, j , k )
815 + u(i , j , k-1) + u(i+1, j , k-1) );
816 const Real uy = 0.25 * ( v(i , j , k ) + v(i , j+1, k )
817 + v(i , j , k-1) + v(i , j+1, k-1) );
818 const Real uz = w(i, j, k);
819 const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
820 const Real t_blank = 0.5 * (t_blank_arr(i, j, k) + t_blank_arr(i, j, k-1));
821 const Real CdM = std::min(drag_coefficient / (windspeed + tiny), 1000.0);
823 zmom_src_arr(i, j, k) -= t_blank * rho_zface * CdM * uz * windspeed;
830 if ((solverChoice.
buildings_type == BuildingsType::ImmersedForcing ) &&
831 ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast)))
834 const Real* dx_arr = geom.CellSize();
835 const Real dx_x = dx_arr[0];
836 const Real dx_y = dx_arr[1];
840 const Real min_t_blank = 0.005;
842 ParallelFor(tbx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
844 const Real ux = u(i, j, k);
845 const Real uy = 0.25 * ( v(i, j , k ) + v(i-1, j , k )
846 + v(i, j+1, k ) + v(i-1, j+1, k ) );
847 const Real uz = 0.25 * ( w(i, j , k ) + w(i-1, j , k )
848 + w(i, j , k+1) + w(i-1, j , k+1) );
849 const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
851 Real t_blank = 0.5 * (t_blank_arr(i, j, k) + t_blank_arr(i-1, j, k));
852 if (t_blank < min_t_blank) { t_blank = 0.0; }
853 const Real dx_z = (z_cc_arr) ? (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-1)) : dx_arr[2];
854 const Real drag_coefficient = alpha_m / std::pow(dx_x*dx_y*dx_z, 1./3.);
855 const Real CdM = std::min(drag_coefficient / (windspeed + tiny), 1000.0);
857 xmom_src_arr(i, j, k) -= t_blank * rho_xface * CdM * ux * windspeed;
859 ParallelFor(tby, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
861 const Real ux = 0.25 * ( u(i , j , k ) + u(i , j-1, k )
862 + u(i+1, j , k ) + u(i+1, j-1, k ) );
863 const Real uy = v(i, j, k);
864 const Real uz = 0.25 * ( w(i , j , k ) + w(i , j-1, k )
865 + w(i , j , k+1) + w(i , j-1, k+1) );
866 const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
868 Real t_blank = 0.5 * (t_blank_arr(i, j, k) + t_blank_arr(i, j-1, k));
869 if (t_blank < min_t_blank) { t_blank = 0.0; }
870 const Real dx_z = (z_cc_arr) ? (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-1)) : dx_arr[2];
871 const Real drag_coefficient = alpha_m / std::pow(dx_x*dx_y*dx_z, 1./3.);
872 const Real CdM = std::min(drag_coefficient / (windspeed + tiny), 1000.0);
874 ymom_src_arr(i, j, k) -= t_blank * rho_yface * CdM * uy * windspeed;
876 ParallelFor(tbz, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
878 const Real ux = 0.25 * ( u(i , j , k ) + u(i+1, j , k )
879 + u(i , j , k-1) + u(i+1, j , k-1) );
880 const Real uy = 0.25 * ( v(i , j , k ) + v(i , j+1, k )
881 + v(i , j , k-1) + v(i , j+1, k-1) );
882 const Real uz = w(i, j, k);
883 const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
885 Real t_blank = 0.5 * (t_blank_arr(i, j, k) + t_blank_arr(i, j, k-1));
886 if (t_blank < min_t_blank) { t_blank = 0.0; }
887 const Real dx_z = (z_nd_arr) ? (z_nd_arr(i,j,k) - z_nd_arr(i,j,k-1)) : dx_arr[2];
888 const Real drag_coefficient = alpha_m / std::pow(dx_x*dx_y*dx_z, 1./3.);
889 const Real CdM = std::min(drag_coefficient / (windspeed + tiny), 1000.0);
891 zmom_src_arr(i, j, k) -= t_blank * rho_zface * CdM * uz * windspeed;
898 if (is_slow_step && (enforce_massflux_x || enforce_massflux_y)) {
901 ParallelFor(tbx, tby,
902 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
903 xmom_src_arr(i, j, k) += tau_inv * (rhoUA_target - rhoUA);
905 [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept {
906 ymom_src_arr(i, j, k) += tau_inv * (rhoVA_target - rhoVA);
void ApplyBndryForcing_Forecast(const SolverChoice &solverChoice, const Geometry geom, const Box &tbx, const Box &tby, const Box &tbz, const Array4< const Real > &z_phys_nd, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< Real > &rho_w_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &rho_w, const Array4< const Real > &rho_u_initial_state, const Array4< const Real > &rho_v_initial_state, const Array4< const Real > &rho_w_initial_state, const Array4< const Real > &cons_initial_state)
Definition: ERF_ApplyBndryForcing_Forecast.cpp:8
void ApplySpongeZoneBCsForMom(const SpongeChoice &spongeChoice, const Geometry geom, const Box &tbx, const Box &tby, const Box &tbz, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< Real > &rho_w_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &rho_w, const Array4< const Real > &r0, const Array4< const Real > &z_phys_nd, const Array4< const Real > &z_phys_cc)
Definition: ERF_ApplySpongeZoneBCs.cpp:118
void ApplySpongeZoneBCsForMom_ReadFromFile(const SpongeChoice &spongeChoice, const Geometry geom, const Box &tbx, const Box &tby, const Array4< const Real > &cell_data, const Array4< const Real > &z_phys_cc, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Vector< Real * > d_sponge_ptrs_at_lev)
Definition: ERF_ApplySpongeZoneBCs_ReadFromFile.cpp:8
void ApplySurfaceTreatment_BulkCoeff_Mom(const Box &tbx, const Box &tby, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &cons_state, const Array4< const Real > &z_phys_nd, const Array4< const Real > &surface_state_arr)
Definition: ERF_ApplySurfaceTreatment_BulkCoeff.cpp:8
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:20
constexpr amrex::Real PI
Definition: ERF_Constants.H:6
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
@ ubar
Definition: ERF_DataStruct.H:97
@ wbar
Definition: ERF_DataStruct.H:97
@ vbar
Definition: ERF_DataStruct.H:97
DirectionSelector< 2 > ZDir
Definition: ERF_DirectionSelector.H:38
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
void NumericalDiffusion_Ymom(const Box &bx, const Real dt, const Real num_diff_coeff, const Array4< const Real > &prim_data, const Array4< const Real > &cell_data, const Array4< Real > &rhs, const Array4< const Real > &mfx_arr, const Array4< const Real > &mfy_arr)
Definition: ERF_NumericalDiffusion.cpp:151
void NumericalDiffusion_Xmom(const Box &bx, const Real dt, const Real num_diff_coeff, const Array4< const Real > &prim_data, const Array4< const Real > &cell_data, const Array4< Real > &rhs, const Array4< const Real > &mfx_arr, const Array4< const Real > &mfy_arr)
Definition: ERF_NumericalDiffusion.cpp:85
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_PlaneAverage.H:14
@ r0_comp
Definition: ERF_IndexDefines.H:63
@ ymom
Definition: ERF_IndexDefines.H:160
@ cons
Definition: ERF_IndexDefines.H:158
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159
@ nr
Definition: ERF_Morrison.H:45
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ yvel
Definition: ERF_IndexDefines.H:142
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
bool rayleigh_damp_V
Definition: ERF_DampingStruct.H:85
amrex::Real rayleigh_dampcoef
Definition: ERF_DampingStruct.H:88
bool rayleigh_damp_W
Definition: ERF_DampingStruct.H:86
RayleighDampingType rayleigh_damping_type
Definition: ERF_DampingStruct.H:101
bool rayleigh_damp_U
Definition: ERF_DampingStruct.H:84
bool do_mom_advection
Definition: ERF_DataStruct.H:1141
amrex::Real coriolis_factor
Definition: ERF_DataStruct.H:1132
static MeshType mesh_type
Definition: ERF_DataStruct.H:1048
amrex::Real if_surf_temp_flux
Definition: ERF_DataStruct.H:1112
bool if_use_most
Definition: ERF_DataStruct.H:1116
DampingChoice dampingChoice
Definition: ERF_DataStruct.H:1058
amrex::Real const_massflux_v
Definition: ERF_DataStruct.H:1206
amrex::Real if_z0
Definition: ERF_DataStruct.H:1111
amrex::Real cosphi
Definition: ERF_DataStruct.H:1133
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > abl_geo_forcing
Definition: ERF_DataStruct.H:1179
bool hindcast_lateral_forcing
Definition: ERF_DataStruct.H:1215
int massflux_klo
Definition: ERF_DataStruct.H:1210
bool custom_w_subsidence
Definition: ERF_DataStruct.H:1139
bool nudging_from_input_sounding
Definition: ERF_DataStruct.H:1149
bool immersed_forcing_substep
Definition: ERF_DataStruct.H:1104
amrex::Real sinphi
Definition: ERF_DataStruct.H:1134
bool have_geo_wind_profile
Definition: ERF_DataStruct.H:1181
amrex::Real const_massflux_u
Definition: ERF_DataStruct.H:1205
amrex::Real if_Olen_in
Definition: ERF_DataStruct.H:1115
bool use_coriolis
Definition: ERF_DataStruct.H:1098
amrex::Real num_diff_coeff
Definition: ERF_DataStruct.H:1165
bool variable_coriolis
Definition: ERF_DataStruct.H:1184
amrex::Real if_Cd_momentum
Definition: ERF_DataStruct.H:1109
bool custom_forcing_prim_vars
Definition: ERF_DataStruct.H:1143
static BuildingsType buildings_type
Definition: ERF_DataStruct.H:1039
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1036
SpongeChoice spongeChoice
Definition: ERF_DataStruct.H:1059
static InitType init_type
Definition: ERF_DataStruct.H:1030
bool hindcast_surface_bcs
Definition: ERF_DataStruct.H:1216
bool has_lat_lon
Definition: ERF_DataStruct.H:1183
bool do_forest_drag
Definition: ERF_DataStruct.H:1202
amrex::Real const_massflux_tau
Definition: ERF_DataStruct.H:1207
int massflux_khi
Definition: ERF_DataStruct.H:1211
bool forest_substep
Definition: ERF_DataStruct.H:1105
int ave_plane
Definition: ERF_DataStruct.H:1186
std::string sponge_type
Definition: ERF_SpongeStruct.H:58
Definition: ERF_MOSTStress.H:40
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_m(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:90
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_h(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:104