Report Typos and Errors
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
particle_methods
pic_opt
pic_opt_particles
sll_m_particle_group_6d.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
18
module
sll_m_particle_group_6d
19
!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
20
#include "sll_memory.h"
21
#include "sll_working_precision.h"
22
23
use
sll_m_cartesian_meshes
,
only
: &
24
sll_t_cartesian_mesh_3d
25
26
use
sll_m_particle_representations
,
only
: &
27
sll_t_particle_6d
28
29
!#ifdef _OPENMP
30
! use omp_lib, only: &
31
! omp_get_num_threads, &
32
! omp_get_thread_num
33
!
34
!#endif
35
36
implicit none
37
38
public
:: &
39
sll_s_particle_6d_group_init
, &
40
sll_s_particle_6d_group_free
, &
41
sll_t_particle_group_6d
42
43
private
44
45
type
::
sll_t_particle_group_6d
46
sll_int32 :: number_particles
47
sll_int32 :: active_particles
48
sll_int32 :: guard_list_size
49
sll_int32,
allocatable
:: num_postprocess_particles(:)
50
sll_real64 :: qoverm
51
type
(
sll_t_cartesian_mesh_3d
),
pointer
:: mesh
52
type
(
sll_t_particle_6d
),
pointer
:: p_list(:)
53
end type
sll_t_particle_group_6d
54
55
contains
56
57
subroutine
sll_s_particle_6d_group_init
( &
58
res, &
59
num_particles, &
60
particle_array_size, &
61
guard_list_size, &
62
qoverm, &
63
mesh)
64
65
type
(
sll_t_particle_group_6d
) :: res
66
sll_int32,
intent(in)
:: num_particles
67
sll_int32,
intent(in)
:: particle_array_size
68
sll_int32,
intent(in)
:: guard_list_size
69
sll_real64,
intent(in)
:: qoverm
70
type
(
sll_t_cartesian_mesh_3d
),
target
:: mesh
71
72
sll_int32 :: ierr
73
sll_int32 :: n_thread
74
sll_int32 :: thread_id
75
76
if
(num_particles > particle_array_size)
then
77
print *,
'sll_f_new_particle_6d_group(): ERROR, num_particles should not '
, &
78
'be greater than the memory size requested, particle_array_size.'
79
stop
80
end if
81
82
res%number_particles = num_particles
83
res%active_particles = num_particles
84
res%guard_list_size = guard_list_size
85
res%qoverm = qoverm
86
87
sll_allocate(res%p_list(particle_array_size), ierr)
88
89
n_thread = 1
90
thread_id = 0
91
!
92
! !$omp parallel PRIVATE(thread_id)
93
!#ifdef _OPENMP
94
! thread_id = OMP_GET_THREAD_NUM()
95
! if (thread_id ==0) then
96
! n_thread = OMP_GET_NUM_THREADS()
97
! endif
98
!#endif
99
! !$omp end parallel
100
!
101
! nn = guard_list_size/n_thread
102
! SLL_ALLOCATE( res%p_guard(1:n_thread), ierr)
103
sll_allocate(res%num_postprocess_particles(1:n_thread), ierr)
104
!
105
! !$omp parallel PRIVATE(thread_id)
106
!#ifdef _OPENMP
107
! thread_id = OMP_GET_THREAD_NUM()
108
!#endif
109
! SLL_ALLOCATE( res%p_guard(thread_id+1)%g_list(1:nn),ierr)
110
! !$omp end parallel
111
112
res%mesh => mesh
113
114
end subroutine
sll_s_particle_6d_group_init
115
116
subroutine
sll_s_particle_6d_group_free
(p_group)
117
118
type
(
sll_t_particle_group_6d
) :: p_group
119
120
! sll_int32 :: ierr
121
! sll_int32 :: thread_id
122
123
! if(.not. associated(p_group) ) then
124
! print *, 'delete_particle_group_2d(): ERROR, passed group was not ', &
125
! 'associated.'
126
! end if
127
128
! thread_id = 0
129
!#ifdef _OPENMP
130
! thread_id = OMP_GET_THREAD_NUM()
131
!#endif
132
deallocate
(p_group%num_postprocess_particles)
133
! SLL_DEALLOCATE(p_group%p_guard(thread_id+1)%g_list, ierr)
134
deallocate
(p_group%p_list)
135
! SLL_DEALLOCATE(p_group%p_guard, ierr)
136
! SLL_DEALLOCATE(p_group, ierr)
137
138
end subroutine
sll_s_particle_6d_group_free
139
140
end module
sll_m_particle_group_6d
sll_m_cartesian_meshes
Cartesian mesh basic types.
Definition:
sll_m_cartesian_meshes.F90:20
sll_m_particle_group_6d
Definition:
sll_m_particle_group_6d.F90:18
sll_m_particle_group_6d::sll_s_particle_6d_group_free
subroutine, public sll_s_particle_6d_group_free(p_group)
Definition:
sll_m_particle_group_6d.F90:117
sll_m_particle_group_6d::sll_s_particle_6d_group_init
subroutine, public sll_s_particle_6d_group_init(res, num_particles, particle_array_size, guard_list_size, qoverm, mesh)
Definition:
sll_m_particle_group_6d.F90:64
sll_m_particle_representations
Definition:
sll_m_particle_representations.F90:18
sll_m_cartesian_meshes::sll_t_cartesian_mesh_3d
3D cartesian mesh
Definition:
sll_m_cartesian_meshes.F90:107
sll_m_particle_group_6d::sll_t_particle_group_6d
Definition:
sll_m_particle_group_6d.F90:45
sll_m_particle_representations::sll_t_particle_6d
Type for a single particle in a 3d volume.
Definition:
sll_m_particle_representations.F90:38
Report Typos and Errors
Generated on Mon Oct 23 2023 19:15:41 for Semi-Lagrangian Library by
1.9.1