ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_TerrainMetrics.cpp File Reference
#include <ERF_TerrainMetrics.H>
#include <ERF_Utils.H>
#include <AMReX_ParmParse.H>
#include <ERF_Constants.H>
#include <ERF_Interpolation_1D.H>
#include <cmath>
Include dependency graph for ERF_TerrainMetrics.cpp:

Functions

void init_default_zphys (int, const Geometry &geom, MultiFab &z_phys_nd, MultiFab &z_phys_cc)
 
void make_terrain_fitted_coords (int lev, const Geometry &geom, MultiFab &z_phys_nd, Vector< Real > const &z_levels_h, GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &phys_bc_type)
 
void init_which_terrain_grid (int lev, Geometry const &geom, MultiFab &z_phys_nd, Vector< Real > const &z_levels_h)
 
void make_J (const Geometry &geom, MultiFab &z_phys_nd, MultiFab &detJ_cc)
 
void make_areas (const Geometry &geom, MultiFab &z_phys_nd, MultiFab &ax, MultiFab &ay, MultiFab &az)
 
void make_zcc (const Geometry &geom, MultiFab &z_phys_nd, MultiFab &z_phys_cc)
 
Real get_dzmin_terrain (MultiFab &z_phys_nd)
 

Function Documentation

◆ get_dzmin_terrain()

Real get_dzmin_terrain ( MultiFab &  z_phys_nd)

Computation min dz at cell-center

653 {
654  auto const& ma_z_nd_arr = z_phys_nd.const_arrays();
655  GpuTuple<Real> min = ParReduce(TypeList<ReduceOpMin>{},
656  TypeList<Real>{},
657  z_phys_nd, IntVect(0),
658  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
659  -> GpuTuple<Real>
660  {
661  amrex::Real dz = Compute_Z_AtWFace(i,j,k+1,ma_z_nd_arr[box_no]) -
662  Compute_Z_AtWFace(i,j,k ,ma_z_nd_arr[box_no]);
663  return { dz };
664  });
665  Real r = (get<0>(min) + std::numeric_limits<amrex::Real>::epsilon());
666  ParallelDescriptor::ReduceRealMin(r);
667  return r;
668 }
const Real dz
Definition: ERF_InitCustomPert_Bubble.H:25
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtWFace(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:374
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12

Referenced by ERF::init_zphys(), ERF::ReadCheckpointFile(), and ERF::remake_zphys().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_default_zphys()

void init_default_zphys ( int  ,
const Geometry &  geom,
MultiFab &  z_phys_nd,
MultiFab &  z_phys_cc 
)

Define a default z_phys so we have it even if a completely regular grid This will be over-written if we use z_levels, or grid stretching, or terrain-fitted grids

16 {
17  const auto& dx = geom.CellSize();
18  Real dz = dx[2];
19 
20  for (MFIter mfi(z_phys_nd,true); mfi.isValid(); ++mfi)
21  {
22  const Box& bx = mfi.growntilebox();
23  const Array4< Real> z_nd_arr = z_phys_nd.array(mfi);
24  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
25  {
26  z_nd_arr(i,j,k) = k * dz;
27  });
28  }
29 
30  for (MFIter mfi(z_phys_cc,true); mfi.isValid(); ++mfi)
31  {
32  const Box& bx = mfi.growntilebox();
33  const Array4< Real> z_cc_arr = z_phys_cc.array(mfi);
34  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
35  {
36  z_cc_arr(i,j,k) = (k + 0.5) * dz;
37  });
38  }
39 }
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+0.5) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=1.0) { Real dT=T_pert *(std::cos(PI *L)+1.0)/2.0;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })

Referenced by ERF::init_stuff().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_which_terrain_grid()

