Manual
Grids
ParticleInCell.AbstractGrid
— TypeParent type for all grid objects, which are used to define the simulation domain, and to convert between coordinate systems. There are three different numbering systems that can refer to a location in the simulation domain:
- The 'physical coordinates' of a point are the real (dimensionalful) coordinates associated with that point. This value can range from the lower bounds to the upper bounds of the simulation. This value will typically take the form
Vector{T}
orNTuple{N, T}
whereT <: Real
. - The 'cell coordinates' of a point is the non-dimensional location of the point in units of cell lengths. This value can range from 0 to num_cells - 1, or outside this range if guard cells are included. The value will typically have the type
NTuple{N, Int}
orNTuple{N, T}
withT <: Real
. - The 'cell index' of a point is the
CartesianIndex
that can be used to index into field arrays at that point. This value must strictly be confined toaxes(field.values)
, which, for any given dimension, will typically range from 1 to numcells + 2*numguard_cells + 1.
The first two types of indexing, physcoords and cellcoords, are independent of the number of guard cells in a given field, and depend only on grid quantities. Thus utilities for converting between these systems require only a reference to a grid object. On the other hand, the utilities for cell_index are specific the field being used, and thus those must be provided an AbstractField
to do the coordinate conversion.
In general, physical coordinates are useful when considering the location of a particle, while the cell index is used to interpolate to and from the particle locations. The cell coordinates are useful for some interpolation algorithms, especially those that are defined for non-uniform grids.
ParticleInCell.UniformCartesianGrid
— TypeUniformCartesianGrid(lower_bounds, upper_bounds, num_cells, periodic)
The simplest grid type, which represents a set of equally spaced rectangular cells. The grid can optionally be periodic in one or more dimensions.
Grid utility functions
ParticleInCell.cell_lengths
— Functioncell_lengths(grid, [cell_coords])
Returns the length of the cell located at cell_coords
. For uniform grids, the cell_coords
argument is optional.
ParticleInCell.sim_lengths
— Functionsim_lengths(grid)
Returns a tuple of the length of the simulation domain in each dimension.
ParticleInCell.cell_coords_to_phys_coords
— Functioncell_coords_to_phys_coords(grid, idxs)
Converts the cell coordinates idxs
to a physical coordinate using the geometry specified in grid
.
For more information on the different types of coordinate systems, see AbstractGrid
.
ParticleInCell.phys_coords_to_cell_coords
— Functionphys_coords_to_cell_coords(grid, xs)
Converts the physical coordinate xs
to a grid coordinate using the geometry specified in grid
.
For more information on the different types of coordinate systems, see AbstractGrid
.
Species
ParticleInCell.AbstractSpecies
— TypeSubtypes of AbstractSpecies
represent a group of macroparticles with a single physical charge and mass. Each particle of a species can optionally represent a different number of physical particles.
ParticleInCell.VariableWeightSpecies
— TypeVariableWeightSpecies(positions::Matrix, momentums::Matrix, weight::Vector charge, mass)
Stores the positions and momentums of particles that share a common charge and mass. Each particle can have a different weight, and the momentum space can have a larger dimension than the position space.
ParticleInCell.electrons
— Functionelectrons(positions::Matrix, momentums::Matrix, weights::Vector)
electrons(positions::Matrix, momentums::Matrix, weight)
Creates an AbsractSpecies
with the given positions
, momentums
, and weight
, and the charge and mass of a single electron.
Species utility functions
ParticleInCell.particle_charge
— Functionparticle_charge(species, idx)
Returns the charge of the macroparticle of species
with index idx
. See also physical_charge
.
ParticleInCell.physical_charge
— Functionphysical_charge(species, idx)
Returns the charge of a physical particle represented by the macroparticle of species
with index idx
. See also particle_charge
.
ParticleInCell.particle_mass
— Functionparticle_mass(species, idx)
Returns the mass of the macroparticle of species
with index idx
. See also physical_mass
.
ParticleInCell.physical_mass
— Functionphysical_mass(species, idx)
Returns the mass of a physical particle represented by the macroparticle of species
with index idx
. See also particle_charge
.
ParticleInCell.particle_weight
— Functionparticle_weight(species, idx)
Returns the number of physical particles represented by the macroparticle of species
with index idx
.
ParticleInCell.particle_position
— Functionparticle_position(species, idx)
Returns the position of the macroparticle of species
with index idx
. See also particle_momentum
and particle_position!
.
ParticleInCell.particle_position!
— Functionparticle_position!(species, idx, value)
Sets the position of the macroparticle of species
with index idx
to value
. See also particle_momentum!
and particle_position
.
ParticleInCell.particle_momentum
— Functionparticle_momentum(species, idx)
Returns the momentum of the macroparticle of species
with index idx
. See also particle_momentum!
, particle_velocity
, and physical_momentum
.
ParticleInCell.particle_momentum!
— Functionparticle_momentum!(species, idx, value)
Sets the momentum of the macroparticle of species
with index idx
to value
. See also particle_momentum
, particle_velocity
, and physical_momentum
.
ParticleInCell.particle_velocity
— Functionparticle_velocity(species, idx)
Returns the velocity of the macroparticle of species
with index idx
. See also particle_momentum
, particle_momentum!
, and physical_momentum
.
ParticleInCell.physical_momentum
— Functionphysical_momentum(species, idx)
Returns the momentum of a physical particle represented by the macroparticle of species
with index idx
. See also particle_momentum
, particle_momentum!
, and particle_velocity
.
Base.eachindex
— Functioneachindex(species)
Returns an iterator to the index of each of the macroparticles contained in species
.
Fields
ParticleInCell.Field
— TypeField(grid, offset, num_elements,[ lower_guard_cells = 0,
[upper_guard_cells = lower_guard_cells + 1]])
Stores the values of a field.
Update steps
ParticleInCell.AbstractUpdateStep
step!
ParticleInCell.ElectrostaticParticlePush
— TypeElectrostaticParticlePush(species, E, timestep, [interpolation_order=1])
An update step that moves and accelerates particles of species
according to the electric field E
. The field will be interpolated to the particle positions using b-splines of interpolation_order
.
ParticleInCell.BorisParticlePush
— TypeBorisParticlePush(species, E, B, timestep)
This update step moves the particles of species
subject to both an electric field, E
, and a magnetic field, B
. The method is frequently referred to by the shorthand accelerate-rotate-accelerate because the acceleration from the electric field is split in half, and applied before and after the magnetic field rotation.
An applied magnetic field forces charged particles to travel in circular orbits perpendicular to the magnetic field. Thus, a simulation using the Boris method only makes sense when there are at least two velocity components, and in this case, the applied magnetic field must be strictly perpendicular to the velocity components. So for example, one could have a 1d2v simulation with a spatial grid along the x-axis, and velocity components in the x and y directions. In this case, the magnetic field would be forced to point along the z axis. If all three velocity components are included in the simulation, then the magnetic field can point in any direction.
For more details on the method, see Sections 4.3 and 4.4 of Birdsall and Langdon, or the Boris' original conference proceedings.
ParticleInCell.LinearSolveStep
— TypeSolves the linear equation A x = b
where x
and b
are fields.
Misc
ParticleInCell.compute_knots
— Functioncompute_knots(degree)
Returns a vector of the knot locations for a polynomial b-spline with degree
.
ParticleInCell.compute_bspline_coeffs
— Functioncompute_bspline_coeffs(degree, [T])
Returns a vector of vectors of the coefficients for the b-spline with polynomial degree
, and unit-spaced knots, centered at zero.