Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_particle_group_3d3v.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_3d3v
36  procedure :: get_v => get_v_3d3v
37  procedure :: get_charge => get_charge_3d3v
38  procedure :: get_mass => get_mass_3d3v
39  procedure :: get_weights => get_weights_3d3v
40  procedure :: get_common_weight => get_common_weight_3d3v
41 
42  ! Setters
43  procedure :: set_x => set_x_3d3v
44  procedure :: set_v => set_v_3d3v
45  procedure :: set_weights => set_weights_3d3v
46  procedure :: set_common_weight => set_common_weight_3d3v
47 
48  ! Initializer
49  procedure :: init => initialize_particle_group_3d3v
50  procedure :: free => delete_particle_group_3d3v
51 
52  procedure :: print => print_particle_group_3d3v
53  procedure :: read => read_particle_group_3d3v
54 
56 
57 contains
58 
59  !----------------------------------------------------------------------!
61  subroutine delete_particle_group_3d3v(self)
62  class( sll_t_particle_group_3d3v ), intent( inout ) :: self
63 
64  deallocate(self%particle_array)
65 
66  end subroutine delete_particle_group_3d3v
67 
68  !----------------------------------------------------------------------!
71  self, &
72  n_particles, &
73  n_total_particles, &
74  charge, &
75  mass, &
76  n_weights)
77  class( sll_t_particle_group_3d3v ), intent( inout ) :: self
78  sll_int32 , intent( in ) :: n_particles
79  sll_int32 , intent( in ) :: n_total_particles
80  sll_real64 , intent( in ) :: charge
81  sll_real64 , intent( in ) :: mass
82  sll_int32 , intent( in ) :: n_weights
83 
84  sll_int32 :: ierr
85 
86  self%n_particles = n_particles
87  self%n_total_particles = n_total_particles
88 
89  sll_allocate(self%particle_array(6+n_weights, self%n_particles), ierr)
90 
91  allocate(self%species, stat=ierr)
92  sll_assert( ierr == 0)
93  call self%species%init( charge, mass)
94 
95  self%n_weights = n_weights
96 
97  end subroutine initialize_particle_group_3d3v
98 
99 
100  !----------------------------------------------------------------------!
103  particle_group, &
104  n_particles, &
105  n_total_particles, &
106  charge, &
107  mass, &
108  n_weights)
109  class( sll_c_particle_group_base ), pointer, intent( out ) :: particle_group
110  sll_int32 , intent( in ) :: n_particles
111  sll_int32 , intent( in ) :: n_total_particles
112  sll_real64 , intent( in ) :: charge
113  sll_real64 , intent( in ) :: mass
114  sll_int32 , intent(in) :: n_weights
115 
116  sll_int32 :: ierr
117 
118  sll_allocate( sll_t_particle_group_3d3v :: particle_group, ierr)
119 
120  select type( particle_group )
121  type is ( sll_t_particle_group_3d3v )
122  call particle_group%init(n_particles, n_total_particles, charge, mass, n_weights)
123  end select
124 
125  end subroutine sll_s_new_particle_group_3d3v_ptr
126 
127 
128  !----------------------------------------------------------------------!
131  particle_group, &
132  n_particles, &
133  n_total_particles, &
134  charge, &
135  mass, &
136  n_weights)
137  class( sll_c_particle_group_base ), allocatable, intent( out ) :: particle_group
138  sll_int32 , intent( in ) :: n_particles
139  sll_int32 , intent( in ) :: n_total_particles
140  sll_real64 , intent( in ) :: charge
141  sll_real64 , intent( in ) :: mass
142  sll_int32 , intent(in) :: n_weights
143 
144  sll_int32 :: ierr
145 
146  sll_allocate( sll_t_particle_group_3d3v :: particle_group, ierr)
147 
148  select type( particle_group )
149  type is ( sll_t_particle_group_3d3v )
150  call particle_group%init(n_particles, n_total_particles, charge, mass, n_weights)
151  end select
152 
153  end subroutine sll_s_new_particle_group_3d3v
154 
155  !----------------------------------------------------------------------!
157  pure function get_x_3d3v( self, i ) result( r )
158  class( sll_t_particle_group_3d3v ), intent( in ) :: self
159  sll_int32 , intent( in ) :: i
160  sll_real64 :: r(3)
161 
162  r = self%particle_array(1:3, i)
163 
164  end function get_x_3d3v
165 
166  !----------------------------------------------------------------------!
168  pure function get_v_3d3v( self, i ) result( r )
169  class( sll_t_particle_group_3d3v ), intent( in ) :: self
170  sll_int32 , intent( in ) :: i
171  sll_real64 :: r(3)
172 
173  r = self%particle_array(4:6, i)
174 
175  end function get_v_3d3v
176 
177 
178  !----------------------------------------------------------------------!
180  pure function get_charge_3d3v( self, i , i_weight) result (r)
181  class( sll_t_particle_group_3d3v ), intent( in ) :: self
182  sll_int32 , intent( in ) :: i
183  sll_int32, optional , intent( in ) :: i_weight
184  sll_real64 :: r
185 
186  sll_int32 :: i_wi
187 
188  i_wi = 1
189  if(present(i_weight)) i_wi = i_weight
190  r = self%species%q * self%particle_array(6+i_wi, i) * self%common_weight
191 
192  end function get_charge_3d3v
193 
194 
195  !----------------------------------------------------------------------!
197  pure function get_mass_3d3v( self, i, i_weight) result (r)
198  class( sll_t_particle_group_3d3v ), intent( in ) :: self
199  sll_int32 , intent( in ) :: i
200  sll_int32, optional , intent( in ) :: i_weight
201  sll_real64 :: r
202 
203  sll_int32 :: i_wi
204 
205  i_wi = 1
206  if(present(i_weight)) i_wi = i_weight
207  r = self%species%m * self%particle_array(6+i_wi, i) * self%common_weight
208 
209  end function get_mass_3d3v
210 
211 
212  !----------------------------------------------------------------------!
214  pure function get_weights_3d3v( self, i) result (r)
215  class( sll_t_particle_group_3d3v ), intent( in ) :: self
216  sll_int32 , intent( in ) :: i
217  sll_real64 :: r(self%n_weights)
218 
219  r = self%particle_array(7:6+self%n_weights, i)
220 
221  end function get_weights_3d3v
222 
223  !----------------------------------------------------------------------!
225  pure function get_common_weight_3d3v( self )result(r)
226  class( sll_t_particle_group_3d3v ), intent( in ) :: self
227  sll_real64 :: r
228 
229  r = self%common_weight
230 
231  end function get_common_weight_3d3v
232 
233  !----------------------------------------------------------------------!
235  subroutine set_x_3d3v( self, i, x )
236  class( sll_t_particle_group_3d3v ), intent( inout ) :: self
237  sll_int32 , intent( in ) :: i
238  sll_real64 , intent( in) :: x(3)
239 
240  self%particle_array(1:3, i) = x
241 
242  end subroutine set_x_3d3v
243 
244  !----------------------------------------------------------------------!
246  subroutine set_v_3d3v( self, i, x )
247  class( sll_t_particle_group_3d3v ), intent( inout ) :: self
248  sll_int32 , intent( in ) :: i
249  sll_real64 , intent( in) :: x(3)
250 
251  self%particle_array(4:6, i) = x
252 
253  end subroutine set_v_3d3v
254 
255  !----------------------------------------------------------------------!
257  subroutine set_weights_3d3v( self, i, x )
258  class( sll_t_particle_group_3d3v ), intent( inout ) :: self
259  sll_int32 , intent( in ) :: i
260  sll_real64 , intent( in) :: x(self%n_weights)
261 
262  self%particle_array(7:6+self%n_weights, i) = x
263 
264  end subroutine set_weights_3d3v
265 
266  !----------------------------------------------------------------------!
268  subroutine set_common_weight_3d3v( self, x )
269  class( sll_t_particle_group_3d3v ), intent( inout ) :: self
270  sll_real64 , intent( in) :: x
271 
272  self%common_weight = x
273 
274  end subroutine set_common_weight_3d3v
275 
276  !----------------------------------------------------------------------!
278  subroutine print_particle_group_3d3v(self, filename)
279  class( sll_t_particle_group_3d3v ), intent(in) :: self
280  character(len=*), intent(in) :: filename
281  sll_int32 :: file_id
282  sll_int32 :: i
283 
284  open(newunit=file_id,file=filename//'_weight')
285  write(file_id,*) self%common_weight
286  close(file_id)
287 
288  open(newunit=file_id,file=filename)
289  write(file_id,*) self%common_weight
290  do i=1, self%n_particles
291  write(file_id,*) self%particle_array(:,i)
292  end do
293  close(file_id)
294 
295  end subroutine print_particle_group_3d3v
296 
297  !----------------------------------------------------------------------!
299  subroutine read_particle_group_3d3v(self, filename)
300  class( sll_t_particle_group_3d3v ), intent(inout) :: self
301  character(len=*), intent(in) :: filename
302  sll_int32 :: file_id
303  sll_int32 :: i, ierr
304 
305  open(newunit=file_id,file=filename//'_weight', status='old', iostat=ierr)
306  if( ierr .ne. 0 ) then
307  sll_error('read_particle_group_3d3v', 'File not found: '//filename//'_weight')
308  end if
309  read(file_id,*) self%common_weight
310  close(file_id)
311 
312  open(newunit=file_id,file=filename, status='old', iostat=ierr)
313  if( ierr .ne. 0 ) then
314  sll_error('read_particle_group_3d3v', 'File not found: '//filename)
315  end if
316  read(file_id,*) self%common_weight
317  do i = 1, self%n_particles
318  read(file_id,*) self%particle_array(:,i)
319  !print*, 'a', i, self%n_particles
320  end do
321  close(file_id)
322 
323  end subroutine read_particle_group_3d3v
324 
325 end module sll_m_particle_group_3d3v
Simple particle group group for 3d3v.
subroutine set_weights_3d3v(self, i, x)
Set weights of particle i.
subroutine print_particle_group_3d3v(self, filename)
Print particle array.
pure real(kind=f64) function get_common_weight_3d3v(self)
Set the common weight.
subroutine set_v_3d3v(self, i, x)
Set velocities of particle i.
subroutine set_common_weight_3d3v(self, x)
Set the common weight.
subroutine set_x_3d3v(self, i, x)
Set positions of particle i.
subroutine, public sll_s_new_particle_group_3d3v_ptr(particle_group, n_particles, n_total_particles, charge, mass, n_weights)
Constructor for abstract type.
pure real(kind=f64) function, dimension(3) get_x_3d3v(self, i)
Get positions of particle i.
subroutine initialize_particle_group_3d3v(self, n_particles, n_total_particles, charge, mass, n_weights)
Initialize particle group.
subroutine delete_particle_group_3d3v(self)
Destructor.
subroutine, public sll_s_new_particle_group_3d3v(particle_group, n_particles, n_total_particles, charge, mass, n_weights)
Constructor for abstract type.
pure real(kind=f64) function, dimension(3) get_v_3d3v(self, i)
Get velocities of particle i.
pure real(kind=f64) function get_mass_3d3v(self, i, i_weight)
Get mass of particle i ( m * particle_weight)
subroutine read_particle_group_3d3v(self, filename)
Read particle array from file.
pure real(kind=f64) function, dimension(self%n_weights) get_weights_3d3v(self, i)
Get weights of particle i.
pure real(kind=f64) function get_charge_3d3v(self, i, i_weight)
Get charge of particle i ( q * particle_weight)
Simple version of a PIC particle group in 2d2v.
    Report Typos and Errors