void init_which_terrain_grid ( int  lev,
Geometry const &  geom,
MultiFab &  z_phys_nd,
Vector< Real > const &  z_levels_h 
)
163 {
164  // User-selected method from inputs file (BTF default)
165  ParmParse pp("erf");
166  int terrain_smoothing = 0;
167  pp.query("terrain_smoothing", terrain_smoothing);
168 
169  if (lev > 0 && terrain_smoothing != 0) {
170  Abort("Must use terrain_smoothing = 0 when doing multilevel");
171  }
172 
173  // Number of ghost cells
174  int ngrow = z_phys_nd.nGrow();
175  IntVect ngrowVect = z_phys_nd.nGrowVect();
176 
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;
181 
182  int imin = domlo_x; // if (geom.isPeriodic(0)) imin -= z_phys_nd.nGrowVect()[0];
183  int jmin = domlo_y; // if (geom.isPeriodic(1)) jmin -= z_phys_nd.nGrowVect()[1];
184 
185  int imax = domhi_x; // if (geom.isPeriodic(0)) imax += z_phys_nd.nGrowVect()[0];
186  int jmax = domhi_y; // if (geom.isPeriodic(1)) jmax += z_phys_nd.nGrowVect()[1];
187 
188  int nz = z_levels_h.size();
189  Real z_top = z_levels_h[nz-1];
190 
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());
194 
195  switch(terrain_smoothing) {
196  case 0: // BTF Method
197  {
198  for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
199  {
200  // Note that this box is nodal because it is based on z_phys_nd
201  const Box& bx = mfi.validbox();
202 
203  int k0 = bx.smallEnd()[2];
204 
205  // Grown box with corrected ghost cells at top
206  Box gbx = mfi.growntilebox(ngrowVect);
207 
208  if (bx.smallEnd(2) == domlo_z) {
209  gbx.setSmall(2,domlo_z);
210  } else {
211  gbx.growLo(2,-1);
212  }
213 
214  // Note that we don't overwrite the values at the high end of the box
215  // regardless of whether the box reaches the top of the domain or not.
216  // In the case of lev > 0, this ensures that the fine nodes at the top
217  // of a fine box are those that are interpolated from the coarse grid
218  if (bx.bigEnd(2) == domhi_z) {
219  gbx.setBig(2,domhi_z);
220  } else {
221  gbx.growHi(2,-1);
222  if (gbx.bigEnd(2) > domhi_z) {
223  gbx.setBig(2,domhi_z);
224  }
225  }
226 
227  Array4<Real> const& z_arr = z_phys_nd.array(mfi);
228  auto const& z_lev = z_levels_d.data();
229 
230  //
231  // Vertical grid stretching using BTF model from p2163 of Klemp2011
232  // z_levels are only defined from k = dom_lo to dom_hi (nodal)
233  //
234  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
235  {
236  int ii = amrex::max(amrex::min(i,imax),imin);
237  int jj = amrex::max(amrex::min(j,jmax),jmin);
238 
239  //
240  // Start with flat z_lev set either with uniform cell size or specified z_levels
241  // If k0 = 0 then z_arr at k0 has already been filled from the terrain data
242  // If k0 > 0 then z_arr at k0 has already been filled from interpolation
243  //
244  Real z = z_lev[k];
245  Real z_sfc = z_arr(ii,jj,k0);
246  Real z_lev_sfc = z_lev[k0];
247 
248  z_arr(i,j,k) = ( (z_sfc - z_lev_sfc) * z_top +
249  (z_top - z_sfc ) * z ) / (z_top - z_lev_sfc);
250  });
251  } // mfi
252 
253  z_phys_nd.FillBoundary(geom.periodicity());
254 
255  for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
256  {
257  // Note that this box is nodal because it is based on z_phys_nd
258  const Box& bx = mfi.validbox();
259  Box gbx = mfi.growntilebox(ngrowVect);
260 
261  int k0 = bx.smallEnd()[2];
262 
263  if (k0 == 0)
264  {
265  Array4<Real> const& z_arr = z_phys_nd.array(mfi);
266 
267  // Fill lateral boundaries below the bottom surface
268  ParallelFor(makeSlab(gbx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int)
269  {
270  z_arr(i,j,-1) = 2.0*z_arr(i,j,0) - z_arr(i,j,1);
271  });
272  }
273  } // mfi
274  break;
275  } // case 0
276 
277  case 1: // STF Method
278  {
279  // Get MultiFab spanning domain with 1 level of ghost cells
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);
282 
283  // Save max height for smoothing
284  Real h_m;
285 
286  // Create 2D MF without allocation
287  MultiFab mf2d;
288  {
289  BoxList bl2d = h_mf.boxArray().boxList();
290  for (auto& b : bl2d) { b.setRange(2,b.smallEnd(2)); }
291  BoxArray ba2d(std::move(bl2d));
292  mf2d = MultiFab(ba2d, h_mf.DistributionMap(), 1, ngrow, MFInfo().SetAlloc(false));
293  }
294 
295  // Get MultiArray4s from the multifabs
296  MultiArray4<Real> const& ma_h_s = h_mf.arrays();
297  MultiArray4<Real> const& ma_h_s_old = h_mf_old.arrays();
298  MultiArray4<Real> const& ma_z_phys = z_phys_nd.arrays();
299 
300  // Bottom boundary
301  int k0 = domlo_z;
302 
303  // Get max value
304  h_m = ParReduce(TypeList<ReduceOpMax>{}, TypeList<Real>{}, mf2d, IntVect(ngrow,ngrow,0),
305  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int) noexcept
306  -> GpuTuple<Real>
307  {
308  // Get Array4s
309  const auto & h = ma_h_s[box_no];
310  const auto & z_arr = ma_z_phys[box_no];
311 
312  int ii = amrex::max(amrex::min(i,imax),imin);
313  int jj = amrex::max(amrex::min(j,jmax),jmin);
314 
315  // Fill the lateral boundaries
316  z_arr(i,j,k0) = z_arr(ii,jj,k0);
317 
318  // Populate h with terrain
319  h(i,j,k0) = z_arr(i,j,k0);
320 
321  // Return height for max
322  return { z_arr(i,j,k0) };
323  });
324  amrex::ParallelDescriptor::ReduceRealMax(h_m);
325 
326  if (h_m < std::numeric_limits<Real>::epsilon()) h_m = 1e-16;
327 
328  // Fill ghost cells (neglects domain boundary if not periodic)
329  h_mf.FillBoundary(geom.periodicity());
330 
331  // Make h_mf copy for old values
332  MultiFab::Copy(h_mf_old, h_mf,0,0,1,h_mf_old.nGrow());
333 
334  // Minimum allowed fractional grid spacing
335  Real gamma_m = 0.5;
336  pp.query("terrain_gamma_m", gamma_m);
337  Real z_H = 2.44*h_m/(1-gamma_m); // Klemp2011 Eqn. 11
338 
339  // Populate h_mf at k>0 with h_s, solving in ordered 2D slices
340  for (int k = domlo_z+1; k <= domhi_z; k++) // skip terrain level
341  {
342  auto const& z_lev_h = z_levels_h.data();
343 
344  Real zz = z_lev_h[k];
345  Real zz_minus = z_lev_h[k-1];
346 
347  // Hybrid attenuation profile, Klemp2011 Eqn. 9
348  Real A;
349  Real foo = cos((PI/2)*(zz/z_H));
350  if(zz < z_H) { A = foo*foo*foo*foo*foo*foo; } // A controls rate of return to atm
351  else { A = 0; }
352  Real foo_minus = cos((PI/2)*(zz_minus/z_H));
353  Real A_minus;
354  if(zz_minus < z_H) { A_minus = foo_minus*foo_minus*foo_minus*foo_minus*foo_minus*foo_minus; } // A controls rate of return to atm
355  else { A_minus = 0; }
356 
357  unsigned maxIter = 50; // M_k in paper
358  unsigned iter = 0;
359  Real threshold = gamma_m;
360  Real diff = 1.e20;
361  while (iter < maxIter && diff > threshold)
362  {
363 
364  diff = ParReduce(TypeList<ReduceOpMin>{}, TypeList<Real>{}, mf2d, IntVect(ngrow,ngrow,0),
365  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int) noexcept
366  -> GpuTuple<Real>
367  {
368  const auto & h_s = ma_h_s[box_no];
369  const auto & h_s_old = ma_h_s_old[box_no];
370 
371  Real beta_k = 0.2*std::min(zz/(2*h_m),1.0); //smoothing coefficient (Eqn. 8)
372 
373  // Clip indices for ghost-cells
374  int ii = amrex::min(amrex::max(i,domlo_x),domhi_x);
375  int jj = amrex::min(amrex::max(j,domlo_y),domhi_y);
376 
377  if (iter == 0) {
378  h_s(i,j,k) = h_s_old(i,j,k-1) + beta_k*(h_s_old(ii+1,jj ,k-1)
379  + h_s_old(ii-1,jj ,k-1)
380  + h_s_old(ii ,jj+1,k-1)
381  + h_s_old(ii ,jj-1,k-1) - 4*h_s_old(ii,jj,k-1));
382  }
383  else {
384  h_s(i,j,k) = h_s_old(i,j,k ) + beta_k*(h_s_old(ii+1,jj ,k )
385  + h_s_old(ii-1,jj ,k )
386  + h_s_old(ii ,jj+1,k )
387  + h_s_old(ii ,jj-1,k ) - 4*h_s_old(ii,jj,k ));
388  }
389 
390  // Minimum vertical grid spacing condition (Klemp2011 Eqn. 7)
391  return { (zz + A * h_s(i,j,k) - (zz_minus + A_minus * h_s(i,j,k-1))) / (zz - zz_minus) };
392 
393  }); //ParReduce
394 
395  MultiFab::Copy(h_mf_old, h_mf,0,0,1,h_mf_old.nGrow());
396 
397  ParallelDescriptor::ReduceRealMin(diff);
398 
399  iter++;
400 
401  //fill ghost points
402  h_mf_old.FillBoundary(geom.periodicity());
403 
404  } //while
405 
406  auto const& z_lev_d = z_levels_d.data();
407 
408  //Populate z_phys_nd by solving z_arr(i,j,k) = z + A*h_s(i,j,k)
409  for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
410  {
411  // Grown box with no z range
412  Box xybx = mfi.growntilebox(ngrow);
413  xybx.setRange(2,0);
414 
415  Array4<Real> const& h_s = h_mf_old.array(mfi);
416  Array4<Real> const& z_arr = z_phys_nd.array(mfi);
417 
418  ParallelFor(xybx, [=] AMREX_GPU_DEVICE (int i, int j, int) {
419 
420  // Location of nodes
421  Real z = z_lev_d[k];
422 
423  // STF model from p2164 of Klemp2011 (Eqn. 4)
424  z_arr(i,j,k) = z + A*h_s(i,j,k);
425 
426  // Fill below the bottom surface
427  if (k == 1) {
428  z_arr(i,j,k0-1) = 2.0*z_arr(i,j,k0) - z_arr(i,j,k);
429  }
430  });
431  } // mfi
432  } // k
433 
434  Gpu::streamSynchronize();
435 
436  break;
437  } // case 1
438 
439  case 2: // Sullivan TF Method
440  {
441  int k0 = 0;
442 
443  for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
444  {
445  // Grown box with corrected ghost cells at top
446  Box gbx = mfi.growntilebox(ngrow);
447  gbx.setRange(2,domlo_z,domhi_z+1);
448 
449  Array4<Real> const& z_arr = z_phys_nd.array(mfi);
450  auto const& z_lev = z_levels_d.data();
451 
452  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
453  {
454  // Vertical grid stretching
455  Real z = z_lev[k];
456 
457  int ii = amrex::max(amrex::min(i,imax),imin);
458  int jj = amrex::max(amrex::min(j,jmax),jmin);
459 
460  // Fill levels using model from Sullivan et. al. 2014
461  int omega = 3; //Used to adjust how rapidly grid lines level out. omega=1 is BTF!
462  z_arr(i,j,k) = z + (std::pow((1. - (z/z_top)),omega) * z_arr(ii,jj,k0));
463 
464  // Fill lateral boundaries and below the bottom surface
465  if (k == k0) {
466  z_arr(i,j,k0 ) = z_arr(ii,jj,k0);
467  }
468  if (k == 1) {
469  z_arr(i,j,k0-1) = 2.0*z_arr(ii,jj,k0) - z_arr(i,j,k);
470  }
471  });
472  } // mfi
473  break;
474  } // case 2
475 
476  case 3: // Debugging Test Method -- applies Sullivan TF starting at k = 1 so that domain does not change size
477  {
478  int k0 = 0;
479 
480  for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
481  {
482  // Grown box with corrected ghost cells at top
483  Box gbx = mfi.growntilebox(ngrow);
484  gbx.setRange(2,domlo_z,domhi_z+1);
485 
486  Array4<Real> const& z_arr = z_phys_nd.array(mfi);
487  auto const& z_lev = z_levels_d.data();
488 
489  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
490  {
491  // Vertical grid stretching
492  Real z = z_lev[k];
493 
494  int ii = amrex::max(amrex::min(i,imax),imin);
495  int jj = amrex::max(amrex::min(j,jmax),jmin);
496 
497  // Fill values outside the lateral boundaries and below the bottom surface (necessary if init_type = "WRFInput")
498  if (k == k0+1) {
499  z_arr(i,j,k) = z + z_arr(ii,jj,k0);
500  } else {
501  // Fill levels using model from Sullivan et. al. 2014
502  int omega = 3; //Used to adjust how rapidly grid lines level out. omega=1 is BTF!
503  z_arr(i,j,k) = z + (std::pow((1. - (z/z_top)),omega) * z_arr(ii,jj,k0));
504  }
505  });
506  gbx.setBig(2,0);
507  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
508  {
509  z_arr(i,j,k ) = 0.0;
510  z_arr(i,j,k-1) = -z_arr(i,j,k+1);
511  });
512  } // mfi
513  break;
514  } // case 3
515  } //switch
516 }
constexpr amrex::Real PI
Definition: ERF_Constants.H:6
ParmParse pp("prob")
@ omega
Definition: ERF_Morrison.H:53

