ARTEMIS
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 
362  void SetRhofieldFromPEC (const amrex::IntVect & dom_lo,
363  const amrex::IntVect & dom_hi,
364  const amrex::IntVect & nguard,
365  const amrex::IntVect & ijk_vec,
366  const int n,
367  amrex::Array4<amrex::Real> const& rho,
368  const amrex::IntVect & is_nodal,
369  amrex::GpuArray<int, 3> const& fbndry_lo,
370  amrex::GpuArray<int, 3> const& fbndry_hi,
373  {
374  amrex::IntVect ijk_mirror;
375  // Loop over all dimensions
376  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
377  ijk_mirror = ijk_vec;
378  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
379  for (int iside = 0; iside < 2; ++iside) {
380  // check if the boundary in this direction is a PEC, otherwise
381  // continue to the next direction
382  if (!( (iside == 0 )
383  ? is_boundary_PEC(fbndry_lo, idim)
384  : is_boundary_PEC(fbndry_hi, idim) )) continue;
385 
386  // grid point ijk_vec is ig number of points before (hence the 1-)
387  // the domain boundary in direction, idim
388  // Note: For nodal fields ig=0 is on the boundary and will be
389  // handled in `SetRhoFieldOnPEC` below. For cell-centered
390  // fields ig=0 is the cell center next to the PEC (inside the
391  // domain) and so it mirrors the first ghost cell.
392  int ig = -1 * get_cell_count_to_boundary(
393  dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
394 
395  if (ig >= is_nodal[idim] && ig < nguard[idim] + is_nodal[idim]) {
396 
397  // Mirror location outside the domain by "ig" number of cells
398  // across PEC boundary in direction, idim, and side, iside
399  ijk_mirror[idim] = ( (iside == 0)
400  ? (dom_lo[idim] - (1 - is_nodal[idim]) - ig)
401  : (dom_hi[idim] + 1 + ig));
402 
403  const bool isReflectingBoundary = ( (iside == 0 )
404  ? is_boundary_reflecting(pbndry_lo, idim)
405  : is_boundary_reflecting(pbndry_hi, idim) );
406  if (isReflectingBoundary) {
407  // reflect density if the particles are reflected
408  rho(ijk_vec,n) += rho(ijk_mirror,n);
409  }
410  else {
411  // add image charge density (opposite sign) if
412  // particles are absorbed
413  rho(ijk_vec,n) -= rho(ijk_mirror,n);
414  }
415  }
416  } // loop over sides
417  } // loop of dimensions
418  }
419 
438  void SetRhofieldOnPEC (const amrex::IntVect & dom_lo,
439  const amrex::IntVect & dom_hi,
440  const amrex::IntVect & ijk_vec, const int n,
441  amrex::Array4<amrex::Real> const& rho,
442  const amrex::IntVect & is_nodal,
443  amrex::GpuArray<int, 3> const& fbndry_lo,
444  amrex::GpuArray<int, 3> const& fbndry_hi)
445  {
446  amrex::IntVect ijk_mirror = ijk_vec;
447  bool OnPECBoundary = false;
448  bool GuardCell = false;
449  // Loop over all the dimensions
450  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
451  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
452  for (int iside = 0; iside < 2; ++iside) {
453  // check if the boundary in this direction is a PEC, otherwise
454  // continue to the next direction
455  if (!( (iside == 0 )
456  ? is_boundary_PEC(fbndry_lo, idim)
457  : is_boundary_PEC(fbndry_hi, idim) )) continue;
458 
459  // grid point ijk_vec is ig number of points pass
460  // the domain boundary in direction, idim
462  dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
463 
464  if (ig == 0 && is_nodal[idim] == 1) OnPECBoundary = true;
465  else if (ig > 0) {
466  // Find mirror location across PEC boundary
467  ijk_mirror[idim] = ( (iside == 0)
468  ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
469  : (dom_hi[idim] + 1 - ig));
470  GuardCell = true;
471  }
472  } // loop over iside
473  } // loop over dimensions
474  if (OnPECBoundary == true) {
475  // if ijk_vec is on a PEC boundary, set rho to 0.
476  rho(ijk_vec,n) = 0._rt;
477  } else if (GuardCell == true) {
478  rho(ijk_vec,n) = -rho(ijk_mirror,n);
479  }
480  }
481 
506  void SetJfieldFromPEC (const int icomp,
507  const amrex::IntVect & dom_lo,
508  const amrex::IntVect & dom_hi,
509  const amrex::IntVect & nguard,
510  const amrex::IntVect & ijk_vec,
511  amrex::Array4<amrex::Real> const& Jfield,
512  const amrex::IntVect & is_nodal,
513  amrex::GpuArray<int, 3> const& fbndry_lo,
514  amrex::GpuArray<int, 3> const& fbndry_hi,
517  {
518  amrex::IntVect ijk_mirror;
519  // Loop over all dimensions
520  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
521  ijk_mirror = ijk_vec;
522  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
523  for (int iside = 0; iside < 2; ++iside) {
524  // check if the boundary in this direction is a PEC, otherwise
525  // continue to the next direction
526  if (!( (iside == 0 )
527  ? is_boundary_PEC(fbndry_lo, idim)
528  : is_boundary_PEC(fbndry_hi, idim) )) continue;
529 
530 #if (defined WARPX_DIM_1D_Z)
531  // For 1D : icomp=0 and icomp=1 (Ex and Ey are tangential to the z boundary)
532  // The logic below ensures that the flags are set right for 1D
533  const bool is_tangent_to_PEC = ( ( icomp == idim+2) ? false : true );
534 #elif (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
535  // For 2D : for icomp==1, (Ey in XZ, Etheta in RZ),
536  // icomp=1 is tangential to both x and z boundaries
537  // The logic below ensures that the flags are set right for 2D
538  const bool is_tangent_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
539  ? false : true );
540 #else
541  const bool is_tangent_to_PEC = ( ( icomp == idim) ? false : true );
542 #endif
543 
544  // grid point ijk_vec is ig number of points "before" (hence the -1)
545  // the domain boundary in direction, idim
546  int ig = -1 * get_cell_count_to_boundary(
547  dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
548 
549  if (ig >= is_nodal[idim] && ig < nguard[idim] + is_nodal[idim]) {
550 
551  // Mirror location outside the domain by "ig" number of cells
552  // across PEC boundary in direction, idim, and side, iside
553  ijk_mirror[idim] = ( (iside == 0)
554  ? (dom_lo[idim] - (1 - is_nodal[idim]) - ig)
555  : (dom_hi[idim] + 1 + ig));
556 
557  const bool isReflectingBoundary = ( (iside == 0 )
558  ? is_boundary_reflecting(pbndry_lo, idim)
559  : is_boundary_reflecting(pbndry_hi, idim) );
560  if (isReflectingBoundary) {
561  // reflect density if the particles are reflected
562  if (is_tangent_to_PEC) Jfield(ijk_vec) += Jfield(ijk_mirror);
563  else Jfield(ijk_vec) -= Jfield(ijk_mirror);
564  }
565  else {
566  // add image charge density (opposite sign) if
567  // particles are absorbed
568  if (is_tangent_to_PEC) Jfield(ijk_vec) -= Jfield(ijk_mirror);
569  else Jfield(ijk_vec) += Jfield(ijk_mirror);
570  }
571  }
572  } // loop over sides
573  } // loop of dimensions
574  }
575 
576 
599  void SetJfieldOnPEC (const int icomp,
600  const amrex::IntVect & dom_lo,
601  const amrex::IntVect & dom_hi,
602  const amrex::IntVect & ijk_vec,
603  amrex::Array4<amrex::Real> const& Jfield,
604  const amrex::IntVect & is_nodal,
605  amrex::GpuArray<int, 3> const& fbndry_lo,
606  amrex::GpuArray<int, 3> const& fbndry_hi )
607  {
608  amrex::IntVect ijk_mirror = ijk_vec;
609  bool OnPECBoundary = false;
610  bool GuardCell = false;
611  amrex::Real sign = 1._rt;
612 
613  // Loop over all dimensions
614  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
615  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
616  for (int iside = 0; iside < 2; ++iside) {
617  // check if the boundary in this direction is a PEC, otherwise
618  // continue to the next direction
619  if (!( (iside == 0 )
620  ? is_boundary_PEC(fbndry_lo, idim)
621  : is_boundary_PEC(fbndry_hi, idim) )) continue;
622 
623 #if (defined WARPX_DIM_1D_Z)
624  // For 1D : icomp=0 and icomp=1 (Ex and Ey are tangential to the z boundary)
625  // The logic below ensures that the flags are set right for 1D
626  const bool is_tangent_to_PEC = ( ( icomp == idim+2) ? false : true );
627 #elif (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
628  // For 2D : for icomp==1, (Ey in XZ, Etheta in RZ),
629  // icomp=1 is tangential to both x and z boundaries
630  // The logic below ensures that the flags are set right for 2D
631  const bool is_tangent_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
632  ? false : true );
633 #else
634  const bool is_tangent_to_PEC = ( ( icomp == idim) ? false : true );
635 #endif
636 
637  // grid point ijk_vec is ig number of points pass
638  // the domain boundary in direction, idim
639  // Note: For nodal fields ig=0 is on the boundary.
641  dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
642 
643  if (ig == 0 && is_nodal[idim] == 1) OnPECBoundary = true;
644  else if (ig > 0) {
645  // Find mirror location across PEC boundary
646  ijk_mirror[idim] = ( (iside == 0)
647  ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
648  : (dom_hi[idim] + 1 - ig));
649  GuardCell = true;
650  // Sign of the normal component in guard cell is inverted
651  if (is_tangent_to_PEC == true) sign *= -1._rt;
652  }
653  } // loop over sides
654  } // loop of dimensions
655 
656  if (OnPECBoundary == true) {
657  // if ijk_vec is on a PEC boundary in any direction, set Bnormal to 0.
658  Jfield(ijk_vec) = 0._rt;
659  } else if (GuardCell == true) {
660  // tangential current is set opposite to the value in the mirror
661  // location while the normal current is set equal.
662  Jfield(ijk_vec) = sign * Jfield(ijk_mirror);
663  }
664  }
665 
666 
668  bool isAnyBoundaryPEC();
680  void ApplyPECtoEfield ( std::array<amrex::MultiFab*, 3> Efield,
681  const int lev, PatchType patch_type,
682  const bool split_pml_field = false);
692  void ApplyPECtoBfield ( std::array<amrex::MultiFab*, 3> Bfield,
693  const int lev, PatchType patch_type);
694 
703  void ApplyPECtoRhofield(amrex::MultiFab* rho, const int lev,
704  PatchType patch_type);
705 
715  amrex::MultiFab* Jz, const int lev,
716  PatchType patch_type);
717 }
718 
719 #endif // WarpX_PEC_KERNELS_H_
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
PatchType
Definition: WarpX.H:74
@ Reflecting
particles are reflected
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:214
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void SetJfieldOnPEC(const int icomp, const amrex::IntVect &dom_lo, const amrex::IntVect &dom_hi, const amrex::IntVect &ijk_vec, amrex::Array4< amrex::Real > const &Jfield, const amrex::IntVect &is_nodal, amrex::GpuArray< int, 3 > const &fbndry_lo, amrex::GpuArray< int, 3 > const &fbndry_hi)
Sets the J field value in guard cells next to a PEC boundary. The current density on the PEC boundary...
Definition: WarpX_PEC.H:599
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void SetJfieldFromPEC(const int icomp, const amrex::IntVect &dom_lo, const amrex::IntVect &dom_hi, const amrex::IntVect &nguard, const amrex::IntVect &ijk_vec, amrex::Array4< amrex::Real > const &Jfield, const amrex::IntVect &is_nodal, amrex::GpuArray< int, 3 > const &fbndry_lo, amrex::GpuArray< int, 3 > const &fbndry_hi, amrex::GpuArray< ParticleBoundaryType, 3 > const &pbndry_lo, amrex::GpuArray< ParticleBoundaryType, 3 > const &pbndry_hi)
Sets the J field value in cells close to a PEC boundary. The current density deposited in the guard c...
Definition: WarpX_PEC.H:506
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void SetRhofieldOnPEC(const amrex::IntVect &dom_lo, const amrex::IntVect &dom_hi, const amrex::IntVect &ijk_vec, const int n, amrex::Array4< amrex::Real > const &rho, const amrex::IntVect &is_nodal, amrex::GpuArray< int, 3 > const &fbndry_lo, amrex::GpuArray< int, 3 > const &fbndry_hi)
Sets the rho field value in guard cells next to a PEC boundary. The charge density on the PEC boundar...
Definition: WarpX_PEC.H:438
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:123
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 SetRhofieldFromPEC(const amrex::IntVect &dom_lo, const amrex::IntVect &dom_hi, const amrex::IntVect &nguard, const amrex::IntVect &ijk_vec, const int n, amrex::Array4< amrex::Real > const &rho, const amrex::IntVect &is_nodal, amrex::GpuArray< int, 3 > const &fbndry_lo, amrex::GpuArray< int, 3 > const &fbndry_hi, amrex::GpuArray< ParticleBoundaryType, 3 > const &pbndry_lo, amrex::GpuArray< ParticleBoundaryType, 3 > const &pbndry_hi)
Sets the rho field value in cells close to a PEC boundary. The charge density deposited in the guard ...
Definition: WarpX_PEC.H:362
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:300
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:18
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:27
int n
Definition: run_libensemble_on_warpx.py:67
@ PEC
perfect electric conductor (PEC) with E_tangential=0
Definition: WarpXAlgorithmSelection.H:148