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