Referenced by make_terrain_fitted_coords().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ make_areas()

void make_areas ( const Geometry &  geom,
MultiFab &  z_phys_nd,
MultiFab &  ax,
MultiFab &  ay,
MultiFab &  az 
)

Computation of area fractions on faces

565 {
566  const auto* dx = geom.CellSize();
567  Real dzInv = 1.0/dx[2];
568 
569  // Domain valid box (z_nd is nodal)
570  const Box& domain = geom.Domain();
571  int domlo_z = domain.smallEnd(2);
572 
573  // The z-faces are always full when using terrain-fitted coordinates
574  az.setVal(1.0);
575 
576  //
577  // x-areas
578  //
579 #ifdef _OPENMP
580 #pragma omp parallel if (Gpu::notInLaunchRegion())
581 #endif
582  for ( MFIter mfi(ax, TilingIfNotGPU()); mfi.isValid(); ++mfi )
583  {
584  Box gbx = mfi.growntilebox(ax.nGrow());
585  if (gbx.smallEnd(2) < domlo_z) {
586  gbx.setSmall(2,domlo_z);
587  }
588 
589  Array4<Real const> z_nd = z_phys_nd.const_array(mfi);
590  Array4<Real > ax_arr = ax.array(mfi);
591  ParallelFor(gbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
592  ax_arr(i, j, k) = .5 * dzInv * (
593  z_nd(i,j,k+1) + z_nd(i,j+1,k+1) - z_nd(i,j,k) - z_nd(i,j+1,k));
594  });
595  }
596 
597  //
598  // y-areas
599  //
600 #ifdef _OPENMP
601 #pragma omp parallel if (Gpu::notInLaunchRegion())
602 #endif
603  for ( MFIter mfi(ay, TilingIfNotGPU()); mfi.isValid(); ++mfi )
604  {
605  Box gbx = mfi.growntilebox(ay.nGrow());
606  if (gbx.smallEnd(2) < domlo_z) {
607  gbx.setSmall(2,domlo_z);
608  }
609 
610  Array4<Real const> z_nd = z_phys_nd.const_array(mfi);
611  Array4<Real > ay_arr = ay.array(mfi);
612  ParallelFor(gbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
613  ay_arr(i, j, k) = .5 * dzInv * (
614  z_nd(i,j,k+1) + z_nd(i+1,j,k+1) - z_nd(i,j,k) - z_nd(i+1,j,k));
615  });
616  }
617 
618  ax.FillBoundary(geom.periodicity());
619  ay.FillBoundary(geom.periodicity());
620  az.FillBoundary(geom.periodicity());
621 }

