MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
Utility_module.f90
Go to the documentation of this file.
1! Copyright 2019
2!
3! For a comprehensive list of the developers that contributed to these codes
4! see the UK-AMOR website.
5!
6! This file is part of UKRmol-in (UKRmol+ suite).
7!
8! UKRmol-in is free software: you can redistribute it and/or modify
9! it under the terms of the GNU General Public License as published by
10! the Free Software Foundation, either version 3 of the License, or
11! (at your option) any later version.
12!
13! UKRmol-in is distributed in the hope that it will be useful,
14! but WITHOUT ANY WARRANTY; without even the implied warranty of
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16! GNU General Public License for more details.
17!
18! You should have received a copy of the GNU General Public License
19! along with UKRmol-in (in source/COPYING). Alternatively, you can also visit
20! <https://www.gnu.org/licenses/>.
21
22!> \brief Utility module
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
27!>
28module utility_module
29
30 use consts_mpi_ci, only: nidx
31 use precisn, only: longint, wp
32 use mpi_gbl, only: mpi_mod_wtime, mpi_xermsg
33
34 implicit none
35
36 public string_hash, get_real_time, get_cpu_time, compute_total_triangular, triangular_index_to_ij
37 public compute_total_box, box_index_to_ij
38
39contains
40
41 !> \brief Calculate triangular area
42 !> \authors A Al-Refaie
43 !> \date 2017
44 !>
45 !> Calculate n*(n+1)/2.
46 !>
47 integer(longint) function compute_total_triangular (n)
48 integer, intent(in) :: n
49
50 compute_total_triangular = n * (n + 1_longint) / 2
51
52 end function compute_total_triangular
53
54
55 !> \brief Calculate rectangular product
56 !> \authors A Al-Refaie
57 !> \date 2017
58 !>
59 !> Calculate width*height.
60 !>
61 integer function compute_total_box (width, height)
62 integer, intent(in) :: width, height
63
64 compute_total_box = width * height
65
66 end function compute_total_box
67
68
69 !> \brief Extract indices from rectangular multi-index
70 !> \authors A Al-Refaie
71 !> \date 2017
72 !>
73 subroutine box_index_to_ij (idx, height, i, j)
74 integer, intent(in) :: idx, height
75 integer, intent(out) :: i, j
76
77 i = mod(idx - 1, height) + 1
78 j = (idx - 1) / height + 1
79
80 end subroutine box_index_to_ij
81
82
83 !> \brief Extract indices from triangular multi-index
84 !> \authors A Al-Refaie
85 !> \date 2017
86 !>
87 subroutine triangular_index_to_ij (idx_f, N, row, column)
88 integer, intent(in) :: idx_f, n
89 integer, intent(out) :: row, column
90 integer :: idx, ii, k, jj
91
92 idx = idx_f - 1
93 ii = n * (n + 1) / 2 - 1 - idx
94 k = int((sqrt(8.0_wp*real(ii,wp) + 1.0_wp) - 1.0_wp) / 2.0_wp, longint)
95 jj = ii - k * (k + 1) / 2
96
97 row = n - 1 - k
98 column = n - 1 - jj
99
100 end subroutine triangular_index_to_ij
101
102
103 !> \brief Calculate a string hash
104 !> \authors A Al-Refaie
105 !> \date 2017
106 !>
107 !> This hash assumes at least 29-bit integers. It is supposedly
108 !> documented in Aho, Sethi, and Ullman, pp. 434-438
109 !>
110 integer function string_hash (str, table_size) result(h)
111 character(len=*), intent(in) :: str
112 integer, intent(in) :: table_size
113 integer :: i, chr, g, mask = int(z"1FFFFFF")
114
115 h = 0
116 do i = 1, len_trim(str)
117 chr = ichar(str(i:i))
118 h = ishft(h, 4) + chr
119 g = ishft(h,-24)
120 h = iand(ieor(h,g), mask)
121 end do
122 h = 1 + modulo(h, table_size)
123
124 end function string_hash
125
126
127 !> \brief Get current (real) time
128 !> \authors A Al-Refaie
129 !> \date 2017
130 !>
131 !> Uses a function from GBTOlib.
132 !>
133 function get_real_time () result(t)
134 real(wp) :: t
135
136 t = mpi_mod_wtime()
137
138 end function get_real_time
139
140
141 !> \brief Get current (CPU) time
142 !> \authors A Al-Refaie
143 !> \date 2017
144 !>
145 function get_cpu_time () result(t)
146 real(wp) :: t
147
148 call cpu_time(t)
149
150 end function get_cpu_time
151
152
153 !> \brief Store 8 integers in given number of integers
154 !> \authors J Benda
155 !> \date 2023
156 !>
157 !> Copy or pack 8 integers to the given packed storage.
158 !>
159 subroutine pack_ints (i1, i2, i3, i4, i5, i6, i7, i8, z)
160
161 integer, intent(in) :: i1, i2, i3, i4, i5, i6, i7, i8
162 integer(longint), intent(out) :: z(:)
163
164 integer(longint), parameter :: max16 = 2_longint**16 - 1
165 integer(longint), parameter :: max32 = 2_longint**32 - 1
166
167 select case (nidx)
168
169 case (2)
170
171 if (i1 > max16 .or. i2 > max16 .or. i3 > max16 .or. i4 > max16 .or. &
172 i5 > max16 .or. i6 > max16 .or. i7 > max16 .or. i8 > max16) then
173 call mpi_xermsg('Utility_module', 'pack_ints', &
174 'Cannot pack integers, values too large. Try increasing NIDX.', 2, 1)
175 end if
176
177 z(1) = i1
178 z(1) = ior(ishft(z(1), 16_longint), int(i2, longint))
179 z(1) = ior(ishft(z(1), 16_longint), int(i3, longint))
180 z(1) = ior(ishft(z(1), 16_longint), int(i4, longint))
181 z(2) = i5
182 z(2) = ior(ishft(z(2), 16_longint), int(i6, longint))
183 z(2) = ior(ishft(z(2), 16_longint), int(i7, longint))
184 z(2) = ior(ishft(z(2), 16_longint), int(i8, longint))
185
186 case (4)
187
188 if (i1 > max32 .or. i2 > max32 .or. i3 > max32 .or. i4 > max32 .or. &
189 i5 > max32 .or. i6 > max32 .or. i7 > max32 .or. i8 > max32) then
190 call mpi_xermsg('Utility_module', 'pack_ints', &
191 'Cannot pack integers, values too large. Try increasing NIDX.', 4, 1)
192 end if
193
194 z(1) = i1
195 z(1) = ior(ishft(z(1), 32_longint), int(i2, longint))
196 z(2) = i3
197 z(2) = ior(ishft(z(2), 32_longint), int(i4, longint))
198 z(3) = i5
199 z(3) = ior(ishft(z(3), 32_longint), int(i6, longint))
200 z(4) = i7
201 z(4) = ior(ishft(z(4), 32_longint), int(i8, longint))
202
203 case (8)
204
205 z(1) = i1
206 z(2) = i2
207 z(3) = i3
208 z(4) = i4
209 z(5) = i5
210 z(6) = i6
211 z(7) = i7
212 z(8) = i8
213
214 case default
215
216 call mpi_xermsg('Utility_module', 'pack_ints', 'Packing implemented only for NIDX = 2, 4 or 8.', 1, 1)
217
218 end select
219
220 end subroutine pack_ints
221
222
223 !> \brief Retrieve 8 integers from given number of integers
224 !> \authors J Benda
225 !> \date 2023
226 !>
227 !> Copy or unpack 8 default integers from the provided packed storage.
228 !>
229 subroutine unpack_ints (z, u)
230
231 integer(longint), intent(in) :: z(:)
232 integer, intent(out) :: u(8)
233
234 integer(longint), parameter :: mask16 = 2_longint**16 - 1
235 integer(longint), parameter :: mask32 = 2_longint**32 - 1
236
237 integer(longint) :: n
238
239 select case (nidx)
240
241 case (2)
242
243 n = z(1); u(4) = int(iand(n, mask16))
244 n = ishft(n, -16_longint); u(3) = int(iand(n, mask16))
245 n = ishft(n, -16_longint); u(2) = int(iand(n, mask16))
246 n = ishft(n, -16_longint); u(1) = int(iand(n, mask16))
247 n = z(2); u(8) = int(iand(n, mask16))
248 n = ishft(n, -16_longint); u(7) = int(iand(n, mask16))
249 n = ishft(n, -16_longint); u(6) = int(iand(n, mask16))
250 n = ishft(n, -16_longint); u(5) = int(iand(n, mask16))
251
252 case (4)
253
254 n = z(1); u(2) = int(iand(n, mask32))
255 n = ishft(n, -32_longint); u(1) = int(iand(n, mask32))
256 n = z(2); u(4) = int(iand(n, mask32))
257 n = ishft(n, -32_longint); u(3) = int(iand(n, mask32))
258 n = z(3); u(6) = int(iand(n, mask32))
259 n = ishft(n, -32_longint); u(5) = int(iand(n, mask32))
260 n = z(4); u(8) = int(iand(n, mask32))
261 n = ishft(n, -32_longint); u(7) = int(iand(n, mask32))
262
263 case (8)
264
265 u(1) = int(z(1))
266 u(2) = int(z(2))
267 u(3) = int(z(3))
268 u(4) = int(z(4))
269 u(5) = int(z(5))
270 u(6) = int(z(6))
271 u(7) = int(z(7))
272 u(8) = int(z(8))
273
274 case default
275
276 call mpi_xermsg('Utility_module', 'unpack_ints', 'Unpacking implemented only for NIDX = 2, 4 or 8.', 1, 1)
277
278 end select
279
280 end subroutine unpack_ints
281
282
283 !> \brief Compare two arrays
284 !> \authors J Benda
285 !> \date 2023
286 !>
287 !> Lexicographically compare two integer arrays of equal length NIDX. Return 0 if the arrays are indentical, -1 if the first
288 !> unequal element is smaller in the first array that in the second array, and +1 otherwise.
289 !>
290 integer function lexicographical_compare (a, b) result (verdict)
291
292 integer(longint), intent(in) :: a(nidx), b(nidx)
293
294 integer :: i
295
296 verdict = 0
297
298 do i = 1, nidx
299 if (a(i) < b(i)) then
300 verdict = -1
301 exit
302 end if
303 if (a(i) > b(i)) then
304 verdict = +1
305 exit
306 end if
307 end do
308
309 end function lexicographical_compare
310
311end module utility_module
MPI SCATCI Constants module.
integer, parameter nidx
Number of long integers used to store up to 8 (packed) integral indices.