Tracer particles
Tracer particles are particles that move along a prescribed velocity field
where \(\mathbf{X}\) is the particle position and \(\mathbf{V}\) is the particle velocity. The velocity is interpolated from a mesh-based field as
where \(\mathbf{v}\) is a velocity field defined on the mesh. Such particles are useful, for example, for numerical integration along field lines.
Tip
The chombo-discharge tracer particle functionality resides in $DISCHARGE_HOME/Source/TracerParticles.
TracerParticleSolver
The tracer particle solver is templated as
/*!
@brief Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity field.
@details The user can set the velocity field through public member functions. This class is templated so the user can
switch tracer particle implementations.
This is a single-phase solver -- i.e. the particles only live on one of the phases. Extensions to multiphase
is certainly possible.
The template requirements on the particle type P are:
1) It must contain a function RealVect& position() (derived from BinItem)
2) It must contain a function const Real& weight() const.
3) It must contain a function RealVect& velocity().
*/
template <typename P>
class TracerParticleSolver
where P is the particle type used for the solver.
The template constraints on P are
It must contain a function
RealVect& position()It must contain a function
const Real& weight() constIt must contain a function
RealVect& velocity().
Users are free to provide their own particle type provided that it meets these template constraints. However, we also define a plug-and-play particle class that meets these requirements, see TracerParticle.
Note
The TracerParticleSolver<P> API is available at https://chombo-discharge.github.io/chombo-discharge/doxygen/html/classTracerParticleSolver.html.
TracerParticle
The TracerParticle type inherits from GenericParticle particle class and is templated as
/*!
@brief A tracer particle class. This is templated for holding extra storage (useful for kernels).
@details The template parameters M and N determine extra storage allocated to the particle. M determines
the number of allocated scalars (Reals) and N determines the number of allocated vectors (RealVects). These
quantities are communicated when remapping particles.
*/
template <size_t M, size_t N>
class TracerParticle : public GenericParticle<M, N>
The class also defines two more members; the weight and a particle velocity. These are accessible as
/*!
@brief Get the particle "weight"
@return m_weight
*/
inline Real&
weight();
/*!
@brief Get the particle "weight"
@return m_weight
*/
inline const Real&
weight() const;
/*!
@brief Get the particle velocity.
@return m_velocity
*/
inline RealVect&
velocity();
/*!
@brief Get the particle velocity.
@return m_velocity
*/
inline const RealVect&
velocity() const;
Note that, just as for GenericParticle, the template arguments M and N indicates the number of scalars and vectors allocated to the particle, see GenericParticle.
These data fields can be used by applications for, e.g., storing integration variables (such as intermediate positions in a Runge-Kutta code).
Initialization
To initialize the solver, one can use the full constructor
/*!
@brief Full contructor
@param[in] a_amr Handle to AmrMesh.
@param[in] a_compGeom Computational geometry.
*/
TracerParticleSolver(const RefCountedPtr<AmrMesh>& a_amr, const RefCountedPtr<ComputationalGeometry> a_compGeom);
Getting the particles
To obtain the solver particles, simply call
/*!
@brief Get all particles.
@return m_particles
*/
virtual ParticleContainer<P>&
getParticles();
/*!
@brief Get all particles. Const version.
@return m_particles
*/
virtual const ParticleContainer<P>&
getParticles() const;
This returns the ParticleContainer<P> holding the particles, see ParticleContainer.
Setting \(\mathbf{v}\)
To set the velocity field on the mesh, use
/*!
@brief Set the tracer particle velocity field.
@param[in] a_velocityField Velocity field.
*/
virtual void
setVelocity(const EBAMRCellData& a_velocityField);
This will associate the input velocity a_velocityField with \(\mathbf{v}\).
Interpolating velocities
To compute \(\mathbf{V} = \mathbf{v}\left(\mathbf{X}\right)\) for all particles that reside in the solver, use
/*!
@brief Interpolate particles velocities.
*/
virtual void
interpolateVelocities();
This will interpolate the velocities to the particle positions using the user-defined interpolation method (see Input options).
If desirable, one can also interpolate a scalar field defined on the mesh onto the particle weight by calling
/*!
@brief Interpolate a scalar field onto the particle weight
*/
virtual void
interpolateWeight(const EBAMRCellData& a_scalar) noexcept;
The interpolation function is set by the user, see Input options. See Particle-mesh for further details.
Deposit particles
To deposit the particles, call
/*!
@brief Deposit particle weight on mesh.
@param[out] a_phi Deposited weight.
*/
virtual void
deposit(EBAMRCellData& a_phi) const noexcept;
This will deposit the particle weights onto the input data holder.
The deposition function is set by the user, see Input options. Complete details regarding how the deposition functions work is available in Particle-mesh.
Input options
Available input options for the tracer particle solver are given in the listing below.
TracerParticleSolver<P>.
All options are run-time configurable.# ====================================================================================================
# TracerParticleSolver class options
# ====================================================================================================
TracerParticleSolver.verbosity = -1 ## Solver verbosity level.
TracerParticleSolver.deposition = cic ## Deposition method. Must be 'ngp' or 'cic'
TracerParticleSolver.interpolation = cic ## Interpolation method. Must be 'ngp' or 'cic'
TracerParticleSolver.deposition_cf = transition ## 'interp', 'halo', 'halo_ngp', 'transition'.
TracerParticleSolver.plot_weight = true ## Turn on/off plotting of the particle weight.
TracerParticleSolver.plot_velocity = true ## Turn on/off plotting of the particle velocities.
TracerParticleSolver.volume_scale = false ## If true, depositions yield density * volume instead of just volume