8 #include "sll_assert.h"
9 #include "sll_memory.h"
10 #include "sll_working_precision.h"
49 logical :: electrostatic = .false.
68 sll_real64,
intent(in) :: dt
69 sll_int32,
intent(in) :: number_steps
72 if(self%electrostatic)
then
73 do i_step = 1, number_steps
74 call self%helper%advect_x( dt*0.5_f64 )
75 call self%helper%advect_vb( dt*0.5_f64 )
76 call self%helper%advect_e( dt )
77 call self%helper%advect_vb( dt*0.5_f64 )
78 call self%helper%advect_x( dt*0.5_f64 )
81 do i_step = 1, number_steps
82 call self%helper%advect_eb( dt*0.5_f64 )
83 call self%helper%advect_x( dt*0.5_f64 )
84 call self%helper%advect_vb( dt*0.5_f64 )
85 call self%helper%advect_e( dt )
86 call self%helper%advect_vb( dt*0.5_f64 )
87 call self%helper%advect_x( dt*0.5_f64 )
88 call self%helper%advect_eb( dt*0.5_f64 )
98 sll_real64,
intent(in) :: dt
99 sll_int32,
intent(in) :: number_steps
103 if(self%electrostatic)
then
104 do i_step = 1,number_steps
105 call self%helper%advect_x( dt )
106 call self%helper%advect_vb( dt )
107 call self%helper%advect_e( dt )
110 do i_step = 1,number_steps
111 call self%helper%advect_x( dt )
112 call self%helper%advect_vb( dt )
113 call self%helper%advect_eb( dt )
114 call self%helper%advect_e( dt )
124 sll_real64,
intent(in) :: dt
125 sll_int32,
intent(in) :: number_steps
129 if(self%electrostatic)
then
130 do i_step = 1,number_steps
131 call self%helper%advect_e( dt )
132 call self%helper%advect_vb( dt )
133 call self%helper%advect_x( dt )
136 do i_step = 1,number_steps
137 call self%helper%advect_e( dt )
138 call self%helper%advect_vb( dt )
139 call self%helper%advect_x( dt )
140 call self%helper%advect_eb( dt )
152 particle_mesh_coupling, &
160 boundary_particles, &
172 sll_real64,
target,
intent( in ) :: phi_dofs(:)
173 sll_real64,
target,
intent( in ) :: efield_dofs(:)
174 sll_real64,
target,
intent( in ) :: bfield_dofs(:)
175 sll_real64,
intent( in ) :: x_min(3)
176 sll_real64,
intent( in ) :: lx(3)
178 sll_int32,
optional,
intent( in ) :: boundary_particles
179 sll_real64,
optional,
intent( in ) :: solver_tolerance
180 sll_real64,
optional,
intent(in) :: betar(2)
181 sll_real64,
optional,
intent(in) :: force_sign
182 logical,
optional,
intent( in ) :: electrostatic
183 sll_real64,
optional,
target,
intent( in ) :: rhob(:)
185 logical,
optional,
intent(in) :: jmean
187 sll_real64 :: solver_tolerance_set
188 sll_real64 :: betar_set(2), force_sign_set
189 sll_int32 :: boundary_particles_set
192 if(
present(boundary_particles) )
then
193 boundary_particles_set = boundary_particles
195 boundary_particles_set = 100
198 if (
present(solver_tolerance) )
then
199 solver_tolerance_set = solver_tolerance
201 solver_tolerance_set = 1d-12
204 if (
present(betar) )
then
210 if(
present(force_sign) )
then
211 force_sign_set = force_sign
213 force_sign_set = 1._f64
216 if(
present(electrostatic) )
then
217 self%electrostatic = electrostatic
220 if (
present(jmean))
then
226 if(
present( control_variate) )
then
227 call self%helper%init(maxwell_solver, &
228 particle_mesh_coupling, &
236 boundary_particles=boundary_particles_set, &
237 solver_tolerance=solver_tolerance_set, &
239 force_sign=force_sign_set, &
240 control_variate = control_variate,&
243 call self%helper%init(maxwell_solver, &
244 particle_mesh_coupling, &
252 boundary_particles=boundary_particles_set, &
253 solver_tolerance=solver_tolerance_set, &
255 force_sign=force_sign_set, &
266 particle_mesh_coupling, &
275 boundary_particles, &
286 sll_real64,
target,
intent( in ) :: phi_dofs(:)
287 sll_real64,
target,
intent( in ) :: efield_dofs(:)
288 sll_real64,
target,
intent( in ) :: bfield_dofs(:)
289 sll_real64,
intent( in ) :: x_min(3)
290 sll_real64,
intent( in ) :: lx(3)
292 character(len=*),
intent( in ) :: filename
293 sll_int32,
optional,
intent( in ) :: boundary_particles
294 sll_real64,
optional,
intent( in ) :: betar(2)
295 sll_real64,
optional,
intent(in) :: force_sign
296 logical,
optional,
intent( in ) :: electrostatic
297 sll_real64,
optional,
target,
intent( in ) :: rhob(:)
299 logical,
optional,
intent(in) :: jmean
301 sll_real64 :: betar_set(2), force_sign_set
302 sll_int32 :: boundary_particles_set
305 if(
present(boundary_particles) )
then
306 boundary_particles_set = boundary_particles
308 boundary_particles_set = 100
311 if(
present(force_sign) )
then
312 force_sign_set = force_sign
314 force_sign_set = 1._f64
317 if(
present(electrostatic) )
then
318 self%electrostatic = electrostatic
321 if(
present(betar) )
then
327 if (
present(jmean))
then
333 if(
present( control_variate) )
then
334 call self%helper%init_from_file( maxwell_solver, &
335 particle_mesh_coupling, &
344 boundary_particles = boundary_particles_set, &
346 force_sign=force_sign_set, &
348 control_variate = control_variate, &
351 call self%helper%init_from_file( maxwell_solver, &
352 particle_mesh_coupling, &
361 boundary_particles = boundary_particles_set, &
363 force_sign=force_sign_set, &
374 call self%helper%free()
Module interface to solve Maxwell's equations in 3D.
Base class for kernel smoothers for accumulation and field evaluation in PIC.
Base class for Hamiltonian splittings.
Particle pusher based on antisymmetric splitting with AVF for 3d3v Vlasov-Maxwell with coordinate tra...
Particle pusher based on antisymmetric splitting with disgradE for 3d3v Vlasov-Maxwell.
subroutine initialize_pic_vm_3d3v(self, maxwell_solver, particle_mesh_coupling, particle_group, phi_dofs, efield_dofs, bfield_dofs, x_min, Lx, filter, boundary_particles, solver_tolerance, betar, force_sign, electrostatic, rhob, control_variate, jmean)
Constructor.
subroutine strang_splitting_pic_vm_3d3v(self, dt, number_steps)
Finalization.
subroutine lie_splitting_pic_vm_3d3v(self, dt, number_steps)
Lie splitting.
subroutine delete_pic_vm_3d3v(self)
Destructor.
subroutine initialize_file_pic_vm_3d3v(self, maxwell_solver, particle_mesh_coupling, particle_group, phi_dofs, efield_dofs, bfield_dofs, x_min, Lx, filter, filename, boundary_particles, betar, force_sign, electrostatic, rhob, control_variate, jmean)
Constructor.
subroutine lie_splitting_back_pic_vm_3d3v(self, dt, number_steps)
Lie splitting (oposite ordering)
Particle pusher based on antisymmetric splitting with AVF for 3d3v Vlasov-Maxwell.
Basic type of a kernel smoother used for PIC simulations.
Type for Hamiltonian splittings.
Hamiltonian splitting type for Vlasov-Maxwell 3d3v.
Time propagator for Vlasov-Maxwell 3d3v.
Helper for implicit time propagators for Vlasov-Maxwell 3d3v.