fokker_planck
moment_kinetics.fokker_planck — Module
Module for including the Full-F Fokker-Planck Collision Operator.
The functions in this module are split into two groups.
The first set of functions implement the weak-form Collision operator using the Rosenbluth-MacDonald-Judd formulation in a divergence form. The Green's functions for the Rosenbluth potentials are used to obtain the Rosenbluth potentials at the boundaries. To find the potentials everywhere else elliptic solves of the PDEs for the Rosenbluth potentials are performed with Dirichlet boundary conditions. These routines provide the default collision operator used in the code.
The second set of functions are used to set up the necessary arrays to compute the Rosenbluth potentials everywhere in vpa, vperp by direct integration of the Green's functions. These functions are supported for the purposes of testing and debugging.
Lower-level routines are provided by functions from moment_kinetics.fokker_planck_calculus.
Parallelisation of the collision operator uses a special 'anysv' region type, see Collision operator and anysv region.
moment_kinetics.fokker_planck.allocate_fokkerplanck_arrays_direct_integration — Method
Function that allocates the required ancilliary arrays for direct integration routines.
sourcemoment_kinetics.fokker_planck.calculate_entropy_production! — Method
Function to calculate entropy production, in place.
sourcemoment_kinetics.fokker_planck.conserving_corrections! — Method
Function that applies numerical-error correcting terms to ensure numerical conservation of the moments density, upar, pressure in the self-collision operator. Modifies the collision operator such that the operator becomes
\[C_{ss} = C^\ast_{ss}[F_s,F_{s}] - \left(x_0 + x_1(v_{\|}-u_{\|})+ x_2(v_\perp^2 +(v_{\|}-u_{\|})^2)\right)F_s\]
where $C^\ast_{ss}[F_s,F_{s}]$ is the weak-form self-collision operator computed using the finite-element implementation, $u_{\|}$ is the parallel velocity of $F_s$, and $x_0,x_1,x_2$ are parameters that are chosen so that $C_{ss}$ conserves density, parallel velocity and pressure of $F_s$.
sourcemoment_kinetics.fokker_planck.density_conserving_correction! — Method
Function that applies a numerical-error correcting term to ensure numerical conservation of the density in the collision operator.
\[C_{ss^\prime} = C^\ast_{ss}[F_s,F_{s^\prime}] - x_0 F_s\]
where $C^\ast_{ss}[F_s,F_{s^\prime}]$ is the weak-form collision operator computed using the finite-element implementation.
sourcemoment_kinetics.fokker_planck.explicit_fokker_planck_collisions_weak_form! — Method
Function for advancing with the explicit, weak-form, self-collision operator.
sourcemoment_kinetics.fokker_planck.explicit_fp_collisions_weak_form_Maxwellian_cross_species! — Method
Function for advancing with the explicit, weak-form, self-collision operator using the existing method for computing the Rosenbluth potentials, with the addition of cross-species collisions against fixed Maxwellian distribution functions where the Rosenbluth potentials are specified using analytical results.
sourcemoment_kinetics.fokker_planck.fokker_planck_collision_operator_weak_form! — Method
Function for evaluating $C_{ss'} = C_{ss'}[F_s,F_{s'}]$
The result is stored in the array fkpl_arrays.CC.
The normalised collision frequency for collisions between species s and s' is defined by
\[\tilde{\nu}_{ss'} = \frac{L_{\mathrm{ref}}}{c_{\mathrm{ref}}}\frac{\gamma_{ss'} n_\mathrm{ref}}{m_s^2 c_\mathrm{ref}^3}\]
with $\gamma_{ss'} = 2 \pi (Z_s Z_{s'})^2 e^4 \ln \Lambda_{ss'} / (4 \pi \epsilon_0)^2$. The input parameter to this code is
\[\tilde{\nu}_{ii} = \frac{L_{\mathrm{ref}}}{c_{\mathrm{ref}}}\frac{\gamma_\mathrm{ref} n_\mathrm{ref}}{m_\mathrm{ref}^2 c_\mathrm{ref}^3}\]
with $\gamma_\mathrm{ref} = 2 \pi e^4 \ln \Lambda_{ii} / (4 \pi \epsilon_0)^2$. This means that $\tilde{\nu}_{ss'} = (Z_s Z_{s'})^2\tilde{\nu}_\mathrm{ref}$ and this conversion is handled explicitly in the code with the charge number input provided by the user.
sourcemoment_kinetics.fokker_planck.fokker_planck_collision_operator_weak_form_Maxwellian_Fsp! — Method
Function for computing the collision operator
\[\sum_{s^\prime} C[F_{s},F_{s^\prime}]\]
when $F_{s^\prime}$ is an analytically specified Maxwellian distribution and the corresponding Rosenbluth potentials are specified using analytical results.
sourcemoment_kinetics.fokker_planck.init_fokker_planck_collisions_direct_integration — Method
Function that initialises the arrays needed to calculate the Rosenbluth potentials by direct integration. As this function is only supported to keep the testing of the direct integration method, the struct 'fka' created here does not contain all of the arrays necessary to compute the weak-form operator. This functionality could be ported if necessary.
sourcemoment_kinetics.fokker_planck.init_fokker_planck_collisions_weak_form — Method
Function that initialises the arrays needed for Fokker Planck collisions using numerical integration to compute the Rosenbluth potentials only at the boundary and using an elliptic solve to obtain the potentials in the rest of the velocity space domain.
sourcemoment_kinetics.fokker_planck.moment_kinetic_collision_frequency_prefactor — Method
Function to account for the normalisation of the moment kinetic normalised distribution function by multiplying the input collision frequency by the normalisation factors.
Only applicable for self collisions.
sourcemoment_kinetics.fokker_planck.setup_fkpl_collisions_input — Method
Function for reading Fokker Planck collision operator input parameters. Structure the namelist as follows.
[fokker_planck_collisions]
use_fokker_planck = true
nuii = 1.0
frequency_option = "manual"sourcemoment_kinetics.fokker_planck.setup_fp_nl_solve — Method
Function to setup nonlinear_solver struct for implicit Fokker-Planck collisions. An input namelist of form
[fokker_planck_collisions_nonlinear_solver]
atol = 1.0e-10
rtol = 0.0
...is used, with the same defaults as the main [nonlinear_solver] namelist, apart from atol and rtol, which are set to their own defaults here.
moment_kinetics.fokker_planck.symmetric_matrix_inverse — Method
Function that solves A x = b for a matrix of the form
\[\begin{array}{ccc} A_{00} & 0 & A_{02} \\ 0 & A_{11} & A_{12} \\ A_{02} & A_{12} & A_{22} \\ \end{array}\]
appropriate for the moment numerical conserving terms used in the Fokker-Planck collision operator.
sourcemoment_kinetics.fokker_planck.symmetric_matrix_inverse — Method
Function that solves A x = b for a matrix of the form
\[\begin{array}{ccc} A_{00} & A_{01} & A_{02} \\ A_{01} & A_{11} & A_{12} \\ A_{02} & A_{12} & A_{22} \\ \end{array}\]
appropriate for moment numerical conserving terms.
source