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

Description

Particle pusher based on Hamiltonian splitting using in Crouseilles, Einkemmer, Faou for 3d3v Vlasov-Maxwell.

Author
Katharina Kormann, IPP

MPI parallelization by domain cloning. Periodic boundaries. Spline DoFs numerated by the point the spline starts. Reference: Kraus, Kormann, Sonnendrücker, Morrison: GEMPIC: Geometric ElectroMagnetic Particle-In-Cell Methods

Derived types and interfaces

type  sll_t_time_propagator_pic_vm_3d3v_cef
 Hamiltonian splitting type for Vlasov-Maxwell 3d3v. More...
 

Functions/Subroutines

subroutine strang_splitting_pic_vm_3d3v (self, dt, number_steps)
 Finalization. More...
 
subroutine lie_splitting_pic_vm_3d3v (self, dt, number_steps)
 Lie splitting. More...
 
subroutine lie_splitting_back_pic_vm_3d3v (self, dt, number_steps)
 Lie splitting (oposite ordering) More...
 
subroutine operatorhp_pic_vm_3d3v (self, dt)
 Push H_p: Equations to be solved $X^{n+1}=X^n+\Delta t V^n$ $M_1 e^n= M_1 e^n- \int_{t^n}^{t^{n+1}} \mathbb{\Lambda}^1(X(\tau))^\top d\tau \mathbb{W}_q V^n$. More...
 
subroutine compute_particle_boundary (self, xold, xnew, vi, wi)
 Helper function for operatorHp. More...
 
subroutine operatorhb_pic_vm_3d3v (self, dt)
 Push H_B: Equations to be solved $(\mathbb{I}-\Delta \frac{\Delta t q}{2 m} \mathbb{B}(X^n,b^n) V^{n+1}=(\mathbb{I}+ \frac{\Delta t q}{2 m} \mathbb{B}(X^n,b^n) ) V^n$ $M_1 e^{n+1}=M_1 e^n+\Delta t C^\top M_2 b^n$. More...
 
subroutine operatorhe_pic_vm_3d3v (self, dt)
 Push H_E: Equations to be solved $V^{n+1}=V^n+\Delta t\mathbb{W}_{\frac{q}{m}} \mathbb{Lambda}^1(X^n) e^n$ $b^{n+1}=b^n-\Delta t C e^n$. More...
 
subroutine initialize_pic_vm_3d3v (self, maxwell_solver, particle_mesh_coupling, particle_group, phi_dofs, efield_dofs, bfield_dofs, x_min, Lx, boundary_particles, betar, electrostatic, rhob, control_variate)
 Constructor. More...
 
subroutine delete_pic_vm_3d3v (self)
 Destructor. More...
 

Function/Subroutine Documentation

◆ compute_particle_boundary()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::compute_particle_boundary ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(inout)  self,
real(kind=f64), dimension(3), intent(in)  xold,
real(kind=f64), dimension(3), intent(inout)  xnew,
real(kind=f64), dimension(3), intent(inout)  vi,
real(kind=f64), dimension(1), intent(in)  wi 
)
private

Helper function for operatorHp.

Parameters
[in,out]selftime propagator object

Definition at line 213 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

Here is the caller graph for this function:

◆ delete_pic_vm_3d3v()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::delete_pic_vm_3d3v ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(inout)  self)
private

Destructor.

Parameters
[in,out]selftime propagator object

Definition at line 482 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

◆ initialize_pic_vm_3d3v()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::initialize_pic_vm_3d3v ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(out)  self,
class(sll_c_maxwell_3d_base), intent(in), target  maxwell_solver,
class(sll_c_particle_mesh_coupling_3d), intent(in), target  particle_mesh_coupling,
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), dimension(3), intent(in)  x_min,
real(kind=f64), dimension(3), intent(in)  Lx,
integer(kind=i32), intent(in), optional  boundary_particles,
real(kind=f64), dimension(2), intent(in), optional  betar,
logical, intent(in), optional  electrostatic,
real(kind=f64), dimension(:), intent(in), optional, target  rhob,
class(sll_t_control_variates), intent(in), optional, target  control_variate 
)
private

Constructor.

Parameters
[out]selftime propagator object
[in]maxwell_solverMaxwell solver
[in]particle_mesh_couplingParticle mesh coupling
[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]betarreciprocal plasma beta
[in]electrostatictrue for electrostatic simulation
[in]rhobcharge at the boundary
[in]control_variateControl variate (if delta f)

Definition at line 398 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

◆ lie_splitting_back_pic_vm_3d3v()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::lie_splitting_back_pic_vm_3d3v ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(inout)  self,
real(kind=f64), intent(in)  dt,
integer(kind=i32), intent(in)  number_steps 
)
private

Lie splitting (oposite ordering)

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

Definition at line 146 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

◆ lie_splitting_pic_vm_3d3v()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::lie_splitting_pic_vm_3d3v ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(inout)  self,
real(kind=f64), intent(in)  dt,
integer(kind=i32), intent(in)  number_steps 
)
private

Lie splitting.

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

Definition at line 128 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

◆ operatorhb_pic_vm_3d3v()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::operatorhb_pic_vm_3d3v ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(inout)  self,
real(kind=f64), intent(in)  dt 
)
private

Push H_B: Equations to be solved $(\mathbb{I}-\Delta \frac{\Delta t q}{2 m} \mathbb{B}(X^n,b^n) V^{n+1}=(\mathbb{I}+ \frac{\Delta t q}{2 m} \mathbb{B}(X^n,b^n) ) V^n$ $M_1 e^{n+1}=M_1 e^n+\Delta t C^\top M_2 b^n$.

Parameters
[in,out]selftime propagator object
[in]dttime step

Definition at line 268 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

◆ operatorhe_pic_vm_3d3v()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::operatorhe_pic_vm_3d3v ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(inout)  self,
real(kind=f64), intent(in)  dt 
)
private

Push H_E: Equations to be solved $V^{n+1}=V^n+\Delta t\mathbb{W}_{\frac{q}{m}} \mathbb{Lambda}^1(X^n) e^n$ $b^{n+1}=b^n-\Delta t C e^n$.

Parameters
[in,out]selftime propagator object
[in]dttime step

Definition at line 327 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

◆ operatorhp_pic_vm_3d3v()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::operatorhp_pic_vm_3d3v ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(inout)  self,
real(kind=f64), intent(in)  dt 
)
private

Push H_p: Equations to be solved $X^{n+1}=X^n+\Delta t V^n$ $M_1 e^n= M_1 e^n- \int_{t^n}^{t^{n+1}} \mathbb{\Lambda}^1(X(\tau))^\top d\tau \mathbb{W}_q V^n$.

Parameters
[in,out]selftime propagator object
[in]dttime step

Definition at line 166 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

Here is the call graph for this function:

◆ strang_splitting_pic_vm_3d3v()

subroutine sll_m_time_propagator_pic_vm_3d3v_cef::strang_splitting_pic_vm_3d3v ( class(sll_t_time_propagator_pic_vm_3d3v_cef), intent(inout)  self,
real(kind=f64), intent(in)  dt,
integer(kind=i32), intent(in)  number_steps 
)
private

Finalization.

Strang splitting

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

Definition at line 109 of file sll_m_time_propagator_pic_vm_3d3v_cef.F90.

    Report Typos and Errors