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

Description

Particle pusher based on Hamiltonian splitting 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_hs
 Hamiltonian splitting type for Vlasov-Maxwell 3d3v. More...
 

Functions/Subroutines

subroutine strang_splitting_pic_vm_3d3v_hs (self, dt, number_steps)
 Finalization. More...
 
subroutine lie_splitting_pic_vm_3d3v_hs (self, dt, number_steps)
 Lie splitting. More...
 
subroutine lie_splitting_back_pic_vm_3d3v_hs (self, dt, number_steps)
 Lie splitting (oposite ordering) More...
 
subroutine operatorhp1_pic_vm_3d3v_hs (self, dt)
 Push Hp1: Equations to solve are \partial_t f + v_1 \partial_{x_1} f = 0 -> xnew = xold + dt V_1 V_new,2 = V_old,2 + \int_0 h V_old,1 B_old \partial_t E_1 = - \int v_1 f(t,x_1, v) dv -> E_{1,new} = E_{1,old} - \int \int v_1 f(t,x_1+s v_1,v) dv ds \partial_t E_2 = 0 -> E_{2,new} = E_{2,old} \partial_t B = 0 => B_new = B_old. More...
 
subroutine compute_particle_boundary_current (self, xold, xnew, vi, wi, qoverm)
 Helper function for operatorHp1. More...
 
subroutine operatorhp2_pic_vm_3d3v_hs (self, dt)
 Push Hp2: Equations to solve are. More...
 
subroutine operatorhp3_pic_vm_3d3v_hs (self, dt)
 Push Hp3: Equations to solve are. More...
 
subroutine operatorhe_pic_vm_3d3v_hs (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 operatorhb_pic_vm_3d3v_hs (self, dt)
 Push H_B: Equations to be solved $M_1 e^{n+1}=M_1 e^n+\Delta t C^\top M_2 b^n$. More...
 
subroutine initialize_pic_vm_3d3v_hs (self, maxwell_solver, particle_mesh_coupling, particle_group, phi_dofs, efield_dofs, bfield_dofs, x_min, Lx, filter, boundary_particles, betar, force_sign, electrostatic, rhob, control_variate, jmean)
 Constructor. More...
 
subroutine delete_pic_vm_3d3v_hs (self)
 Destructor. More...
 

Function/Subroutine Documentation

◆ compute_particle_boundary_current()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::compute_particle_boundary_current ( class(sll_t_time_propagator_pic_vm_3d3v_hs), intent(inout)  self,
real(kind=f64), dimension(3), intent(inout)  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,
real(kind=f64), intent(in)  qoverm 
)
private

Helper function for operatorHp1.

Parameters
[in,out]selftime propagator object

Definition at line 300 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

Here is the caller graph for this function:

◆ delete_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::delete_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), intent(inout)  self)
private

Destructor.

Parameters
[in,out]selftime propagator object

Definition at line 706 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

◆ initialize_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::initialize_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), 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,
class(sll_c_filter_base_3d), target  filter,
integer(kind=i32), intent(in), optional  boundary_particles,
real(kind=f64), dimension(2), intent(in), optional  betar,
real(kind=f64), intent(in), optional  force_sign,
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,
logical, intent(in), optional  jmean 
)
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]force_signsign of particle force
[in]electrostatictrue for electrostatic simulation
[in]rhobcharge at the boundary
[in]control_variateControl variate (if delta f)
[in]jmeanlogical for mean value of current

Definition at line 601 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

◆ lie_splitting_back_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::lie_splitting_back_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), 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 185 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

◆ lie_splitting_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::lie_splitting_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), 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 157 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

◆ operatorhb_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::operatorhb_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), intent(inout)  self,
real(kind=f64), intent(in)  dt 
)
private

Push H_B: Equations to be solved $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 588 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

◆ operatorhe_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::operatorhe_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), 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 507 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

◆ operatorhp1_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::operatorhp1_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), intent(inout)  self,
real(kind=f64), intent(in)  dt 
)
private

Push Hp1: Equations to solve are \partial_t f + v_1 \partial_{x_1} f = 0 -> xnew = xold + dt V_1 V_new,2 = V_old,2 + \int_0 h V_old,1 B_old \partial_t E_1 = - \int v_1 f(t,x_1, v) dv -> E_{1,new} = E_{1,old} - \int \int v_1 f(t,x_1+s v_1,v) dv ds \partial_t E_2 = 0 -> E_{2,new} = E_{2,old} \partial_t B = 0 => B_new = B_old.

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

Definition at line 219 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

Here is the call graph for this function:

◆ operatorhp2_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::operatorhp2_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), intent(inout)  self,
real(kind=f64), intent(in)  dt 
)
private

Push Hp2: Equations to solve are.

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

Definition at line 354 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

◆ operatorhp3_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::operatorhp3_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), intent(inout)  self,
real(kind=f64), intent(in)  dt 
)
private

Push Hp3: Equations to solve are.

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

Definition at line 430 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

◆ strang_splitting_pic_vm_3d3v_hs()

subroutine sll_m_time_propagator_pic_vm_3d3v_hs::strang_splitting_pic_vm_3d3v_hs ( class(sll_t_time_propagator_pic_vm_3d3v_hs), 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 122 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.

    Report Typos and Errors