Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_particle_group_2d.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! CALVI project team
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
19 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
20 #include "sll_memory.h"
21 #include "sll_working_precision.h"
22 
23  use sll_m_cartesian_meshes, only: &
25 
29 
30 #ifdef _OPENMP
31  use omp_lib, only: &
32  omp_get_num_threads, &
33  omp_get_thread_num
34 
35 #endif
36  implicit none
37 
38  public :: &
42  sll_o_delete, &
44 
45  private
46 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 
49  sll_int32 :: number_particles
50  sll_int32 :: active_particles
51  sll_int32 :: guard_list_size
52  ! an array indexed by the thread number, of the number of particles
53  ! to post-process after the main loop
54  sll_int32, dimension(:), pointer :: num_postprocess_particles
55  sll_real64 :: qoverm
56  type(sll_t_cartesian_mesh_2d), pointer :: mesh
57  type(sll_t_particle_2d), dimension(:), pointer :: p_list
58  type(sll_t_particle_2d_guard_ptr), dimension(:), pointer :: p_guard
60 
61  interface sll_o_delete
62  module procedure sll_s_particle_2d_group_free
63  end interface sll_o_delete
64 
65 contains
66 
68  res, &
69  num_particles, &
70  particle_array_size, &
71  guard_list_size, &
72  qoverm, &
73  mesh)
74 
75  type(sll_t_particle_group_2d) :: res
76  sll_int32, intent(in) :: num_particles
77  sll_int32, intent(in) :: particle_array_size
78  sll_int32, intent(in) :: guard_list_size
79  sll_real64, intent(in) :: qoverm
80  type(sll_t_cartesian_mesh_2d), target :: mesh
81  sll_int32 :: ierr
82  sll_int32 :: n_thread
83  sll_int32 :: thread_id
84  sll_int32 :: nn
85 
86  if (num_particles > particle_array_size) then
87  print *, 'sll_f_new_particle_2d_group(): ERROR, num_particles should not ', &
88  'be greater than the memory size requested, particle_array_size.'
89  stop
90  end if
91 
92  res%number_particles = num_particles
93  res%active_particles = num_particles
94  res%guard_list_size = guard_list_size
95  res%qoverm = qoverm
96 
97  sll_allocate(res%p_list(particle_array_size), ierr)
98 
99  n_thread = 1
100  thread_id = 0
101 
102 !$omp parallel PRIVATE(thread_id)
103 #ifdef _OPENMP
104  thread_id = omp_get_thread_num()
105  if (thread_id == 0) then
106  n_thread = omp_get_num_threads()
107  end if
108 #endif
109 !$omp end parallel
110 
111  nn = guard_list_size/n_thread
112  sll_allocate(res%p_guard(1:n_thread), ierr)
113  sll_allocate(res%num_postprocess_particles(1:n_thread), ierr)
114 
115 !$omp parallel PRIVATE(thread_id)
116 #ifdef _OPENMP
117  thread_id = omp_get_thread_num()
118 #endif
119  sll_allocate(res%p_guard(thread_id + 1)%g_list(1:nn), ierr)
120 !$omp end parallel
121 
122  res%mesh => mesh
123 
124  end subroutine sll_s_particle_2d_group_init
125 
127  num_particles, &
128  particle_array_size, &
129  guard_list_size, &
130  qoverm, &
131  mesh) result(res)
132 
133  type(sll_t_particle_group_2d), pointer :: res
134  sll_int32, intent(in) :: num_particles
135  sll_int32, intent(in) :: particle_array_size
136  sll_int32, intent(in) :: guard_list_size
137  sll_real64, intent(in) :: qoverm
138  type(sll_t_cartesian_mesh_2d), pointer :: mesh
139  sll_int32 :: ierr
140 
141  sll_allocate(res, ierr)
143  res, &
144  num_particles, &
145  particle_array_size, &
146  guard_list_size, &
147  qoverm, &
148  mesh)
149 
150  end function sll_f_new_particle_2d_group
151 
152  subroutine sll_s_particle_2d_group_free(p_group)
153  type(sll_t_particle_group_2d), pointer :: p_group
154  sll_int32 :: ierr
155  sll_int32 :: thread_id
156 
157  if (.not. associated(p_group)) then
158  print *, 'delete_particle_group_2d(): ERROR, passed group was not ', &
159  'associated.'
160  end if
161 
162  thread_id = 0
163 #ifdef _OPENMP
164  thread_id = omp_get_thread_num()
165 #endif
166  sll_deallocate(p_group%num_postprocess_particles, ierr)
167  sll_deallocate(p_group%p_guard(thread_id + 1)%g_list, ierr)
168  sll_deallocate(p_group%p_list, ierr)
169  sll_deallocate(p_group%p_guard, ierr)
170  sll_deallocate(p_group, ierr)
171  end subroutine sll_s_particle_2d_group_free
172 
173 end module sll_m_particle_group_2d
Cartesian mesh basic types.
subroutine, public sll_s_particle_2d_group_init(res, num_particles, particle_array_size, guard_list_size, qoverm, mesh)
subroutine, public sll_s_particle_2d_group_free(p_group)
type(sll_t_particle_group_2d) function, pointer, public sll_f_new_particle_2d_group(num_particles, particle_array_size, guard_list_size, qoverm, mesh)
    Report Typos and Errors