Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_particle_group_1d2v.F90
Go to the documentation of this file.
1 
6 
7  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_assert.h"
9 #include "sll_errors.h"
10 #include "sll_memory.h"
11 #include "sll_working_precision.h"
12 
13  use sll_m_particle_group_base, only: &
16 
17  implicit none
18 
19  public :: &
23 
24  private
25  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
26 
29  !sll_int32 :: n_particles !< number of particle
30  sll_real64, allocatable :: particle_array(:,:)
31  sll_real64 :: common_weight = 1.0_f64
32 
33  contains
34  ! Getters
35  procedure :: get_x => get_x_1d2v
36  procedure :: get_v => get_v_1d2v
37  procedure :: get_charge => get_charge_1d2v
38  procedure :: get_mass => get_mass_1d2v
39  procedure :: get_weights => get_weights_1d2v
40  procedure :: get_common_weight => get_common_weight_1d2v
41 
42  ! Setters
43  procedure :: set_x => set_x_1d2v
44  procedure :: set_v => set_v_1d2v
45  procedure :: set_weights => set_weight_1d2v
46  procedure :: set_common_weight => set_common_weight_1d2v
47 
48  ! Initializer
49  procedure :: init => initialize_particle_group_1d2v
50  procedure :: free => delete_particle_group_1d2v
51 
52  procedure :: print => print_particle_group_1d2v
53  procedure :: read => read_particle_group_1d2v
54 
56 
57 
58 contains
59 
60  !----------------------------------------------------------------------!
62  subroutine delete_particle_group_1d2v(self)
63  class( sll_t_particle_group_1d2v ), intent(inout) :: self
64 
65  deallocate(self%particle_array)
66 
67  end subroutine delete_particle_group_1d2v
68 
69  !----------------------------------------------------------------------!
72  self, &
73  n_particles, &
74  n_total_particles, &
75  charge, &
76  mass, &
77  n_weights)
78  class( sll_t_particle_group_1d2v ), intent( inout) :: self
79  sll_int32 , intent( in ) :: n_particles
80  sll_int32 , intent( in ) :: n_total_particles
81  sll_real64 , intent( in ) :: charge
82  sll_real64 , intent( in ) :: mass
83  sll_int32 , intent( in ) :: n_weights
84 
85  sll_int32 :: ierr
86 
87  self%n_particles = n_particles
88  self%n_total_particles = n_total_particles
89 
90  sll_allocate(self%particle_array(3+n_weights, n_particles), ierr)
91 
92  allocate(self%species, stat=ierr)
93  sll_assert( ierr == 0)
94  call self%species%init( charge, mass )
95 
96  self%n_weights = n_weights
97 
98  end subroutine initialize_particle_group_1d2v
99 
100 
101 
102  !----------------------------------------------------------------------!
105  particle_group, &
106  n_particles, &
107  n_total_particles, &
108  charge, &
109  mass, &
110  n_weights)
111  class( sll_c_particle_group_base ), pointer, intent( out ) :: particle_group
112  sll_int32 , intent( in ) :: n_particles
113  sll_int32 , intent( in ) :: n_total_particles
114  sll_real64 , intent( in ) :: charge
115  sll_real64 , intent( in ) :: mass
116  sll_int32 , intent( in ) :: n_weights
117 
118  sll_int32 :: ierr
119 
120  sll_allocate( sll_t_particle_group_1d2v :: particle_group, ierr)
121 
122  select type (particle_group)
123  type is (sll_t_particle_group_1d2v)
124  call particle_group%init( n_particles, n_total_particles, charge, mass, n_weights )
125  end select
126 
127  end subroutine sll_s_new_particle_group_1d2v_ptr
128 
129 
130  !----------------------------------------------------------------------!
133  particle_group, &
134  n_particles, &
135  n_total_particles, &
136  charge, &
137  mass, &
138  n_weights)
139  class( sll_c_particle_group_base ), allocatable, intent( out ) :: particle_group
140  sll_int32 , intent( in ) :: n_particles
141  sll_int32 , intent( in ) :: n_total_particles
142  sll_real64 , intent( in ) :: charge
143  sll_real64 , intent( in ) :: mass
144  sll_int32 , intent( in ) :: n_weights
145 
146 
147  allocate( sll_t_particle_group_1d2v :: particle_group)
148 
149  select type (particle_group)
150  type is (sll_t_particle_group_1d2v)
151  call particle_group%init( n_particles, n_total_particles, charge, mass, n_weights )
152  end select
153 
154  end subroutine sll_s_new_particle_group_1d2v
155 
156 
157 
158  !----------------------------------------------------------------------!
160  pure function get_x_1d2v( self, i ) result( r )
161  class( sll_t_particle_group_1d2v ), intent( in ) :: self
162  sll_int32 , intent( in ) :: i
163  sll_real64 :: r(3)
164 
165  r = 1.0_f64
166  r(1) = self%particle_array(1, i)
167 
168  end function get_x_1d2v
169 
170  !----------------------------------------------------------------------!
172  pure function get_v_1d2v( self, i ) result( r )
173  class( sll_t_particle_group_1d2v ), intent( in ) :: self
174  sll_int32 , intent( in ) :: i
175  sll_real64 :: r(3)
176 
177  r = 1.0_f64
178  r(1:2) = self%particle_array( 2:3, i )
179 
180  end function get_v_1d2v
181 
182 
183  !----------------------------------------------------------------------!
185  pure function get_charge_1d2v( self, i , i_weight) result (r)
186  class( sll_t_particle_group_1d2v ), intent( in ) :: self
187  sll_int32 , intent( in ) :: i
188  sll_int32, optional , intent( in ) :: i_weight
189  sll_real64 :: r
190 
191  sll_int32 :: i_wi
192 
193  i_wi = 1
194  if(present(i_weight)) i_wi = i_weight
195  r = self%species%q * self%particle_array(3+i_wi, i) * self%common_weight
196 
197  end function get_charge_1d2v
198 
199 
200  !----------------------------------------------------------------------!
202  pure function get_mass_1d2v( self, i, i_weight) result (r)
203  class( sll_t_particle_group_1d2v ), intent( in ) :: self
204  sll_int32 , intent( in ) :: i
205  sll_int32, optional , intent( in ) :: i_weight
206  sll_real64 :: r
207 
208  sll_int32 :: i_wi
209 
210  i_wi = 1
211  if(present(i_weight)) i_wi = i_weight
212  r = self%species%m * self%particle_array( 3+i_wi, i) * self%common_weight
213 
214  end function get_mass_1d2v
215 
216  !----------------------------------------------------------------------!
218  pure function get_weights_1d2v( self, i) result (r)
219  class( sll_t_particle_group_1d2v ), intent( in ) :: self
220  sll_int32 , intent( in ) :: i
221  sll_real64 :: r(self%n_weights)
222 
223  r = self%particle_array(4:3+self%n_weights, i)
224 
225  end function get_weights_1d2v
226 
227  !----------------------------------------------------------------------!
229  pure function get_common_weight_1d2v( self )result(r)
230  class( sll_t_particle_group_1d2v ), intent( in ) :: self
231  sll_real64 :: r
232 
233  r = self%common_weight
234 
235  end function get_common_weight_1d2v
236 
237  !----------------------------------------------------------------------!
239  subroutine set_x_1d2v( self, i, x )
240  class( sll_t_particle_group_1d2v ), intent( inout ) :: self
241  sll_int32 , intent( in ) :: i
242  sll_real64 , intent( in) :: x(3)
243 
244  self%particle_array(1, i) = x(1)
245 
246  end subroutine set_x_1d2v
247 
248  !----------------------------------------------------------------------!
250  subroutine set_v_1d2v( self, i, x )
251  class( sll_t_particle_group_1d2v ), intent( inout ) :: self
252  sll_int32 , intent( in ) :: i
253  sll_real64 , intent( in) :: x(3)
254 
255  self%particle_array(2:3, i) = x(1:2)
256 
257  end subroutine set_v_1d2v
258 
259  !----------------------------------------------------------------------!
261  subroutine set_weight_1d2v( self, i, x )
262  class( sll_t_particle_group_1d2v ), intent( inout ) :: self
263  sll_int32 , intent( in ) :: i
264  sll_real64 , intent( in) :: x(self%n_weights)
265 
266  self%particle_array(4:3+self%n_weights, i) = x
267 
268  end subroutine set_weight_1d2v
269 
270 
271  !----------------------------------------------------------------------!
273  subroutine set_common_weight_1d2v( self, x )
274  class( sll_t_particle_group_1d2v ), intent( inout ) :: self
275  sll_real64 , intent( in):: x
276 
277  self%common_weight = x
278 
279  end subroutine set_common_weight_1d2v
280 
281 
282  !----------------------------------------------------------------------!
284  subroutine print_particle_group_1d2v(self, filename)
285  class( sll_t_particle_group_1d2v ), intent(in) :: self
286  character(len=*), intent(in) :: filename
287  sll_int32 :: file_id
288  sll_int32 :: i
289  character(len=256) :: fmt
290 
291  open(newunit=file_id,file=filename//'_weight')
292  if (self%n_weights == 1) then
293  fmt = '(2g24.16,2g24.16,2g24.16,2g24.16)'
294  elseif (self%n_weights == 3) then
295  fmt = '(2g24.16,2g24.16,2g24.16,2g24.16,2g24.16,2g24.16)'
296  else
297  sll_error('read_particle_group_1d2v', 'Unknown number of weights.' )
298  end if
299 
300  write(file_id,trim(fmt)) self%common_weight
301  close(file_id)
302 
303  open(newunit=file_id,file=filename)
304  if (self%n_weights == 1) then
305  fmt = '(2g24.16,2g24.16,2g24.16,2g24.16)'
306  elseif (self%n_weights == 3) then
307  fmt = '(2g24.16,2g24.16,2g24.16,2g24.16,2g24.16,2g24.16)'
308  else
309  sll_error('read_particle_group_1d2v', 'Unknown number of weights.' )
310  end if
311 
312  do i=1, self%n_particles
313  write(file_id,trim(fmt)) self%particle_array(:,i)
314  end do
315  close(file_id)
316 
317 
318 
319  end subroutine print_particle_group_1d2v
320 
321 
322  !----------------------------------------------------------------------!
324  subroutine read_particle_group_1d2v(self, filename)
325  class( sll_t_particle_group_1d2v ), intent(inout) :: self
326  character(len=*), intent(in) :: filename
327  sll_int32 :: file_id
328  sll_int32 :: i, ierr
329  character(len=256) :: fmt
330 
331 
332  open(newunit=file_id,file=filename//'_weight', status='old', iostat=ierr)
333  if( ierr .ne. 0 ) then
334  sll_error('read_particle_group_1d2v', 'File not found: '//filename//'_weight')
335  end if
336  if (self%n_weights == 1) then
337  fmt = '(2g24.16,2g24.16,2g24.16,2g24.16)'
338  elseif (self%n_weights == 3) then
339  fmt = '(2g24.16,2g24.16,2g24.16,2g24.16,2g24.16,2g24.16)'
340  else
341  sll_error('read_particle_group_1d2v', 'Unknown number of weights.' )
342  end if
343 
344  read(file_id,trim(fmt)) self%common_weight
345  close(file_id)
346 
347 
348  open(newunit=file_id,file=filename, status='old', iostat=ierr)
349  if( ierr .ne. 0 ) then
350  sll_error('read_particle_group_1d2v', 'File not found: '//filename)
351  end if
352  if (self%n_weights == 1) then
353  fmt = '(2g24.16,2g24.16,2g24.16,2g24.16)'
354  elseif (self%n_weights == 3) then
355  fmt = '(2g24.16,2g24.16,2g24.16,2g24.16,2g24.16,2g24.16)'
356  else
357  sll_error('read_particle_group_1d2v', 'Unknown number of weights.' )
358  end if
359 
360 
361  do i = 1, self%n_particles
362  read(file_id,trim(fmt)) self%particle_array(:,i)
363  !print*, 'a', i, self%n_particles
364  end do
365  close(file_id)
366 
367 
368 
369  end subroutine read_particle_group_1d2v
370 
371 
372 
373 
374 end module sll_m_particle_group_1d2v
Simple particle group type for 1d2v.
subroutine set_v_1d2v(self, i, x)
Set velocity of particle i.
pure real(kind=f64) function get_mass_1d2v(self, i, i_weight)
Get mass of particle (m * particle_weight)
subroutine initialize_particle_group_1d2v(self, n_particles, n_total_particles, charge, mass, n_weights)
Initialization of the particle group.
pure real(kind=f64) function, dimension(3) get_v_1d2v(self, i)
Get velocities.
subroutine delete_particle_group_1d2v(self)
Destructor.
pure real(kind=f64) function get_common_weight_1d2v(self)
Set the common weight.
subroutine read_particle_group_1d2v(self, filename)
Read particle array from file.
pure real(kind=f64) function, dimension(3) get_x_1d2v(self, i)
Get position.
subroutine print_particle_group_1d2v(self, filename)
Print particle array to file.
subroutine set_x_1d2v(self, i, x)
Set position of particle i.
subroutine, public sll_s_new_particle_group_1d2v(particle_group, n_particles, n_total_particles, charge, mass, n_weights)
Constructor for allocatable.
subroutine set_weight_1d2v(self, i, x)
Set weights of particle i.
pure real(kind=f64) function, dimension(self%n_weights) get_weights_1d2v(self, i)
Get particle weights.
pure real(kind=f64) function get_charge_1d2v(self, i, i_weight)
Get charge of particle (q * particle_weight)
subroutine, public sll_s_new_particle_group_1d2v_ptr(particle_group, n_particles, n_total_particles, charge, mass, n_weights)
Constructor for pointer.
subroutine set_common_weight_1d2v(self, x)
Set the common weight.
Simple version of a PIC particle group in 1d2v.
    Report Typos and Errors