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, Real z_offset)
 
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

654 {
655  auto const& ma_z_nd_arr = z_phys_nd.const_arrays();
656  GpuTuple<Real> min = ParReduce(TypeList<ReduceOpMin>{},
657  TypeList<Real>{},
658  z_phys_nd, IntVect(0),
659  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
660  -> GpuTuple<Real>
661  {
662  amrex::Real dz = Compute_Z_AtWFace(i,j,k+1,ma_z_nd_arr[box_no]) -
663  Compute_Z_AtWFace(i,j,k ,ma_z_nd_arr[box_no]);
664  return { dz };
665  });
666  Real r = (get<0>(min) + std::numeric_limits<amrex::Real>::epsilon());
667  ParallelDescriptor::ReduceRealMin(r);
668  return r;
669 }
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:376
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,
Real  z_offset 
)

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

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

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

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

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

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

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

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

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

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: