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

Description

Particle pusher based on energy and charge-conserving discrete gradient method, implicit.

Author
Katharina Kormann, IPP

MPI parallelization by domain cloning. Periodic boundaries. Spline DoFs numerated by the point the spline starts. Reference: Kormann, Sonnendrücker, Energy-conserving time propagation for a structure-preserving particle-in-cell Vlasov–Maxwell solver, Journal of Computational Physics 425, 109890, 2021

Derived types and interfaces

type  sll_t_time_propagator_pic_vm_1d2v_disgradec
 

Functions/Subroutines

subroutine strang_splitting_pic_vm_1d2v_disgradec (self, dt, number_steps)
 Strang splitting. More...
 
subroutine lie_splitting_pic_vm_1d2v_disgradec (self, dt, number_steps)
 Lie splitting. More...
 
subroutine lie_splitting_back_pic_vm_1d2v_disgradec (self, dt, number_steps)
 Lie splitting. More...
 
subroutine initialize_pic_vm_1d2v_disgradec (self, maxwell_solver, kernel_smoother_0, kernel_smoother_1, particle_group, phi_dofs, efield_dofs, bfield_dofs, x_min, Lx, filter, boundary_particles, solver_tolerance, iter_tolerance, max_iter, force_sign, betar, electrostatic, jmean)
 Constructor. More...
 
subroutine initialize_file_pic_vm_1d2v_disgradec (self, maxwell_solver, kernel_smoother_0, kernel_smoother_1, particle_group, phi_dofs, efield_dofs, bfield_dofs, x_min, Lx, filter, filename, boundary_particles, force_sign, betar, electrostatic, jmean)
 Constructor. More...
 
subroutine delete_pic_vm_1d2v_disgradec (self)
 
subroutine reinit_fields_disgradec (self)
 

Function/Subroutine Documentation

◆ delete_pic_vm_1d2v_disgradec()

subroutine sll_m_time_propagator_pic_vm_1d2v_disgradec::delete_pic_vm_1d2v_disgradec ( class(sll_t_time_propagator_pic_vm_1d2v_disgradec), intent(inout)  self)
private
Parameters
[in,out]selftime splitting object

Definition at line 339 of file sll_m_time_propagator_pic_vm_1d2v_disgradEC.F90.

◆ initialize_file_pic_vm_1d2v_disgradec()

subroutine sll_m_time_propagator_pic_vm_1d2v_disgradec::initialize_file_pic_vm_1d2v_disgradec ( class(sll_t_time_propagator_pic_vm_1d2v_disgradec), intent(out)  self,
class(sll_c_maxwell_1d_base), intent(in), target  maxwell_solver,
class(sll_c_particle_mesh_coupling_1d), intent(in), target  kernel_smoother_0,
class(sll_c_particle_mesh_coupling_1d), intent(in), target  kernel_smoother_1,
class(sll_t_particle_array), intent(in), target  particle_group,
real(kind=f64), dimension(:), intent(in), target  phi_dofs,
real(kind=f64), dimension(:,:), intent(in), target  efield_dofs,
real(kind=f64), dimension(:), intent(in), target  bfield_dofs,
real(kind=f64), intent(in)  x_min,
real(kind=f64), intent(in)  Lx,
class(sll_c_filter_base_1d), intent(in), target  filter,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in), optional  boundary_particles,
real(kind=f64), intent(in), optional  force_sign,
real(kind=f64), dimension(2), intent(in), optional  betar,
logical, intent(in), optional  electrostatic,
logical, intent(in), optional  jmean 
)
private

Constructor.

Parameters
[out]selftime propagator object
[in]maxwell_solverMaxwell solver
[in]kernel_smoother_0Kernel smoother
[in]kernel_smoother_1Kernel smoother
[in]particle_groupParticle group
[in]phi_dofsarray for the coefficients of phi
[in]efield_dofsarray for the coefficients of the efields
[in]bfield_dofsarray for the coefficients of the bfield
[in]x_minLower bound of x domain
[in]lxLength of the domain in x direction.
[in]boundary_particlesparticle boundary conditions
[in]force_signsign of particle force
[in]betarreciprocal plasma beta
[in]electrostatictrue for electrostatic simulation
[in]jmeanlogical for mean value of current

Definition at line 248 of file sll_m_time_propagator_pic_vm_1d2v_disgradEC.F90.

◆ initialize_pic_vm_1d2v_disgradec()

