ARTEMIS
Functions
PEC Namespace Reference

Functions

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. More...
 
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, is reflecting. More...
 
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.e. a value of +1 means the current cell is outside of the simulation domain by 1 cell. Note that the high side domain boundary is between cell dom_hi and dom_hi+1 for cell centered grids and on cell dom_hi+1 for nodal grid. This is why (dom_hi[idim] + is_nodal[idim]) is used below. More...
 
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 components in the guard cells outside the domain boundary are set equal and opposite to the field in the valid cells at their mirrored locations. The normal Efield components in the guard cells are set equal to the field in the valid cells at their mirrored locations. The number or depth of guard cells updated is equal to the shape factor of particles in each dimension. For corner cells with mixed boundaries, the mirror location could be outside valid region, while still ensuring PEC condition is maintained across the PEC boundary, and the necessary sign change is accounted for depending on if the component, icomp, is tangential or normal to the PEC boundary. More...
 
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 value of the guard cells outside the domain boundary are set equal (and opposite) to the respective field components in the valid cells at their mirrored locations. The number or depth of guard cells updated is equal to the shape factor of particles in each dimension. More...
 
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 cells are either reflected back into the simulation domain (if a reflecting particle boundary is used), or the opposite charge density is deposited back in the domain to capture the effect of an image charge. More...
 
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 boundary is set to 0 while values in the guard cells are set equal and opposite to their mirror location inside the domain - representing image charges. More...
 
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 cells are either reflected back into the simulation domain (if a reflecting particle boundary is used), or the opposite current density is deposited back in the domain to capture the effect of an image charge. More...
 
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 is set to 0 while values in the guard cells are set equal (and opposite) to their mirror location inside the domain - representing image charges - in the normal (tangential) direction. More...
 
bool isAnyBoundaryPEC ()
 
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 and opposite to the valid cell field value at the respective mirror locations. More...
 
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 are set equal and opposite to the valid cell field value at the respective mirror locations. More...
 
void ApplyPECtoRhofield (amrex::MultiFab *rho, const int lev, PatchType patch_type)
 Reflects charge density deposited over the PEC boundary back into the simulation domain. More...
 
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. More...
 

Function Documentation

◆ ApplyPECtoBfield()

void PEC::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 are set equal and opposite to the valid cell field value at the respective mirror locations.

Parameters
[in,out]BfieldBoundary values of normal Bfield are set to zero.
[in]levlevel of the Multifab
[in]patch_typecoarse or fine

◆ ApplyPECtoEfield()

void PEC::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 and opposite to the valid cell field value at the respective mirror locations.

Parameters
[in,out]EfieldBoundary values of tangential Efield are set to zero
[in]levlevel of the Multifab
[in]patch_typecoarse or fine
[in]split_pml_fieldwhether pml the multifab is the regular Efield or split pml field

◆ ApplyPECtoJfield()

void PEC::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.

Parameters
[in,out]Jx,Jy,JzMultifabs containing the current density
[in]levlevel of the Multifab
[in]patch_typecoarse or fine

◆ ApplyPECtoRhofield()

void PEC::ApplyPECtoRhofield ( amrex::MultiFab rho,
const int  lev,
PatchType  patch_type 
)

Reflects charge density deposited over the PEC boundary back into the simulation domain.

Parameters
[in,out]rhoMultifab containing the charge density
[in]levlevel of the Multifab
[in]patch_typecoarse or fine

◆ get_cell_count_to_boundary()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE int PEC::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.e. a value of +1 means the current cell is outside of the simulation domain by 1 cell. Note that the high side domain boundary is between cell dom_hi and dom_hi+1 for cell centered grids and on cell dom_hi+1 for nodal grid. This is why (dom_hi[idim] + is_nodal[idim]) is used below.

Parameters
[in]dom_lo,dom_hiDomain boundaries
[in]ijk_vecCell coordinates
[in]is_nodalWhether the field of interest is nodal
[in]idimDimension of interest
[in]iside0 for low and 1 for high
Returns
number of grid points to the boundary

◆ is_boundary_PEC()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool PEC::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.

Parameters
[in]fboundaryValue containing boundary type
[in]dirdirection
Returns
1 if the boundary type is PEC else 0

◆ is_boundary_reflecting()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool PEC::is_boundary_reflecting ( amrex::GpuArray< ParticleBoundaryType, 3 > const &  pboundary,
int  dir 
)

Determines if the particle boundary condition stored in pboundary in direction, dir, is reflecting.

Parameters
[in]pboundaryValue containing boundary type
[in]dirdirection
Returns
1 if the boundary type is reflecting else 0

