Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_particle_group_1d1v.F90
Go to the documentation of this file.
1 
6 
7  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_assert.h"
9 #include "sll_memory.h"
10 #include "sll_working_precision.h"
11 
12  use sll_m_particle_group_base, only: &
15 
16  implicit none
17 
18  public :: &
22 
23 
24  private
25  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
28  !sll_int32 :: n_particles !< number of particle
29  sll_real64, allocatable :: particle_array(:,:)
30  sll_real64 :: common_weight = 1.0_f64
31 
32  contains
33  ! Getters
34  procedure :: get_x => get_x_1d1v
35  procedure :: get_v => get_v_1d1v
36  procedure :: get_charge => get_charge_1d1v
37  procedure :: get_mass => get_mass_1d1v
38  procedure :: get_weights => get_weights_1d1v
39  procedure :: get_common_weight => get_common_weight_1d1v
40 
41  ! Setters
42  procedure :: set_x => set_x_1d1v
43  procedure :: set_v => set_v_1d1v
44  procedure :: set_weights => set_weight_1d1v
45  procedure :: set_common_weight => set_common_weight_1d1v
46 
47  ! Initializer
48  procedure :: init => initialize_particle_group_1d1v
49  procedure :: free => delete_particle_group_1d1v
50 
51  procedure :: print => print_particle_group_1d1v
52 
54 
55 contains
56 
57  !----------------------------------------------------------------------!
59  subroutine delete_particle_group_1d1v(self)
60  class( sll_t_particle_group_1d1v ), intent(inout) :: self
61 
62  deallocate(self%particle_array)
63 
64  end subroutine delete_particle_group_1d1v
65 
66  !----------------------------------------------------------------------!
69  self, &
70  n_particles, &
71  n_total_particles, &
72  charge, &
73  mass, &
74  n_weights)
75  class( sll_t_particle_group_1d1v ), intent( inout) :: self
76  sll_int32 , intent( in ) :: n_particles
77  sll_int32 , intent( in ) :: n_total_particles
78  sll_real64 , intent( in ) :: charge
79  sll_real64 , intent( in ) :: mass
80  sll_int32 , intent( in ) :: n_weights
81 
82  sll_int32 :: ierr
83 
84  self%n_particles = n_particles
85  self%n_total_particles = n_total_particles
86 
87  sll_allocate(self%particle_array(2+n_weights, n_particles), ierr)
88 
89  allocate(self%species, stat=ierr)
90  sll_assert( ierr == 0)
91  call self%species%init( charge, mass )
92 
93  self%n_weights = n_weights
94 
95  end subroutine initialize_particle_group_1d1v
96 
97 
98 
99  !----------------------------------------------------------------------!
102  particle_group, &
103  n_particles, &
104  n_total_particles, &
105  charge, &
106  mass, &
107  n_weights)
108  class( sll_c_particle_group_base ), pointer, intent( out ) :: particle_group
109  sll_int32 , intent( in ) :: n_particles
110  sll_int32 , intent( in ) :: n_total_particles
111  sll_real64 , intent( in ) :: charge
112  sll_real64 , intent( in ) :: mass
113  sll_int32 , intent( in ) :: n_weights
114 
115  sll_int32 :: ierr
116 
117  sll_allocate( sll_t_particle_group_1d1v :: particle_group, ierr)
118 
119  select type (particle_group)
120  type is (sll_t_particle_group_1d1v)
121  call particle_group%init( n_particles, n_total_particles, charge, mass, n_weights )
122  end select
123 
124  end subroutine sll_s_new_particle_group_1d1v_ptr
125 
126  !----------------------------------------------------------------------!
129  particle_group, &
130  n_particles, &
131  n_total_particles, &
132  charge, &
133  mass, &
134  n_weights)
135  class( sll_c_particle_group_base ), allocatable, intent( out ) :: particle_group
136  sll_int32 , intent( in ) :: n_particles
137  sll_int32 , intent( in ) :: n_total_particles
138  sll_real64 , intent( in ) :: charge
139  sll_real64 , intent( in ) :: mass
140  sll_int32 , intent( in ) :: n_weights
141 
142 
143  allocate( sll_t_particle_group_1d1v :: particle_group)
144 
145  select type (particle_group)
146  type is (sll_t_particle_group_1d1v)
147  call particle_group%init( n_particles, n_total_particles, charge, mass, n_weights )
148  end select
149 
150  end subroutine sll_s_new_particle_group_1d1v
151 
152 
153  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
155  pure function get_x_1d1v(self,i) result(r)
156  class( sll_t_particle_group_1d1v ), intent( in ) :: self
157  sll_int32 , intent( in ) :: i
158  sll_real64 :: r(3)
159 
160  r = 1.0_f64
161  r(1) = self%particle_array(1, i)
162 
163  end function get_x_1d1v
164 
165  !----------------------------------------------------------------------!
167  pure function get_v_1d1v( self, i ) result( r )
168  class( sll_t_particle_group_1d1v ), intent( in ) :: self
169  sll_int32 , intent( in ) :: i
170  sll_real64 :: r(3)
171 
172  r = 1.0_f64
173  r(1) = self%particle_array( 2, i )
174 
175  end function get_v_1d1v
176 
177  !----------------------------------------------------------------------!
179  pure function get_charge_1d1v( self, i , i_weight) result (r)
180  class( sll_t_particle_group_1d1v ), intent( in ) :: self
181  sll_int32 , intent( in ) :: i
182  sll_int32, optional , intent( in ) :: i_weight
183  sll_real64 :: r
184 
185  sll_int32 :: i_wi
186 
187  i_wi = 1
188  if(present(i_weight)) i_wi = i_weight
189  r = self%species%q * self%particle_array(2+i_wi, i) * self%common_weight
190 
191  end function get_charge_1d1v
192 
193  !----------------------------------------------------------------------!
195  pure function get_mass_1d1v( self, i, i_weight) result (r)
196  class( sll_t_particle_group_1d1v ), intent( in ) :: self
197  sll_int32 , intent( in ) :: i
198  sll_int32, optional , intent( in ) :: i_weight
199  sll_real64 :: r
200 
201  sll_int32 :: i_wi
202 
203  i_wi = 1
204  if(present(i_weight)) i_wi = i_weight
205  r = self%species%m * self%particle_array( 2+i_wi, i) * self%common_weight
206 
207  end function get_mass_1d1v
208 
209  !----------------------------------------------------------------------!
211  pure function get_weights_1d1v( self, i) result (r)
212  class( sll_t_particle_group_1d1v ), intent( in ) :: self
213  sll_int32 , intent( in ) :: i
214  sll_real64 :: r(self%n_weights)
215 
216  r = self%particle_array(3:2+self%n_weights, i)
217 
218  end function get_weights_1d1v
219 
220  !----------------------------------------------------------------------!
222  pure function get_common_weight_1d1v( self ) result(r)
223  class( sll_t_particle_group_1d1v ), intent( in ) :: self
224  sll_real64 :: r
225 
226  r = self%common_weight
227 
228  end function get_common_weight_1d1v
229 
230  !----------------------------------------------------------------------!
232  subroutine set_x_1d1v( self, i, x )
233  class( sll_t_particle_group_1d1v ), intent( inout ) :: self
234  sll_int32 , intent( in ) :: i
235  sll_real64 , intent( in) :: x(3)
236 
237  self%particle_array(1, i) = x(1)
238 
239  end subroutine set_x_1d1v
240 
241  !----------------------------------------------------------------------!
243  subroutine set_v_1d1v( self, i, x )
244  class( sll_t_particle_group_1d1v ), intent( inout ) :: self
245  sll_int32 , intent( in ) :: i
246  sll_real64 , intent( in) :: x(3)
247 
248  self%particle_array(2, i) = x(1)
249 
250  end subroutine set_v_1d1v
251 
252  !----------------------------------------------------------------------!
254  subroutine set_weight_1d1v( self, i, x )
255  class( sll_t_particle_group_1d1v ), intent( inout ) :: self
256  sll_int32 , intent( in ) :: i
257  sll_real64 , intent( in) :: x(self%n_weights)
258 
259  self%particle_array(3:2+self%n_weights, i) = x
260 
261  end subroutine set_weight_1d1v
262 
263 
264  !----------------------------------------------------------------------!
266  subroutine set_common_weight_1d1v( self, x )
267  class( sll_t_particle_group_1d1v ), intent( inout ) :: self
268  sll_real64 , intent( in):: x
269 
270  self%common_weight = x
271 
272  end subroutine set_common_weight_1d1v
273 
274 
275  !----------------------------------------------------------------------!
277  subroutine print_particle_group_1d1v(self, filename)
278  class( sll_t_particle_group_1d1v ), intent(in) :: self
279  character(len=*), intent(in) :: filename
280  sll_int32 :: file_id
281 
282  open(newunit=file_id,file=filename)
283  write(file_id,*) self%particle_array
284  close(file_id)
285 
286  end subroutine print_particle_group_1d1v
287 
288 
289 
290 end module sll_m_particle_group_1d1v
Simple particle group type for 1d1v.
pure real(kind=f64) function, dimension(3) get_x_1d1v(self, i)
Get position.
pure real(kind=f64) function, dimension(self%n_weights) get_weights_1d1v(self, i)
Get particle weights.
subroutine set_x_1d1v(self, i, x)
Set position of particle i.
pure real(kind=f64) function get_charge_1d1v(self, i, i_weight)
Get charge of particle (q * particle_weight)
subroutine, public sll_s_new_particle_group_1d1v_ptr(particle_group, n_particles, n_total_particles, charge, mass, n_weights)
Constructor for pointer.
pure real(kind=f64) function get_common_weight_1d1v(self)
Set the common weight.
pure real(kind=f64) function, dimension(3) get_v_1d1v(self, i)
Get velocity.
subroutine initialize_particle_group_1d1v(self, n_particles, n_total_particles, charge, mass, n_weights)
Initialization of the particle group.
subroutine print_particle_group_1d1v(self, filename)
Print particle array to file.
subroutine delete_particle_group_1d1v(self)
print particle array to file
subroutine set_v_1d1v(self, i, x)
Set velocity of particle i.
subroutine set_weight_1d1v(self, i, x)
Set weights of particle i.
subroutine set_common_weight_1d1v(self, x)
Set the common weight.
subroutine, public sll_s_new_particle_group_1d1v(particle_group, n_particles, n_total_particles, charge, mass, n_weights)
Constructor for allocatable.
pure real(kind=f64) function get_mass_1d1v(self, i, i_weight)
Get mass of particle (m * particle_weight)
Simple version of a PIC particle group in 1d1v.
    Report Typos and Errors