Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_particle_group_base.F90
Go to the documentation of this file.
2 
3  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
4 #include "sll_memory.h"
5 #include "sll_working_precision.h"
6 
7  implicit none
8 
9  public :: &
12  sll_t_species, &
17 
18  private
19  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
20 
21 
22 
23 
24  !============================================================================
25  ! Particle species
26  !============================================================================
27  type :: sll_t_species
28  character(len=64) :: name
29  sll_real64 :: q
30  sll_real64 :: m
31 
32  contains
33  procedure :: q_over_m
34  procedure :: init => initialize_species
35 
36  end type sll_t_species
37 
38 
39  !============================================================================
40  ! Particle group
41  !============================================================================
42  type, abstract :: sll_c_particle_group_base
43 
44  class( sll_t_species ), allocatable :: species
45  sll_int32 :: id
46  sll_int32 :: n_particles
47  sll_int32 :: n_total_particles
48  sll_int32 :: n_weights
49 
50 
51  contains
52  ! Getters ( for single particles )
53  procedure( i_get_coords ), deferred :: get_x
54  procedure( i_get_coords ), deferred :: get_v
55  procedure( i_get_scalar ), deferred :: get_charge
56  procedure( i_get_scalar ), deferred :: get_mass
57  procedure( i_get_array ), deferred :: get_weights
58  procedure( get_scalar ), deferred :: get_common_weight
59  procedure :: get_box
60  procedure :: get_boxnd
61  procedure :: get_xbox
62  procedure :: get_patch
63  procedure :: get_xpatch
64 
65  ! Getters ( for chunk of particles )
66  ! All chunk-getters are calling the single versions by default
67  ! We use pointers to avoid copying. Note that no setters are necessary since we use pointers to access the values.
68  ! For chunks, only one chget_x is provided this gives the coordinates in the format that the positions are stored.
69  procedure :: chget_x
70  procedure :: chget_v
71  procedure :: chget_weights
72  procedure :: chget_box
73  procedure :: chget_boxnd
74  procedure :: chget_patch
75 
76  ! Setters ( for single particles )
77  procedure( i_set_coords ), deferred :: set_x
78  procedure( i_set_coords ), deferred :: set_v
79  procedure( i_set_array ), deferred :: set_weights
80  procedure( set_scalar ), deferred :: set_common_weight
81 
82  procedure :: set_box
83  procedure :: set_boxnd
84  procedure :: set_patch
85 
86  procedure( empty ), deferred :: free
87 
88  ! ! Getters for the whole group
89  ! procedure( get_all_coords), deferred :: get_all_x
90  ! procedure( get_all_coords), deferred :: get_all_v
91 
92  procedure :: print
93  procedure :: read
94 
96 
97  !============================================================================
98  ! Particle array
99  !============================================================================
101 
102  sll_int32 :: n_species
103  class(sll_c_particle_group_base), allocatable :: group(:)
104 
105  end type sll_t_particle_array
106 
107  !----------------------------------------------------------------------------
108  abstract interface
109  pure function i_get_int( self, i ) result( r )
112  class( sll_c_particle_group_base ), intent( in ) :: self
113  sll_int32 , intent( in ) :: i
114  sll_int32 :: r
115  end function i_get_int
116  end interface
117 
118  !----------------------------------------------------------------------------
119  abstract interface
120  pure function i_get_intnd( self, i ) result( r )
123  class( sll_c_particle_group_base ), intent( in ) :: self
124  sll_int32 , intent( in ) :: i
125  sll_int32 :: r(3)
126  end function i_get_intnd
127  end interface
128 
129  !----------------------------------------------------------------------------
130  abstract interface
131  pure function i_get_scalar( self, i , i_weight) result( r )
134  class( sll_c_particle_group_base ), intent( in ) :: self
135  sll_int32 , intent( in ) :: i
136  sll_int32, optional , intent( in ) :: i_weight
137  sll_real64 :: r
138  end function i_get_scalar
139  end interface
140  !----------------------------------------------------------------------------
141  abstract interface
142  pure function get_scalar( self ) result(r)
145  class( sll_c_particle_group_base ), intent( in ) :: self
146  sll_real64 :: r
147  end function get_scalar
148  end interface
149 
150 
151  !----------------------------------------------------------------------------
152  abstract interface
153  pure function i_get_coords( self, i ) result( r )
156  class( sll_c_particle_group_base ), intent( in ) :: self
157  sll_int32 , intent( in ) :: i
158  sll_real64 :: r(3)
159  end function i_get_coords
160  end interface
161 
162  !----------------------------------------------------------------------------
163  abstract interface
164  pure function i_get_array( self, i ) result( r )
167  class( sll_c_particle_group_base ), intent( in ) :: self
168  sll_int32 , intent( in ) :: i
169  sll_real64 :: r(self%n_weights)
170  end function i_get_array
171  end interface
172 
173  !----------------------------------------------------------------------------
174  abstract interface
175  subroutine i_set_coords( self, i, x )
178  class( sll_c_particle_group_base ), intent( inout ) :: self
179  sll_int32 , intent( in ) :: i
180  sll_real64 , intent( in ) :: x(3)
181  end subroutine i_set_coords
182  end interface
183 
184  !----------------------------------------------------------------------------
185  abstract interface
186  subroutine i_set_scalar( self, i, x )
189  class( sll_c_particle_group_base ), intent( inout ) :: self
190  sll_int32 , intent( in ) :: i
191  sll_real64 , intent( in ) :: x
192  end subroutine i_set_scalar
193  end interface
194 
195  !----------------------------------------------------------------------------
196  abstract interface
197  subroutine set_scalar( self, x )
200  class( sll_c_particle_group_base ), intent( inout ) :: self
201  sll_real64 , intent( in ) :: x
202  end subroutine set_scalar
203  end interface
204 
205  !----------------------------------------------------------------------------
206  abstract interface
207  subroutine i_set_array( self, i, x )
210  class( sll_c_particle_group_base ), intent( inout ) :: self
211  sll_int32 , intent( in ) :: i
212  sll_real64 , intent( in ) :: x(self%n_weights)
213  end subroutine i_set_array
214  end interface
215 
216  !---------------------------------------------------------------------------!
217  abstract interface
218  subroutine empty(self)
220  class( sll_c_particle_group_base ), intent( inout ) :: self
221 
222  end subroutine empty
223  end interface
224 
226 contains
228 
229  !----------------------------------------------------------------------------
230  pure function q_over_m( self ) result( r )
231  class( sll_t_species ), intent( in ) :: self
232  sll_real64 :: r
233 
234  r = self%q / self%m
235  end function q_over_m
236 
237 
238  !----------------------------------------------------------------------------
239  subroutine initialize_species(&
240  self, &
241  species_charge, &
242  species_mass &
243  )
244  class(sll_t_species), intent ( out ) :: self
245  sll_real64, intent ( in ) :: species_charge
246  sll_real64, intent ( in ) :: species_mass
247 
248 
249  self%q = species_charge
250  self%m = species_mass
251 
252  end subroutine initialize_species
253 
254 
255 
256  !----------------------------------------------------------------------------!
257  ! Single particle gets and setters for boxes and patches
258  !----------------------------------------------------------------------------!
259 
260 
261  pure function get_box( self, i ) result( r )
262  class(sll_c_particle_group_base), intent ( in ) :: self
263  sll_int32, intent ( in ) :: i
264 
265  sll_int32 :: r
266 
267  r = 1
268 
269  end function get_box
270 
271 
272  pure function get_boxnd( self, i ) result( r )
273  class(sll_c_particle_group_base), intent ( in ) :: self
274  sll_int32, intent ( in ) :: i
275 
276  sll_int32 :: r(3)
277 
278  r = 1
279 
280  end function get_boxnd
281 
282 
283  pure function get_xbox( self, i ) result( r )
284  class(sll_c_particle_group_base), intent ( in ) :: self
285  sll_int32, intent ( in ) :: i
286 
287  sll_real64 :: r(3)
288 
289  r = self%get_x( i )
290 
291  end function get_xbox
292 
294  pure function get_patch( self, i ) result( r )
295  class(sll_c_particle_group_base), intent ( in ) :: self
296  sll_int32, intent ( in ) :: i
297 
298  sll_int32 :: r
299 
300  r = 1
301 
302  end function get_patch
303 
305  pure function get_xpatch( self, i ) result( r )
306  class(sll_c_particle_group_base), intent ( in ) :: self
307  sll_int32, intent ( in ) :: i
308 
309  sll_real64 :: r(3)
310 
311  r = self%get_x( i )
312 
313  end function get_xpatch
314 
315 
317  subroutine set_box( self, i , box)
318  class(sll_c_particle_group_base), intent ( in ) :: self
319  sll_int32, intent ( in ) :: i
320  sll_int32, intent ( out) :: box
321 
322 
323  end subroutine set_box
324 
326  subroutine set_boxnd( self, i , box)
327  class(sll_c_particle_group_base), intent ( in ) :: self
328  sll_int32, intent ( in ) :: i
329  sll_int32, intent ( out) :: box(:)
330 
331 
332  end subroutine set_boxnd
333 
334 
336  subroutine set_patch( self, i , box)
337  class(sll_c_particle_group_base), intent ( in ) :: self
338  sll_int32, intent ( in ) :: i
339  sll_int32, intent ( out) :: box
340 
341 
342  end subroutine set_patch
343 
344  !----------------------------------------------------------------------------!
345  ! Getters for chunks of particles
346  !----------------------------------------------------------------------------!
347 
348 
349  subroutine chget_x( self, chunk, x )
350  class(sll_c_particle_group_base), intent ( in ) :: self
351  sll_int32, intent ( in ) :: chunk(1:2)
352  sll_real64, pointer, intent( out ) :: x(:,:)
353 
354  sll_int32 :: i
355 
356  do i=chunk(1),chunk(2)
357  x(:,i-chunk(1)+1) = self%get_x( i )
358  end do
359 
360  end subroutine chget_x
361 
362 
363  subroutine chget_v( self, chunk, v )
364  class(sll_c_particle_group_base), intent ( in ) :: self
365  sll_int32, intent ( in ) :: chunk(1:2)
366  sll_real64, pointer, intent( out ) :: v(:,:)
367 
368  sll_int32 :: i
369 
370  do i=chunk(1),chunk(2)
371  v(:,i-chunk(1)+1) = self%get_v( i )
372  end do
373 
374  end subroutine chget_v
375 
376 
377  subroutine chget_weights( self, chunk, weights )
378  class(sll_c_particle_group_base), intent ( in ) :: self
379  sll_int32, intent ( in ) :: chunk(1:2)
380  sll_real64, pointer, intent( out ) :: weights(:,:)
381 
382  sll_int32 :: i
383 
384  do i=chunk(1),chunk(2)
385  weights(:,i-chunk(1)+1) = self%get_weights( i )
386  end do
387 
388  end subroutine chget_weights
389 
390 
391  subroutine chget_box( self, chunk, box )
392  class(sll_c_particle_group_base), intent ( in ) :: self
393  sll_int32, intent ( in ) :: chunk(1:2)
394  sll_int32, pointer, intent( out ) :: box(:)
395 
396  sll_int32 :: i
397 
398  do i=chunk(1),chunk(2)
399  box(i-chunk(1)+1) = self%get_box( i )
400  end do
401 
402  end subroutine chget_box
403 
404 
405  subroutine chget_boxnd( self, chunk, boxnd )
406  class(sll_c_particle_group_base), intent ( in ) :: self
407  sll_int32, intent ( in ) :: chunk(1:2)
408  sll_int32, pointer, intent( out ) :: boxnd(:,:)
409 
410  sll_int32 :: i
411 
412  do i=chunk(1),chunk(2)
413  boxnd(:,i-chunk(1)+1) = self%get_boxnd( i )
414  end do
415 
416  end subroutine chget_boxnd
417 
418 
419  subroutine chget_patch( self, chunk, patch )
420  class(sll_c_particle_group_base), intent ( in ) :: self
421  sll_int32, intent ( in ) :: chunk(1:2)
422  sll_int32, pointer, intent( out ) :: patch(:)
423 
424  sll_int32 :: i
425 
426  do i=chunk(1),chunk(2)
427  patch(i-chunk(1)+1) = self%get_patch( i )
428  end do
429 
430  end subroutine chget_patch
431 
433  subroutine print( self, filename )
434  class(sll_c_particle_group_base), intent ( in ) :: self
435  character(len=*), intent(in) :: filename
436 
437  end subroutine print
438 
440  subroutine read( self, filename )
441  class(sll_c_particle_group_base), intent ( inout ) :: self
442  character(len=*), intent(in) :: filename
443 
444  end subroutine read
445 
446  !----------------------------------------------------------------------------!
447  ! Index conversion helper functions
448  !----------------------------------------------------------------------------!
449 
451  pure function sll_f_index_1dto2d( num_pts, ind1d) result( ind2d )
452  sll_int32, intent( in ) :: num_pts(2)
453  sll_int32, intent( in ) :: ind1d
454  sll_int32 :: ind2d(2)
455 
456  ind2d(2) = ( ind1d/num_pts(1) )
457  ind2d(1) = ind1d - ind2d(2)*num_pts(1)
458  ind2d(2) = ind2d(2)+1
459 
460  end function sll_f_index_1dto2d
461 
463  pure function sll_f_index_2dto1d( num_pts, ind2d) result( ind1d )
464  sll_int32, intent( in ) :: num_pts(2)
465  sll_int32, intent( in ) :: ind2d(2)
466  sll_int32 :: ind1d
467 
468 
469  ind1d = ind2d(1)+ (ind2d(2)-1) * num_pts(1)
470 
471  end function sll_f_index_2dto1d
472 
473 
475  pure function sll_f_index_1dto3d( num_pts, ind1d) result( ind3d )
476  sll_int32, intent( in ) :: num_pts(3)
477  sll_int32, intent( in ) :: ind1d
478  sll_int32 :: ind3d(3)
479 
480  sll_int32 :: ind
481 
482  ind3d(3) = ( ind1d/(num_pts(1)*num_pts(2)) )
483  ind = ind1d - ind3d(3)*(num_pts(1)*num_pts(2))
484  ind3d(3) = ind3d(3) +1
485  ind3d(2) = ( ind/num_pts(1) )
486  ind3d(1) = ind - ind3d(2)*num_pts(1)
487  ind3d(2) = ind3d(2)+1
488 
489  end function sll_f_index_1dto3d
490 
492  pure function sll_f_index_3dto1d( num_pts, ind3d) result( ind1d )
493  sll_int32, intent( in ) :: num_pts(3)
494  sll_int32, intent( in ) :: ind3d(3)
495  sll_int32 :: ind1d
496 
497 
498  ind1d = ind3d(1) + (ind3d(2)-1) * num_pts(1) + (ind3d(3)-1)* num_pts(1)*num_pts(2)
499 
500  end function sll_f_index_3dto1d
501 
502 
503 end module sll_m_particle_group_base
subroutine set_boxnd(self, i, box)
Set the number of the box from nd index (by default nothing is done)
subroutine chget_x(self, chunk, x)
subroutine set_patch(self, i, box)
Set the number of the patch (by default nothing is done)
subroutine read(self, filename)
Dummy read function ( can be overwritten to read for restart )
pure integer(kind=i32) function get_box(self, i)
pure integer(kind=i32) function, dimension(2), public sll_f_index_1dto2d(num_pts, ind1d)
Helper function to compute the 2d index for a tensor product grid from the 1d index.
subroutine print(self, filename)
Dummy print function ( can be overwritten to print for debugging )
pure real(kind=f64) function, dimension(3) get_xpatch(self, i)
Get the coordinate in the patch-local system. Default: get_x.
pure real(kind=f64) function, dimension(3) get_xbox(self, i)
pure integer(kind=i32) function, dimension(3) get_boxnd(self, i)
subroutine chget_box(self, chunk, box)
subroutine initialize_species(self, species_charge, species_mass)
pure integer(kind=i32) function, public sll_f_index_3dto1d(num_pts, ind3d)
Helper function to compute the 1d index for a tensor product grid from the 2d index.
subroutine chget_v(self, chunk, v)
pure real(kind=f64) function q_over_m(self)
pure integer(kind=i32) function get_patch(self, i)
Get the number of the patch (per default we only have one patch and therefore always return a 1)
subroutine chget_boxnd(self, chunk, boxnd)
subroutine set_box(self, i, box)
Set the number of the box (by default nothing is done)
subroutine chget_weights(self, chunk, weights)
pure integer(kind=i32) function, dimension(3), public sll_f_index_1dto3d(num_pts, ind1d)
Helper function to compute the 2d index for a tensor product grid from the 1d index.
pure integer(kind=i32) function, public sll_f_index_2dto1d(num_pts, ind2d)
Helper function to compute the 1d index for a tensor product grid from the 2d index.
subroutine chget_patch(self, chunk, patch)
Module to select the kind parameter.
    Report Typos and Errors