List of JustPIC functions
Here an overview of all functions in JustPIC.jl, , for a complete list see here:
JustPIC.Particles Type
Particles{Backend,N,I,T1,T2} <: AbstractParticlesA struct representing a collection of particles.
Arguments
backend: The backend used for particle computations (CPUBackend, CUDABackend, AMDGPUBackend).coords: Coordinates of the particlesindex: Helper array flaggin active particlesnxcell: Initial number of particles per cellmax_xcell: Maximum number of particles per cellmin_xcell: Minimum number of particles per cellnp: Number of particles
JustPIC.nphases Method
nphases(x::PhaseRatios)Return the number of phases in x::PhaseRatios.
JustPIC._2D.advection! Method
advection!(particles::Particles, method::AbstractAdvectionIntegrator, V, grid_vi::NTuple{N,NTuple{N,T}}, dt)Advects the particles using the advection scheme defined by method.
Arguments
particles: Particles object to be advected.method: Time integration method (EulerorRungeKutta2).V: Tuple containingVx,Vy; andVzin 3D.grid_vi: Tuple containing the grids corresponding toVx,Vy; andVzin 3D.dt: Time step.
JustPIC._2D.advection_LinP! Method
advection!(particles::Particles, method::AbstractAdvectionIntegrator, V, grid_vi::NTuple{N,NTuple{N,T}}, dt)Advects the particles using the advection scheme defined by method.
Arguments
particles: Particles object to be advected.method: Time integration method (EulerorRungeKutta2).V: Tuple containingVx,Vy; andVzin 3D.grid_vi: Tuple containing the grids corresponding toVx,Vy; andVzin 3D.dt: Time step.
JustPIC._2D.advection_MQS! Method
advection!(particles::Particles, method::AbstractAdvectionIntegrator, V, grid_vi::NTuple{N,NTuple{N,T}}, dt)Advects the particles using the advection scheme defined by method.
Arguments
particles: Particles object to be advected.method: Time integration method (EulerorRungeKutta2).V: Tuple containingVx,Vy; andVzin 3D.grid_vi: Tuple containing the grids corresponding toVx,Vy; andVzin 3D.dt: Time step.
JustPIC._2D.centroid2particle! Method
centroid2particle!(Fp, xci, F, particles)Interpolates properties F that are defined on a mesh at center points with location xci to particles Fp.
JustPIC._2D.checkpointing_particles Method
checkpointing_particles(dst, particles;phases=nothing, phase_ratios=nothing, chain=nothing, t=nothing, dt=nothing, particle_args=nothing)Save the state of particles and related data to a checkpoint file in a jld2 format. The name of the checkpoint file is particles_checkpoint.jld2.
Arguments
dst: The destination directory where the checkpoint file will be saved.particles: The array of particles to be saved.
Keyword Arguments
phases: The array of phases associated with the particles. If nothing is stated, the default isnothing.phase_ratios: The array of phase ratios. If nothing is stated, the default isnothing.chain: The chain data to be saved. If nothing is stated, the default isnothing.t: The current time to be saved. If nothing is stated, the default isnothing.dt: The timestep to be saved. If nothing is stated, the default isnothing.particle_args: Additional particle arguments to be saved. If nothing is stated, the default isnothing.kwargs: Additional keyword arguments to be saved in the checkpoint file.
JustPIC._2D.fill_chain_from_chain! Method
fill_chain!(chain::MarkerChain, topo_x, topo_y)Fill the given chain of markers with topographical data.
Arguments
chain::MarkerChain: The chain of markers to be filled.topo_x: The x-coordinates of the topography.topo_y: The y-coordinates of the topography.
Description
This function populates the chain with markers based on the provided topographical data (topo_x and topo_y). The function modifies the chain in place.
JustPIC._2D.force_injection! Method
force_injection!(particles::Particles{Backend}, p_new) where {Backend}Forces the injection of new particles into the existing particles collection. This function modifies the particles in place by adding the new particles specified in p_new.
Arguments
particles::Particles{Backend}: The existing collection of particles to which new particles will be added. The type of backend is specified by theBackendparameter.p_new: The new particles to be injected into the existing collection.
JustPIC._2D.force_injection! Method
force_injection!(particles::Particles{Backend}, p_new, fields::NTuple{N, Any}, values::NTuple{N, Any}) where {Backend, N}Forcefully injects new particles into the particles object. This function modifies the particles object in place.
Arguments
particles::Particles{Backend}: The particles object to be modified.p_new: The new particles to be injected.fields::NTuple{N, Any}: A tuple containing the fields to be updated.values::NTuple{N, Any}: A tuple containing the values corresponding to the fields.
Returns
- Nothing. This function modifies the
particlesobject in place.
JustPIC._2D.init_particles Method
init_particles( backend, nxcell, max_xcell, min_xcell, coords::NTuple{N,AbstractArray}, dxᵢ::NTuple{N,T}, nᵢ::NTuple{N,I})Initialize the particles object.
Arguments
backend: Device backendnxcell: Initial number of particles per cellmax_xcell: Maximum number of particles per cellmin_xcell: Minimum number of particles per cellxvi: Grid cells vertices
JustPIC._2D.inject_particles! Method
inject_particles!(particles::Particles, args, grid)Injects particles if the number of particles in a given cell is such that n < particles.min_xcell.
Arguments
particles: The particles object.args:CellArrayss containing particle fields.grid: The grid cell vertices.
JustPIC._2D.interpolate_velocity_to_markerchain! Method
interpolate_velocity_to_markerchain!(chain::MarkerChain, chain_V::NTuple{N, CellArray}, V, grid_vi::NTuple{N, NTuple{N, T}}) where {N, T}Interpolates the velocity field to the positions of the marker chain.
Arguments
chain::MarkerChain: The marker chain object containing the particle coordinates and indices.chain_V::NTuple{N, CellArray}: The output velocity field at the marker chain positions.V: The velocity field to be interpolated.grid_vi::NTuple{N, NTuple{N, T}}: The grid information for each dimension.
JustPIC._2D.lerp Method
lerp(v, t::NTuple{nD,T}) where {nD,T}Linearly interpolates the value v between the elements of the tuple t. This function is specialized for tuples of length nD.
Arguments
v: The value to be interpolated.t: The tuple of values to interpolate between.
JustPIC._2D.move_particles! Method
move_particles!(particles::AbstractParticles, grid, args)Move particles in the given particles container according to the provided grid and particles fields in args.
Arguments
particles: The container of particles to be moved.grid: The grid used for particle movement.args:CellArrayss containing particle fields.
JustPIC._2D.particle2centroid! Method
particle2centroid!(F, Fp, xci::NTuple, particles::Particles)Interpolates properties Fp from particles to the grid F at center points that are defined by 1D coordinate arrays in xci
JustPIC._2D.semilagrangian_advection! Method
semilagrangian_advection!F, F0, integrator, V, grid_vi, grid, dt)Performs semi-Lagrangian advection by backtracking particle positions in a velocity field. This function updates the positions and/or properties of particles according to the semi-Lagrangian scheme.
Arguments
F: The new state of the grid field (e.g., density, temperature).F0: The current state of the grid field (used for interpolation).integrator: The numerical integrator to use for advection (e.g., Euler, Rk2, RK4).V: The velocity field at the particle positions.grid_vi: The grid cell indices for the velocity field.grid: The spatial grid information.dt: The time step for the advection.
JustPIC._2D.semilagrangian_advection_LinP! Method
semilagrangian_advection_LinP!, F0, integrator, V, grid_vi, grid, dt)Performs semi-Lagrangian advection by backtracking particle positions in a velocity field. This function updates the positions and/or properties of particles according to the semi-Lagrangian scheme.
Arguments
F: The new state of the grid field (e.g., density, temperature).F0: The current state of the grid field (used for interpolation).integrator: The numerical integrator to use for advection (e.g., Euler, Rk2, RK4).V: The velocity field at the particle positions.grid_vi: The grid cell indices for the velocity field.grid: The spatial grid information.dt: The time step for the advection.
JustPIC._2D.semilagrangian_advection_MQS! Method
semilagrangian_advection_MQS!, F0, integrator, V, grid_vi, grid, dt)Performs semi-Lagrangian advection by backtracking particle positions in a velocity field. This function updates the positions and/or properties of particles according to the semi-Lagrangian scheme.
Arguments
F: The new state of the grid field (e.g., density, temperature).F0: The current state of the grid field (used for interpolation).integrator: The numerical integrator to use for advection (e.g., Euler, Rk2, RK4).V: The velocity field at the particle positions.grid_vi: The grid cell indices for the velocity field.grid: The spatial grid information.dt: The time step for the advection.
JustPIC._2D.@idx Macro
@idx(args...)Make a linear range from 1 to args[i], with i ∈ [1, ..., n]
JustPIC._3D.advection! Method
advection!(particles::Particles, method::AbstractAdvectionIntegrator, V, grid_vi::NTuple{N,NTuple{N,T}}, dt)Advects the particles using the advection scheme defined by method.
Arguments
particles: Particles object to be advected.method: Time integration method (EulerorRungeKutta2).V: Tuple containingVx,Vy; andVzin 3D.grid_vi: Tuple containing the grids corresponding toVx,Vy; andVzin 3D.dt: Time step.
JustPIC._3D.advection_LinP! Method
advection!(particles::Particles, method::AbstractAdvectionIntegrator, V, grid_vi::NTuple{N,NTuple{N,T}}, dt)Advects the particles using the advection scheme defined by method.
Arguments
particles: Particles object to be advected.method: Time integration method (EulerorRungeKutta2).V: Tuple containingVx,Vy; andVzin 3D.grid_vi: Tuple containing the grids corresponding toVx,Vy; andVzin 3D.dt: Time step.
JustPIC._3D.advection_MQS! Method
advection!(particles::Particles, method::AbstractAdvectionIntegrator, V, grid_vi::NTuple{N,NTuple{N,T}}, dt)Advects the particles using the advection scheme defined by method.
Arguments
particles: Particles object to be advected.method: Time integration method (EulerorRungeKutta2).V: Tuple containingVx,Vy; andVzin 3D.grid_vi: Tuple containing the grids corresponding toVx,Vy; andVzin 3D.dt: Time step.
JustPIC._3D.centroid2particle! Method
centroid2particle!(Fp, xci, F, particles)Interpolates properties F that are defined on a mesh at center points with location xci to particles Fp.
JustPIC._3D.checkpointing_particles Method
checkpointing_particles(dst, particles;phases=nothing, phase_ratios=nothing, chain=nothing, t=nothing, dt=nothing, particle_args=nothing)Save the state of particles and related data to a checkpoint file in a jld2 format. The name of the checkpoint file is particles_checkpoint.jld2.
Arguments
dst: The destination directory where the checkpoint file will be saved.particles: The array of particles to be saved.
Keyword Arguments
phases: The array of phases associated with the particles. If nothing is stated, the default isnothing.phase_ratios: The array of phase ratios. If nothing is stated, the default isnothing.chain: The chain data to be saved. If nothing is stated, the default isnothing.t: The current time to be saved. If nothing is stated, the default isnothing.dt: The timestep to be saved. If nothing is stated, the default isnothing.particle_args: Additional particle arguments to be saved. If nothing is stated, the default isnothing.kwargs: Additional keyword arguments to be saved in the checkpoint file.
JustPIC._3D.fill_chain_from_chain! Method
fill_chain!(chain::MarkerChain, topo_x, topo_y)Fill the given chain of markers with topographical data.
Arguments
chain::MarkerChain: The chain of markers to be filled.topo_x: The x-coordinates of the topography.topo_y: The y-coordinates of the topography.
Description
This function populates the chain with markers based on the provided topographical data (topo_x and topo_y). The function modifies the chain in place.
JustPIC._3D.force_injection! Method
force_injection!(particles::Particles{Backend}, p_new) where {Backend}Forces the injection of new particles into the existing particles collection. This function modifies the particles in place by adding the new particles specified in p_new.
Arguments
particles::Particles{Backend}: The existing collection of particles to which new particles will be added. The type of backend is specified by theBackendparameter.p_new: The new particles to be injected into the existing collection.
JustPIC._3D.force_injection! Method
force_injection!(particles::Particles{Backend}, p_new, fields::NTuple{N, Any}, values::NTuple{N, Any}) where {Backend, N}Forcefully injects new particles into the particles object. This function modifies the particles object in place.
Arguments
particles::Particles{Backend}: The particles object to be modified.p_new: The new particles to be injected.fields::NTuple{N, Any}: A tuple containing the fields to be updated.values::NTuple{N, Any}: A tuple containing the values corresponding to the fields.
Returns
- Nothing. This function modifies the
particlesobject in place.
JustPIC._3D.init_particles Method
init_particles( backend, nxcell, max_xcell, min_xcell, coords::NTuple{N,AbstractArray}, dxᵢ::NTuple{N,T}, nᵢ::NTuple{N,I})Initialize the particles object.
Arguments
backend: Device backendnxcell: Initial number of particles per cellmax_xcell: Maximum number of particles per cellmin_xcell: Minimum number of particles per cellxvi: Grid cells vertices
JustPIC._3D.inject_particles! Method
inject_particles!(particles::Particles, args, grid)Injects particles if the number of particles in a given cell is such that n < particles.min_xcell.
Arguments
particles: The particles object.args:CellArrayss containing particle fields.grid: The grid cell vertices.
JustPIC._3D.interpolate_velocity_to_markerchain! Method
interpolate_velocity_to_markerchain!(chain::MarkerChain, chain_V::NTuple{N, CellArray}, V, grid_vi::NTuple{N, NTuple{N, T}}) where {N, T}Interpolates the velocity field to the positions of the marker chain.
Arguments
chain::MarkerChain: The marker chain object containing the particle coordinates and indices.chain_V::NTuple{N, CellArray}: The output velocity field at the marker chain positions.V: The velocity field to be interpolated.grid_vi::NTuple{N, NTuple{N, T}}: The grid information for each dimension.
JustPIC._3D.lerp Method
lerp(v, t::NTuple{nD,T}) where {nD,T}Linearly interpolates the value v between the elements of the tuple t. This function is specialized for tuples of length nD.
Arguments
v: The value to be interpolated.t: The tuple of values to interpolate between.
JustPIC._3D.move_particles! Method
move_particles!(particles::AbstractParticles, grid, args)Move particles in the given particles container according to the provided grid and particles fields in args.
Arguments
particles: The container of particles to be moved.grid: The grid used for particle movement.args:CellArrayss containing particle fields.
JustPIC._3D.particle2centroid! Method
particle2centroid!(F, Fp, xci::NTuple, particles::Particles)Interpolates properties Fp from particles to the grid F at center points that are defined by 1D coordinate arrays in xci
JustPIC._3D.semilagrangian_advection! Method
semilagrangian_advection!F, F0, integrator, V, grid_vi, grid, dt)Performs semi-Lagrangian advection by backtracking particle positions in a velocity field. This function updates the positions and/or properties of particles according to the semi-Lagrangian scheme.
Arguments
F: The new state of the grid field (e.g., density, temperature).F0: The current state of the grid field (used for interpolation).integrator: The numerical integrator to use for advection (e.g., Euler, Rk2, RK4).V: The velocity field at the particle positions.grid_vi: The grid cell indices for the velocity field.grid: The spatial grid information.dt: The time step for the advection.
JustPIC._3D.semilagrangian_advection_LinP! Method
semilagrangian_advection_LinP!, F0, integrator, V, grid_vi, grid, dt)Performs semi-Lagrangian advection by backtracking particle positions in a velocity field. This function updates the positions and/or properties of particles according to the semi-Lagrangian scheme.
Arguments
F: The new state of the grid field (e.g., density, temperature).F0: The current state of the grid field (used for interpolation).integrator: The numerical integrator to use for advection (e.g., Euler, Rk2, RK4).V: The velocity field at the particle positions.grid_vi: The grid cell indices for the velocity field.grid: The spatial grid information.dt: The time step for the advection.
JustPIC._3D.semilagrangian_advection_MQS! Method
semilagrangian_advection_MQS!, F0, integrator, V, grid_vi, grid, dt)Performs semi-Lagrangian advection by backtracking particle positions in a velocity field. This function updates the positions and/or properties of particles according to the semi-Lagrangian scheme.
Arguments
F: The new state of the grid field (e.g., density, temperature).F0: The current state of the grid field (used for interpolation).integrator: The numerical integrator to use for advection (e.g., Euler, Rk2, RK4).V: The velocity field at the particle positions.grid_vi: The grid cell indices for the velocity field.grid: The spatial grid information.dt: The time step for the advection.
JustPIC._3D.@idx Macro
@idx(args...)Make a linear range from 1 to args[i], with i ∈ [1, ..., n]