167 int terrain_smoothing = 0;
168 pp.query(
"terrain_smoothing", terrain_smoothing);
170 if (lev > 0 && terrain_smoothing != 0) {
171 Abort(
"Must use terrain_smoothing = 0 when doing multilevel");
175 int ngrow = z_phys_nd.nGrow();
176 IntVect ngrowVect = z_phys_nd.nGrowVect();
178 const Box& domain = geom.Domain();
179 int domlo_x = domain.smallEnd(0);
int domhi_x = domain.bigEnd(0) + 1;
180 int domlo_y = domain.smallEnd(1);
int domhi_y = domain.bigEnd(1) + 1;
181 int domlo_z = domain.smallEnd(2);
int domhi_z = domain.bigEnd(2) + 1;
189 int nz = z_levels_h.size();
190 Real z_top = z_levels_h[nz-1];
192 Gpu::DeviceVector<Real> z_levels_d;
193 z_levels_d.resize(nz);
194 Gpu::copy(Gpu::hostToDevice, z_levels_h.begin(), z_levels_h.end(), z_levels_d.begin());
196 switch(terrain_smoothing) {
199 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
202 const Box& bx = mfi.validbox();
204 int k0 = bx.smallEnd()[2];
207 Box gbx = mfi.growntilebox(ngrowVect);
209 if (bx.smallEnd(2) == domlo_z) {
210 gbx.setSmall(2,domlo_z);
219 if (bx.bigEnd(2) == domhi_z) {
220 gbx.setBig(2,domhi_z);
223 if (gbx.bigEnd(2) > domhi_z) {
224 gbx.setBig(2,domhi_z);
228 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
229 auto const& z_lev = z_levels_d.data();
235 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
237 int ii = amrex::max(amrex::min(i,imax),imin);
238 int jj = amrex::max(amrex::min(j,jmax),jmin);
246 Real z_sfc = z_arr(ii,jj,k0);
247 Real z_lev_sfc = z_lev[k0];
249 z_arr(i,j,k) = ( (z_sfc - z_lev_sfc) * z_top +
250 (z_top - z_sfc ) *
z ) / (z_top - z_lev_sfc);
254 z_phys_nd.FillBoundary(geom.periodicity());
256 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
259 const Box& bx = mfi.validbox();
260 Box gbx = mfi.growntilebox(ngrowVect);
262 int k0 = bx.smallEnd()[2];
266 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
269 ParallelFor(makeSlab(gbx,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
271 z_arr(i,j,-1) =
two*z_arr(i,j,0) - z_arr(i,j,1);
281 MultiFab h_mf( z_phys_nd.boxArray(), z_phys_nd.DistributionMap(), 1, ngrow+1);
282 MultiFab h_mf_old(z_phys_nd.boxArray(), z_phys_nd.DistributionMap(), 1, ngrow+1);
290 BoxList bl2d = h_mf.boxArray().boxList();
291 for (
auto& b : bl2d) { b.setRange(2,b.smallEnd(2)); }
292 BoxArray ba2d(std::move(bl2d));
293 mf2d = MultiFab(ba2d, h_mf.DistributionMap(), 1, ngrow, MFInfo().SetAlloc(
false));
297 MultiArray4<Real>
const& ma_h_s = h_mf.arrays();
298 MultiArray4<Real>
const& ma_h_s_old = h_mf_old.arrays();
299 MultiArray4<Real>
const& ma_z_phys = z_phys_nd.arrays();
305 h_m = ParReduce(TypeList<ReduceOpMax>{}, TypeList<Real>{}, mf2d, IntVect(ngrow,ngrow,0),
306 [=] AMREX_GPU_DEVICE (
int box_no,
int i,
int j,
int) noexcept
310 const auto & h = ma_h_s[box_no];
311 const auto & z_arr = ma_z_phys[box_no];
313 int ii = amrex::max(amrex::min(i,imax),imin);
314 int jj = amrex::max(amrex::min(j,jmax),jmin);
317 z_arr(i,j,k0) = z_arr(ii,jj,k0);
320 h(i,j,k0) = z_arr(i,j,k0);
323 return { z_arr(i,j,k0) };
325 amrex::ParallelDescriptor::ReduceRealMax(h_m);
330 h_mf.FillBoundary(geom.periodicity());
333 MultiFab::Copy(h_mf_old, h_mf,0,0,1,h_mf_old.nGrow());
337 pp.query(
"terrain_gamma_m", gamma_m);
338 Real z_H =
Real(2.44)*h_m/(1-gamma_m);
341 for (
int k = domlo_z+1; k <= domhi_z; k++)
343 auto const& z_lev_h = z_levels_h.data();
345 Real zz = z_lev_h[k];
346 Real zz_minus = z_lev_h[k-1];
350 Real foo = cos((
PI/2)*(zz/z_H));
351 if(zz < z_H) { A = foo*foo*foo*foo*foo*foo; }
353 Real foo_minus = cos((
PI/2)*(zz_minus/z_H));
355 if(zz_minus < z_H) { A_minus = foo_minus*foo_minus*foo_minus*foo_minus*foo_minus*foo_minus; }
356 else { A_minus = 0; }
358 unsigned maxIter = 50;
360 Real threshold = gamma_m;
362 while (iter < maxIter && diff > threshold)
365 diff = ParReduce(TypeList<ReduceOpMin>{}, TypeList<Real>{}, mf2d, IntVect(ngrow,ngrow,0),
366 [=] AMREX_GPU_DEVICE (
int box_no,
int i,
int j,
int) noexcept
369 const auto & h_s = ma_h_s[box_no];
370 const auto & h_s_old = ma_h_s_old[box_no];
375 int ii = amrex::min(amrex::max(i,domlo_x),domhi_x);
376 int jj = amrex::min(amrex::max(j,domlo_y),domhi_y);
379 h_s(i,j,k) = h_s_old(i,j,k-1) + beta_k*(h_s_old(ii+1,jj ,k-1)
380 + h_s_old(ii-1,jj ,k-1)
381 + h_s_old(ii ,jj+1,k-1)
382 + h_s_old(ii ,jj-1,k-1) - 4*h_s_old(ii,jj,k-1));
385 h_s(i,j,k) = h_s_old(i,j,k ) + beta_k*(h_s_old(ii+1,jj ,k )
386 + h_s_old(ii-1,jj ,k )
387 + h_s_old(ii ,jj+1,k )
388 + h_s_old(ii ,jj-1,k ) - 4*h_s_old(ii,jj,k ));
392 return { (zz + A * h_s(i,j,k) - (zz_minus + A_minus * h_s(i,j,k-1))) / (zz - zz_minus) };
396 MultiFab::Copy(h_mf_old, h_mf,0,0,1,h_mf_old.nGrow());
398 ParallelDescriptor::ReduceRealMin(diff);
403 h_mf_old.FillBoundary(geom.periodicity());
407 auto const& z_lev_d = z_levels_d.data();
410 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
413 Box xybx = mfi.growntilebox(ngrow);
416 Array4<Real>
const& h_s = h_mf_old.array(mfi);
417 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
419 ParallelFor(xybx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int) {
425 z_arr(i,j,k) =
z + A*h_s(i,j,k);
429 z_arr(i,j,k0-1) =
two*z_arr(i,j,k0) - z_arr(i,j,k);
435 Gpu::streamSynchronize();
444 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
447 Box gbx = mfi.growntilebox(ngrow);
448 gbx.setRange(2,domlo_z,domhi_z+1);
450 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
451 auto const& z_lev = z_levels_d.data();
453 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
458 int ii = amrex::max(amrex::min(i,imax),imin);
459 int jj = amrex::max(amrex::min(j,jmax),jmin);
463 z_arr(i,j,k) =
z + (std::pow((
one - (z/z_top)),
omega) * z_arr(ii,jj,k0));
467 z_arr(i,j,k0 ) = z_arr(ii,jj,k0);
470 z_arr(i,j,k0-1) =
two*z_arr(ii,jj,k0) - z_arr(i,j,k);
481 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
484 Box gbx = mfi.growntilebox(ngrow);
485 gbx.setRange(2,domlo_z,domhi_z+1);
487 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
488 auto const& z_lev = z_levels_d.data();
490 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
495 int ii = amrex::max(amrex::min(i,imax),imin);
496 int jj = amrex::max(amrex::min(j,jmax),jmin);
500 z_arr(i,j,k) =
z + z_arr(ii,jj,k0);
504 z_arr(i,j,k) =
z + (std::pow((
one - (z/z_top)),
omega) * z_arr(ii,jj,k0));
508 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
510 z_arr(i,j,k ) =
zero;
511 z_arr(i,j,k-1) = -z_arr(i,j,k+1);
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real PI
Definition: ERF_Constants.H:16
@ omega
Definition: ERF_Morrison.H:53