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 
94  void SetEfieldOnPEC (const int icomp, const amrex::IntVect & dom_lo,
95  const amrex::IntVect &dom_hi,
96  const amrex::IntVect &ijk_vec, const int n,
97  amrex::Array4<amrex::Real> const& Efield,
98  const amrex::IntVect& is_nodal,
99  amrex::GpuArray<int, 3> const& fbndry_lo,
100  amrex::GpuArray<int, 3> const& fbndry_hi )
101  {
102  // Tangential Efield componentes in guard cells set equal and opposite to cells
103  // in the mirror locations across the PEC boundary, whereas normal E-field
104  // components are set equal to values in the mirror locations across the PEC
105  // boundary. Here we just initialize it.
106  amrex::IntVect ijk_mirror = ijk_vec;
107  bool OnPECBoundary = false;
108  bool GuardCell = false;
109  amrex::Real sign = 1._rt;
110  // Loop over all the dimensions
111  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
112  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
113  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
114  for (int iside = 0; iside < 2; ++iside) {
115  const bool isPECBoundary = ( (iside == 0)
116  ? is_boundary_PEC(fbndry_lo, idim)
117  : is_boundary_PEC(fbndry_hi, idim) );
118 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
119  // For 2D : for icomp==1, (Ey in XZ, Etheta in RZ),
120  // icomp=1 is tangential to both x and z boundaries
121  // The logic below ensures that the flags are set right for 2D
122  const bool is_tangent_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
123  ? false : true );
124 #else
125  const bool is_tangent_to_PEC = ( ( icomp == idim) ? false : true );
126 #endif
127  if (isPECBoundary == true) {
128  // guard cell outside the domain by "ig" number cells, in direction, idim
129  const int ig = ( (iside == 0)
130  ? (dom_lo[idim] - ijk_vec[idim])
131  : (ijk_vec[idim] - (dom_hi[idim] + is_nodal[idim]) ) );
132  if (ig == 0) {
133  if (is_tangent_to_PEC == true and is_nodal[idim] == 1) {
134  OnPECBoundary = true;
135  }
136  } else if (ig > 0) {
137  // Find mirror location across PEC boundary
138  ijk_mirror[idim] = ( ( iside == 0)
139  ? (dom_lo[idim] + ig)
140  : (dom_hi[idim] + is_nodal[idim] - ig));
141  GuardCell = true;
142  // tangential components are inverted across PEC boundary
143  if (is_tangent_to_PEC == true) sign *= -1._rt;
144  }
145  } // is PEC boundary
146  } // loop over iside
147  } // loop over dimensions
148  if (OnPECBoundary == true) {
149  // if ijk_vec is on a PEC boundary in any direction, set Etangential to 0.
150  Efield(ijk_vec,n) = 0._rt;
151  } else if (GuardCell == true) {
152  Efield(ijk_vec,n) = sign * Efield(ijk_mirror,n);
153  }
154  }
155 
210  void SetBfieldOnPEC (const int icomp, const amrex::IntVect & dom_lo,
211  const amrex::IntVect & dom_hi,
212  const amrex::IntVect & ijk_vec, const int n,
213  amrex::Array4<amrex::Real> const& Bfield,
214  const amrex::IntVect & is_nodal,
215  amrex::GpuArray<int, 3> const& fbndry_lo,
216  amrex::GpuArray<int, 3> const& fbndry_hi )
217  {
218  amrex::IntVect ijk_mirror = ijk_vec;
219  bool OnPECBoundary = false;
220  bool GuardCell = false;
221  amrex::Real sign = 1._rt;
222  // Loop over all dimensions
223  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
224  // Loop over sides, iside = 0 (lo), iside = 1 (hi)
225  for (int iside = 0; iside < 2; ++iside) {
226  const bool isPECBoundary = ( (iside == 0 )
227  ? is_boundary_PEC(fbndry_lo, idim)
228  : is_boundary_PEC(fbndry_hi, idim) );
229  if (isPECBoundary == true) {
230 #if (defined WARPX_DIM_XZ) || (defined WARPX_DIM_RZ)
231  // For 2D : for icomp==1, (By in XZ, Btheta in RZ),
232  // icomp=1 is not normal to x or z boundary
233  // The logic below ensures that the flags are set right for 2D
234  const bool is_normal_to_PEC = ( (icomp == AMREX_SPACEDIM*idim)
235  ? true : false );
236 #else
237  const bool is_normal_to_PEC = ( ( icomp == idim) ? true : false );
238 #endif
239  // guard cell outside the domain by "ig" number cells, in direction, idim
240  const int ig = ( (iside == 0)
241  ? (dom_lo[idim] - ijk_vec[idim])
242  : (ijk_vec[idim] - (dom_hi[idim] + is_nodal[idim]) ) );
243  if (ig == 0) {
244  // Only normal component is set to 0
245  if (is_normal_to_PEC == true and is_nodal[idim]==1) {
246  OnPECBoundary = true;
247  }
248  } else if ( ig > 0) {
249  // Mirror location inside the domain by "ig" number of cells
250  // across PEC boundary in direction, idim, and side, iside
251  ijk_mirror[idim] = ( (iside == 0)
252  ? (dom_lo[idim] + ig)
253  : (dom_hi[idim] + is_nodal[idim] - ig));
254  GuardCell = true;
255  // Sign of the normal component in guard cell is inverted
256  if (is_normal_to_PEC == true) sign *= -1._rt;
257  }
258  } // if PEC Boundary
259  } // loop over sides
260  } // loop of dimensions
261 
262  if (OnPECBoundary == true) {
263  // if ijk_vec is on a PEC boundary in any direction, set Bnormal to 0.
264  Bfield(ijk_vec,n) = 0._rt;
265  } else if (GuardCell == true) {
266  // Bnormal and Btangential is set opposite and equal to the value
267  // in the mirror location, respectively.
268  Bfield(ijk_vec,n) = sign * Bfield(ijk_mirror,n);
269  }
270  }
271 
273  bool isAnyBoundaryPEC();
285  void ApplyPECtoEfield ( std::array<amrex::MultiFab*, 3> Efield,
286  const int lev, PatchType patch_type,
287  const bool split_pml_field = false);
297  void ApplyPECtoBfield ( std::array<amrex::MultiFab*, 3> Bfield,
298  const int lev, PatchType patch_type);
299 }
300 
301 #endif // WarpX_PEC_KERNELS_H_
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:210
bool isAnyBoundaryPEC()
Definition: WarpX_PEC.cpp:18
perfect electric conductor (PEC) with E_tangential=0
Definition: WarpXAlgorithmSelection.H:117
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
int n
Definition: run_libensemble_on_warpx.py:67
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:124
Definition: WarpX_PEC.H:20
PatchType
Definition: WarpX.H:69
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:94
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
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 is PEC in direction, dir, is PEC.
Definition: WarpX_PEC.H:32