Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_particle_group_2d2v.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  private
24  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
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_2d2v
35  procedure :: get_v => get_v_2d2v
36  procedure :: get_charge => get_charge_2d2v
37  procedure :: get_mass => get_mass_2d2v
38  procedure :: get_weights => get_weights_2d2v
39  procedure :: get_common_weight => get_common_weight_2d2v
40 
41  ! Setters
42  procedure :: set_x => set_x_2d2v
43  procedure :: set_v => set_v_2d2v
44  procedure :: set_weights => set_weights_2d2v
45  procedure :: set_common_weight => set_common_weight_2d2v
46 
47  ! Initializer
48  procedure :: init => initialize_particle_group_2d2v
49  procedure :: free => delete_particle_group_2d2v
50 
51 
52  procedure :: print => print_particle_group_2d2v
53 
55 
56 contains
57 
58  !----------------------------------------------------------------------!
60  subroutine delete_particle_group_2d2v(self)
61  class( sll_t_particle_group_2d2v ), intent( inout ) :: self
62 
63  deallocate(self%particle_array)
64 
65  end subroutine delete_particle_group_2d2v
66 
67  !----------------------------------------------------------------------!
70  self, &
71  n_particles, &
72  n_total_particles, &
73  charge, &
74  mass, &
75  n_weights)
76  class( sll_t_particle_group_2d2v ), intent( inout ) :: self
77  sll_int32 , intent( in ) :: n_particles
78  sll_int32 , intent( in ) :: n_total_particles
79  sll_real64 , intent( in ) :: charge
80  sll_real64 , intent( in ) :: mass
81  sll_int32 , intent( in ) :: n_weights
82 
83  sll_int32 :: ierr
84 
85  self%n_particles = n_particles
86  self%n_total_particles = n_total_particles
87 
88  sll_allocate(self%particle_array(4+n_weights, self%n_particles), ierr)
89 
90  allocate(self%species, stat=ierr)
91  sll_assert( ierr == 0)
92  call self%species%init( charge, mass)
93 
94  self%n_weights = n_weights
95 
96  end subroutine initialize_particle_group_2d2v
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_2d2v :: particle_group, ierr)
118 
119  select type( particle_group )
120  type is ( sll_t_particle_group_2d2v )
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_2d2v_ptr
125 
126 
127  !----------------------------------------------------------------------!
130  particle_group, &
131  n_particles, &
132  n_total_particles, &
133  charge, &
134  mass, &
135  n_weights)
136  class( sll_c_particle_group_base ), allocatable, intent( out ) :: particle_group
137  sll_int32 , intent( in ) :: n_particles
138  sll_int32 , intent( in ) :: n_total_particles
139  sll_real64 , intent( in ) :: charge
140  sll_real64 , intent( in ) :: mass
141  sll_int32 , intent(in) :: n_weights
142 
143  sll_int32 :: ierr
144 
145  sll_allocate( sll_t_particle_group_2d2v :: particle_group, ierr)
146 
147  select type( particle_group )
148  type is ( sll_t_particle_group_2d2v )
149  call particle_group%init(n_particles, n_total_particles, charge, mass, n_weights)
150  end select
151 
152  end subroutine sll_s_new_particle_group_2d2v
153 
154  !----------------------------------------------------------------------!
156  pure function get_x_2d2v( self, i ) result( r )
157  class( sll_t_particle_group_2d2v ), intent( in ) :: self
158  sll_int32 , intent( in ) :: i
159  sll_real64 :: r(3)
160 
161  r = 1.0_f64
162  r(1:2) = self%particle_array(1:2, i)
163 
164  end function get_x_2d2v
165 
166  !----------------------------------------------------------------------!
168  pure function get_v_2d2v( self, i ) result( r )
169  class( sll_t_particle_group_2d2v ), intent( in ) :: self
170  sll_int32 , intent( in ) :: i
171  sll_real64 :: r(3)
172 
173  r = 1.0_f64
174  r(1:2) = self%particle_array(3:4, i)
175 
176  end function get_v_2d2v
177 
178 
179  !----------------------------------------------------------------------!
181  pure function get_charge_2d2v( self, i , i_weight) result (r)
182  class( sll_t_particle_group_2d2v ), intent( in ) :: self
183  sll_int32 , intent( in ) :: i
184  sll_int32, optional , intent( in ) :: i_weight
185  sll_real64 :: r
186 
187  sll_int32 :: i_wi
188 
189  i_wi = 1
190  if(present(i_weight)) i_wi = i_weight
191  r = self%species%q * self%particle_array(4+i_wi, i) * self%common_weight
192 
193  end function get_charge_2d2v
194 
195 
196  !----------------------------------------------------------------------!
198  pure function get_mass_2d2v( self, i, i_weight) result (r)
199  class( sll_t_particle_group_2d2v ), intent( in ) :: self
200  sll_int32 , intent( in ) :: i
201  sll_int32, optional , intent( in ) :: i_weight
202  sll_real64 :: r
203 
204  sll_int32 :: i_wi
205 
206  i_wi = 1
207  if(present(i_weight)) i_wi = i_weight
208  r = self%species%m * self%particle_array(4+i_wi, i) * self%common_weight
209 
210  end function get_mass_2d2v
211 
212 
213  !----------------------------------------------------------------------!
215  pure function get_weights_2d2v( self, i) result (r)
216  class( sll_t_particle_group_2d2v ), intent( in ) :: self
217  sll_int32 , intent( in ) :: i
218  sll_real64 :: r(self%n_weights)
219 
220  r = self%particle_array(5:4+self%n_weights, i)
221 
222  end function get_weights_2d2v
223 
224  !----------------------------------------------------------------------!
226  pure function get_common_weight_2d2v( self )result(r)
227  class( sll_t_particle_group_2d2v ), intent( in ) :: self
228  sll_real64 :: r
229 
230  r = self%common_weight
231 
232  end function get_common_weight_2d2v
233 
234  !----------------------------------------------------------------------!
236  subroutine set_x_2d2v( self, i, x )
237  class( sll_t_particle_group_2d2v ), intent( inout ) :: self
238  sll_int32 , intent( in ) :: i
239  sll_real64 , intent( in) :: x(3)
240 
241  self%particle_array(1:2, i) = x(1:2)
242 
243  end subroutine set_x_2d2v
244 
245  !----------------------------------------------------------------------!
247  subroutine set_v_2d2v( self, i, x )
248  class( sll_t_particle_group_2d2v ), intent( inout ) :: self
249  sll_int32 , intent( in ) :: i
250  sll_real64 , intent( in) :: x(3)
251 
252  self%particle_array(3:4, i) = x(1:2)
253 
254  end subroutine set_v_2d2v
255 
256  !----------------------------------------------------------------------!
258  subroutine set_weights_2d2v( self, i, x )
259  class( sll_t_particle_group_2d2v ), intent( inout ) :: self
260  sll_int32 , intent( in ) :: i
261  sll_real64 , intent( in) :: x(self%n_weights)
262 
263  self%particle_array(5:4+self%n_weights, i) = x
264 
265  end subroutine set_weights_2d2v
266 
267  !----------------------------------------------------------------------!
269  subroutine set_common_weight_2d2v( self, x )
270  class( sll_t_particle_group_2d2v ), intent( inout ) :: self
271  sll_real64 , intent( in) :: x
272 
273  self%common_weight = x
274 
275  end subroutine set_common_weight_2d2v
276 
277 
278  !----------------------------------------------------------------------!
280  subroutine print_particle_group_2d2v(self, filename)
281  class( sll_t_particle_group_2d2v ), intent(in) :: self
282  character(len=*), intent(in) :: filename
283  sll_int32 :: file_id
284 
285  open(newunit=file_id,file=filename)
286  write(file_id,*) self%particle_array
287  close(file_id)
288 
289  end subroutine print_particle_group_2d2v
290 
291 
292 end module sll_m_particle_group_2d2v
Simple particle group group for 2d2v.
subroutine print_particle_group_2d2v(self, filename)
Print particle array to file.
subroutine, public sll_s_new_particle_group_2d2v_ptr(particle_group, n_particles, n_total_particles, charge, mass, n_weights)
Constructor for abstract type.
pure real(kind=f64) function, dimension(3) get_v_2d2v(self, i)
Get velocities of particle i.
subroutine initialize_particle_group_2d2v(self, n_particles, n_total_particles, charge, mass, n_weights)
Initialize particle group.
pure real(kind=f64) function get_common_weight_2d2v(self)
Set the common weight.
subroutine set_v_2d2v(self, i, x)
Set velocities of particle i.
subroutine set_common_weight_2d2v(self, x)
Set the common weight.
pure real(kind=f64) function get_mass_2d2v(self, i, i_weight)
Get mass of particle i ( m * particle_weight)
subroutine delete_particle_group_2d2v(self)
print particle array to file
subroutine set_weights_2d2v(self, i, x)
Set weights of particle i.
pure real(kind=f64) function get_charge_2d2v(self, i, i_weight)
Get charge of particle i ( q * particle_weight)
pure real(kind=f64) function, dimension(3) get_x_2d2v(self, i)
Get positions of particle i.
subroutine set_x_2d2v(self, i, x)
Set positions of particle i.
subroutine, public sll_s_new_particle_group_2d2v(particle_group, n_particles, n_total_particles, charge, mass, n_weights)
Constructor for abstract type.
pure real(kind=f64) function, dimension(self%n_weights) get_weights_2d2v(self, i)
Get weights of particle i.
Simple version of a PIC particle group in 2d2v.
    Report Typos and Errors