version 3.11-dev
Dumux::Particles::FokkerPlanck< Scalar, GridGeometry > Class Template Reference

A particle solver for the Fokker-Planck equation. More...

#include <dumux/particles/fokkerplanck.hh>

Description

template<class Scalar, class GridGeometry>
class Dumux::Particles::FokkerPlanck< Scalar, GridGeometry >

A cloud of particles, each particle associated with a given mass \(m_i\), is evolved such that the probability density function \(\rho(\mathbf{x}, t)\) of the particle position \( \mathbf{X}_t \) approximates the solution of the Fokker-Planck equation:

\[ \frac{\partial \rho}{\partial t} + \nabla \cdot (\boldsymbol{\mu} \rho + D \nabla \rho) = 0 \]

Particles are evolved with the Dumux::Particles::FokkerPlanck::run method, which updates each particles position \(\mathbf{X}_t\) with the Euler scheme:

\[ \Delta\mathbf{X}_t = \boldsymbol{\mu}(\mathbf{X}_t, t) \Delta t + \sqrt{2 D} \Delta \mathbf{W}_t \]

The particle density on the grid can be computed with the Dumux::Particles::FokkerPlanck::computeDensity method as

\[ \rho_E = \frac{1}{|E|} \sum_{i=1}^{N_E} m_i, \]

where \(|E|\) is the volume of the grid element \(E\) and \(m_i\) is the mass of the particle \(i\) in element \(E\).

The Fokker-Planck equation for particle positions is equivalent to the advection-diffusion equation for the particle density.

Classes

struct  ParticleData
 data associated with each particle id More...
 

Public Member Functions

 FokkerPlanck (std::shared_ptr< const GridGeometry > gg, std::shared_ptr< Cloud > cloud)
 
template<class Func >
void init (const Func &func)
 For each particle apply a function to initialize it's state. More...
 
void run (const Scalar dt)
 Evolve particle state with Euler scheme. More...
 
void computeDensity (std::vector< Scalar > &density) const
 computes the particle density in the grid More...
 
void setVelocityData (const std::vector< Velocity > &velocities)
 set the velocity data for each grid element More...
 
std::shared_ptr< const CloudparticleCloud () const
 access the underlying particle cloud More...
 

Constructor & Destructor Documentation

◆ FokkerPlanck()

template<class Scalar , class GridGeometry >
Dumux::Particles::FokkerPlanck< Scalar, GridGeometry >::FokkerPlanck ( std::shared_ptr< const GridGeometry >  gg,
std::shared_ptr< Cloud cloud 
)
inline

Member Function Documentation

◆ computeDensity()

template<class Scalar , class GridGeometry >
void Dumux::Particles::FokkerPlanck< Scalar, GridGeometry >::computeDensity ( std::vector< Scalar > &  density) const
inline
Parameters
densitythe density data to fill for each grid element

Computes a discrete piece-wise constant approximation of the density with cell values given as

\[ \rho_E = \frac{1}{|E|} \sum_{i=1}^{N_E} m_i, \]

where \(|E|\) is the volume of the grid element \(E\) and \(m_i\) is the mass of the particle \(i\) in element \(E\).

◆ init()

template<class Scalar , class GridGeometry >
template<class Func >
void Dumux::Particles::FokkerPlanck< Scalar, GridGeometry >::init ( const Func &  func)
inline
Note
The function may modify the particle and its associated data

◆ particleCloud()

template<class Scalar , class GridGeometry >
std::shared_ptr< const Cloud > Dumux::Particles::FokkerPlanck< Scalar, GridGeometry >::particleCloud ( ) const
inline

◆ run()

template<class Scalar , class GridGeometry >
void Dumux::Particles::FokkerPlanck< Scalar, GridGeometry >::run ( const Scalar  dt)
inline
Parameters
dtthe time step size

Updates the particle position \(\mathbf{X}_t\) by

\[ \Delta\mathbf{X}_t = \boldsymbol{\mu}(\mathbf{X}_t, t) \Delta t + \sqrt{2 D} \Delta \mathbf{W}_t \]

where \(\boldsymbol{\mu}(\mathbf{X}_t, t)\) is the velocity field, \(D\) is the diffusion coefficient, and \(\Delta \mathbf{W}_t\) is a Wiener increment. Note that \(\Delta \mathbf{W}_t\) is normal-distributed with mean zero and variance \(\Delta t\): \(\Delta \mathbf{W}_t \sim \mathcal{N}(0, \Delta t)\).

◆ setVelocityData()

template<class Scalar , class GridGeometry >
void Dumux::Particles::FokkerPlanck< Scalar, GridGeometry >::setVelocityData ( const std::vector< Velocity > &  velocities)
inline
Parameters
velocitiesthe velocity data
Note
to set a constant velocity field, use a vector of size 1

The documentation for this class was generated from the following file: