WarpX
WarpX_PEC.H
Go to the documentation of this file.
1 #ifndef WARPX_PEC_KERNELS_H_
2 #define WARPX_PEC_KERNELS_H_
3 
4 #include "WarpX.H"
6 
7 #include <AMReX_Array.H>
8 #include <AMReX_Array4.H>
9 #include <AMReX_Config.H>
10 #include <AMReX_Extension.H>
11 #include <AMReX_GpuQualifiers.H>
12 #include <AMReX_IntVect.H>
13 #include <AMReX_REAL.H>
14 
15 #include <AMReX_BaseFwd.H>
16 
17 #include <array>
18 #include <memory>
19 
20 namespace PEC {
21 using namespace amrex;
32  bool is_boundary_PEC (amrex::GpuArray<int, 3> const& fboundary, int dir) {
33  return ( fboundary[dir] == FieldBoundaryType::PEC );
34  }
35 
47  return ( pboundary[dir] == ParticleBoundaryType::Reflecting );
48  }
49 
68  const amrex::IntVect& dom_hi, const amrex::IntVect& ijk_vec,
69  const amrex::IntVect& is_nodal, const int idim, const int iside)
70  {
71  return ((iside == 0) ? (dom_lo[idim] - ijk_vec[idim])
72  : (ijk_vec[idim] - (dom_hi[idim] + is_nodal[idim])));
73  }
74 
142  void SetEfieldOnPEC (const int icomp, const amrex::IntVect & dom_lo,
143  const amrex::IntVect &dom_hi,
144  const amrex::IntVect &ijk_vec, const int n,
145  amrex::Array4<amrex::Real> const& Efield,
146  const amrex::IntVect& is_nodal,
147  amrex::GpuArray<int, 3> const& fbndry_lo,
148  amrex::GpuArray<int, 3> const& fbndry_hi )
149  {
150  // Tangential Efield componentes in guard cells set equal and opposite to cells
151  // in the mirror locations across the PEC boundary, whereas normal E-field
152  // components are set equal to values in the mirror locations across the PEC
153  // boundary. Here we just initialize it.
154  amrex::IntVect ijk_mirror = ijk_vec;
155  bool OnPECBoundary = false;
156  bool GuardCell = false;
157  amrex::Real sign = 1._rt;
158  // Loop over all the dimensions
159  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
160  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
161  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
162  for (int iside = 0; iside < 2; ++iside) {
163  const bool isPECBoundary = ( (iside == 0)
164  ? is_boundary_PEC(fbndry_lo, idim)
165  : is_boundary_PEC(fbndry_hi, idim) );
166 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
167  // For 2D : for icomp==1, (Ey in XZ, Etheta in RZ),
168  // icomp=1 is tangential to both x and z boundaries
169  // The logic below ensures that the flags are set right for 2D
170  const bool is_tangent_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
171  ? false : true );
172 #elif (defined WARPX_DIM_1D_Z)
173  // For 1D : icomp=0 and icomp=1 (Ex and Ey are tangential to the z boundary)
174  // The logic below ensures that the flags are set right for 1D
175  const bool is_tangent_to_PEC = ( ( icomp == idim+2) ? false : true );
176 #else
177  const bool is_tangent_to_PEC = ( ( icomp == idim) ? false : true );
178 #endif
179  if (isPECBoundary == true) {
180  // grid point ijk_vec is ig number of points pass the
181  // domain boundary in direction, idim
183  dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
184 
185  if (ig == 0) {
186  if (is_tangent_to_PEC == true and is_nodal[idim] == 1) {
187  OnPECBoundary = true;
188  }
189  } else if (ig > 0) {
190  // Find mirror location across PEC boundary
191  ijk_mirror[idim] = ( ( iside == 0)
192  ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
193  : (dom_hi[idim] + 1 - ig));
194  GuardCell = true;
195  // tangential components are inverted across PEC boundary
196  if (is_tangent_to_PEC == true) sign *= -1._rt;
197  }
198  } // is PEC boundary
199  } // loop over iside
200  } // loop over dimensions
201  if (OnPECBoundary == true) {
202  // if ijk_vec is on a PEC boundary in any direction, set Etangential to 0.
203  Efield(ijk_vec,n) = 0._rt;
204  } else if (GuardCell == true) {
205  Efield(ijk_vec,n) = sign * Efield(ijk_mirror,n);
206  }
207  }
208 
272  void SetBfieldOnPEC (const int icomp, const amrex::IntVect & dom_lo,
273  const amrex::IntVect & dom_hi,
274  const amrex::IntVect & ijk_vec, const int n,
275  amrex::Array4<amrex::Real> const& Bfield,
276  const amrex::IntVect & is_nodal,
277  amrex::GpuArray<int, 3> const& fbndry_lo,
278  amrex::GpuArray<int, 3> const& fbndry_hi )
279  {
280  amrex::IntVect ijk_mirror = ijk_vec;
281  bool OnPECBoundary = false;
282  bool GuardCell = false;
283  amrex::Real sign = 1._rt;
284  // Loop over all dimensions
285  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
286  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
287  for (int iside = 0; iside < 2; ++iside) {
288  const bool isPECBoundary = ( (iside == 0 )
289  ? is_boundary_PEC(fbndry_lo, idim)
290  : is_boundary_PEC(fbndry_hi, idim) );
291  if (isPECBoundary == true) {
292 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
293  // For 2D : for icomp==1, (By in XZ, Btheta in RZ),
294  // icomp=1 is not normal to x or z boundary
295  // The logic below ensures that the flags are set right for 2D
296  const bool is_normal_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
297  ? true : false );
298 #elif (defined WARPX_DIM_1D_Z)
299  // For 1D : icomp=0 and icomp=1 (Bx and By are not normal to the z boundary)
300  // The logic below ensures that the flags are set right for 1D
301  const bool is_normal_to_PEC = ( ( icomp == idim+2) ? true : false );
302 #else
303  const bool is_normal_to_PEC = ( ( icomp == idim) ? true : false );
304 #endif
305 
306  // grid point ijk_vec is ig number of points pass the
307  // domain boundary in direction, idim
309  dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
310 
311  if (ig == 0) {
312  // Only normal component is set to 0
313  if (is_normal_to_PEC == true and is_nodal[idim]==1) {
314  OnPECBoundary = true;
315  }
316  } else if ( ig > 0) {
317  // Mirror location inside the domain by "ig" number of cells
318  // across PEC boundary in direction, idim, and side, iside
319  ijk_mirror[idim] = ( (iside == 0)
320  ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
321  : (dom_hi[idim] + 1 - ig));
322  GuardCell = true;
323  // Sign of the normal component in guard cell is inverted
324  if (is_normal_to_PEC == true) sign *= -1._rt;
325  }
326  } // if PEC Boundary
327  } // loop over sides
328  } // loop of dimensions
329 
330  if (OnPECBoundary == true) {
331  // if ijk_vec is on a PEC boundary in any direction, set Bnormal to 0.
332  Bfield(ijk_vec,n) = 0._rt;
333  } else if (GuardCell == true) {
334  // Bnormal and Btangential is set opposite and equal to the value
335  // in the mirror location, respectively.
336  Bfield(ijk_vec,n) = sign * Bfield(ijk_mirror,n);
337  }
338  }
339 
340 
362  void SetRhoOrJfieldFromPEC (const int n,
363  const amrex::IntVect & ijk_vec,
365  amrex::GpuArray<GpuArray<int, 2>, AMREX_SPACEDIM> const& mirrorfac,
366  amrex::GpuArray<GpuArray<amrex::Real, 2>, AMREX_SPACEDIM> const& psign,
367  amrex::GpuArray<GpuArray<bool, 2>, AMREX_SPACEDIM> const& is_pec,
368  amrex::GpuArray<bool, AMREX_SPACEDIM> const& tangent_to_bndy,
369  amrex::Box const& fabbox)
370  {
371  // The boundary is handled in 2 steps:
372  // 1) The cells internal to the domain are updated using the
373  // current deposited in the guard cells
374  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
375  {
376  for (int iside = 0; iside < 2; ++iside)
377  {
378  if (!is_pec[idim][iside]) continue;
379 
380  // Get the mirror guard cell index
381  amrex::IntVect iv_mirror = ijk_vec;
382  iv_mirror[idim] = mirrorfac[idim][iside] - ijk_vec[idim];
383 
384  // On the PEC boundary the current density is set to 0
385  if (ijk_vec == iv_mirror) field(ijk_vec, n) = 0._rt;
386  // otherwise update the internal cell if the mirror guard cell exists
387  else if (fabbox.contains(iv_mirror))
388  {
389  field(ijk_vec,n) += psign[idim][iside] * field(iv_mirror,n);
390  }
391  }
392  }
393  // 2) The guard cells are updated with the appropriate image
394  // charge current based on the current in the valid cells
395  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
396  {
397  for (int iside = 0; iside < 2; ++iside)
398  {
399  if (!is_pec[idim][iside]) continue;
400 
401  amrex::IntVect iv_mirror = ijk_vec;
402  iv_mirror[idim] = mirrorfac[idim][iside] - ijk_vec[idim];
403  if (ijk_vec != iv_mirror && fabbox.contains(iv_mirror))
404  {
405  if (tangent_to_bndy[idim])
406  field(iv_mirror, n) = -field(ijk_vec, n);
407  else
408  field(iv_mirror, n) = field(ijk_vec, n);
409  }
410  }
411  }
412  }
413 
414 
416  bool isAnyBoundaryPEC();
428  void ApplyPECtoEfield ( std::array<amrex::MultiFab*, 3> Efield,
429  const int lev, PatchType patch_type,
430  const bool split_pml_field = false);
440  void ApplyPECtoBfield ( std::array<amrex::MultiFab*, 3> Bfield,
441  const int lev, PatchType patch_type);
442 
451  void ApplyPECtoRhofield(amrex::MultiFab* rho, const int lev,
452  PatchType patch_type);
453 
463  amrex::MultiFab* Jz, const int lev,
464  PatchType patch_type);
465 }
466 
467 #endif // WarpX_PEC_KERNELS_H_
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
PatchType
Definition: WarpX.H:72
@ Reflecting
particles are reflected
AMREX_GPU_HOST_DEVICE bool contains(const IntVect &p) const noexcept
Definition: WarpX_PEC.H:20
void ApplyPECtoRhofield(amrex::MultiFab *rho, const int lev, PatchType patch_type)
Reflects charge density deposited over the PEC boundary back into the simulation domain.
Definition: WarpX_PEC.cpp:225
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void SetRhoOrJfieldFromPEC(const int n, const amrex::IntVect &ijk_vec, amrex::Array4< amrex::Real > const &field, amrex::GpuArray< GpuArray< int, 2 >, AMREX_SPACEDIM > const &mirrorfac, amrex::GpuArray< GpuArray< amrex::Real, 2 >, AMREX_SPACEDIM > const &psign, amrex::GpuArray< GpuArray< bool, 2 >, AMREX_SPACEDIM > const &is_pec, amrex::GpuArray< bool, AMREX_SPACEDIM > const &tangent_to_bndy, amrex::Box const &fabbox)
Sets the rho or J field value in cells close to and on a PEC boundary. The charge/current density dep...
Definition: WarpX_PEC.H:362
void ApplyPECtoBfield(std::array< amrex::MultiFab *, 3 > Bfield, const int lev, PatchType patch_type)
Sets the normal component of the magnetic field at the PEC boundary to zero. The guard cell values ar...
Definition: WarpX_PEC.cpp:124
AMREX_GPU_DEVICE AMREX_FORCE_INLINE int get_cell_count_to_boundary(const amrex::IntVect &dom_lo, const amrex::IntVect &dom_hi, const amrex::IntVect &ijk_vec, const amrex::IntVect &is_nodal, const int idim, const int iside)
Calculates the number of grid points the given index is pass the domain boundary i....
Definition: WarpX_PEC.H:67
void ApplyPECtoJfield(amrex::MultiFab *Jx, amrex::MultiFab *Jy, amrex::MultiFab *Jz, const int lev, PatchType patch_type)
Reflects current density deposited over the PEC boundary back into the simulation domain.
Definition: WarpX_PEC.cpp:315
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void SetEfieldOnPEC(const int icomp, const amrex::IntVect &dom_lo, const amrex::IntVect &dom_hi, const amrex::IntVect &ijk_vec, const int n, amrex::Array4< amrex::Real > const &Efield, const amrex::IntVect &is_nodal, amrex::GpuArray< int, 3 > const &fbndry_lo, amrex::GpuArray< int, 3 > const &fbndry_hi)
Sets the electric field value tangential to the PEC boundary to zero. The tangential Efield component...
Definition: WarpX_PEC.H:142
AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool is_boundary_PEC(amrex::GpuArray< int, 3 > const &fboundary, int dir)
Determines if the field boundary condition stored in fboundary in direction, dir, is PEC.
Definition: WarpX_PEC.H:32
AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool is_boundary_reflecting(amrex::GpuArray< ParticleBoundaryType, 3 > const &pboundary, int dir)
Determines if the particle boundary condition stored in pboundary in direction, dir,...
Definition: WarpX_PEC.H:46
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void SetBfieldOnPEC(const int icomp, const amrex::IntVect &dom_lo, const amrex::IntVect &dom_hi, const amrex::IntVect &ijk_vec, const int n, amrex::Array4< amrex::Real > const &Bfield, const amrex::IntVect &is_nodal, amrex::GpuArray< int, 3 > const &fbndry_lo, amrex::GpuArray< int, 3 > const &fbndry_hi)
Sets the magnetic field value normal to the PEC boundary to zero. The tangential (and normal) field v...
Definition: WarpX_PEC.H:272
bool isAnyBoundaryPEC()
Definition: WarpX_PEC.cpp:19
void ApplyPECtoEfield(std::array< amrex::MultiFab *, 3 > Efield, const int lev, PatchType patch_type, const bool split_pml_field=false)
Sets the tangential electric field at the PEC boundary to zero. The guard cell values are set equal a...
Definition: WarpX_PEC.cpp:28
int n
Definition: run_libensemble_on_warpx.py:67
string field
Definition: video_yt.py:31
@ PEC
perfect electric conductor (PEC) with E_tangential=0
Definition: WarpXAlgorithmSelection.H:134