Referenced by ERF::update_terrain_arrays().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ make_J()

void make_J ( const Geometry &  geom,
MultiFab &  z_phys_nd,
MultiFab &  detJ_cc 
)

Computation of detJ at cell-center

525 {
526  const auto *dx = geom.CellSize();
527  Real dzInv = 1.0/dx[2];
528 
529  // Domain valid box (z_nd is nodal)
530  const Box& domain = geom.Domain();
531  int domlo_z = domain.smallEnd(2);
532 
533  // Number of ghost cells
534  int ngrow= detJ_cc.nGrow();
535 
536 #ifdef _OPENMP
537 #pragma omp parallel if (Gpu::notInLaunchRegion())
538 #endif
539  for ( MFIter mfi(detJ_cc, TilingIfNotGPU()); mfi.isValid(); ++mfi )
540  {
541  Box gbx = mfi.growntilebox(ngrow);
542  if (gbx.smallEnd(2) < domlo_z) {
543  gbx.setSmall(2,domlo_z);
544  }
545  Array4<Real const> z_nd = z_phys_nd.const_array(mfi);
546  Array4<Real > detJ = detJ_cc.array(mfi);
547  ParallelFor(gbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
548  detJ(i, j, k) = .25 * dzInv * (
549  z_nd(i,j,k+1) + z_nd(i+1,j,k+1) + z_nd(i,j+1,k+1) + z_nd(i+1,j+1,k+1)
550  -z_nd(i,j,k ) - z_nd(i+1,j,k ) - z_nd(i,j+1,k ) - z_nd(i+1,j+1,k ) );
551  });
552  }
553  detJ_cc.FillBoundary(geom.periodicity());
554 }

