#include <memory>#include <vector>#include <array>#include <random>#include <cmath>#include <iostream>#include <dune/common/fvector.hh>#include <dune/common/exceptions.hh>#include <dumux/parallel/parallel_for.hh>#include <dumux/geometry/intersectingentities.hh>#include <dumux/geometry/geometryintersection.hh>#include <dumux/particles/particle.hh>#include <dumux/particles/simpleparticlecloud.hh>Go to the source code of this file.
A particle solver for the Fokker-Planck equation.
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 | |
| class | Dumux::Particles::FokkerPlanck< Scalar, GridGeometry > |
| A particle solver for the Fokker-Planck equation. More... | |
| struct | Dumux::Particles::FokkerPlanck< Scalar, GridGeometry >::ParticleData |
| data associated with each particle id More... | |
Namespaces | |
| namespace | Dumux |
| namespace | Dumux::Particles |