ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERFPCParticleToMesh.H
Go to the documentation of this file.
1 /*! \file ERFPCParticleToMesh.H
2  * \brief Template implementation of ERFPC::ERFPCParticleToMesh
3  *
4  * Separated from ERFPC.H to keep the header lightweight.
5  * Include this in .cpp files that call ERFPCParticleToMesh().
6  */
7 
8 #ifndef ERFPC_PARTICLE_TO_MESH_H_
9 #define ERFPC_PARTICLE_TO_MESH_H_
10 
11 #include <ERFPC.H>
12 #include <ERF_TerrainConversion.H>
13 
14 template<typename ValueFunc>
15 void ERFPC::ERFPCParticleToMesh(amrex::MultiFab& a_mf,
16  const amrex::MultiFab& a_z_phys_nd,
17  int a_lev, int a_comp,
18  ValueFunc&& value_func) const
19 {
20  using namespace amrex;
21 
22  AMREX_ASSERT(OK());
23  AMREX_ASSERT(numParticlesOutOfRange(*this, 0) == 0);
24 
25  const auto& geom = Geom(a_lev);
26  const auto plo = geom.ProbLoArray();
27  const auto dxi = geom.InvCellSizeArray();
28  const Real dx = geom.CellSize(0);
29  const Real dy = geom.CellSize(1);
30  const Box domain = geom.Domain();
31  const int k_max = domain.bigEnd(AMREX_SPACEDIM-1);
32 
33  a_mf.setVal(0.0);
34 
35 #ifdef AMREX_USE_OMP
36 #pragma omp parallel if (Gpu::notInLaunchRegion())
37 #endif
38  for (ParConstIterType pti(*this, a_lev); pti.isValid(); ++pti)
39  {
40  const int grid = pti.index();
41  const auto& ptile = ParticlesAt(a_lev, pti);
42  const auto& aos = ptile.GetArrayOfStructs();
43  const int np = aos.numParticles();
44  if (np == 0) { continue; }
45 
46  auto rho = a_mf[grid].array();
47  auto zheight = a_z_phys_nd[grid].array();
48 
49  const auto ptd = ptile.getConstParticleTileData();
50 
51  ParallelFor(np, [=] AMREX_GPU_DEVICE (int i)
52  {
53  const auto& p = ptd.m_aos[i];
54  if (p.id() <= 0) { return; }
55 
56  const Real pval = static_cast<Real>(value_func(ptd, i));
57 
58  const Real lx = static_cast<Real>((p.pos(0) - plo[0]) * dxi[0]) + Real(0.5);
59  const int ix = static_cast<int>(Math::floor(lx)) - 1;
60  const Real wx1 = lx - static_cast<Real>(ix + 1);
61  const Real wx0 = Real(1.0) - wx1;
62 
63  const Real ly = static_cast<Real>((p.pos(1) - plo[1]) * dxi[1]) + Real(0.5);
64  const int iy = static_cast<int>(Math::floor(ly)) - 1;
65  const Real wy1 = ly - static_cast<Real>(iy + 1);
66  const Real wy0 = Real(1.0) - wy1;
67 
68  const int kz = static_cast<int>(amrex::Math::floor((p.pos(AMREX_SPACEDIM-1) - plo[AMREX_SPACEDIM-1])
69  * dxi[AMREX_SPACEDIM-1]));
70  const int ix_p = static_cast<int>(Math::floor((p.pos(0) - plo[0]) * dxi[0]));
71  const int iy_p = static_cast<int>(Math::floor((p.pos(1) - plo[1]) * dxi[1]));
72  const Real fx = static_cast<Real>((p.pos(0) - plo[0]) * dxi[0]) - static_cast<Real>(ix_p);
73  const Real fy = static_cast<Real>((p.pos(1) - plo[1]) * dxi[1]) - static_cast<Real>(iy_p);
74 
75  auto zface = [=] AMREX_GPU_DEVICE (int k_face) -> Real {
76  return zheight(ix_p , iy_p , k_face) * (Real(1.0)-fx) * (Real(1.0)-fy)
77  + zheight(ix_p+1, iy_p , k_face) * fx * (Real(1.0)-fy)
78  + zheight(ix_p , iy_p+1, k_face) * (Real(1.0)-fx) * fy
79  + zheight(ix_p+1, iy_p+1, k_face) * fx * fy;
80  };
81 
82  // Effective k bounds for safe zface access: must respect both the
83  // geometry domain (k_max) and the local fab's z extent. With AMR
84  // partial-z fabs (especially with many small boxes from large rank
85  // counts), a particle at the top/bottom of its fab can have kz
86  // such that kz+2 / kz-1 read outside the fab. Use the fallback
87  // (one-sided) span formula in those cases.
88  const int kz_hi_safe = amrex::min(k_max, zheight.end[2] - 2);
89  const int kz_lo_safe = amrex::max(0, zheight.begin[2] + 1);
90 
91  const Real z_lo_k = zface(kz);
92  const Real z_hi_k = zface(kz + 1);
93  const Real z_ctr_k = Real(0.5) * (z_lo_k + z_hi_k);
94 
95  // Particle pos(2) is computational zeta; convert to physical z so
96  // the cell-center comparison and weight span are in the same units
97  // as zface() (which is built from the height array).
98  const Real p_z = static_cast<Real>(ERF::ParticlePos::z_from_zeta(
99  p.pos(0), p.pos(1), p.pos(AMREX_SPACEDIM-1),
100  plo, dxi, zheight));
101 
102  int kz_lo, kz_hi;
103  Real wz_lo, wz_hi;
104  if (p_z >= z_ctr_k) {
105  kz_lo = kz;
106  kz_hi = kz + 1;
107  Real span;
108  if (kz < kz_hi_safe) {
109  const Real z_hi2 = zface(kz + 2);
110  const Real z_ctr_hi = Real(0.5) * (z_hi_k + z_hi2);
111  span = z_ctr_hi - z_ctr_k;
112  } else {
113  span = z_hi_k - z_lo_k;
114  }
115  wz_hi = (p_z - z_ctr_k) / span;
116  wz_lo = Real(1.0) - wz_hi;
117  } else {
118  kz_hi = kz;
119  kz_lo = kz - 1;
120  Real span;
121  if (kz > kz_lo_safe) {
122  const Real z_lo2 = zface(kz - 1);
123  const Real z_ctr_lo = Real(0.5) * (z_lo2 + z_lo_k);
124  span = z_ctr_k - z_ctr_lo;
125  } else {
126  span = z_hi_k - z_lo_k;
127  }
128  wz_lo = (z_ctr_k - p_z) / span;
129  wz_hi = Real(1.0) - wz_lo;
130  }
131 
132  // Per-corner node-averaged cell thickness for volume normalization.
133  auto dz_cell = [=] AMREX_GPU_DEVICE (int ic, int jc, int kc) -> Real {
134  const Real dz = Real(0.25) * (
135  zheight(ic , jc , kc+1) - zheight(ic , jc , kc)
136  + zheight(ic+1, jc , kc+1) - zheight(ic+1, jc , kc)
137  + zheight(ic , jc+1, kc+1) - zheight(ic , jc+1, kc)
138  + zheight(ic+1, jc+1, kc+1) - zheight(ic+1, jc+1, kc) );
139  return dz;
140  };
141 
142  const Real inv_dxdy = Real(1.0) / (dx * dy);
143 
144  auto deposit = [=] AMREX_GPU_DEVICE (int ic, int jc, int kc, Real w_horiz, Real w_z)
145  {
146  // Clamp kc only for dz lookup; still write to ghost rho cell.
147  const int kc_dz = amrex::max(0, amrex::min(k_max, kc));
148  const Real dz = dz_cell(ic, jc, kc_dz);
149  const Real contrib = pval * w_horiz * w_z * inv_dxdy / dz;
150  Gpu::Atomic::AddNoRet(&rho(ic, jc, kc, a_comp), contrib);
151  };
152 
153  deposit(ix , iy , kz_lo, wx0*wy0, wz_lo);
154  deposit(ix+1, iy , kz_lo, wx1*wy0, wz_lo);
155  deposit(ix , iy+1, kz_lo, wx0*wy1, wz_lo);
156  deposit(ix+1, iy+1, kz_lo, wx1*wy1, wz_lo);
157  deposit(ix , iy , kz_hi, wx0*wy0, wz_hi);
158  deposit(ix+1, iy , kz_hi, wx1*wy0, wz_hi);
159  deposit(ix , iy+1, kz_hi, wx0*wy1, wz_hi);
160  deposit(ix+1, iy+1, kz_hi, wx1*wy1, wz_hi);
161  });
162  }
163 
164  a_mf.SumBoundary(geom.periodicity());
165 }
166 
167 #endif
const Real dy
Definition: ERF_InitCustomPert_ABL.H:24
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
rho
Definition: ERF_InitCustomPert_Bubble.H:106
const Real dz
Definition: ERF_InitCustomPert_Bubble.H:25
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+myhalf) *dx[0];const Real y=(j+myhalf) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - one)/(two *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(one+deltaT, inv_gm1);const Real T=(one+deltaT) *T_inf;const Real p=std::pow(rho_norm, Gamma)/Gamma *rho_0 *a_inf *a_inf;const Real rho_theta=rho_0 *rho_norm *(T *std::pow(p_0/p, rdOcp));state_pert(i, j, k, RhoTheta_comp)=rho_theta - getRhoThetagivenP(p_hse(i, j, k));const Real r2d_xy=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc));state_pert(i, j, k, RhoScalar_comp)=fourth *(one+std::cos(PI *std::min(r2d_xy, R)/R));})
Real * p
Definition: ERF_InitCustomPert_SquallLine.H:61
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real z_from_zeta(amrex::Real x, amrex::Real y, amrex::Real zeta, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &plo, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &dxi, amrex::Array4< amrex::Real const > const &height_arr) noexcept
Definition: ERF_TerrainConversion.H:51
Definition: ERF_ConsoleIO.cpp:12