1 #ifndef WARPX_PEC_KERNELS_H_
2 #define WARPX_PEC_KERNELS_H_
9 #include <AMReX_Config.H>
21 using namespace amrex;
71 return ((iside == 0) ? (dom_lo[idim] - ijk_vec[idim])
72 : (ijk_vec[idim] - (dom_hi[idim] + is_nodal[idim])));
155 bool OnPECBoundary =
false;
156 bool GuardCell =
false;
157 amrex::Real sign = 1._rt;
159 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
161 for (
int iside = 0; iside < 2; ++iside) {
162 const bool isPECBoundary = ( (iside == 0)
165 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
169 const bool is_tangent_to_PEC = (icomp != AMREX_SPACEDIM*idim);
170 #elif (defined WARPX_DIM_1D_Z)
173 const bool is_tangent_to_PEC = (icomp != idim+2);
175 const bool is_tangent_to_PEC = (icomp != idim);
181 dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
184 if (is_tangent_to_PEC && is_nodal[idim] == 1) {
185 OnPECBoundary =
true;
189 ijk_mirror[idim] = ( ( iside == 0)
190 ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
191 : (dom_hi[idim] + 1 - ig));
194 if (is_tangent_to_PEC) sign *= -1._rt;
201 Efield(ijk_vec,
n) = 0._rt;
202 }
else if (GuardCell) {
203 Efield(ijk_vec,
n) = sign * Efield(ijk_mirror,
n);
280 bool OnPECBoundary =
false;
281 bool GuardCell =
false;
282 amrex::Real sign = 1._rt;
284 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
286 for (
int iside = 0; iside < 2; ++iside) {
287 const bool isPECBoundary = ( (iside == 0 )
291 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
295 const bool is_normal_to_PEC = (icomp == (AMREX_SPACEDIM*idim));
296 #elif (defined WARPX_DIM_1D_Z)
299 const bool is_normal_to_PEC = (icomp == (idim+2));
301 const bool is_normal_to_PEC = (icomp == idim);
307 dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
311 if (is_normal_to_PEC && is_nodal[idim]==1) {
312 OnPECBoundary =
true;
314 }
else if ( ig > 0) {
317 ijk_mirror[idim] = ( (iside == 0)
318 ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
319 : (dom_hi[idim] + 1 - ig));
322 if (is_normal_to_PEC) sign *= -1._rt;
330 Bfield(ijk_vec,
n) = 0._rt;
331 }
else if (GuardCell) {
334 Bfield(ijk_vec,
n) = sign * Bfield(ijk_mirror,
n);
372 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim)
374 for (
int iside = 0; iside < 2; ++iside)
376 if (!is_pec[idim][iside])
continue;
380 iv_mirror[idim] = mirrorfac[idim][iside] - ijk_vec[idim];
383 if (ijk_vec == iv_mirror)
field(ijk_vec,
n) = 0._rt;
385 else if (fabbox.
contains(iv_mirror))
387 field(ijk_vec,
n) += psign[idim][iside] *
field(iv_mirror,
n);
393 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim)
395 for (
int iside = 0; iside < 2; ++iside)
397 if (!is_pec[idim][iside])
continue;
400 iv_mirror[idim] = mirrorfac[idim][iside] - ijk_vec[idim];
401 if (ijk_vec != iv_mirror && fabbox.
contains(iv_mirror))
403 if (tangent_to_bndy[idim])
435 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim)
437 for (
int iside = 0; iside < 2; ++iside)
439 if (!is_pec[idim][iside])
continue;
443 iv_mirror[idim] = mirrorfac[idim][iside] - ijk_vec[idim];
447 if (ijk_vec == iv_mirror) {
448 iv_mirror[idim] += (iside == 0) ? 1 : -1;
452 else if (fabbox.
contains(iv_mirror))
476 bool split_pml_field =
false);
#define AMREX_FORCE_INLINE
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