166 int terrain_smoothing = 0;
167 pp.query(
"terrain_smoothing", terrain_smoothing);
169 if (lev > 0 && terrain_smoothing != 0) {
170 Abort(
"Must use terrain_smoothing = 0 when doing multilevel");
174 int ngrow = z_phys_nd.nGrow();
175 IntVect ngrowVect = z_phys_nd.nGrowVect();
177 const Box& domain = geom.Domain();
178 int domlo_x = domain.smallEnd(0);
int domhi_x = domain.bigEnd(0) + 1;
179 int domlo_y = domain.smallEnd(1);
int domhi_y = domain.bigEnd(1) + 1;
180 int domlo_z = domain.smallEnd(2);
int domhi_z = domain.bigEnd(2) + 1;
188 int nz = z_levels_h.size();
189 Real z_top = z_levels_h[nz-1];
191 Gpu::DeviceVector<Real> z_levels_d;
192 z_levels_d.resize(nz);
193 Gpu::copy(Gpu::hostToDevice, z_levels_h.begin(), z_levels_h.end(), z_levels_d.begin());
195 switch(terrain_smoothing) {
198 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
201 const Box& bx = mfi.validbox();
203 int k0 = bx.smallEnd()[2];
206 Box gbx = mfi.growntilebox(ngrowVect);
208 if (bx.smallEnd(2) == domlo_z) {
209 gbx.setSmall(2,domlo_z);
218 if (bx.bigEnd(2) == domhi_z) {
219 gbx.setBig(2,domhi_z);
222 if (gbx.bigEnd(2) > domhi_z) {
223 gbx.setBig(2,domhi_z);
227 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
228 auto const& z_lev = z_levels_d.data();
234 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
236 int ii = amrex::max(amrex::min(i,imax),imin);
237 int jj = amrex::max(amrex::min(j,jmax),jmin);
245 Real z_sfc = z_arr(ii,jj,k0);
246 Real z_lev_sfc = z_lev[k0];
248 z_arr(i,j,k) = ( (z_sfc - z_lev_sfc) * z_top +
249 (z_top - z_sfc ) *
z ) / (z_top - z_lev_sfc);
253 z_phys_nd.FillBoundary(geom.periodicity());
255 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
258 const Box& bx = mfi.validbox();
259 Box gbx = mfi.growntilebox(ngrowVect);
261 int k0 = bx.smallEnd()[2];
265 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
268 ParallelFor(makeSlab(gbx,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
270 z_arr(i,j,-1) = 2.0*z_arr(i,j,0) - z_arr(i,j,1);
280 MultiFab h_mf( z_phys_nd.boxArray(), z_phys_nd.DistributionMap(), 1, ngrow+1);
281 MultiFab h_mf_old(z_phys_nd.boxArray(), z_phys_nd.DistributionMap(), 1, ngrow+1);
289 BoxList bl2d = h_mf.boxArray().boxList();
290 for (
auto& b : bl2d) {
293 BoxArray ba2d(std::move(bl2d));
294 mf2d = MultiFab(ba2d, h_mf.DistributionMap(), 1, ngrow, MFInfo().SetAlloc(
false));
298 MultiArray4<Real>
const& ma_h_s = h_mf.arrays();
299 MultiArray4<Real>
const& ma_h_s_old = h_mf_old.arrays();
300 MultiArray4<Real>
const& ma_z_phys = z_phys_nd.arrays();
306 h_m = ParReduce(TypeList<ReduceOpMax>{}, TypeList<Real>{}, mf2d, IntVect(ngrow,ngrow,0),
307 [=] AMREX_GPU_DEVICE (
int box_no,
int i,
int j,
int) noexcept
311 const auto & h = ma_h_s[box_no];
312 const auto & z_arr = ma_z_phys[box_no];
314 int ii = amrex::max(amrex::min(i,imax),imin);
315 int jj = amrex::max(amrex::min(j,jmax),jmin);
318 z_arr(i,j,k0) = z_arr(ii,jj,k0);
321 h(i,j,k0) = z_arr(i,j,k0);
324 return { z_arr(i,j,k0) };
326 amrex::ParallelDescriptor::ReduceRealMax(h_m);
331 h_mf.FillBoundary(geom.periodicity());
334 MultiFab::Copy(h_mf_old, h_mf,0,0,1,h_mf_old.nGrow());
338 pp.query(
"terrain_gamma_m", gamma_m);
339 Real z_H = 2.44*h_m/(1-gamma_m);
342 for (
int k = domlo_z+1; k <= domhi_z; k++)
344 auto const& z_lev_h = z_levels_h.data();
346 Real zz = z_lev_h[k];
347 Real zz_minus = z_lev_h[k-1];
351 Real foo = cos((
PI/2)*(zz/z_H));
352 if(zz < z_H) { A = foo*foo*foo*foo*foo*foo; }
354 Real foo_minus = cos((
PI/2)*(zz_minus/z_H));
356 if(zz_minus < z_H) { A_minus = foo_minus*foo_minus*foo_minus*foo_minus*foo_minus*foo_minus; }
357 else { A_minus = 0; }
359 unsigned maxIter = 50;
361 Real threshold = gamma_m;
363 while (iter < maxIter && diff > threshold)
366 diff = ParReduce(TypeList<ReduceOpMin>{}, TypeList<Real>{}, mf2d, IntVect(ngrow,ngrow,0),
367 [=] AMREX_GPU_DEVICE (
int box_no,
int i,
int j,
int) noexcept
370 const auto & h_s = ma_h_s[box_no];
371 const auto & h_s_old = ma_h_s_old[box_no];
373 Real beta_k = 0.2*std::min(zz/(2*h_m),1.0);
376 int ii = amrex::min(amrex::max(i,domlo_x),domhi_x);
377 int jj = amrex::min(amrex::max(j,domlo_y),domhi_y);
380 h_s(i,j,k) = h_s_old(i,j,k-1) + beta_k*(h_s_old(ii+1,jj ,k-1)
381 + h_s_old(ii-1,jj ,k-1)
382 + h_s_old(ii ,jj+1,k-1)
383 + h_s_old(ii ,jj-1,k-1) - 4*h_s_old(ii,jj,k-1));
386 h_s(i,j,k) = h_s_old(i,j,k ) + beta_k*(h_s_old(ii+1,jj ,k )
387 + h_s_old(ii-1,jj ,k )
388 + h_s_old(ii ,jj+1,k )
389 + h_s_old(ii ,jj-1,k ) - 4*h_s_old(ii,jj,k ));
393 return { (zz + A * h_s(i,j,k) - (zz_minus + A_minus * h_s(i,j,k-1))) / (zz - zz_minus) };
397 MultiFab::Copy(h_mf_old, h_mf,0,0,1,h_mf_old.nGrow());
399 ParallelDescriptor::ReduceRealMin(diff);
404 h_mf_old.FillBoundary(geom.periodicity());
408 auto const& z_lev_d = z_levels_d.data();
411 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
414 Box xybx = mfi.growntilebox(ngrow);
417 Array4<Real>
const& h_s = h_mf_old.array(mfi);
418 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
420 ParallelFor(xybx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int) {
426 z_arr(i,j,k) =
z + A*h_s(i,j,k);
430 z_arr(i,j,k0-1) = 2.0*z_arr(i,j,k0) - z_arr(i,j,k);
436 Gpu::streamSynchronize();
445 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
448 Box gbx = mfi.growntilebox(ngrow);
449 gbx.setRange(2,domlo_z,domhi_z+1);
451 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
452 auto const& z_lev = z_levels_d.data();
454 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
459 int ii = amrex::max(amrex::min(i,imax),imin);
460 int jj = amrex::max(amrex::min(j,jmax),jmin);
464 z_arr(i,j,k) =
z + (std::pow((1. - (z/z_top)),
omega) * z_arr(ii,jj,k0));
468 z_arr(i,j,k0 ) = z_arr(ii,jj,k0);
471 z_arr(i,j,k0-1) = 2.0*z_arr(ii,jj,k0) - z_arr(i,j,k);
482 for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
485 Box gbx = mfi.growntilebox(ngrow);
486 gbx.setRange(2,domlo_z,domhi_z+1);
488 Array4<Real>
const& z_arr = z_phys_nd.array(mfi);
489 auto const& z_lev = z_levels_d.data();
491 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
496 int ii = amrex::max(amrex::min(i,imax),imin);
497 int jj = amrex::max(amrex::min(j,jmax),jmin);
501 z_arr(i,j,k) =
z + z_arr(ii,jj,k0);
505 z_arr(i,j,k) =
z + (std::pow((1. - (z/z_top)),
omega) * z_arr(ii,jj,k0));
509 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
512 z_arr(i,j,k-1) = -z_arr(i,j,k+1);
constexpr amrex::Real PI
Definition: ERF_Constants.H:6
@ omega
Definition: ERF_Morrison.H:53