◆ isAnyBoundaryPEC()

bool PEC::isAnyBoundaryPEC ( )

Returns 1 if any domain boundary is set to PEC, else returns 0.

◆ SetBfieldOnPEC()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void PEC::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 value of the guard cells outside the domain boundary are set equal (and opposite) to the respective field components in the valid cells at their mirrored locations. The number or depth of guard cells updated is equal to the shape factor of particles in each dimension.

For 3D : x component is tangential to the y-boundary and z-boundary y component is tangential to the x-boundary and z-boundary z component is tangential to the x-boundary and y-boundary x component is normal to the x-boundary y component is normal to the y-boundary z component is normal to the z-boundary where, x-boundary is the yz-plane at x=xmin and x=xmax y-boundary is the xz-plane at y=ymin and y=ymax z-boundary is the xy-plane at z=zmin and z=zmax

For 2D : WarpX uses X-Z as the two dimensions x component is tangential to the z-boundary y component is tangential to the x-boundary and z-boundary z component is tangential to the x-boundary x component is normal to the x-boundary y component is not normal to any boundary (Only xz dimensions in 2D) z component is normal to the z-boundary where, x-boundary is along the line z at x=xmin and x=xmax z-boundary is along the line x at z=zmin and z=zmax

For 1D : WarpX uses Z as the only dimension x component is tangential to the z-boundary y component is tangential to the z-boundary z component is not tangential to the z-boundary x component is not normal to any boundary (Only z dimension in 1D) y component is not normal to any boundary (Only z dimension in 1D) z component is normal to the z-boundary where, z-boundary is a point at z=zmin and z=zmax

For RZ : WarpX uses R-Z as the two dimensions r component is tangential to the z-boundary theta_component is tangential to the r-boundary and z-boundary z component is tangential to the r-boundary r component is normal to the r-boundary theta_component is not normal to any boundary (on RZ dimensions are modeled) z component is normal to the z-boundary where, r-boundary is along the line z at r=rmin and r=rmax z-boundary is along the line r at z=zmin and z=zmax

Parameters
[in]icompcomponent of the Bfield being updated (0=x, 1=y, 2=z in Cartesian) (0=r, 1=theta, 2=z in RZ)
[in]dom_loindex value of the lower domain boundary (cell-centered)
[in]dom_hiindex value of the higher domain boundary (cell-centered)
[in]ijk_vecindices along the x(i), y(j), z(k) of Efield Array4
[in]nindex of the MultiFab component being updated
[in]Bfieldfield data to be updated if (ijk) is at the boundary or a guard cell
[in]is_nodalstaggering of the field data being updated.
[in]fbndry_loField boundary type at the lower boundaries
[in]fbndry_hiField boundary type at the upper boundaries

◆ SetEfieldOnPEC()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void PEC::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 components in the guard cells outside the domain boundary are set equal and opposite to the field in the valid cells at their mirrored locations. The normal Efield components in the guard cells are set equal to the field in the valid cells at their mirrored locations. The number or depth of guard cells updated is equal to the shape factor of particles in each dimension. For corner cells with mixed boundaries, the mirror location could be outside valid region, while still ensuring PEC condition is maintained across the PEC boundary, and the necessary sign change is accounted for depending on if the component, icomp, is tangential or normal to the PEC boundary.

For 3D : x component is tangential to the y-boundary and z-boundary y component is tangential to the x-boundary and z-boundary z component is tangential to the x-boundary and y-boundary x component is normal to the x-boundary y component is normal to the y-boundary z component is normal to the z-boundary where, x-boundary is the yz-plane at x=xmin and x=xmax y-boundary is the xz-plane at y=ymin and y=ymax z-boundary is the xy-plane at z=zmin and z=zmax

For 2D : WarpX uses X-Z as the two dimensions x component is tangential to the z-boundary y component is tangential to the x-boundary and z-boundary z component is tangential to the x-boundary x component is normal to the x-boundary y component is not normal to any boundary (Only xz dimensions in 2D) z component is normal to the z-boundary where, x-boundary is along the line z at x=xmin and x=xmax z-boundary is along the line x at z=zmin and z=zmax

For 1D : WarpX uses Z as the only dimension x component is tangential to the z-boundary y component is tangential to the z-boundary z component is not tangential to the z-boundary x component is not normal to any boundary (Only z dimension in 1D) y component is not normal to any boundary (Only z dimension in 1D) z component is normal to the z-boundary where, z-boundary is a point at z=zmin and z=zmax

