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_pic_poisson_base Module Reference

Description

Base class for Poisson solver for particle methods.

Author
Katharina Kormann, IPP

This base class gives an abstract interface to the basic functions for accumulation of charge as well as the evaluation of a function at particle positions.

Todo:

Try to integrate also general elliptic and quasi neutral solvers.

Allow for more than one weights

Functions to compute derivatives

Derived types and interfaces

type  sll_c_pic_poisson
 Basic type of Poisson solver for PIC simulations. More...
 
interface  add_single
 
interface  eval_single
 
interface  eval_component_single
 
interface  empty
 
interface  compute_energy
 
interface  update_dofs_function
 
interface  linear_combination
 

Functions/Subroutines

subroutine add_charge_vector (self, n_part, position, marker_charge)
 
subroutine evaluate_rho_vector (self, position, n_part, func_value)
 
subroutine evaluate_phi_vector (self, position, n_part, func_value)
 
subroutine evaluate_field_vector (self, position, n_part, components, func_value)
 
subroutine solve_phi (self)
 Just compute phi, by default it calls the combined procedure solve. More...
 
subroutine solve_fields (self)
 Just compute the efield, by default it calls the combined procedure solve. More...
 

Function/Subroutine Documentation

◆ add_charge_vector()

subroutine sll_m_pic_poisson_base::add_charge_vector ( class (sll_c_pic_poisson), intent(inout)  self,
integer(kind=i32), intent(in)  n_part,
real(kind=f64), dimension(self%dim, n_part), intent(in)  position,
real(kind=f64), dimension(n_part), intent(in)  marker_charge 
)
Parameters
[in,out]selfPic Poisson solver object
[in]n_partNumber of particles whos positions are given
[in]positionPosition of the particle
[in]marker_chargeParticle weights times charge

Definition at line 137 of file sll_m_pic_poisson_base.F90.

Here is the caller graph for this function:

◆ evaluate_field_vector()

subroutine sll_m_pic_poisson_base::evaluate_field_vector ( class (sll_c_pic_poisson), intent(inout)  self,
real(kind=f64), dimension(:,:), intent(in)  position,
integer(kind=i32), intent(in)  n_part,
integer(kind=i32), dimension(:), intent(in)  components,
real(kind=f64), dimension(:,:), intent(out)  func_value 
)
private
Parameters
[in,out]selfPic Poisson solver object
[in]positionPosition of the particle (size (selfdim, n_part))
[in]n_partNumber of particles whos positions are given
[in]componentsComponents of the field to be computed
[out]func_valueValue(s) of the electric fields at given positions

Definition at line 193 of file sll_m_pic_poisson_base.F90.

Here is the caller graph for this function:

◆ evaluate_phi_vector()

subroutine sll_m_pic_poisson_base::evaluate_phi_vector ( class (sll_c_pic_poisson), intent(inout)  self,
real(kind=f64), dimension(:,:), intent(in)  position,
integer(kind=i32), intent(in)  n_part,
real(kind=f64), dimension(n_part), intent(out)  func_value 
)
private
Parameters
[in,out]selfPic Poisson solver object
[in]positionPosition of the particles(size (selfdim, n_part))
[in]n_partNumber of particles whos positions are given
[out]func_valueValue(s) of the electric fields at given position

Definition at line 174 of file sll_m_pic_poisson_base.F90.

Here is the caller graph for this function:

◆ evaluate_rho_vector()

subroutine sll_m_pic_poisson_base::evaluate_rho_vector ( class (sll_c_pic_poisson), intent(inout)  self,
real(kind=f64), dimension(:, :), intent(in)  position,
integer(kind=i32), intent(in)  n_part,
real(kind=f64), dimension(n_part), intent(out)  func_value 
)
private
Parameters
[in,out]selfPic Poisson solver object
[in]positionPosition of the particles (size (selfdim, n_part))
[in]n_partNumber of particles whos positions are given
[out]func_valueValue(s) of the electric fields at given position

Definition at line 155 of file sll_m_pic_poisson_base.F90.

Here is the caller graph for this function:

◆ solve_fields()

subroutine sll_m_pic_poisson_base::solve_fields ( class (sll_c_pic_poisson), intent(inout)  self)
private

Just compute the efield, by default it calls the combined procedure solve.

Parameters
[in,out]selfPic Poisson solver object

Definition at line 220 of file sll_m_pic_poisson_base.F90.

◆ solve_phi()

subroutine sll_m_pic_poisson_base::solve_phi ( class (sll_c_pic_poisson), intent(inout)  self)
private

Just compute phi, by default it calls the combined procedure solve.

Parameters
[in,out]selfPic Poisson solver object

Definition at line 212 of file sll_m_pic_poisson_base.F90.

    Report Typos and Errors