ARTEMIS
ParticleScraper.H
Go to the documentation of this file.
1 /* Copyright 2021 Andrew Myers
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef PARTICLESCRAPER_H_
8 #define PARTICLESCRAPER_H_
9 
12 
14 
15 #include <AMReX.H>
16 #include <AMReX_MultiFab.H>
17 #include <AMReX_Particle.H>
18 #include <AMReX_RandomEngine.H>
19 #include <AMReX_REAL.H>
20 #include <AMReX_TypeTraits.H>
21 #include <AMReX_Vector.H>
22 
23 #include <type_traits>
24 #include <utility>
25 
26 
27 
66 void
67 scrapeParticles (PC& pc, const amrex::Vector<const amrex::MultiFab*>& distance_to_eb, int lev, F&& f)
68 {
69  scrapeParticles(pc, distance_to_eb, lev, lev, std::forward<F>(f));
70 }
71 
109 void
110 scrapeParticles (PC& pc, const amrex::Vector<const amrex::MultiFab*>& distance_to_eb, F&& f)
111 {
112  scrapeParticles(pc, distance_to_eb, 0, pc.finestLevel(), std::forward<F>(f));
113 }
114 
154 void
156  int lev_min, int lev_max, F&& f)
157 {
158  BL_PROFILE("scrapeParticles");
159 
160  for (int lev = lev_min; lev <= lev_max; ++lev)
161  {
162  const auto plo = pc.Geom(lev).ProbLoArray();
163  const auto dxi = pc.Geom(lev).InvCellSizeArray();
164  for(WarpXParIter pti(pc, lev); pti.isValid(); ++pti)
165  {
166  const auto getPosition = GetParticlePosition(pti);
167  auto& tile = pti.GetParticleTile();
168  auto ptd = tile.getParticleTileData();
169  const auto np = tile.numParticles();
170  amrex::Particle<0,0> * const particles = tile.GetArrayOfStructs()().data();
171  auto phi = (*distance_to_eb[lev])[pti].array(); // signed distance function
173  [=] AMREX_GPU_DEVICE (const int ip, amrex::RandomEngine const& engine) noexcept
174  {
175  // skip particles that are already flagged for removal
176  if (particles[ip].id() < 0) return;
177 
178  amrex::ParticleReal xp, yp, zp;
179  getPosition(ip, xp, yp, zp);
180 
181  int i, j, k;
182  amrex::Real W[AMREX_SPACEDIM][2];
183  ablastr::particles::compute_weights_nodal(xp, yp, zp, plo, dxi, i, j, k, W);
184 
185  amrex::Real phi_value = ablastr::particles::interp_field_nodal(i, j, k, W, phi);
186 
187  if (phi_value < 0.0)
188  {
189  amrex::RealVect normal = DistanceToEB::interp_normal(i, j, k, W, phi, dxi);
190 
191  // the closest point on the surface to pos is pos - grad phi(pos) * phi(pos)
192  amrex::RealVect pos;
193 #if (defined WARPX_DIM_3D)
194  pos[0] = xp - normal[0]*phi_value;
195  pos[1] = yp - normal[1]*phi_value;
196  pos[2] = zp - normal[2]*phi_value;
197 #elif (defined WARPX_DIM_XZ)
198  pos[0] = xp - normal[0]*phi_value;
199  pos[1] = zp - normal[1]*phi_value;
200 #elif (defined WARPX_DIM_RZ)
201  pos[0] = std::sqrt(xp*xp + yp*yp) - normal[0]*phi_value;
202  pos[1] = zp - normal[1]*phi_value;
203 #elif (defined WARPX_DIM_1D_Z)
204  pos[0] = zp - normal[0]*phi_value;
205 #endif
206  DistanceToEB::normalize(normal);
207 
208  f(ptd, ip, pos, normal, engine);
209  }
210  });
211  }
212  }
213 }
214 
215 #endif
#define BL_PROFILE(a)
#define AMREX_GPU_DEVICE
void scrapeParticles(PC &pc, const amrex::Vector< const amrex::MultiFab * > &distance_to_eb, int lev, F &&f)
Interact particles with the embedded boundary walls.
Definition: ParticleScraper.H:67
Definition: WarpXParticleContainer.H:52
bool isValid() const noexcept
list data
open text file and read data #####
Definition: Excitation_Flag_Generator.py:24
AMREX_GPU_HOST_DEVICE AMREX_INLINE amrex::Real interp_field_nodal(int i, int j, int k, const amrex::Real W[AMREX_SPACEDIM][2], amrex::Array4< const amrex::Real > const &scalar_field) noexcept
Interpolate nodal field value based on surrounding indices and weights.
Definition: NodalFieldGather.H:92
AMREX_GPU_HOST_DEVICE AMREX_INLINE void compute_weights_nodal(const amrex::ParticleReal xp, const amrex::ParticleReal yp, const amrex::ParticleReal zp, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &plo, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &dxi, int &i, int &j, int &k, amrex::Real W[AMREX_SPACEDIM][2]) noexcept
Compute weight of each surrounding node in interpolating a nodal field to the given coordinates.
Definition: NodalFieldGather.H:32
AMREX_ATTRIBUTE_FLATTEN_FOR void ParallelForRNG(T n, L const &f) noexcept
i
Definition: check_interp_points_and_weights.py:174
value
Definition: updateAMReX.py:141
f
Definition: write_atomic_data_cpp.py:88
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition: GetAndSetPosition.H:53