Referenced by ERF::update_terrain_arrays().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ make_terrain_fitted_coords()

void make_terrain_fitted_coords ( int  lev,
const Geometry &  geom,
MultiFab &  z_phys_nd,
Vector< Real > const &  z_levels_h,
GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &  phys_bc_type 
)

Computation of the terrain grid from BTF, STF, or Sullivan TF model

49 {
50  const Box& domain = geom.Domain();
51 
52  int domlo_z = domain.smallEnd(2);
53  int domhi_z = domain.bigEnd(2) + 1;
54 
55  // Just in case ...
56  z_phys_nd.setDomainBndry(1.234e20,0,1,geom);
57 
58  // ****************************************************************************
59 
60  if (lev == 0) {
61  BoxArray ba(z_phys_nd.boxArray());
62  bool all_boxes_touch_bottom = true;
63  for (int i = 0; i < ba.size(); i++) {
64  if (ba[i].smallEnd(2) != domlo_z) {
65  all_boxes_touch_bottom = false;
66  }
67  }
68 
69  if (all_boxes_touch_bottom) {
70  init_which_terrain_grid(lev, geom, z_phys_nd, z_levels_h);
71  } else {
72 
73  BoxArray ba_new(domain);
74  ChopGrids2D(ba_new, domain, ParallelDescriptor::NProcs());
75 
76  DistributionMapping dm_new(ba_new);
77  ba_new.surroundingNodes();
78 
79  MultiFab z_phys_nd_new(ba_new, dm_new, 1, z_phys_nd.nGrowVect());
80 
81  z_phys_nd_new.ParallelCopy(z_phys_nd,0,0,1,z_phys_nd.nGrowVect(),z_phys_nd.nGrowVect());
82 
83  init_which_terrain_grid(lev, geom, z_phys_nd_new, z_levels_h);
84 
85  z_phys_nd.ParallelCopy(z_phys_nd_new,0,0,1,z_phys_nd.nGrowVect(),z_phys_nd.nGrowVect());
86  }
87  } else { // lev > 0
88  init_which_terrain_grid(lev, geom, z_phys_nd, z_levels_h);
89  }
90 
91  //
92  // Fill ghost layers and corners (including periodic) -- no matter what level
93  //
94  z_phys_nd.FillBoundary(geom.periodicity());
95 
96  if (phys_bc_type[Orientation(0,Orientation::low )] == ERF_BC::symmetry ||
97  phys_bc_type[Orientation(0,Orientation::high)] == ERF_BC::symmetry ||
98  phys_bc_type[Orientation(1,Orientation::low )] == ERF_BC::symmetry ||
99  phys_bc_type[Orientation(1,Orientation::high)] == ERF_BC::symmetry) {
100 
101  const auto& dom_lo = lbound(convert(domain,IntVect(1,1,1)));
102  const auto& dom_hi = ubound(convert(domain,IntVect(1,1,1)));
103 
104  for (MFIter mfi(z_phys_nd,true); mfi.isValid(); ++mfi) {
105  const Box& bx = mfi.growntilebox();
106  const Array4< Real> z_nd_arr = z_phys_nd.array(mfi);
107  if (phys_bc_type[Orientation(0,Orientation::low)] == ERF_BC::symmetry && bx.smallEnd(0) == dom_lo.x) {
108  ParallelFor(makeSlab(bx,0,1), [=] AMREX_GPU_DEVICE (int , int j, int k)
109  {
110  z_nd_arr(dom_lo.x-1,j,k) = z_nd_arr(dom_lo.x+1,j,k);
111  });
112  }
113  if (phys_bc_type[Orientation(0,Orientation::high)] == ERF_BC::symmetry && bx.bigEnd(0) == dom_hi.x) {
114  ParallelFor(makeSlab(bx,0,1), [=] AMREX_GPU_DEVICE (int , int j, int k)
115  {
116  z_nd_arr(dom_hi.x+1,j,k) = z_nd_arr(dom_hi.x-1,j,k);
117  });
118  }
119  if (phys_bc_type[Orientation(1,Orientation::low)] == ERF_BC::symmetry && bx.smallEnd(1) == dom_lo.y) {
120  ParallelFor(makeSlab(bx,1,1), [=] AMREX_GPU_DEVICE (int i, int , int k)
121  {
122  z_nd_arr(i,dom_lo.y-1,k) = z_nd_arr(i,dom_lo.y+1,k);
123  });
124  }
125  if (phys_bc_type[Orientation(1,Orientation::high)] == ERF_BC::symmetry && bx.bigEnd(1) == dom_hi.y) {
126  ParallelFor(makeSlab(bx,1,1), [=] AMREX_GPU_DEVICE (int i, int , int k)
127  {
128  z_nd_arr(i,dom_hi.y+1,k) = z_nd_arr(i,dom_hi.y-1,k);
129  });
130  }
131  }
132  }
133 
134  //********************************************************************************
135  // Populate domain boundary cells in z-direction
136  //********************************************************************************
137  int ngrow = z_phys_nd.nGrow();
138 
139  for ( MFIter mfi(z_phys_nd, TilingIfNotGPU()); mfi.isValid(); ++mfi )
140  {
141  // Only set values above top of domain if this box reaches that far
142  Box nd_bx = mfi.tilebox();
143 
144  // Note that domhi_z is already nodal in the z-direction
145  if (nd_bx.bigEnd(2) >= domhi_z) {
146  // Grown box with no z range
147  Box bx_zhi = mfi.growntilebox(ngrow);
148  bx_zhi.setSmall(2,domhi_z+1);
149  Array4<Real> const& z_arr = z_phys_nd.array(mfi);
150 
151  // Extrapolate top layer
152  ParallelFor(bx_zhi, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
153  z_arr(i,j,k) = z_arr(i,j,domhi_z)
154  + (k-domhi_z) * (z_arr(i,j,domhi_z) - z_arr(i,j,domhi_z-1));
155  });
156  }
157  }
158 } // make_terrain_fitted_coords
void ChopGrids2D(BoxArray &ba, const Box &domain, int target_size)
Definition: ERF_ChopGrids.cpp:21
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
void init_which_terrain_grid(int lev, Geometry const &geom, MultiFab &z_phys_nd, Vector< Real > const &z_levels_h)
Definition: ERF_TerrainMetrics.cpp:161

