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  for (int iside = 0; iside < 2; ++iside) {
162  const bool isPECBoundary = ( (iside == 0)
163  ? is_boundary_PEC(fbndry_lo, idim)
164  : is_boundary_PEC(fbndry_hi, idim) );
165 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
166  // For 2D : for icomp==1, (Ey in XZ, Etheta in RZ),
167  // icomp=1 is tangential to both x and z boundaries
168  // The logic below ensures that the flags are set right for 2D
169  const bool is_tangent_to_PEC = (icomp != AMREX_SPACEDIM*idim);
170 #elif (defined WARPX_DIM_1D_Z)
171  // For 1D : icomp=0 and icomp=1 (Ex and Ey are tangential to the z boundary)
172  // The logic below ensures that the flags are set right for 1D
173  const bool is_tangent_to_PEC = (icomp != idim+2);
174 #else
175  const bool is_tangent_to_PEC = (icomp != idim);
176 #endif
177  if (isPECBoundary) {
178  // grid point ijk_vec is ig number of points pass the
179  // domain boundary in direction, idim
180  const int ig = get_cell_count_to_boundary(
181  dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
182 
183  if (ig == 0) {
184  if (is_tangent_to_PEC && is_nodal[idim] == 1) {
185  OnPECBoundary = true;
186  }
187  } else if (ig > 0) {
188  // Find mirror location across PEC boundary
189  ijk_mirror[idim] = ( ( iside == 0)
190  ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
191  : (dom_hi[idim] + 1 - ig));
192  GuardCell = true;
193  // tangential components are inverted across PEC boundary
194  if (is_tangent_to_PEC) sign *= -1._rt;
195  }
196  } // is PEC boundary
197  } // loop over iside
198  } // loop over dimensions
199  if (OnPECBoundary) {
200  // if ijk_vec is on a PEC boundary in any direction, set Etangential to 0.
201  Efield(ijk_vec,n) = 0._rt;
202  } else if (GuardCell) {
203  Efield(ijk_vec,n) = sign * Efield(ijk_mirror,n);
204  }
205  }
206 
207 
271  void SetBfieldOnPEC (const int icomp, const amrex::IntVect & dom_lo,
272  const amrex::IntVect & dom_hi,
273  const amrex::IntVect & ijk_vec, const int n,
274  amrex::Array4<amrex::Real> const& Bfield,
275  const amrex::IntVect & is_nodal,
276  amrex::GpuArray<int, 3> const& fbndry_lo,
277  amrex::GpuArray<int, 3> const& fbndry_hi )
278  {
279  amrex::IntVect ijk_mirror = ijk_vec;
280  bool OnPECBoundary = false;
281  bool GuardCell = false;
282  amrex::Real sign = 1._rt;
283  // Loop over all dimensions
284  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
285  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
286  for (int iside = 0; iside < 2; ++iside) {
287  const bool isPECBoundary = ( (iside == 0 )
288  ? is_boundary_PEC(fbndry_lo, idim)
289  : is_boundary_PEC(fbndry_hi, idim) );
290  if (isPECBoundary) {
291 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
292  // For 2D : for icomp==1, (By in XZ, Btheta in RZ),
293  // icomp=1 is not normal to x or z boundary
294  // The logic below ensures that the flags are set right for 2D
295  const bool is_normal_to_PEC = (icomp == (AMREX_SPACEDIM*idim));
296 #elif (defined WARPX_DIM_1D_Z)
297  // For 1D : icomp=0 and icomp=1 (Bx and By are not normal to the z boundary)
298  // The logic below ensures that the flags are set right for 1D
299  const bool is_normal_to_PEC = (icomp == (idim+2));
300 #else
301  const bool is_normal_to_PEC = (icomp == idim);
302 #endif
303 
304  // grid point ijk_vec is ig number of points pass the
305  // domain boundary in direction, idim
306  const int ig = get_cell_count_to_boundary(
307  dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
308 
309  if (ig == 0) {
310  // Only normal component is set to 0
311  if (is_normal_to_PEC && is_nodal[idim]==1) {
312  OnPECBoundary = true;
313  }
314  } else if ( ig > 0) {
315  // Mirror location inside the domain by "ig" number of cells
316  // across PEC boundary in direction, idim, and side, iside
317  ijk_mirror[idim] = ( (iside == 0)
318  ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
319  : (dom_hi[idim] + 1 - ig));
320  GuardCell = true;
321  // Sign of the normal component in guard cell is inverted
322  if (is_normal_to_PEC) sign *= -1._rt;
323  }
324  } // if PEC Boundary
325  } // loop over sides
326  } // loop of dimensions
327 
328  if (OnPECBoundary) {
329  // if ijk_vec is on a PEC boundary in any direction, set Bnormal to 0.
330  Bfield(ijk_vec,n) = 0._rt;
331  } else if (GuardCell) {
332  // Bnormal and Btangential is set opposite and equal to the value
333  // in the mirror location, respectively.
334  Bfield(ijk_vec,n) = sign * Bfield(ijk_mirror,n);
335  }
336  }
337 
338 
360  void SetRhoOrJfieldFromPEC (const int n,
361  const amrex::IntVect & ijk_vec,
363  amrex::GpuArray<GpuArray<int, 2>, AMREX_SPACEDIM> const& mirrorfac,
364  amrex::GpuArray<GpuArray<amrex::Real, 2>, AMREX_SPACEDIM> const& psign,
365  amrex::GpuArray<GpuArray<bool, 2>, AMREX_SPACEDIM> const& is_pec,
366  amrex::GpuArray<bool, AMREX_SPACEDIM> const& tangent_to_bndy,
367  amrex::Box const& fabbox)
368  {
369  // The boundary is handled in 2 steps:
370  // 1) The cells internal to the domain are updated using the
371  // current deposited in the guard cells
372  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
373  {
374  for (int iside = 0; iside < 2; ++iside)
375  {
376  if (!is_pec[idim][iside]) continue;
377 
378  // Get the mirror guard cell index
379  amrex::IntVect iv_mirror = ijk_vec;
380  iv_mirror[idim] = mirrorfac[idim][iside] - ijk_vec[idim];
381 
382  // On the PEC boundary the charge/current density is set to 0
383  if (ijk_vec == iv_mirror) field(ijk_vec, n) = 0._rt;
384  // otherwise update the internal cell if the mirror guard cell exists
385  else if (fabbox.contains(iv_mirror))
386  {
387  field(ijk_vec,n) += psign[idim][iside] * field(iv_mirror,n);
388  }
389  }
390  }
391  // 2) The guard cells are updated with the appropriate image
392  // charge based on the charge/current in the valid cells
393  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
394  {
395  for (int iside = 0; iside < 2; ++iside)
396  {
397  if (!is_pec[idim][iside]) continue;
398 
399  amrex::IntVect iv_mirror = ijk_vec;
400  iv_mirror[idim] = mirrorfac[idim][iside] - ijk_vec[idim];
401  if (ijk_vec != iv_mirror && fabbox.contains(iv_mirror))
402  {
403  if (tangent_to_bndy[idim])
404  field(iv_mirror, n) = -field(ijk_vec, n);
405  else
406  field(iv_mirror, n) = field(ijk_vec, n);
407  }
408  }
409  }
410  }
411 
412 
428  void SetNeumannOnPEC (const int n,
429  const amrex::IntVect & ijk_vec,
431  amrex::GpuArray<GpuArray<int, 2>, AMREX_SPACEDIM> const& mirrorfac,
432  amrex::GpuArray<GpuArray<bool, 2>, AMREX_SPACEDIM> const& is_pec,
433  amrex::Box const& fabbox )
434  {
435  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
436  {
437  for (int iside = 0; iside < 2; ++iside)
438  {
439  if (!is_pec[idim][iside]) continue;
440 
441  // Get the mirror guard cell index
442  amrex::IntVect iv_mirror = ijk_vec;
443  iv_mirror[idim] = mirrorfac[idim][iside] - ijk_vec[idim];
444 
445  // On the PEC boundary the field value is set equal to the
446  // first value in the domain (nodal fields)
447  if (ijk_vec == iv_mirror) {
448  iv_mirror[idim] += (iside == 0) ? 1 : -1;
449  if (fabbox.contains(iv_mirror)) field(ijk_vec, n) = field(iv_mirror, n);
450  }
451  // otherwise set the mirror guard cell equal to the internal cell value
452  else if (fabbox.contains(iv_mirror))
453  {
454  field(iv_mirror, n) = field(ijk_vec, n);
455  }
456  }
457  }
458  }
459 
460 
462  bool isAnyBoundaryPEC();
474  void ApplyPECtoEfield ( std::array<amrex::MultiFab*, 3> Efield,
475  int lev, PatchType patch_type,
476  bool split_pml_field = false);
486  void ApplyPECtoBfield ( std::array<amrex::MultiFab*, 3> Bfield,
487  int lev, PatchType patch_type);
488 
497  void ApplyPECtoRhofield(amrex::MultiFab* rho, int lev,
498  PatchType patch_type);
499 
509  amrex::MultiFab* Jz, int lev,
510  PatchType patch_type);
511 
520  int lev, PatchType patch_type);
521 }
522 #endif // WarpX_PEC_KERNELS_H_
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
PatchType
Definition: WarpX.H:78
@ 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, 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:360
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void SetNeumannOnPEC(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< bool, 2 >, AMREX_SPACEDIM > const &is_pec, amrex::Box const &fabbox)
This function sets the given field value on a PEC boundary to enforce a Neumann boundary condition (z...
Definition: WarpX_PEC.H:428
void ApplyPECtoJfield(amrex::MultiFab *Jx, amrex::MultiFab *Jy, amrex::MultiFab *Jz, int lev, PatchType patch_type)
Reflects current density deposited over the PEC boundary back into the simulation domain.
Definition: WarpX_PEC.cpp:305
void ApplyPECtoBfield(std::array< amrex::MultiFab *, 3 > Bfield, 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
void ApplyPECtoElectronPressure(amrex::MultiFab *Pefield, int lev, PatchType patch_type)
Apply the PEC boundary to the electron pressure field.
Definition: WarpX_PEC.cpp:496
void ApplyPECtoEfield(std::array< amrex::MultiFab *, 3 > Efield, int lev, PatchType patch_type, 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
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
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:271
bool isAnyBoundaryPEC()
Definition: WarpX_PEC.cpp:19
int n
Definition: run_libensemble_on_warpx.py:70
string field
Definition: video_yt.py:31
@ PEC
perfect electric conductor (PEC) with E_tangential=0
Definition: WarpXAlgorithmSelection.H:135