For RZ : WarpX uses R-Z as the two dimensions r component is tangential to the z-boundary theta_component is tangential to the r-boundary and z-boundary z component is tangential to the r-boundary r component is normal to the r-boundary theta_component is not normal to any boundary (on RZ dimensions are modeled) z component is normal to the z-boundary where, r-boundary is along the line z at r=rmin and r=rmax z-boundary is along the line r at z=zmin and z=zmax

Parameters
[in]icompcomponent of the Efield being updated (0=x, 1=y, 2=z in Cartesian) (0=r, 1=theta, 2=z in RZ)
[in]dom_loindex value of the lower domain boundary (cell-centered)
[in]dom_hiindex value of the higher domain boundary (cell-centered)
[in]ijk_vecindices along the x(i), y(j), z(k) of Efield Array4
[in]nindex of the MultiFab component being updated
[in]Efieldfield data to be updated if (ijk) is at the boundary or a guard cell
[in]is_nodalstaggering of the field data being updated.
[in]fbndry_loField boundary type at the lower boundaries
[in]fbndry_hiField boundary type at the upper boundaries

◆ SetJfieldFromPEC()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void PEC::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 cells are either reflected back into the simulation domain (if a reflecting particle boundary is used), or the opposite current density is deposited back in the domain to capture the effect of an image charge.

Parameters
[in]icompcomponent of the Jfield being updated (0=x, 1=y, 2=z in Cartesian) (0=r, 1=theta, 2=z in RZ)
[in]dom_loindex value of the lower domain boundary (cell-centered)
[in]dom_hiindex value of the higher domain boundary (cell-centered)
[in]nguardnumber of guard cells along each dimension
[in]ijk_vecindices along the x(i), y(j), z(k) of the rho Array4
[in]Jfieldfield data to be updated if (ijk) is at the boundary or a guard cell
[in]is_nodalstaggering of the field data being updated.
[in]fbndry_loField boundary type at the lower boundaries
[in]fbndry_hiField boundary type at the upper boundaries
[in]pbndry_loParticle boundary type at the lower boundaries
[in]pbndry_hiParticle boundary type at the upper boundaries

◆ SetJfieldOnPEC()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void PEC::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 is set to 0 while values in the guard cells are set equal (and opposite) to their mirror location inside the domain - representing image charges - in the normal (tangential) direction.

Parameters
[in]icompcomponent of the Jfield being updated (0=x, 1=y, 2=z in Cartesian) (0=r, 1=theta, 2=z in RZ)
[in]dom_loindex value of the lower domain boundary (cell-centered)
[in]dom_hiindex value of the higher domain boundary (cell-centered)
[in]nguardnumber of guard cells along each dimension
[in]ijk_vecindices along the x(i), y(j), z(k) of the rho Array4
[in]Jfieldfield data to be updated if (ijk) is at the boundary or a guard cell
[in]is_nodalstaggering of the field data being updated.
[in]fbndry_loField boundary type at the lower boundaries
[in]fbndry_hiField boundary type at the upper boundaries

◆ SetRhofieldFromPEC()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void PEC::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 cells are either reflected back into the simulation domain (if a reflecting particle boundary is used), or the opposite charge density is deposited back in the domain to capture the effect of an image charge.

Parameters
[in]dom_loindex value of the lower domain boundary (cell-centered)
[in]dom_hiindex value of the higher domain boundary (cell-centered)
[in]nguardnumber of guard cells along each dimension
[in]ijk_vecindices along the x(i), y(j), z(k) of the rho Array4
[in]nindex of the MultiFab component being updated
[in]rhofield data to be updated if (ijk) is at the boundary or a guard cell
[in]is_nodalstaggering of the field data being updated.
[in]fbndry_loField boundary type at the lower boundaries
[in]fbndry_hiField boundary type at the upper boundaries
[in]pbndry_loParticle boundary type at the lower boundaries
[in]pbndry_hiParticle boundary type at the upper boundaries

◆ SetRhofieldOnPEC()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void PEC::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 boundary is set to 0 while values in the guard cells are set equal and opposite to their mirror location inside the domain - representing image charges.

Parameters
[in]dom_loindex value of the lower domain boundary (cell-centered)
[in]dom_hiindex value of the higher domain boundary (cell-centered)
[in]ijk_vecindices along the x(i), y(j), z(k) of Efield Array4
[in]nindex of the MultiFab component being updated
[in]rhofield data to be updated if (ijk) is at the boundary or a guard cell
[in]is_nodalstaggering of the field data being updated.
[in]fbndry_loField boundary type at the lower boundaries
[in]fbndry_hiField boundary type at the upper boundaries