Referenced by ERF::init_zphys(), and ERF::remake_zphys().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ make_zcc()

void make_zcc ( const Geometry &  geom,
MultiFab &  z_phys_nd,
MultiFab &  z_phys_cc 
)

Computation of z_phys at cell-center

630 {
631 #ifdef _OPENMP
632 #pragma omp parallel if (Gpu::notInLaunchRegion())
633 #endif
634  for ( MFIter mfi(z_phys_cc, TilingIfNotGPU()); mfi.isValid(); ++mfi )
635  {
636  Box gbx = mfi.growntilebox();
637  Array4<Real const> z_nd = z_phys_nd.const_array(mfi);
638  Array4<Real > z_cc = z_phys_cc.array(mfi);
639  ParallelFor(gbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
640  z_cc(i, j, k) = .125 * ( z_nd(i,j,k ) + z_nd(i+1,j,k ) + z_nd(i,j+1,k ) + z_nd(i+1,j+1,k )
641  +z_nd(i,j,k+1) + z_nd(i+1,j,k+1) + z_nd(i,j+1,k+1) + z_nd(i+1,j+1,k+1) );
642  });
643  }
644  z_phys_cc.FillBoundary(geom.periodicity());
645 }

Referenced by ERF::post_timestep(), and ERF::update_terrain_arrays().

Here is the call graph for this function:
Here is the caller graph for this function: