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
64 deallocate(self%particle_array)
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
86 self%n_particles = n_particles
87 self%n_total_particles = n_total_particles
89 sll_allocate(self%particle_array(6+n_weights, self%n_particles), ierr)
91 allocate(self%species, stat=ierr)
92 sll_assert( ierr == 0)
93 call self%species%init( charge, mass)
95 self%n_weights = n_weights
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
120 select type( particle_group )
122 call particle_group%init(n_particles, n_total_particles, charge, mass, n_weights)
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
148 select type( particle_group )
150 call particle_group%init(n_particles, n_total_particles, charge, mass, n_weights)
159 sll_int32 ,
intent( in ) :: i
162 r = self%particle_array(1:3, i)
170 sll_int32 ,
intent( in ) :: i
173 r = self%particle_array(4:6, i)
182 sll_int32 ,
intent( in ) :: i
183 sll_int32,
optional ,
intent( in ) :: i_weight
189 if(
present(i_weight)) i_wi = i_weight
190 r = self%species%q * self%particle_array(6+i_wi, i) * self%common_weight
199 sll_int32 ,
intent( in ) :: i
200 sll_int32,
optional ,
intent( in ) :: i_weight
206 if(
present(i_weight)) i_wi = i_weight
207 r = self%species%m * self%particle_array(6+i_wi, i) * self%common_weight
216 sll_int32 ,
intent( in ) :: i
217 sll_real64 :: r(self%n_weights)
219 r = self%particle_array(7:6+self%n_weights, i)
229 r = self%common_weight
237 sll_int32 ,
intent( in ) :: i
238 sll_real64 ,
intent( in) :: x(3)
240 self%particle_array(1:3, i) = x
248 sll_int32 ,
intent( in ) :: i
249 sll_real64 ,
intent( in) :: x(3)
251 self%particle_array(4:6, i) = x
259 sll_int32 ,
intent( in ) :: i
260 sll_real64 ,
intent( in) :: x(self%n_weights)
262 self%particle_array(7:6+self%n_weights, i) = x
270 sll_real64 ,
intent( in) :: x
272 self%common_weight = x
280 character(len=*),
intent(in) :: filename
284 open(newunit=file_id,file=filename//
'_weight')
285 write(file_id,*) self%common_weight
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)
301 character(len=*),
intent(in) :: filename
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')
309 read(file_id,*) self%common_weight
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)
316 read(file_id,*) self%common_weight
317 do i = 1, self%n_particles
318 read(file_id,*) self%particle_array(:,i)
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.