Afivo  0.3
Data Types | Functions/Subroutines
m_af_particles Module Reference

This module contains routines related to , which can interpolate 'to' the grid and 'from' the grid (useful for e.g. particle simulations). The interpolation for meshes is called prolongation, see m_aX_prolong. More...

Data Types

interface  subr_particle_id
 To get a particle id. More...
 
interface  subr_particle_rw
 To get a particle's coordinates and weight. More...
 

Functions/Subroutines

subroutine, public af_particles_to_grid (tree, iv, n_particles, get_id, get_rw, order, density, fill_gc, iv_tmp, offset_particles)
 Map a list of particles to a density. The order can be zero (map particle to the containing cell) or one (use bi/tri-linear interpolation). Note that ghost cells are automatically filled by this routine. More...
 
subroutine particles_to_grid_0 (tree, iv, get_rw, ids, threads, n_particles, density, p_offset)
 
subroutine particles_to_grid_1 (tree, iv, get_rw, ids, threads, n_particles, density, p_offset)
 Add weights to the cell centers using linear interpolation. More...
 
subroutine tree_add_from_ghostcells (tree, iv)
 
subroutine add_from_ghostcells (boxes, id, iv)
 
subroutine add_as_density (tree, iv_from, iv_to)
 Convert particle weights to densities and add to another variable. More...
 
subroutine add_as_density_box (box, iv_from, iv_to)
 Convert particle weights to densities and add to another variable. More...
 

Detailed Description

This module contains routines related to , which can interpolate 'to' the grid and 'from' the grid (useful for e.g. particle simulations). The interpolation for meshes is called prolongation, see m_aX_prolong.

Note that the particle coordinates are transferred via subroutines. This has two advantages: first, data does not need to be copied, saving memory. Second, the particle id can be cached in the particle's data structured.

Function/Subroutine Documentation

◆ af_particles_to_grid()

subroutine, public m_af_particles::af_particles_to_grid ( type(af_t), intent(inout)  tree,
integer, intent(in)  iv,
integer, intent(in)  n_particles,
procedure(subr_particle_id get_id,
procedure(subr_particle_rw get_rw,
integer, intent(in)  order,
logical, intent(in), optional  density,
logical, intent(in), optional  fill_gc,
integer, intent(in), optional  iv_tmp,
integer, intent(in), optional  offset_particles 
)

Map a list of particles to a density. The order can be zero (map particle to the containing cell) or one (use bi/tri-linear interpolation). Note that ghost cells are automatically filled by this routine.

Parameters
[in]ivVariable to store density
[in]n_particlesThe number of particles
get_idTo get the particle id
get_rwTo get the particle position and weight
[in]orderOrder of interpolation
[in]densityDivide by cell area/volume (default: true)
[in]fill_gcFill ghost cells afterwards (default: true)
[in]iv_tmpUse temporary variable to convert to density. This can be faster, and is slightly more accurate for cylindrical coordinate systems, due to the way ghost cells are exchanged near refinement boundaries
[in]offset_particlesStart offset for indexing the particles (if zero, start at index 1)

Definition at line 39 of file m_af_particles.f90.

◆ particles_to_grid_0()

subroutine m_af_particles::particles_to_grid_0 ( type(af_t), intent(inout)  tree,
integer, intent(in)  iv,
procedure(subr_particle_rw get_rw,
integer, dimension(n_particles), intent(in)  ids,
integer, dimension(n_particles), intent(in)  threads,
integer, intent(in)  n_particles,
logical, intent(in)  density,
integer, intent(in)  p_offset 
)
Parameters
[in]ivVariable to store particle density
[in]p_offsetOffset for particle indexing

Definition at line 184 of file m_af_particles.f90.

◆ particles_to_grid_1()

subroutine m_af_particles::particles_to_grid_1 ( type(af_t), intent(inout)  tree,
integer, intent(in)  iv,
procedure(subr_particle_rw get_rw,
integer, dimension(n_particles), intent(in)  ids,
integer, dimension(n_particles), intent(in)  threads,
integer, intent(in)  n_particles,
logical, intent(in)  density,
integer, intent(in)  p_offset 
)

Add weights to the cell centers using linear interpolation.

Todo:
Support cylindrical coordinates
Parameters
[in]ivVariable to store particle density
[in]densityAdd particle as a density
[in]p_offsetOffset for particle indexing

Definition at line 239 of file m_af_particles.f90.

◆ tree_add_from_ghostcells()

subroutine m_af_particles::tree_add_from_ghostcells ( type(af_t), intent(inout)  tree,
integer, intent(in)  iv 
)
Parameters
[in]ivIndex of variable

Definition at line 322 of file m_af_particles.f90.

◆ add_from_ghostcells()

subroutine m_af_particles::add_from_ghostcells ( type(box_t), dimension(:), intent(inout)  boxes,
integer, intent(in)  id,
integer, intent(in)  iv 
)
private
Parameters
[in]idIndex of box
[in]ivIndex of variable

Definition at line 339 of file m_af_particles.f90.

◆ add_as_density()

subroutine m_af_particles::add_as_density ( type(af_t), intent(inout)  tree,
integer, intent(in)  iv_from,
integer, intent(in)  iv_to 
)
private

Convert particle weights to densities and add to another variable.

Definition at line 405 of file m_af_particles.f90.

◆ add_as_density_box()

subroutine m_af_particles::add_as_density_box ( type(box_t), intent(inout)  box,
integer, intent(in)  iv_from,
integer, intent(in)  iv_to 
)
private

Convert particle weights to densities and add to another variable.

Definition at line 424 of file m_af_particles.f90.