subroutine sll_m_time_propagator_pic_vm_1d2v_disgradec::initialize_pic_vm_1d2v_disgradec ( class(sll_t_time_propagator_pic_vm_1d2v_disgradec), intent(out)  self,
class(sll_c_maxwell_1d_base), intent(in), target  maxwell_solver,
class(sll_c_particle_mesh_coupling_1d), intent(in), target  kernel_smoother_0,
class(sll_c_particle_mesh_coupling_1d), intent(in), target  kernel_smoother_1,
class(sll_t_particle_array), intent(in), target  particle_group,
real(kind=f64), dimension(:), intent(in), target  phi_dofs,
real(kind=f64), dimension(:,:), intent(in), target  efield_dofs,
real(kind=f64), dimension(:), intent(in), target  bfield_dofs,
real(kind=f64), intent(in)  x_min,
real(kind=f64), intent(in)  Lx,
class(sll_c_filter_base_1d), intent(in), target  filter,
integer(kind=i32), intent(in), optional  boundary_particles,
real(kind=f64), intent(in), optional  solver_tolerance,
real(kind=f64), intent(in), optional  iter_tolerance,
integer(kind=i32), intent(in), optional  max_iter,
real(kind=f64), intent(in), optional  force_sign,
real(kind=f64), dimension(2), intent(in), optional  betar,
logical, intent(in), optional  electrostatic,
logical, intent(in), optional  jmean 
)
private

Constructor.

Parameters
[out]selftime splitting object
[in]maxwell_solverMaxwell solver
[in]kernel_smoother_0Kernel smoother
[in]kernel_smoother_1Kernel smoother
[in]particle_groupParticle group
[in]phi_dofsarray for the coefficients of the scalar potential
[in]efield_dofsarray for the coefficients of the efields
[in]bfield_dofsarray for the coefficients of the bfield
[in]x_minLower bound of x domain
[in]lxLength of the domain in x direction.
[in]boundary_particlesparticle boundary conditions
[in]solver_tolerancesolver tolerance
[in]iter_toleranceiteration tolerance
[in]max_itermaximal number of iterations
[in]force_signsign of particle force
[in]betarreciprocal plasma beta
[in]electrostatictrue for electrostatic simulation
[in]jmeanlogical for mean value of current

Definition at line 136 of file sll_m_time_propagator_pic_vm_1d2v_disgradEC.F90.

◆ lie_splitting_back_pic_vm_1d2v_disgradec()

subroutine sll_m_time_propagator_pic_vm_1d2v_disgradec::lie_splitting_back_pic_vm_1d2v_disgradec ( class(sll_t_time_propagator_pic_vm_1d2v_disgradec), intent(inout)  self,
real(kind=f64), intent(in)  dt,
integer(kind=i32), intent(in)  number_steps 
)
private

Lie splitting.

Parameters
[in,out]selftime splitting object
[in]dttime step
[in]number_stepsnumber of time steps

Definition at line 112 of file sll_m_time_propagator_pic_vm_1d2v_disgradEC.F90.

◆ lie_splitting_pic_vm_1d2v_disgradec()

subroutine sll_m_time_propagator_pic_vm_1d2v_disgradec::lie_splitting_pic_vm_1d2v_disgradec ( class(sll_t_time_propagator_pic_vm_1d2v_disgradec), intent(inout)  self,
real(kind=f64), intent(in)  dt,
integer(kind=i32), intent(in)  number_steps 
)
private

Lie splitting.

Parameters
[in,out]selftime splitting object
[in]dttime step
[in]number_stepsnumber of time steps

Definition at line 88 of file sll_m_time_propagator_pic_vm_1d2v_disgradEC.F90.

◆ reinit_fields_disgradec()

subroutine sll_m_time_propagator_pic_vm_1d2v_disgradec::reinit_fields_disgradec ( class(sll_t_time_propagator_pic_vm_1d2v_disgradec), intent(inout)  self)
private
Parameters
[in,out]selftime splitting object

Definition at line 347 of file sll_m_time_propagator_pic_vm_1d2v_disgradEC.F90.

◆ strang_splitting_pic_vm_1d2v_disgradec()

subroutine sll_m_time_propagator_pic_vm_1d2v_disgradec::strang_splitting_pic_vm_1d2v_disgradec ( class(sll_t_time_propagator_pic_vm_1d2v_disgradec), intent(inout)  self,
real(kind=f64), intent(in)  dt,
integer(kind=i32), intent(in)  number_steps 
)
private

Strang splitting.

Parameters
[in,out]selftime splitting object
[in]dttime step
[in]number_stepsnumber of time steps

Definition at line 61 of file sll_m_time_propagator_pic_vm_1d2v_disgradEC.F90.

    Report Typos and Errors