The Elk Code
modomp.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2015 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 module modomp
7 
8 use omp_lib
9 
10 ! maximum number of OpenMP threads available
11 integer maxthd
12 ! maximum number of OpenMP threads for the first nesting level
13 integer maxthd1
14 ! maximum number of threads available to MKL
15 integer maxthdmkl
16 ! maximum OpenMP nesting level
17 integer maxlvl
18 ! number of active OpenMP threads for each nesting level
19 integer, allocatable, protected :: nathd(:)
20 
21 interface
22 
23 integer function mkl_set_num_threads_local(num_threads)
24 integer, intent(in) :: num_threads
25 end function
26 
27 end interface
28 
29 contains
30 
31 subroutine omp_init
32 implicit none
33 ! determine the maximum number of available threads
34 select case(maxthd)
35 case(:-1)
36 ! set the number of threads equal to a fraction of the number of processors
37  maxthd=omp_get_num_procs()/abs(maxthd)
38  if (maxthd < 1) maxthd=1
39  call omp_set_num_threads(maxthd)
40 case(0)
41 ! use the system default number of threads
42  maxthd=omp_get_max_threads()
43 case default
44 ! use the number of threads specified in the input file
45  call omp_set_num_threads(maxthd)
46 end select
47 ! determine the maximum number of threads available at first nesting level
48 select case(maxthd1)
49 case(:-1)
51  if (maxthd1 < 1) maxthd1=1
52 case(0)
54 case default
56 end select
57 ! switch off dynamic allocation of threads
58 call omp_set_dynamic(.false.)
59 ! allow nested parallelism (deprecated in OpenMP version 5)
60 call omp_set_nested(.true.)
61 ! set the maximum nesting level
62 call omp_set_max_active_levels(maxlvl)
63 ! allocate the number of active threads array
64 if (allocated(nathd)) deallocate(nathd)
65 allocate(nathd(0:maxlvl))
66 ! initialise the number of active threads
67 call omp_reset
68 end subroutine
69 
70 subroutine omp_reset
71 implicit none
72 ! number of active threads at each nesting level
73 nathd(0)=1
74 nathd(1:maxlvl)=0
75 end subroutine
76 
77 subroutine holdthd(nloop,nthd)
78 implicit none
79 ! arguments
80 integer, intent(in) :: nloop
81 integer, intent(out) :: nthd
82 ! local variables
83 integer lvl,na,n
84 ! current nesting level
85 lvl=omp_get_level()
86 if ((lvl < 0).or.(lvl >= maxlvl)) then
87  nthd=1
88  return
89 end if
90 ! determine number of active threads at the current nesting level
91 na=nathd(lvl)
92 if (na < 1) na=1
93 if (na > maxthd) na=maxthd
94 ! number of threads allowed for this loop
95 nthd=maxthd/na
96 if (mod(maxthd,na) > 0) nthd=nthd+1
97 if (lvl == 0) then
98  if (nthd > maxthd1) nthd=maxthd1
99 end if
100 if (nthd < 1) nthd=1
101 if (nthd > nloop) nthd=nloop
102 if (nthd > maxthd) nthd=maxthd
103 ! add to number of active threads in next nesting level
104 n=nathd(lvl+1)+nthd
105 if (n < 0) n=0
106 if (n > maxthd) n=maxthd
107 !$OMP ATOMIC WRITE
108 nathd(lvl+1)=n
109 end subroutine
110 
111 subroutine freethd(nthd)
112 implicit none
113 ! arguments
114 integer, intent(in) :: nthd
115 ! local variables
116 integer lvl,n
117 ! current nesting level
118 lvl=omp_get_level()
119 if ((lvl < 0).or.(lvl >= maxlvl)) return
120 ! subtract from the number of active threads in next nesting level
121 n=nathd(lvl+1)-nthd
122 if (n < 0) n=0
123 if (n > maxthd) n=maxthd
124 !$OMP ATOMIC WRITE
125 nathd(lvl+1)=n
126 end subroutine
127 
128 end module
129 
integer maxthd1
Definition: modomp.f90:13
integer maxthd
Definition: modomp.f90:11
Definition: modomp.f90:6
integer maxthdmkl
Definition: modomp.f90:15
subroutine omp_init
Definition: modomp.f90:32
integer maxlvl
Definition: modomp.f90:17
integer, dimension(:), allocatable, protected nathd
Definition: modomp.f90:19
subroutine omp_reset
Definition: modomp.f90:71
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78