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) {
162 for (
int iside = 0; iside < 2; ++iside) {
163 const bool isPECBoundary = ( (iside == 0)
166 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
170 const bool is_tangent_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
172 #elif (defined WARPX_DIM_1D_Z)
175 const bool is_tangent_to_PEC = ( ( icomp == idim+2) ?
false :
true );
177 const bool is_tangent_to_PEC = ( ( icomp == idim) ?
false :
true );
179 if (isPECBoundary ==
true) {
183 dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
186 if (is_tangent_to_PEC ==
true and is_nodal[idim] == 1) {
187 OnPECBoundary =
true;
191 ijk_mirror[idim] = ( ( iside == 0)
192 ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
193 : (dom_hi[idim] + 1 - ig));
196 if (is_tangent_to_PEC ==
true) sign *= -1._rt;
201 if (OnPECBoundary ==
true) {
203 Efield(ijk_vec,
n) = 0._rt;
204 }
else if (GuardCell ==
true) {
205 Efield(ijk_vec,
n) = sign * Efield(ijk_mirror,
n);
281 bool OnPECBoundary =
false;
282 bool GuardCell =
false;
283 amrex::Real sign = 1._rt;
285 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
287 for (
int iside = 0; iside < 2; ++iside) {
288 const bool isPECBoundary = ( (iside == 0 )
291 if (isPECBoundary ==
true) {
292 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
296 const bool is_normal_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
298 #elif (defined WARPX_DIM_1D_Z)
301 const bool is_normal_to_PEC = ( ( icomp == idim+2) ?
true :
false );
303 const bool is_normal_to_PEC = ( ( icomp == idim) ?
true :
false );
309 dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
313 if (is_normal_to_PEC ==
true and is_nodal[idim]==1) {
314 OnPECBoundary =
true;
316 }
else if ( ig > 0) {
319 ijk_mirror[idim] = ( (iside == 0)
320 ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
321 : (dom_hi[idim] + 1 - ig));
324 if (is_normal_to_PEC ==
true) sign *= -1._rt;
330 if (OnPECBoundary ==
true) {
332 Bfield(ijk_vec,
n) = 0._rt;
333 }
else if (GuardCell ==
true) {
336 Bfield(ijk_vec,
n) = sign * Bfield(ijk_mirror,
n);
376 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
377 ijk_mirror = ijk_vec;
379 for (
int iside = 0; iside < 2; ++iside) {
393 dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
395 if (ig >= is_nodal[idim] && ig < nguard[idim] + is_nodal[idim]) {
399 ijk_mirror[idim] = ( (iside == 0)
400 ? (dom_lo[idim] - (1 - is_nodal[idim]) - ig)
401 : (dom_hi[idim] + 1 + ig));
403 const bool isReflectingBoundary = ( (iside == 0 )
406 if (isReflectingBoundary) {
408 rho(ijk_vec,
n) += rho(ijk_mirror,
n);
413 rho(ijk_vec,
n) -= rho(ijk_mirror,
n);
447 bool OnPECBoundary =
false;
448 bool GuardCell =
false;
450 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
452 for (
int iside = 0; iside < 2; ++iside) {
462 dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
464 if (ig == 0 && is_nodal[idim] == 1) OnPECBoundary =
true;
467 ijk_mirror[idim] = ( (iside == 0)
468 ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
469 : (dom_hi[idim] + 1 - ig));
474 if (OnPECBoundary ==
true) {
476 rho(ijk_vec,
n) = 0._rt;
477 }
else if (GuardCell ==
true) {
478 rho(ijk_vec,
n) = -rho(ijk_mirror,
n);
520 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
521 ijk_mirror = ijk_vec;
523 for (
int iside = 0; iside < 2; ++iside) {
530 #if (defined WARPX_DIM_1D_Z)
533 const bool is_tangent_to_PEC = ( ( icomp == idim+2) ?
false :
true );
534 #elif (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
538 const bool is_tangent_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
541 const bool is_tangent_to_PEC = ( ( icomp == idim) ?
false :
true );
547 dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
549 if (ig >= is_nodal[idim] && ig < nguard[idim] + is_nodal[idim]) {
553 ijk_mirror[idim] = ( (iside == 0)
554 ? (dom_lo[idim] - (1 - is_nodal[idim]) - ig)
555 : (dom_hi[idim] + 1 + ig));
557 const bool isReflectingBoundary = ( (iside == 0 )
560 if (isReflectingBoundary) {
562 if (is_tangent_to_PEC) Jfield(ijk_vec) += Jfield(ijk_mirror);
563 else Jfield(ijk_vec) -= Jfield(ijk_mirror);
568 if (is_tangent_to_PEC) Jfield(ijk_vec) -= Jfield(ijk_mirror);
569 else Jfield(ijk_vec) += Jfield(ijk_mirror);
609 bool OnPECBoundary =
false;
610 bool GuardCell =
false;
611 amrex::Real sign = 1._rt;
614 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
616 for (
int iside = 0; iside < 2; ++iside) {
623 #if (defined WARPX_DIM_1D_Z)
626 const bool is_tangent_to_PEC = ( ( icomp == idim+2) ?
false :
true );
627 #elif (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
631 const bool is_tangent_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
634 const bool is_tangent_to_PEC = ( ( icomp == idim) ?
false :
true );
641 dom_lo, dom_hi, ijk_vec, is_nodal, idim, iside);
643 if (ig == 0 && is_nodal[idim] == 1) OnPECBoundary =
true;
646 ijk_mirror[idim] = ( (iside == 0)
647 ? (dom_lo[idim] + ig - (1 - is_nodal[idim]))
648 : (dom_hi[idim] + 1 - ig));
651 if (is_tangent_to_PEC ==
true) sign *= -1._rt;
656 if (OnPECBoundary ==
true) {
658 Jfield(ijk_vec) = 0._rt;
659 }
else if (GuardCell ==
true) {
662 Jfield(ijk_vec) = sign * Jfield(ijk_mirror);
682 const bool split_pml_field =
false);
#define AMREX_FORCE_INLINE
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