Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines
sll_m_poisson_1d_hmf Module Reference

Description

Module to solve Poisson equation for the HMF model.

Author
Pierre Navaro

You could find more information in

"On numerical Landau damping for splitting methods applied to the Vlasov-HMF model"

by Erwan Faou, Romain Horsin and Frédéric Rousset. https://arxiv.org/abs/1510.06555

Derived types and interfaces

type  sll_t_poisson_1d_hmf
 Implementation of the poisson 1d solver for the Vlasov-HMF model. More...
 

Functions/Subroutines

subroutine, public sll_s_poisson_1d_hmf_init (self, eta1_min, eta1_max, nc_eta1, error)
 Initialize the poisson 1d hmf solver. More...
 
subroutine, public sll_s_solve_poisson_1d_hmf (self, field, rhs)
 Solve the 1d equation for the Vlasov-HMF model. More...
 
subroutine compute_phi_from_rho_1d_hmf (poisson, phi, rho)
 solves More...
 
subroutine compute_e_from_rho_1d_hmf (poisson, E, rho)
 solves More...
 
subroutine, public sll_s_free_poisson_1d_hmf (self, error)
 Free the poisson 1d hmf solver. More...
 

Function/Subroutine Documentation

◆ compute_e_from_rho_1d_hmf()

subroutine sll_m_poisson_1d_hmf::compute_e_from_rho_1d_hmf ( class(sll_t_poisson_1d_hmf poisson,
real(kind=f64), dimension(:), intent(out)  E,
real(kind=f64), dimension(:), intent(in)  rho 
)
private

solves

\[ E = -\nabla \phi \]

with

\[ -\Delta \phi = \rho \]

in 2d

Definition at line 133 of file sll_m_poisson_1d_hmf.F90.

Here is the call graph for this function:

◆ compute_phi_from_rho_1d_hmf()

subroutine sll_m_poisson_1d_hmf::compute_phi_from_rho_1d_hmf ( class(sll_t_poisson_1d_hmf poisson,
real(kind=f64), dimension(:), intent(out)  phi,
real(kind=f64), dimension(:), intent(in)  rho 
)
private

solves

\[ -\Delta \phi = \rho \]

in 2d

Definition at line 123 of file sll_m_poisson_1d_hmf.F90.

◆ sll_s_free_poisson_1d_hmf()

subroutine, public sll_m_poisson_1d_hmf::sll_s_free_poisson_1d_hmf ( class(sll_t_poisson_1d_hmf self,
integer(kind=i32), intent(inout)  error 
)

Free the poisson 1d hmf solver.

Parameters
[in,out]errorerror code

Definition at line 143 of file sll_m_poisson_1d_hmf.F90.

◆ sll_s_poisson_1d_hmf_init()

subroutine, public sll_m_poisson_1d_hmf::sll_s_poisson_1d_hmf_init ( class(sll_t_poisson_1d_hmf), intent(out)  self,
real(kind=f64), intent(in)  eta1_min,
real(kind=f64), intent(in)  eta1_max,
integer(kind=i32), intent(in)  nc_eta1,
integer(kind=i32), intent(out)  error 
)

Initialize the poisson 1d hmf solver.

Parameters
[out]selfSolver structure
[in]nc_eta1number of cells
[out]errorerror code
[in]eta1_minleft corner
[in]eta1_maxright corner

Definition at line 76 of file sll_m_poisson_1d_hmf.F90.

Here is the call graph for this function:

◆ sll_s_solve_poisson_1d_hmf()

subroutine, public sll_m_poisson_1d_hmf::sll_s_solve_poisson_1d_hmf ( class(sll_t_poisson_1d_hmf), intent(inout)  self,
real(kind=f64), dimension(:), intent(out)  field,
real(kind=f64), dimension(:), intent(in)  rhs 
)

Solve the 1d equation for the Vlasov-HMF model.

Definition at line 99 of file sll_m_poisson_1d_hmf.F90.

Here is the call graph for this function:
Here is the caller graph for this function:
    Report Typos and Errors