WarpX
Loading...
Searching...
No Matches
PEC Namespace Reference

Functions

void ApplyPECtoEfield (std::array< amrex::MultiFab *, 3 > Efield, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_lo, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_hi, FieldBoundaryType bc_type, const amrex::IntVect &ng_fieldgather, const amrex::Geometry &geom, int lev, PatchType patch_type, const amrex::Vector< amrex::IntVect > &ref_ratios, 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.
 
void ApplyPECtoBfield (std::array< amrex::MultiFab *, 3 > Bfield, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_lo, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_hi, FieldBoundaryType bc_type, const amrex::IntVect &ng_fieldgather, const amrex::Geometry &geom, int lev, PatchType patch_type, const amrex::Vector< amrex::IntVect > &ref_ratios, bool split_pml_field=false)
 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.
 
void ApplyReflectiveBoundarytoRhofield (amrex::MultiFab *rho, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_lo, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_hi, const amrex::Array< ParticleBoundaryType, 3 > &particle_boundary_lo, const amrex::Array< ParticleBoundaryType, 3 > &particle_boundary_hi, const amrex::Geometry &geom, int lev, PatchType patch_type, const amrex::Vector< amrex::IntVect > &ref_ratios)
 Reflects charge density deposited over the PEC boundary back into the simulation domain.
 
void ApplyReflectiveBoundarytoJfield (amrex::MultiFab *Jx, amrex::MultiFab *Jy, amrex::MultiFab *Jz, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_lo, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_hi, const amrex::Array< ParticleBoundaryType, 3 > &particle_boundary_lo, const amrex::Array< ParticleBoundaryType, 3 > &particle_boundary_hi, const amrex::Geometry &geom, int lev, PatchType patch_type, const amrex::Vector< amrex::IntVect > &ref_ratios)
 Reflects current density deposited over the PEC boundary back into the simulation domain.
 
void ApplyPECtoElectronPressure (amrex::MultiFab *Pefield, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_lo, const amrex::Array< FieldBoundaryType, 3 > &field_boundary_hi, const amrex::Geometry &geom, int lev, PatchType patch_type, const amrex::Vector< amrex::IntVect > &ref_ratios)
 Apply the PEC boundary to the electron pressure field.
 

Function Documentation

◆ ApplyPECtoBfield()

void PEC::ApplyPECtoBfield ( std::array< amrex::MultiFab *, 3 > Bfield,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_lo,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_hi,
FieldBoundaryType bc_type,
const amrex::IntVect & ng_fieldgather,
const amrex::Geometry & geom,
int lev,
PatchType patch_type,
const amrex::Vector< amrex::IntVect > & ref_ratios,
bool split_pml_field = false )

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]field_boundary_loBoundary types of the "low" field boundaries
[in]field_boundary_hiBoundary types of the "high" field boundaries
[in]bc_typeBoundary condition type to be applied
[in]ng_fieldgathernumber of guard cells used by field gather
[in]geomgeometry object of level "lev"
[in]levlevel of the Multifab
[in]patch_typecoarse or fine
[in]ref_ratiosvector containing the refinement ratios of the refinement levels
[in]split_pml_fieldwhether the MultiFab is the regular Bfield or the split PML field

◆ ApplyPECtoEfield()

void PEC::ApplyPECtoEfield ( std::array< amrex::MultiFab *, 3 > Efield,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_lo,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_hi,
FieldBoundaryType bc_type,
const amrex::IntVect & ng_fieldgather,
const amrex::Geometry & geom,
int lev,
PatchType patch_type,
const amrex::Vector< amrex::IntVect > & ref_ratios,
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]field_boundary_loBoundary types of the "low" boundaries
[in]field_boundary_hiBoundary types of the "high" boundaries
[in]bc_typeBoundary condition type to be applied
[in]ng_fieldgathernumber of guard cells used by field gather
[in]geomgeometry object of level "lev"
[in]levlevel of the Multifab
[in]patch_typecoarse or fine
[in]ref_ratiosvector containing the refinement ratios of the refinement levels
[in]split_pml_fieldwhether the MultiFab is the regular Efield or the split PML field

◆ ApplyPECtoElectronPressure()

