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_1d2v_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_filter_base_1d, only: &
17 
18  use sll_m_maxwell_1d_base, only: &
20 
21  use sll_m_particle_group_base, only: &
23 
26 
27  use sll_m_time_propagator_base, only: &
29 
32 
33 
34  implicit none
35 
36  public :: &
38 
39  private
40  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
41 
45  logical :: electrostatic = .false.
46 
47  contains
48  procedure :: lie_splitting => lie_splitting_pic_vm_1d2v_disgrade
49  procedure :: lie_splitting_back => lie_splitting_back_pic_vm_1d2v_disgrade
50  procedure :: strang_splitting => strang_splitting_pic_vm_1d2v_disgrade
51 
52  procedure :: reinit_fields
53 
54  procedure :: init => initialize_pic_vm_1d2v_disgrade
55  procedure :: init_from_file => initialize_file_pic_vm_1d2v_disgrade
56  procedure :: free => delete_pic_vm_1d2v_disgrade
57 
59 
60 contains
61  subroutine reinit_fields( self )
62  class(sll_t_time_propagator_pic_vm_1d2v_disgrade), intent(inout) :: self
63 
64  call self%helper%reinit_fields()
65 
66  end subroutine reinit_fields
67 
68 
70  subroutine strang_splitting_pic_vm_1d2v_disgrade(self,dt, number_steps)
71  class(sll_t_time_propagator_pic_vm_1d2v_disgrade), intent(inout) :: self
72  sll_real64, intent(in) :: dt
73  sll_int32, intent(in) :: number_steps
74  !local variables
75  sll_int32 :: i_step
76 
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 )
84  end do
85  else
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 )
94  end do
95  end if
96 
98 
99 
101  subroutine lie_splitting_pic_vm_1d2v_disgrade(self,dt, number_steps)
102  class(sll_t_time_propagator_pic_vm_1d2v_disgrade), intent(inout) :: self
103  sll_real64, intent(in) :: dt
104  sll_int32, intent(in) :: number_steps
105  !local variables
106  sll_int32 :: i_step
107 
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 )
113  end do
114  else
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 )
120  end do
121  end if
122 
124 
125 
127  subroutine lie_splitting_back_pic_vm_1d2v_disgrade(self,dt, number_steps)
128  class(sll_t_time_propagator_pic_vm_1d2v_disgrade), intent(inout) :: self
129  sll_real64, intent(in) :: dt
130  sll_int32, intent(in) :: number_steps
131  !local variables
132  sll_int32 :: i_step
133 
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 )
139  end do
140  else
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 )
146  end do
147  end if
148 
150 
151 
152  !---------------------------------------------------------------------------!
155  self, &
156  maxwell_solver, &
157  kernel_smoother_0, &
158  kernel_smoother_1, &
159  particle_group, &
160  phi_dofs, &
161  efield_dofs, &
162  bfield_dofs, &
163  x_min, &
164  Lx, &
165  filter, &
166  boundary_particles, &
167  solver_tolerance, &
168  force_sign, &
169  control_variate, &
170  i_weight, &
171  betar, &
172  electrostatic, &
173  jmean)
174  class(sll_t_time_propagator_pic_vm_1d2v_disgrade), intent(out) :: self
175  class(sll_c_maxwell_1d_base), target, intent(in) :: maxwell_solver
176  class(sll_c_particle_mesh_coupling_1d), target, intent(in) :: kernel_smoother_0
177  class(sll_c_particle_mesh_coupling_1d), target, intent(in) :: kernel_smoother_1
178  class(sll_t_particle_array), target, intent(in) :: particle_group
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
184  class( sll_c_filter_base_1d ), intent( in ), target :: filter
185  sll_int32, optional, intent( in ) :: boundary_particles
186  sll_real64, intent(in), optional :: solver_tolerance
187  sll_real64, optional, intent(in) :: force_sign
188  class(sll_t_control_variates), optional, target, intent(in) :: control_variate
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
193  !local variables
194  sll_int32 :: boundary_particles_set
195  sll_real64 :: solver_tolerance_set, force_sign_set, betar_set(2)
196  logical :: jmean_set
197 
198  if (present(boundary_particles)) then
199  boundary_particles_set = boundary_particles
200  else
201  boundary_particles_set = 100
202  end if
203 
204  if (present(solver_tolerance) ) then
205  solver_tolerance_set = solver_tolerance
206  else
207  solver_tolerance_set = 1d-12
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(betar) )then
217  betar_set = betar
218  else
219  betar_set = 1._f64
220  end if
221 
222  if( present(electrostatic) )then
223  self%electrostatic = electrostatic
224  end if
225 
226  if (present(jmean)) then
227  jmean_set = jmean
228  else
229  jmean_set = .false.
230  end if
231 
232  if( present( control_variate) ) then
233  call self%helper%init(maxwell_solver, &
234  kernel_smoother_0, &
235  kernel_smoother_1, &
236  particle_group, &
237  phi_dofs, &
238  efield_dofs, &
239  bfield_dofs, &
240  x_min, &
241  lx, &
242  filter, &
243  boundary_particles=boundary_particles_set, &
244  solver_tolerance=solver_tolerance_set, &
245  force_sign=force_sign_set, &
246  control_variate = control_variate,&
247  i_weight=i_weight, &
248  betar=betar_set, &
249  jmean=jmean_set)
250  else
251  call self%helper%init(maxwell_solver, &
252  kernel_smoother_0, &
253  kernel_smoother_1, &
254  particle_group, &
255  phi_dofs, &
256  efield_dofs, &
257  bfield_dofs, &
258  x_min, &
259  lx, &
260  filter, &
261  boundary_particles=boundary_particles_set, &
262  solver_tolerance=solver_tolerance_set, &
263  betar=betar_set, &
264  force_sign=force_sign_set, &
265  jmean=jmean_set)
266  end if
267 
268 
269  end subroutine initialize_pic_vm_1d2v_disgrade
270 
271 
274  self, &
275  maxwell_solver, &
276  kernel_smoother_0, &
277  kernel_smoother_1, &
278  particle_group, &
279  phi_dofs, &
280  efield_dofs, &
281  bfield_dofs, &
282  x_min, &
283  Lx, &
284  filter, &
285  filename, &
286  boundary_particles, &
287  force_sign, &
288  control_variate, &
289  i_weight, &
290  betar, &
291  electrostatic, &
292  jmean)
293  class(sll_t_time_propagator_pic_vm_1d2v_disgrade), intent(out) :: self
294  class(sll_c_maxwell_1d_base), target, intent(in) :: maxwell_solver
295  class(sll_c_particle_mesh_coupling_1d), target, intent(in) :: kernel_smoother_0
296  class(sll_c_particle_mesh_coupling_1d), target, intent(in) :: kernel_smoother_1
297  class(sll_t_particle_array), target, intent(in) :: particle_group
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
303  class( sll_c_filter_base_1d ), intent( in ), target :: filter
304  character(len=*), intent(in) :: filename
305  sll_int32, optional, intent( in ) :: boundary_particles
306  sll_real64, optional, intent(in) :: force_sign
307  class(sll_t_control_variates), optional, target, intent(in) :: control_variate
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
312  !local variables
313  sll_int32 :: boundary_particles_set
314  sll_real64 :: force_sign_set, betar_set(2)
315  logical :: jmean_set
316 
317 
318  if (present(boundary_particles)) then
319  boundary_particles_set = boundary_particles
320  else
321  boundary_particles_set = 100
322  end if
323 
324  if( present(force_sign) )then
325  force_sign_set = force_sign
326  else
327  force_sign_set = 1._f64
328  end if
329 
330  if( present(betar) )then
331  betar_set = betar
332  else
333  betar_set = 1._f64
334  end if
335 
336  if( present(electrostatic) )then
337  self%electrostatic = electrostatic
338  else
339  self%electrostatic = .false.
340  end if
341 
342  if (present(jmean)) then
343  jmean_set = jmean
344  else
345  jmean_set = .false.
346  end if
347 
348  if( present( control_variate) ) then
349  call self%helper%init_from_file( maxwell_solver, &
350  kernel_smoother_0, &
351  kernel_smoother_1, &
352  particle_group, &
353  phi_dofs, &
354  efield_dofs, &
355  bfield_dofs, &
356  x_min, &
357  lx, &
358  filter, &
359  filename, &
360  boundary_particles = boundary_particles_set, &
361  force_sign=force_sign_set, &
362  control_variate = control_variate, &
363  i_weight=i_weight, &
364  betar=betar_set, &
365  jmean=jmean_set)
366  else
367  call self%helper%init_from_file( maxwell_solver, &
368  kernel_smoother_0, &
369  kernel_smoother_1, &
370  particle_group, &
371  phi_dofs, &
372  efield_dofs, &
373  bfield_dofs, &
374  x_min, &
375  lx, &
376  filter, &
377  filename, &
378  boundary_particles = boundary_particles_set, &
379  betar=betar_set, &
380  force_sign=force_sign_set, &
381  jmean=jmean_set)
382  end if
383 
384 
386 
387 
388  !---------------------------------------------------------------------------!
391  class(sll_t_time_propagator_pic_vm_1d2v_disgrade), intent( inout ) :: self
392  call self%helper%free()
393 
394  end subroutine delete_pic_vm_1d2v_disgrade
395 
396 
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 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 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.
    Report Typos and Errors