Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_time_propagator_pic_vm_2d3v_hs.F90
Go to the documentation of this file.
1 
7 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_assert.h"
9 #include "sll_memory.h"
10 #include "sll_working_precision.h"
11 
12  use sll_m_collective, only: &
15 
16  use sll_m_control_variate, only: &
18 
19  use sll_m_time_propagator_base, only: &
21 
24 
25  use sll_m_maxwell_2d_fem_fft, only: &
27 
28  use sll_m_particle_group_base, only: &
31 
32  use mpi, only: &
33  mpi_sum
34 
35  implicit none
36 
37  public :: &
41 
42  private
43 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
44 
47  type(sll_t_maxwell_2d_fem_fft), pointer :: maxwell_solver
48  type(sll_t_particle_mesh_coupling_spline_2d_feec), pointer :: particle_mesh_coupling
49  class(sll_t_particle_array), pointer :: particle_group
50 
51  sll_int32 :: spline_degree
52  sll_real64 :: lx(2)
53  sll_real64 :: x_min(2)
54  !sll_real64 :: delta_x !< Grid spacing
55  sll_int32 :: ntotal
56 
57  sll_real64 :: betar
58 
59  sll_real64, pointer :: efield_dofs(:)
60  sll_real64, pointer :: bfield_dofs(:)
61  sll_real64, allocatable :: j_dofs(:)
62  sll_real64, allocatable :: j_dofs_local(:)
63 
64  ! For control variate
65  class(sll_t_control_variates), pointer :: control_variate => null()
66  sll_int32 :: i_weight = 1
67 
68  contains
69  procedure :: operatorhp1 => operatorhp1_pic_vm_2d3v_hs
70  procedure :: operatorhp2 => operatorhp2_pic_vm_2d3v_hs
71  procedure :: operatorhp3 => operatorhp3_pic_vm_2d3v_hs
72  procedure :: operatorhe => operatorhe_pic_vm_2d3v_hs
73  procedure :: operatorhb => operatorhb_pic_vm_2d3v_hs
74  procedure :: lie_splitting => lie_splitting_pic_vm_2d3v_hs
75  procedure :: lie_splitting_back => lie_splitting_back_pic_vm_2d3v_hs
76  procedure :: strang_splitting => strang_splitting_pic_vm_2d3v_hs
77 
78  procedure :: init => initialize_pic_vm_2d3v_hs
79  procedure :: free => delete_pic_vm_2d3v_hs
80 
82 
83 contains
84 
85 
87  subroutine strang_splitting_pic_vm_2d3v_hs(self,dt, number_steps)
88  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(inout) :: self
89  sll_real64, intent(in) :: dt
90  sll_int32, intent(in) :: number_steps
91 
92  sll_int32 :: i_step
93 
94  do i_step = 1, number_steps
95  call self%operatorHB(0.5_f64*dt)
96  call self%operatorHE(0.5_f64*dt)
97  call self%operatorHp3(0.5_f64*dt)
98  call self%operatorHp2(0.5_f64*dt)
99  call self%operatorHp1(dt)
100  call self%operatorHp2(0.5_f64*dt)
101  call self%operatorHp3(0.5_f64*dt)
102  call self%operatorHE(0.5_f64*dt)
103  call self%operatorHB(0.5_f64*dt)
104 
105  end do
106 
107  end subroutine strang_splitting_pic_vm_2d3v_hs
108 
110  subroutine lie_splitting_pic_vm_2d3v_hs(self,dt, number_steps)
111  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(inout) :: self
112  sll_real64, intent(in) :: dt
113  sll_int32, intent(in) :: number_steps
114 
115  sll_int32 :: i_step
116 
117  do i_step = 1,number_steps
118  call self%operatorHE(dt)
119  call self%operatorHB(dt)
120  call self%operatorHp1(dt)
121  call self%operatorHp2(dt)
122  call self%operatorHp3(dt)
123  end do
124 
125 
126  end subroutine lie_splitting_pic_vm_2d3v_hs
127 
129  subroutine lie_splitting_back_pic_vm_2d3v_hs(self,dt, number_steps)
130  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(inout) :: self
131  sll_real64, intent(in) :: dt
132  sll_int32, intent(in) :: number_steps
133 
134  sll_int32 :: i_step
135 
136  do i_step = 1,number_steps
137  call self%operatorHp3(dt)
138  call self%operatorHp2(dt)
139  call self%operatorHp1(dt)
140  call self%operatorHB(dt)
141  call self%operatorHE(dt)
142  end do
143 
144  end subroutine lie_splitting_back_pic_vm_2d3v_hs
145 
146 
147  !---------------------------------------------------------------------------!
155  subroutine operatorhp1_pic_vm_2d3v_hs(self, dt)
156  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(inout) :: self
157  sll_real64, intent(in) :: dt
158 
159  !local variables
160  sll_int32 :: i_part, i_sp
161  sll_real64 :: x_new(3), vi(3), wi(1), x_old(3), wall(3)
162  sll_int32 :: n_cells
163  sll_real64 :: qoverm
164  sll_int32 :: i_weight
165 
166  i_weight = self%i_weight
167 
168  n_cells = self%maxwell_solver%n_total
169 
170  ! Here we have to accumulate j and integrate over the time interval.
171  ! At each k=1,...,n_grid, we have for s \in [0,dt]:
172  ! j_k(s) = \sum_{i=1,..,N_p} q_i N((x_k+sv_{1,k}-x_i)/h) v_k,
173  ! where h is the grid spacing and N the normalized B-spline
174  ! In order to accumulate the integrated j, we normalize the values of x to the grid spacing, calling them y, we have
175  ! j_k(s) = \sum_{i=1,..,N_p} q_i N(y_k+s/h v_{1,k}-y_i) v_k.
176  ! Now, we want the integral
177  ! \int_{0..dt} j_k(s) d s = \sum_{i=1,..,N_p} q_i v_k \int_{0..dt} N(y_k+s/h v_{1,k}-y_i) ds = \sum_{i=1,..,N_p} q_i v_k \int_{0..dt} N(y_k + w v_{1,k}-y_i) dw
178 
179 
180  self%j_dofs_local = 0.0_f64
181 
182  ! For each particle compute the index of the first DoF on the grid it contributes to and its position (normalized to cell size one). Note: j_dofs(_local) does not hold the values for j itself but for the integrated j.
183  ! Then update particle position: X_new = X_old + dt * V
184  do i_sp = 1, self%particle_group%n_species
185  qoverm = self%particle_group%group(i_sp)%species%q_over_m();
186  do i_part=1,self%particle_group%group(i_sp)%n_particles
187  ! Read out particle position and velocity
188  x_old = self%particle_group%group(i_sp)%get_x(i_part)
189  vi = self%particle_group%group(i_sp)%get_v(i_part)
190 
191  ! Then update particle position: X_new = X_old + dt * V
192  x_new(1) = x_old(1) + dt * vi(1)
193  x_new(2:3) = x_old(2:3)
194 
195  ! Get charge for accumulation of j
196  wi = self%particle_group%group(i_sp)%get_charge(i_part, i_weight)
197 
198  call self%particle_mesh_coupling%add_current_update_v_component1( x_old, x_new(1), wi(1), qoverm, &
199  self%bfield_dofs, vi, self%j_dofs_local)
200 
201  x_new(1) = modulo(x_new(1), self%Lx(1))
202  call self%particle_group%group(i_sp)%set_x(i_part, x_new)
203  call self%particle_group%group(i_sp)%set_v(i_part, vi)
204 
205  ! Update particle weights
206  if (self%particle_group%group(i_sp)%n_weights == 3 ) then
207  wall = self%particle_group%group(i_sp)%get_weights(i_part)
208  wall(3) = self%control_variate%cv(i_sp)%update_df_weight( x_new, vi, 0.0_f64, wall(1), wall(2) )
209  call self%particle_group%group(i_sp)%set_weights( i_part, wall )
210  end if
211  end do
212  end do
213 
214  self%j_dofs = 0.0_f64
215  ! MPI to sum up contributions from each processor
216  call sll_o_collective_allreduce( sll_v_world_collective, self%j_dofs_local, &
217  n_cells, mpi_sum, self%j_dofs)
218 
219  ! Update the electric field.
220  call self%maxwell_solver%compute_E_from_j(self%j_dofs, 1, self%efield_dofs(1:self%ntotal))
221 
222  end subroutine operatorhp1_pic_vm_2d3v_hs
223 
224  !---------------------------------------------------------------------------!
226  subroutine operatorhp2_pic_vm_2d3v_hs(self, dt)
227  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(inout) :: self
228  sll_real64, intent(in) :: dt
229 
230  !local variables
231  sll_int32 :: i_part, i_sp
232  sll_real64 :: x_new(3), vi(3), wi(1), x_old(3), wall(3)
233  sll_int32 :: n_cells
234  sll_real64 :: qoverm
235  sll_int32 :: i_weight
236 
237  i_weight = self%i_weight
238 
239  n_cells = self%maxwell_solver%n_total
240 
241  self%j_dofs_local = 0.0_f64
242 
243  ! For each particle compute the index of the first DoF on the grid it contributes to and its position (normalized to cell size one). Note: j_dofs(_local) does not hold the values for j itself but for the integrated j.
244  ! Then update particle position: X_new = X_old + dt * V
245  do i_sp = 1, self%particle_group%n_species
246  qoverm = self%particle_group%group(i_sp)%species%q_over_m();
247  do i_part=1,self%particle_group%group(i_sp)%n_particles
248  ! Read out particle position and velocity
249  x_old = self%particle_group%group(i_sp)%get_x(i_part)
250  vi = self%particle_group%group(i_sp)%get_v(i_part)
251 
252  ! Then update particle position: X_new = X_old + dt * V
253  x_new(2) = x_old(2) + dt * vi(2)
254  x_new(1) = x_old(1)
255 
256  ! Get charge for accumulation of j
257  wi = self%particle_group%group(i_sp)%get_charge(i_part, i_weight)
258 
259  call self%particle_mesh_coupling%add_current_update_v_component2( x_old, x_new(2), wi(1), qoverm, &
260  self%bfield_dofs, vi, self%j_dofs_local)
261 
262  x_new(2) = modulo(x_new(2), self%Lx(2))
263  call self%particle_group%group(i_sp)%set_x(i_part, x_new)
264  call self%particle_group%group(i_sp)%set_v(i_part, vi)
265 
266 
267  ! Update particle weights
268  if (self%particle_group%group(i_sp)%n_weights == 3 ) then
269  wall = self%particle_group%group(i_sp)%get_weights(i_part)
270  wall(3) = self%control_variate%cv(i_sp)%update_df_weight( x_new, vi, 0.0_f64, wall(1), wall(2) )
271  call self%particle_group%group(i_sp)%set_weights( i_part, wall )
272  end if
273  end do
274 
275  end do
276 
277  self%j_dofs = 0.0_f64
278  ! MPI to sum up contributions from each processor
279  call sll_o_collective_allreduce( sll_v_world_collective, self%j_dofs_local, &
280  n_cells, mpi_sum, self%j_dofs)
281 
282  ! Update the electric field.
283  call self%maxwell_solver%compute_E_from_j(self%j_dofs, 2, self%efield_dofs(self%ntotal+1:2*self%ntotal))
284 
285  end subroutine operatorhp2_pic_vm_2d3v_hs
286 
287  !---------------------------------------------------------------------------!
289  subroutine operatorhp3_pic_vm_2d3v_hs(self, dt)
290  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(inout) :: self
291  sll_real64, intent(in) :: dt
292 
293  !local variables
294  sll_int32 :: i_part, i_sp
295  sll_real64 :: vi(3), wi(1), x_old(3), wall(3)
296  sll_int32 :: n_cells
297  sll_real64 :: qoverm
298  sll_real64 :: bfield(2)
299  sll_int32 :: i_weight
300 
301  i_weight = self%i_weight
302 
303  n_cells = self%maxwell_solver%n_total
304 
305  self%j_dofs_local = 0.0_f64
306 
307  ! For each particle compute the index of the first DoF on the grid it contributes to and its position (normalized to cell size one). Note: j_dofs(_local) does not hold the values for j itself but for the integrated j.
308  ! Then update particle position: X_new = X_old + dt * V
309  do i_sp = 1, self%particle_group%n_species
310  qoverm = self%particle_group%group(i_sp)%species%q_over_m();
311  do i_part=1,self%particle_group%group(i_sp)%n_particles
312  ! Read out particle position and velocity
313  x_old = self%particle_group%group(i_sp)%get_x(i_part)
314  vi = self%particle_group%group(i_sp)%get_v(i_part)
315  ! Get charge for accumulation of j
316  wi = self%particle_group%group(i_sp)%get_charge(i_part, i_weight)
317 
318  call self%particle_mesh_coupling%evaluate &
319  (x_old(1:2), [self%spline_degree, self%spline_degree-1], &
320  self%bfield_dofs(1:self%ntotal), bfield(1))
321  call self%particle_mesh_coupling%evaluate &
322  (x_old(1:2), [self%spline_degree-1, self%spline_degree], &
323  self%bfield_dofs(self%ntotal+1:self%ntotal*2), bfield(2))
324 
325  vi(1) = vi(1) - dt *qoverm*vi(3)*bfield(2)
326  vi(2) = vi(2) + dt *qoverm*vi(3)*bfield(1)
327 
328  ! was pp
329  call self%particle_mesh_coupling%add_charge( x_old(1:2), wi(1)*vi(3), &
330  [self%spline_degree, self%spline_degree], self%j_dofs_local )
331 
332  call self%particle_group%group(i_sp)%set_v(i_part, vi)
333 
334 
335  ! Update particle weights
336  if (self%particle_group%group(i_sp)%n_weights == 3 ) then
337  wall = self%particle_group%group(i_sp)%get_weights(i_part)
338  wall(3) = self%control_variate%cv(i_sp)%update_df_weight( x_old, vi, 0.0_f64, wall(1), wall(2) )
339  call self%particle_group%group(i_sp)%set_weights( i_part, wall )
340  end if
341  end do
342 
343  end do
344 
345  self%j_dofs = 0.0_f64
346  ! MPI to sum up contributions from each processor
347  call sll_o_collective_allreduce( sll_v_world_collective, self%j_dofs_local, &
348  n_cells, mpi_sum, self%j_dofs)
349 
350  ! Update the electric field.
351  self%j_dofs = self%j_dofs*dt
352  call self%maxwell_solver%compute_E_from_j(self%j_dofs, 3, self%efield_dofs(2*self%ntotal+1:3*self%ntotal))
353 
354  end subroutine operatorhp3_pic_vm_2d3v_hs
355 
356 
357  !---------------------------------------------------------------------------!
363  subroutine operatorhe_pic_vm_2d3v_hs(self, dt)
364  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(inout) :: self
365  sll_real64, intent(in) :: dt
366 
367  !local variables
368  sll_int32 :: i_part, i_sp
369  sll_real64 :: v_new(3), xi(3), wall(3)
370  sll_real64 :: efield(3)
371  sll_real64 :: qm
372  sll_int32 :: i_weight
373 
374  i_weight = self%i_weight
375 
376  ! TODO: Optimize using evaluate multiple
377 
378  !write(25,*) self%efield_dofs
379 
380  do i_sp = 1, self%particle_group%n_species
381  qm = self%particle_group%group(i_sp)%species%q_over_m();
382  ! V_new = V_old + dt * E
383  do i_part=1,self%particle_group%group(i_sp)%n_particles
384  v_new = self%particle_group%group(i_sp)%get_v(i_part)
385  ! Evaluate efields at particle position
386  xi = self%particle_group%group(i_sp)%get_x(i_part)
387 
388  !print*, self%efield_dofs
389  call self%particle_mesh_coupling%evaluate &
390  (xi(1:2), [self%spline_degree-1, self%spline_degree], &
391  self%efield_dofs(1:self%ntotal), efield(1))
392  call self%particle_mesh_coupling%evaluate &
393  (xi(1:2), [self%spline_degree, self%spline_degree-1], &
394  self%efield_dofs(self%ntotal+1:2*self%ntotal), efield(2))
395  call self%particle_mesh_coupling%evaluate &
396  (xi(1:2), [self%spline_degree, self%spline_degree], &
397  self%efield_dofs(2*self%ntotal+1:3*self%ntotal), efield(3))
398 
399  v_new = v_new + dt* qm * efield
400  call self%particle_group%group(i_sp)%set_v(i_part, v_new)
401 
402 
403  ! Update particle weights
404  if (self%particle_group%group(i_sp)%n_weights == 3 ) then
405  wall = self%particle_group%group(i_sp)%get_weights(i_part)
406  wall(3) = self%control_variate%cv(i_sp)%update_df_weight( xi, v_new, 0.0_f64, wall(1), wall(2) )
407  call self%particle_group%group(i_sp)%set_weights( i_part, wall )
408  end if
409  end do
410  end do
411 
412  ! Update bfield
413  call self%maxwell_solver%compute_B_from_E( &
414  dt, self%efield_dofs, self%bfield_dofs)
415 
416  end subroutine operatorhe_pic_vm_2d3v_hs
417 
418 
419  !---------------------------------------------------------------------------!
425  subroutine operatorhb_pic_vm_2d3v_hs(self, dt)
426  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(inout) :: self
427  sll_real64, intent(in) :: dt
428 
429  ! Update efield2
430  call self%maxwell_solver%compute_E_from_B(&
431  dt*self%betar, self%bfield_dofs, self%efield_dofs)
432 
433  end subroutine operatorhb_pic_vm_2d3v_hs
434 
435 
436  !---------------------------------------------------------------------------!
439  self, &
440  maxwell_solver, &
441  particle_mesh_coupling, &
442  particle_group, &
443  efield_dofs, &
444  bfield_dofs, &
445  x_min, &
446  Lx, &
447  control_variate, &
448  betar)
449  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent(out) :: self
450  type(sll_t_maxwell_2d_fem_fft), pointer, intent(in) :: maxwell_solver
451  type(sll_t_particle_mesh_coupling_spline_2d_feec), pointer, intent(in) :: particle_mesh_coupling
452  class(sll_t_particle_array), pointer, intent(in) :: particle_group
453  sll_real64, pointer, intent(in) :: efield_dofs(:)
454  sll_real64, pointer, intent(in) :: bfield_dofs(:)
455  sll_real64, intent(in) :: x_min(2)
456  sll_real64, intent(in) :: lx(2)
457  class(sll_t_control_variates), optional, target, intent(in) :: control_variate
458  sll_real64, optional, intent(in) :: betar
459  !local variables
460  sll_int32 :: ierr
461 
462  self%maxwell_solver => maxwell_solver
463  self%particle_mesh_coupling => particle_mesh_coupling
464  self%particle_group => particle_group
465  self%efield_dofs => efield_dofs
466  self%bfield_dofs => bfield_dofs
467 
468  self%ntotal =self%particle_mesh_coupling%n_dofs
469  sll_allocate(self%j_dofs(self%ntotal), ierr)
470  sll_allocate(self%j_dofs_local(self%ntotal), ierr)
471 
472  self%spline_degree = self%particle_mesh_coupling%spline_degree
473  self%x_min = x_min
474  self%Lx = lx
475 
476  if (present(control_variate)) then
477  allocate(self%control_variate )
478  allocate(self%control_variate%cv(particle_group%n_species) )
479  self%control_variate => control_variate
480  self%i_weight = 3
481  end if
482 
483  if (present(betar)) then
484  self%betar = betar!32.89_f64
485  else
486  self%betar = 1.0_f64
487  end if
488 
489  end subroutine initialize_pic_vm_2d3v_hs
490 
491  !---------------------------------------------------------------------------!
493  subroutine delete_pic_vm_2d3v_hs(self)
494  class(sll_t_time_propagator_pic_vm_2d3v_hs), intent( inout ) :: self
495 
496  deallocate(self%j_dofs)
497  deallocate(self%j_dofs_local)
498  self%maxwell_solver => null()
499  self%particle_mesh_coupling => null()
500  self%particle_group => null()
501  self%efield_dofs => null()
502  self%bfield_dofs => null()
503 
504  end subroutine delete_pic_vm_2d3v_hs
505 
506 
507  !---------------------------------------------------------------------------!
510  splitting, &
511  maxwell_solver, &
512  particle_mesh_coupling, &
513  particle_group, &
514  efield_dofs, &
515  bfield_dofs, &
516  x_min, &
517  Lx, &
518  control_variate, &
519  betar)
520  class(sll_c_time_propagator_base), allocatable, intent(out) :: splitting
521  type(sll_t_maxwell_2d_fem_fft), pointer, intent(in) :: maxwell_solver
522  type(sll_t_particle_mesh_coupling_spline_2d_feec), pointer, intent(in) :: particle_mesh_coupling
523  class(sll_t_particle_array),pointer, intent(in) :: particle_group
524  sll_real64, pointer, intent(in) :: efield_dofs(:)
525  sll_real64, pointer, intent(in) :: bfield_dofs(:)
526  sll_real64, intent(in) :: x_min(2)
527  sll_real64, intent(in) :: lx(2)
528  class(sll_t_control_variates), optional, target, intent(in) :: control_variate
529  sll_real64, intent( in ), optional :: betar
530  !local variables
531  sll_int32 :: ierr
532 
533 
534  sll_allocate(sll_t_time_propagator_pic_vm_2d3v_hs :: splitting, ierr)
535 
536  select type (splitting)
538  if (present(control_variate) ) then
539  if (present(betar) ) then
540  call splitting%init(&
541  maxwell_solver, &
542  particle_mesh_coupling, &
543  particle_group, &
544  efield_dofs, &
545  bfield_dofs, &
546  x_min, &
547  lx, &
549  betar=betar)
550  else
551  call splitting%init(&
552  maxwell_solver, &
553  particle_mesh_coupling, &
554  particle_group, &
555  efield_dofs, &
556  bfield_dofs, &
557  x_min, &
558  lx, &
560  end if
561  else
562  if (present(betar) ) then
563  call splitting%init(&
564  maxwell_solver, &
565  particle_mesh_coupling, &
566  particle_group, &
567  efield_dofs, &
568  bfield_dofs, &
569  x_min, &
570  lx, &
571  betar=betar)
572  else
573  call splitting%init(&
574  maxwell_solver, &
575  particle_mesh_coupling, &
576  particle_group, &
577  efield_dofs, &
578  bfield_dofs, &
579  x_min, &
580  lx)
581  end if
582  end if
583 
584  end select
585 
587 
588  !---------------------------------------------------------------------------!
591  splitting, &
592  maxwell_solver, &
593  particle_mesh_coupling, &
594  particle_group, &
595  efield_dofs, &
596  bfield_dofs, &
597  x_min, &
598  Lx)
599  class(sll_c_time_propagator_base), pointer, intent(out) :: splitting
600  type(sll_t_maxwell_2d_fem_fft), pointer, intent(in) :: maxwell_solver
601  type(sll_t_particle_mesh_coupling_spline_2d_feec), pointer, intent(in) :: particle_mesh_coupling
602  class(sll_t_particle_array),pointer, intent(in) :: particle_group
603  sll_real64, pointer, intent(in) :: efield_dofs(:)
604  sll_real64, pointer, intent(in) :: bfield_dofs(:)
605  sll_real64, intent(in) :: x_min(2)
606  sll_real64, intent(in) :: lx(2)
607  !local variables
608  sll_int32 :: ierr
609 
610  sll_allocate(sll_t_time_propagator_pic_vm_2d3v_hs :: splitting, ierr)
611 
612  select type (splitting)
614  call splitting%init(&
615  maxwell_solver, &
616  particle_mesh_coupling, &
617  particle_group, &
618  efield_dofs, &
619  bfield_dofs, &
620  x_min, &
621  lx)
622  end select
623 
625 
626 
Combines values from all processes and distributes the result back to all processes.
Parallelizing facility.
type(sll_t_collective_t), pointer, public sll_v_world_collective
Control of subset assignment. Processes with the same color are in the same new communicator.
Module interface to solve Maxwell's equations in 2D The linear systems are solved based on FFT diagno...
Particle mesh coupling for 3d with splines of arbitrary degree placed on a uniform tensor product mes...
Base class for Hamiltonian splittings.
Particle pusher based on Hamiltonian splitting for 2d3v Vlasov-Maxwell.
subroutine operatorhe_pic_vm_2d3v_hs(self, dt)
Push H_E: Equations to be solved \partial_t f + E_1 \partial_{v_1} f + E_2 \partial_{v_2} f = 0 -> V_...
subroutine, public sll_s_new_time_propagator_pic_vm_2d3v_hs(splitting, maxwell_solver, particle_mesh_coupling, particle_group, efield_dofs, bfield_dofs, x_min, Lx, control_variate, betar)
Constructor for allocatable abstract type.
subroutine lie_splitting_back_pic_vm_2d3v_hs(self, dt, number_steps)
Lie splitting (oposite ordering)
subroutine initialize_pic_vm_2d3v_hs(self, maxwell_solver, particle_mesh_coupling, particle_group, efield_dofs, bfield_dofs, x_min, Lx, control_variate, betar)
Constructor.
subroutine lie_splitting_pic_vm_2d3v_hs(self, dt, number_steps)
Lie splitting.
subroutine, public sll_s_new_time_propagator_pic_vm_2d3v_hs_ptr(splitting, maxwell_solver, particle_mesh_coupling, particle_group, efield_dofs, bfield_dofs, x_min, Lx)
Constructor for pointer abstract type.
subroutine operatorhp2_pic_vm_2d3v_hs(self, dt)
Push Hp2: Equations to solve are.
subroutine operatorhb_pic_vm_2d3v_hs(self, dt)
Push H_B: Equations to be solved V_new = V_old \partial_t E_1 = 0 -> E_{1,new} = E_{1,...
subroutine operatorhp1_pic_vm_2d3v_hs(self, dt)
Push Hp1: Equations to solve are TODO: UPDATE \partial_t f + v_1 \partial_{x_1} f = 0 -> X_new = X_ol...
subroutine operatorhp3_pic_vm_2d3v_hs(self, dt)
Push Hp3: Equations to solve are.
subroutine strang_splitting_pic_vm_2d3v_hs(self, dt, number_steps)
Finalization.
real(kind=f64) function, dimension(size(particle, 2)) control_variate(particle)
Particle mesh coupling in 3d based on (arbitrary degree) spline on a tensor product uniform mesh.
    Report Typos and Errors