void PEC::ApplyPECtoElectronPressure ( amrex::MultiFab * Pefield,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_lo,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_hi,
const amrex::Geometry & geom,
int lev,
PatchType patch_type,
const amrex::Vector< amrex::IntVect > & ref_ratios )

Apply the PEC boundary to the electron pressure field.

Parameters
[in,out]PefieldMultifab containing the electron pressure
[in]field_boundary_loBoundary types of the "low" field boundaries
[in]field_boundary_hiBoundary types of the "high" field boundaries
[in]geomgeometry object of level "lev"
[in]levlevel of the Multifab
[in]patch_typecoarse or fine
[in]ref_ratiosvector containing the refinement ratios of the refinement levels

◆ ApplyReflectiveBoundarytoJfield()

void PEC::ApplyReflectiveBoundarytoJfield ( amrex::MultiFab * Jx,
amrex::MultiFab * Jy,
amrex::MultiFab * Jz,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_lo,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_hi,
const amrex::Array< ParticleBoundaryType, 3 > & particle_boundary_lo,
const amrex::Array< ParticleBoundaryType, 3 > & particle_boundary_hi,
const amrex::Geometry & geom,
int lev,
PatchType patch_type,
const amrex::Vector< amrex::IntVect > & ref_ratios )

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

Step 1: Reflect the J field values deposited to the guard cells to their mirror locations inside the domain at PEC and PMC boundaries. Step 2: Set the J field values in the guard cells consistent with the assumed symmetries associated with PEC and PMC boundaries.

Parameters
[in,out]Jx,Jy,JzMultifabs containing the current density
[in]field_boundary_loBoundary types of the "low" field boundaries
[in]field_boundary_hiBoundary types of the "high" field boundaries
[in]particle_boundary_loBoundary types of the "low" particle boundaries
[in]particle_boundary_hiBoundary types of the "high" particle boundaries
[in]geomgeometry object of level "lev"
[in]levlevel of the Multifab
[in]patch_typecoarse or fine
[in]ref_ratiosvector containing the refinement ratios of the refinement levels

PEC: This is an anti-symmetry boundary. Jparallel/Jperp to a boundary deposited to guard cells is subtracted/added from/to its mirror location inside the domain, which is equivalent to depositing J associated with the image charge of the opposite sign on the other side of the PEC boundary. PMC: This is a symmetry boundary. Jparallel/Jperp to a boundary deposited to guard cells is added/subtracted to/from its mirror location inside the domain, which is equivalent to depositing J associated with the image charge of the opposite sign on the other side of the PEC boundary.

◆ ApplyReflectiveBoundarytoRhofield()

void PEC::ApplyReflectiveBoundarytoRhofield ( amrex::MultiFab * rho,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_lo,
const amrex::Array< FieldBoundaryType, 3 > & field_boundary_hi,
const amrex::Array< ParticleBoundaryType, 3 > & particle_boundary_lo,
const amrex::Array< ParticleBoundaryType, 3 > & particle_boundary_hi,
const amrex::Geometry & geom,
int lev,
PatchType patch_type,
const amrex::Vector< amrex::IntVect > & ref_ratios )

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

Step 1: Reflect the Rho field values deposited to the guard cells to their mirror locations inside the domain at PEC and PMC boundaries. Step 2: Set the Rho field values in the guard cells consistent with the assumed symmetries associated with PEC and PMC boundaries.

Parameters
[in,out]rhoMultifab containing the charge density
[in]field_boundary_loBoundary types of the "low" field boundaries
[in]field_boundary_hiBoundary types of the "high" field boundaries
[in]particle_boundary_loBoundary types of the "low" particle boundaries
[in]particle_boundary_hiBoundary types of the "high" particle boundaries
[in]geomgeometry object of level "lev"
[in]levlevel of the Multifab
[in]patch_typecoarse or fine
[in]ref_ratiosvector containing the refinement ratios of the refinement levels

PEC: This is an anti-symmetry boundary. Rho deposited to guard cells is subtracted from its mirror location inside the domain, which is equivalent to depositing Rho associated with the image charge of the opposite sign on the other side of the PEC boundary. PMC: This is a symmetry boundary. Rho deposited to guard cells is Added to its mirror location inside the domain, which is equivalent to depositing Rho associated with the image charge of the same sign on the other side of the PMC boundary.