8 #include "sll_assert.h"
9 #include "sll_memory.h"
10 #include "sll_working_precision.h"
45 logical :: electrostatic = .false.
64 call self%helper%reinit_fields()
72 sll_real64,
intent(in) :: dt
73 sll_int32,
intent(in) :: number_steps
77 if(self%electrostatic)
then
78 do i_step = 1, number_steps
79 call self%helper%advect_x( dt*0.5_f64 )
80 call self%helper%advect_vb( dt*0.5_f64 )
81 call self%helper%advect_e( dt )
82 call self%helper%advect_vb( dt*0.5_f64 )
83 call self%helper%advect_x( dt*0.5_f64 )
86 do i_step = 1, number_steps
87 call self%helper%advect_eb( dt*0.5_f64 )
88 call self%helper%advect_x( dt*0.5_f64 )
89 call self%helper%advect_vb( dt*0.5_f64 )
90 call self%helper%advect_e( dt )
91 call self%helper%advect_vb( dt*0.5_f64 )
92 call self%helper%advect_x( dt*0.5_f64 )
93 call self%helper%advect_eb( dt*0.5_f64 )
103 sll_real64,
intent(in) :: dt
104 sll_int32,
intent(in) :: number_steps
108 if(self%electrostatic)
then
109 do i_step = 1,number_steps
110 call self%helper%advect_e( dt )
111 call self%helper%advect_vb( dt )
112 call self%helper%advect_x( dt )
115 do i_step = 1,number_steps
116 call self%helper%advect_e( dt )
117 call self%helper%advect_vb( dt )
118 call self%helper%advect_x( dt )
119 call self%helper%advect_eb( dt )
129 sll_real64,
intent(in) :: dt
130 sll_int32,
intent(in) :: number_steps
134 if(self%electrostatic)
then
135 do i_step = 1,number_steps
136 call self%helper%advect_e( dt )
137 call self%helper%advect_vb( dt )
138 call self%helper%advect_x( dt )
141 do i_step = 1,number_steps
142 call self%helper%advect_e( dt )
143 call self%helper%advect_vb( dt )
144 call self%helper%advect_x( dt )
145 call self%helper%advect_eb( dt )
166 boundary_particles, &
179 sll_real64,
target,
intent(in) :: phi_dofs(:)
180 sll_real64,
target,
intent(in) :: efield_dofs(:,:)
181 sll_real64,
target,
intent(in) :: bfield_dofs(:)
182 sll_real64,
intent(in) :: x_min
183 sll_real64,
intent(in) :: lx
185 sll_int32,
optional,
intent( in ) :: boundary_particles
186 sll_real64,
intent(in),
optional :: solver_tolerance
187 sll_real64,
optional,
intent(in) :: force_sign
189 sll_int32,
optional,
intent(in) :: i_weight
190 sll_real64,
optional,
intent(in) :: betar(2)
191 logical,
optional :: electrostatic
192 logical,
optional,
intent(in) :: jmean
194 sll_int32 :: boundary_particles_set
195 sll_real64 :: solver_tolerance_set, force_sign_set, betar_set(2)
198 if (
present(boundary_particles))
then
199 boundary_particles_set = boundary_particles
201 boundary_particles_set = 100
204 if (
present(solver_tolerance) )
then
205 solver_tolerance_set = solver_tolerance
207 solver_tolerance_set = 1d-12
210 if(
present(force_sign) )
then
211 force_sign_set = force_sign
213 force_sign_set = 1._f64
216 if(
present(betar) )
then
222 if(
present(electrostatic) )
then
223 self%electrostatic = electrostatic
226 if (
present(jmean))
then
232 if(
present( control_variate) )
then
233 call self%helper%init(maxwell_solver, &
243 boundary_particles=boundary_particles_set, &
244 solver_tolerance=solver_tolerance_set, &
245 force_sign=force_sign_set, &
246 control_variate = control_variate,&
251 call self%helper%init(maxwell_solver, &
261 boundary_particles=boundary_particles_set, &
262 solver_tolerance=solver_tolerance_set, &
264 force_sign=force_sign_set, &
286 boundary_particles, &
298 sll_real64,
target,
intent(in) :: phi_dofs(:)
299 sll_real64,
target,
intent(in) :: efield_dofs(:,:)
300 sll_real64,
target,
intent(in) :: bfield_dofs(:)
301 sll_real64,
intent(in) :: x_min
302 sll_real64,
intent(in) :: lx
304 character(len=*),
intent(in) :: filename
305 sll_int32,
optional,
intent( in ) :: boundary_particles
306 sll_real64,
optional,
intent(in) :: force_sign
308 sll_int32,
optional,
intent(in) :: i_weight
309 sll_real64,
optional,
intent( in ) :: betar(2)
310 logical,
optional,
intent(in) :: electrostatic
311 logical,
optional,
intent(in) :: jmean
313 sll_int32 :: boundary_particles_set
314 sll_real64 :: force_sign_set, betar_set(2)
318 if (
present(boundary_particles))
then
319 boundary_particles_set = boundary_particles
321 boundary_particles_set = 100
324 if(
present(force_sign) )
then
325 force_sign_set = force_sign
327 force_sign_set = 1._f64
330 if(
present(betar) )
then
336 if(
present(electrostatic) )
then
337 self%electrostatic = electrostatic
339 self%electrostatic = .false.
342 if (
present(jmean))
then
348 if(
present( control_variate) )
then
349 call self%helper%init_from_file( maxwell_solver, &
360 boundary_particles = boundary_particles_set, &
361 force_sign=force_sign_set, &
362 control_variate = control_variate, &
367 call self%helper%init_from_file( maxwell_solver, &
378 boundary_particles = boundary_particles_set, &
380 force_sign=force_sign_set, &
392 call self%helper%free()
Module interface to solve Maxwell's equations in 1D.
Base class for kernel smoothers for accumulation and field evaluation in PIC.
Base class for Hamiltonian splittings.
Particle pusher based on energy (not charge-conserving) discrete gradient method, semi-implicit.
subroutine reinit_fields(self)
Finalization.
subroutine initialize_file_pic_vm_1d2v_disgrade(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, control_variate, i_weight, betar, electrostatic, jmean)
Constructor.
subroutine delete_pic_vm_1d2v_disgrade(self)
Destructor.
subroutine lie_splitting_back_pic_vm_1d2v_disgrade(self, dt, number_steps)
Lie splitting.
subroutine strang_splitting_pic_vm_1d2v_disgrade(self, dt, number_steps)
Strang splitting.
subroutine initialize_pic_vm_1d2v_disgrade(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, force_sign, control_variate, i_weight, betar, electrostatic, jmean)
Constructor.
subroutine lie_splitting_pic_vm_1d2v_disgrade(self, dt, number_steps)
Lie splitting.
Basic type of a kernel smoother used for PIC simulations.
Type for Hamiltonian splittings.
Hamiltonian splitting type for Vlasov-Maxwell 1d2v.
Hamiltonian splitting type for Vlasov-Maxwell 1d2v.