8 #include "sll_assert.h"
9 #include "sll_errors.h"
10 #include "sll_memory.h"
11 #include "sll_working_precision.h"
30 sll_real64,
allocatable :: particle_array(:,:)
31 sll_real64 :: common_weight = 1.0_f64
65 deallocate(self%particle_array)
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
87 self%n_particles = n_particles
88 self%n_total_particles = n_total_particles
90 sll_allocate(self%particle_array(3+n_weights, n_particles), ierr)
92 allocate(self%species, stat=ierr)
93 sll_assert( ierr == 0)
94 call self%species%init( charge, mass )
96 self%n_weights = n_weights
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
122 select type (particle_group)
124 call particle_group%init( n_particles, n_total_particles, charge, mass, n_weights )
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
149 select type (particle_group)
151 call particle_group%init( n_particles, n_total_particles, charge, mass, n_weights )
162 sll_int32 ,
intent( in ) :: i
166 r(1) = self%particle_array(1, i)
174 sll_int32 ,
intent( in ) :: i
178 r(1:2) = self%particle_array( 2:3, i )
187 sll_int32 ,
intent( in ) :: i
188 sll_int32,
optional ,
intent( in ) :: i_weight
194 if(
present(i_weight)) i_wi = i_weight
195 r = self%species%q * self%particle_array(3+i_wi, i) * self%common_weight
204 sll_int32 ,
intent( in ) :: i
205 sll_int32,
optional ,
intent( in ) :: i_weight
211 if(
present(i_weight)) i_wi = i_weight
212 r = self%species%m * self%particle_array( 3+i_wi, i) * self%common_weight
220 sll_int32 ,
intent( in ) :: i
221 sll_real64 :: r(self%n_weights)
223 r = self%particle_array(4:3+self%n_weights, i)
233 r = self%common_weight
241 sll_int32 ,
intent( in ) :: i
242 sll_real64 ,
intent( in) :: x(3)
244 self%particle_array(1, i) = x(1)
252 sll_int32 ,
intent( in ) :: i
253 sll_real64 ,
intent( in) :: x(3)
255 self%particle_array(2:3, i) = x(1:2)
263 sll_int32 ,
intent( in ) :: i
264 sll_real64 ,
intent( in) :: x(self%n_weights)
266 self%particle_array(4:3+self%n_weights, i) = x
275 sll_real64 ,
intent( in):: x
277 self%common_weight = x
286 character(len=*),
intent(in) :: filename
289 character(len=256) :: fmt
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)'
297 sll_error(
'read_particle_group_1d2v',
'Unknown number of weights.' )
300 write(file_id,trim(fmt)) self%common_weight
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)'
309 sll_error(
'read_particle_group_1d2v',
'Unknown number of weights.' )
312 do i=1, self%n_particles
313 write(file_id,trim(fmt)) self%particle_array(:,i)
326 character(len=*),
intent(in) :: filename
329 character(len=256) :: fmt
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')
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)'
341 sll_error(
'read_particle_group_1d2v',
'Unknown number of weights.' )
344 read(file_id,trim(fmt)) self%common_weight
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)
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)'
357 sll_error(
'read_particle_group_1d2v',
'Unknown number of weights.' )
361 do i = 1, self%n_particles
362 read(file_id,trim(fmt)) self%particle_array(:,i)
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.