Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Particle pusher based on Hamiltonian splitting for 3d3v Vlasov-Maxwell.
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... | |
|
private |
Helper function for operatorHp1.
[in,out] | self | time propagator object |
Definition at line 300 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
private |
Destructor.
[in,out] | self | time propagator object |
Definition at line 706 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
private |
Constructor.
[out] | self | time propagator object |
[in] | maxwell_solver | Maxwell solver |
[in] | particle_mesh_coupling | Particle mesh coupling |
[in] | particle_group | Particle group |
[in] | phi_dofs | array for the coefficients of the scalar potential |
[in] | efield_dofs | array for the coefficients of the efields |
[in] | bfield_dofs | array for the coefficients of the bfield |
[in] | x_min | Lower bound of x domain |
[in] | lx | Length of the domain in x direction. |
[in] | boundary_particles | particle boundary conditions |
[in] | betar | reciprocal plasma beta |
[in] | force_sign | sign of particle force |
[in] | electrostatic | true for electrostatic simulation |
[in] | rhob | charge at the boundary |
[in] | control_variate | Control variate (if delta f) |
[in] | jmean | logical for mean value of current |
Definition at line 601 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
private |
Lie splitting (oposite ordering)
[in,out] | self | time propagator object |
[in] | dt | time step |
[in] | number_steps | number of time steps |
Definition at line 185 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
private |
Lie splitting.
[in,out] | self | time propagator object |
[in] | dt | time step |
[in] | number_steps | number of time steps |
Definition at line 157 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
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$.
[in,out] | self | time propagator object |
[in] | dt | time step |
Definition at line 588 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
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$.
[in,out] | self | time propagator object |
[in] | dt | time step |
Definition at line 507 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
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.
[in,out] | self | time propagator object |
[in] | dt | time step |
Definition at line 219 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
private |
Push Hp2: Equations to solve are.
[in,out] | self | time propagator object |
[in] | dt | time step |
Definition at line 354 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
private |
Push Hp3: Equations to solve are.
[in,out] | self | time propagator object |
[in] | dt | time step |
Definition at line 430 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.
|
private |
Finalization.
Strang splitting
[in,out] | self | time propagator object |
[in] | dt | time step |
[in] | number_steps | number of time steps |
Definition at line 122 of file sll_m_time_propagator_pic_vm_3d3v_hs.F90.