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_3d3v_disgradE.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_control_variate, only: &
14 
15  use sll_m_time_propagator_base, only: &
17 
20 
23 
26 
27  use sll_m_maxwell_3d_base, only: &
29 
30  use sll_m_particle_group_base, only: &
33 
34  use sll_m_filter_base_3d, only: &
36 
37  implicit none
38 
39  public :: &
41 
42  private
43  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
44 
47 
49  logical :: electrostatic = .false.
50 
51  contains
52  procedure :: lie_splitting => lie_splitting_pic_vm_3d3v
53  procedure :: lie_splitting_back => lie_splitting_back_pic_vm_3d3v
54  procedure :: strang_splitting => strang_splitting_pic_vm_3d3v
55 
56  procedure :: init => initialize_pic_vm_3d3v
57  procedure :: init_from_file => initialize_file_pic_vm_3d3v
58  procedure :: free => delete_pic_vm_3d3v
59 
61 
62 contains
63 
64 
66  subroutine strang_splitting_pic_vm_3d3v(self,dt, number_steps)
67  class(sll_t_time_propagator_pic_vm_3d3v_disgrade), intent(inout) :: self
68  sll_real64, intent(in) :: dt
69  sll_int32, intent(in) :: number_steps
70  !local variables
71  sll_int32 :: i_step
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 )
79  end do
80  else
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 )
89  end do
90  end if
91 
92  end subroutine strang_splitting_pic_vm_3d3v
93 
94 
96  subroutine lie_splitting_pic_vm_3d3v(self,dt, number_steps)
97  class(sll_t_time_propagator_pic_vm_3d3v_disgrade), intent(inout) :: self
98  sll_real64, intent(in) :: dt
99  sll_int32, intent(in) :: number_steps
100  !local variables
101  sll_int32 :: i_step
102 
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 )
108  end do
109  else
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 )
115  end do
116  end if
117 
118  end subroutine lie_splitting_pic_vm_3d3v
119 
120 
122  subroutine lie_splitting_back_pic_vm_3d3v(self,dt, number_steps)
123  class(sll_t_time_propagator_pic_vm_3d3v_disgrade), intent(inout) :: self
124  sll_real64, intent(in) :: dt
125  sll_int32, intent(in) :: number_steps
126  !local variables
127  sll_int32 :: i_step
128 
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 )
134  end do
135  else
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 )
141  end do
142  end if
143 
144  end subroutine lie_splitting_back_pic_vm_3d3v
145 
146 
147  !---------------------------------------------------------------------------!
150  self, &
151  maxwell_solver, &
152  particle_mesh_coupling, &
153  particle_group, &
154  phi_dofs, &
155  efield_dofs, &
156  bfield_dofs, &
157  x_min, &
158  Lx, &
159  filter, &
160  boundary_particles, &
161  solver_tolerance, &
162  betar, &
163  force_sign, &
164  electrostatic, &
165  rhob, &
166  control_variate, &
167  jmean)
168  class(sll_t_time_propagator_pic_vm_3d3v_disgrade), intent( out ) :: self
169  class(sll_c_maxwell_3d_base), target, intent( in ) :: maxwell_solver
170  class(sll_c_particle_mesh_coupling_3d), target, intent(in) :: particle_mesh_coupling
171  class(sll_t_particle_array), target, intent( in ) :: particle_group
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)
177  class(sll_c_filter_base_3d), target :: filter
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(:)
184  class(sll_t_control_variates), optional, target, intent(in) :: control_variate
185  logical, optional, intent(in) :: jmean
186  !local variables
187  sll_real64 :: solver_tolerance_set
188  sll_real64 :: betar_set(2), force_sign_set
189  sll_int32 :: boundary_particles_set
190  logical :: jmean_set
191 
192  if( present(boundary_particles) )then
193  boundary_particles_set = boundary_particles
194  else
195  boundary_particles_set = 100
196  end if
197 
198  if (present(solver_tolerance) ) then
199  solver_tolerance_set = solver_tolerance
200  else
201  solver_tolerance_set = 1d-12
202  end if
203 
204  if (present(betar) ) then
205  betar_set = betar
206  else
207  betar_set = 1d0
208  end if
209 
210  if( present(force_sign) )then
211  force_sign_set = force_sign
212  else
213  force_sign_set = 1._f64
214  end if
215 
216  if( present(electrostatic) )then
217  self%electrostatic = electrostatic
218  end if
219 
220  if (present(jmean)) then
221  jmean_set = jmean
222  else
223  jmean_set = .false.
224  end if
225 
226  if( present( control_variate) ) then
227  call self%helper%init(maxwell_solver, &
228  particle_mesh_coupling, &
229  particle_group, &
230  phi_dofs, &
231  efield_dofs, &
232  bfield_dofs, &
233  x_min, &
234  lx, &
235  filter, &
236  boundary_particles=boundary_particles_set, &
237  solver_tolerance=solver_tolerance_set, &
238  betar=betar_set, &
239  force_sign=force_sign_set, &
240  control_variate = control_variate,&
241  jmean=jmean_set)
242  else
243  call self%helper%init(maxwell_solver, &
244  particle_mesh_coupling, &
245  particle_group, &
246  phi_dofs, &
247  efield_dofs, &
248  bfield_dofs, &
249  x_min, &
250  lx, &
251  filter, &
252  boundary_particles=boundary_particles_set, &
253  solver_tolerance=solver_tolerance_set, &
254  betar=betar_set, &
255  force_sign=force_sign_set, &
256  jmean=jmean_set)
257  end if
258 
259  end subroutine initialize_pic_vm_3d3v
260 
261 
264  self, &
265  maxwell_solver, &
266  particle_mesh_coupling, &
267  particle_group, &
268  phi_dofs, &
269  efield_dofs, &
270  bfield_dofs, &
271  x_min, &
272  Lx, &
273  filter, &
274  filename, &
275  boundary_particles, &
276  betar, &
277  force_sign, &
278  electrostatic, &
279  rhob, &
280  control_variate,&
281  jmean)
282  class(sll_t_time_propagator_pic_vm_3d3v_disgrade), intent( out ) :: self
283  class(sll_c_maxwell_3d_base), target, intent( in ) :: maxwell_solver
284  class(sll_c_particle_mesh_coupling_3d), target, intent(in) :: particle_mesh_coupling
285  class(sll_t_particle_array), target, intent( in ) :: particle_group
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)
291  class(sll_c_filter_base_3d), target :: filter
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(:)
298  class(sll_t_control_variates), optional, target, intent(in) :: control_variate
299  logical, optional, intent(in) :: jmean
300  !local variables
301  sll_real64 :: betar_set(2), force_sign_set
302  sll_int32 :: boundary_particles_set
303  logical :: jmean_set
304 
305  if( present(boundary_particles) )then
306  boundary_particles_set = boundary_particles
307  else
308  boundary_particles_set = 100
309  end if
310 
311  if( present(force_sign) )then
312  force_sign_set = force_sign
313  else
314  force_sign_set = 1._f64
315  end if
316 
317  if( present(electrostatic) )then
318  self%electrostatic = electrostatic
319  end if
320 
321  if( present(betar) )then
322  betar_set = betar
323  else
324  betar_set = 1d0
325  end if
326 
327  if (present(jmean)) then
328  jmean_set = jmean
329  else
330  jmean_set = .false.
331  end if
332 
333  if( present( control_variate) ) then
334  call self%helper%init_from_file( maxwell_solver, &
335  particle_mesh_coupling, &
336  particle_group, &
337  phi_dofs, &
338  efield_dofs, &
339  bfield_dofs, &
340  x_min, &
341  lx, &
342  filename, &
343  filter, &
344  boundary_particles = boundary_particles_set, &
345  betar=betar_set, &
346  force_sign=force_sign_set, &
347  rhob = rhob, &
348  control_variate = control_variate, &
349  jmean=jmean_set)
350  else
351  call self%helper%init_from_file( maxwell_solver, &
352  particle_mesh_coupling, &
353  particle_group, &
354  phi_dofs, &
355  efield_dofs, &
356  bfield_dofs, &
357  x_min, &
358  lx, &
359  filename, &
360  filter, &
361  boundary_particles = boundary_particles_set, &
362  betar=betar_set, &
363  force_sign=force_sign_set, &
364  rhob = rhob, &
365  jmean=jmean_set)
366  end if
367 
368  end subroutine initialize_file_pic_vm_3d3v
369 
370 
372  subroutine delete_pic_vm_3d3v(self)
373  class(sll_t_time_propagator_pic_vm_3d3v_disgrade), intent( inout ) :: self
374  call self%helper%free()
375 
376  end subroutine delete_pic_vm_3d3v
377 
378 
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 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.
    Report Typos and Errors