Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_particle_group_4d.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_4d), dimension(:), pointer :: p_list
58  type(sll_t_particle_4d_guard_ptr), dimension(:), pointer :: p_guard
60 
61  interface sll_o_delete
62  module procedure sll_s_particle_4d_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_4d) :: 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 
82  sll_int32 :: ierr
83  sll_int32 :: n_thread
84  sll_int32 :: thread_id
85  sll_int32 :: nn
86 
87  if (num_particles > particle_array_size) then
88  print *, 'sll_f_new_particle_4d_group(): ERROR, num_particles should not ', &
89  'be greater than the memory size requested, particle_array_size.'
90  stop
91  end if
92 
93  res%number_particles = num_particles
94  res%active_particles = num_particles
95  res%guard_list_size = guard_list_size
96  res%qoverm = qoverm
97 
98  sll_allocate(res%p_list(particle_array_size), ierr)
99 
100  n_thread = 1
101  thread_id = 0
102 
103 !$omp parallel PRIVATE(thread_id)
104 #ifdef _OPENMP
105  thread_id = omp_get_thread_num()
106  if (thread_id == 0) then
107  n_thread = omp_get_num_threads()
108  end if
109 #endif
110 !$omp end parallel
111 
112  nn = guard_list_size/n_thread
113  sll_allocate(res%p_guard(1:n_thread), ierr)
114  sll_allocate(res%num_postprocess_particles(1:n_thread), ierr)
115 
116 !$omp parallel PRIVATE(thread_id)
117 #ifdef _OPENMP
118  thread_id = omp_get_thread_num()
119 #endif
120  sll_allocate(res%p_guard(thread_id + 1)%g_list(1:nn), ierr)
121 !$omp end parallel
122 
123  res%mesh => mesh
124 
125  end subroutine sll_s_particle_4d_group_init
126 
128  num_particles, &
129  particle_array_size, &
130  guard_list_size, &
131  qoverm, &
132  mesh) result(res)
133 
134  type(sll_t_particle_group_4d), pointer :: res
135  sll_int32, intent(in) :: num_particles
136  sll_int32, intent(in) :: particle_array_size
137  sll_int32, intent(in) :: guard_list_size
138  sll_real64, intent(in) :: qoverm
139  type(sll_t_cartesian_mesh_2d), pointer :: mesh
140  sll_int32 :: ierr
141 
142  sll_allocate(res, ierr)
143 
145  res, &
146  num_particles, &
147  particle_array_size, &
148  guard_list_size, &
149  qoverm, &
150  mesh)
151 
152  end function sll_f_new_particle_4d_group
153 
154  subroutine sll_s_particle_4d_group_free(p_group)
155  type(sll_t_particle_group_4d), pointer :: p_group
156  sll_int32 :: ierr
157  sll_int32 :: thread_id
158 
159  if (.not. associated(p_group)) then
160  print *, 'delete_particle_group_2d(): ERROR, passed group was not ', &
161  'associated.'
162  end if
163 
164  thread_id = 0
165 #ifdef _OPENMP
166  thread_id = omp_get_thread_num()
167 #endif
168  sll_deallocate(p_group%num_postprocess_particles, ierr)
169  sll_deallocate(p_group%p_guard(thread_id + 1)%g_list, ierr)
170  sll_deallocate(p_group%p_list, ierr)
171  sll_deallocate(p_group%p_guard, ierr)
172  sll_deallocate(p_group, ierr)
173  end subroutine sll_s_particle_4d_group_free
174 
175 end module sll_m_particle_group_4d
Cartesian mesh basic types.
type(sll_t_particle_group_4d) function, pointer, public sll_f_new_particle_4d_group(num_particles, particle_array_size, guard_list_size, qoverm, mesh)
subroutine, public sll_s_particle_4d_group_free(p_group)
subroutine, public sll_s_particle_4d_group_init(res, num_particles, particle_array_size, guard_list_size, qoverm, mesh)
    Report Typos and Errors