libSUFR
a LIBrary of Some Useful Fortran Routines
All Classes Namespaces Files Functions Variables Pages
sorting.f90
Go to the documentation of this file.
1!> \file sorting.f90 Procedures to sort arrays
2
3
4! Copyright (c) 2002-2025 Marc van der Sluys - Nikhef/Utrecht University - marc.vandersluys.nl
5!
6! This file is part of the libSUFR package,
7! see: http://libsufr.sourceforge.net/
8!
9! This is free software: you can redistribute it and/or modify it under the terms of the European Union
10! Public Licence 1.2 (EUPL 1.2). This software is distributed in the hope that it will be useful, but
11! WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12! PURPOSE. See the EU Public Licence for more details. You should have received a copy of the European
13! Union Public Licence along with this code. If not, see <https://www.eupl.eu/1.2/en/>.
14!
15!
16
17
18
19
20!***********************************************************************************************************************************
21!> \brief Procedures for sorting
22
24 implicit none
25 save
26
27contains
28
29
30 !*********************************************************************************************************************************
31 !> \brief Return a list of indices index_list that sorts the members of array to ascending value, using a Quicksort algorithm.
32 !!
33 !! \param array Array of size n with values that must be sorted - use dble(array) for other variable types
34 !!
35 !! \param index_list List with indices of array values, sorted to ascending value, same dimension and size as array. (output)
36 !! array(index_list) gives the sorted array.
37 !!
38 !! \param mask Mask to apply to array. If present, index_list will have zeroes after the last meaningful entry (optional)
39 !! \param index_n Number of meaningful elements in index_list, after applying mask (output, optional)
40 !!
41 !! \note This routine does not need to be called directly, but is implicitly called by sort_array().
42 !! \see Numerical Recipes in Fortran 77, Sect. 8.2, 8.4.
43
44 subroutine sorted_index_list(array, index_list, mask, index_n)
45 use sufr_kinds, only: double
47
48 implicit none
49 real(double), intent(in) :: array(:)
50 integer, intent(out) :: index_list(:)
51 logical, intent(in), optional :: mask(:)
52 integer, intent(out),optional :: index_n
53
54 integer, parameter :: mm=7, nstack=50
55 integer :: ii,index_i,ir,itemp,istack(nstack), jj,jstack, kk,ll, loc_list(size(index_list)), loc_index_n
56 real(double) :: arr_i
57 logical :: loc_mask(size(array))
58
59
60 if(size(array).ne.size(index_list)) &
61 call quit_program_error('libSUFR sorted_index_list(): array and index_list must have the same size',0)
62
63 loc_mask = .true.
64 if(present(mask)) loc_mask = mask
65
66 index_list = 0
67 loc_list = 0
68 do jj=1,size(array)
69 loc_list(jj) = jj
70 end do
71
72
73 ! Starting values:
74 jstack = 0
75 ll = 1
76 ir = size(array)
77
78 mainloop: do
79 if(ir-ll.lt.mm) then
80
81 do jj=ll+1,ir
82 index_i = loc_list(jj)
83 arr_i = array(index_i)
84 subloop1: do ii=jj-1,ll,-1
85 if(array(loc_list(ii)).le.arr_i) exit subloop1
86 loc_list(ii+1) = loc_list(ii)
87 end do subloop1
88 if(ii.eq.0) ii = ll-1 ! if subloop1 completes without match for arr_i
89
90 loc_list(ii+1) = index_i
91 end do
92
93
94 !*************************************************
95 ! Done - apply mask and return to caller routine
96 if(jstack.eq.0) then
97
98 loc_index_n = 0
99 do ii=1,size(array)
100 if(loc_mask(loc_list(ii))) then
101 loc_index_n = loc_index_n + 1
102 index_list(loc_index_n) = loc_list(ii)
103 if(present(index_n)) index_n = loc_index_n
104 end if
105 end do
106
107 return
108
109 end if
110 !*************************************************
111
112
113 ir = istack(jstack)
114 ll = istack(jstack-1)
115 jstack = jstack - 2
116
117 else
118
119 kk = (ll+ir)/2
120 itemp = loc_list(kk)
121 loc_list(kk) = loc_list(ll+1)
122 loc_list(ll+1) = itemp
123
124 if(array(loc_list(ll)).gt.array(loc_list(ir))) then
125 itemp = loc_list(ll)
126 loc_list(ll) = loc_list(ir)
127 loc_list(ir) = itemp
128 end if
129
130 if(array(loc_list(ll+1)).gt.array(loc_list(ir))) then
131 itemp = loc_list(ll+1)
132 loc_list(ll+1) = loc_list(ir)
133 loc_list(ir) = itemp
134 end if
135
136 if(array(loc_list(ll)).gt.array(loc_list(ll+1))) then
137 itemp = loc_list(ll)
138 loc_list(ll) = loc_list(ll+1)
139 loc_list(ll+1) = itemp
140 end if
141
142 ii = ll+1
143 jj = ir
144 index_i = loc_list(ll+1)
145 arr_i = array(index_i)
146
147
148 subloop2: do
149 ii = ii+1
150 if(array(loc_list(ii)).lt.arr_i) cycle subloop2
151
152
153 subloop3: do
154 jj = jj-1
155 if(array(loc_list(jj)).le.arr_i) exit
156 end do subloop3
157
158
159 if(jj.ge.ii) then
160 itemp = loc_list(ii)
161 loc_list(ii) = loc_list(jj)
162 loc_list(jj) = itemp
163 cycle subloop2
164 end if
165
166 exit subloop2
167
168 end do subloop2
169
170
171 loc_list(ll+1) = loc_list(jj)
172 loc_list(jj) = index_i
173 jstack = jstack + 2
174
175 if(jstack.gt.nstack) write(0,'(A)')' sorted_index_list(): nstack is too small'
176
177 if(ir-ii+1.ge.jj-ll) then
178 istack(jstack) = ir
179 istack(jstack-1) = ii
180 ir = jj-1
181 else
182 istack(jstack) = jj-1
183 istack(jstack-1) = ll
184 ll = ii
185 end if
186
187 end if
188
189 end do mainloop
190
191 end subroutine sorted_index_list
192 !*********************************************************************************************************************************
193
194
195
196 !*********************************************************************************************************************************
197 !> \brief Sort an array to ascending value.
198 !!
199 !! \param array Input: array to be sorted, output: sorted array (double)
200 !!
201 !! \note Uses sorted_index_list()
202
203 subroutine sort_array(array)
204 use sufr_kinds, only: double
205
206 implicit none
207 real(double), intent(inout) :: array(:)
208
209 real(double) :: array1(size(array))
210 integer :: index_list(size(array))
211
212 array1 = array
213 call sorted_index_list(array1,index_list)
214 array = array1(index_list)
215
216 end subroutine sort_array
217 !*********************************************************************************************************************************
218
219
220
221 !*********************************************************************************************************************************
222 !> \brief Sort an array to ascending value - single-precision wrapper for sort_array
223 !!
224 !! \param array Input: array to be sorted, output: sorted array (double)
225 !!
226 !! \todo Check whether array = array(index_list) works
227
228 subroutine sort_array_sp(array)
229 implicit none
230 real, intent(inout) :: array(:)
231
232 real :: array1(size(array))
233 integer :: index_list(size(array))
234
235 call sorted_index_list(dble(array), index_list)
236 array1 = array
237 array = array1(index_list) ! CHECK: would array = array(index_list) work?
238
239 end subroutine sort_array_sp
240 !*********************************************************************************************************************************
241
242
243
244 !*********************************************************************************************************************************
245 !> \brief Sort strings alphabetically
246 !!
247 !! \param narr Number of strings in array str
248 !! \param strlen Length of the strings in array str
249 !! \param strarr Array of n strings with length len
250 !! \param index_list Sorting index (output)
251
252 subroutine sort_string_array(narr, strlen, strarr, index_list)
253 implicit none
254 integer, intent(in) :: narr, strlen
255 character, intent(in) :: strarr(narr)*(strlen)
256 integer, intent(out) :: index_list(narr)
257
258 integer :: i,l,ic
259 real :: score(narr)
260
261 do i=1,narr
262 score(i) = 0.
263 do l=1,strlen
264 ic = ichar(strarr(i)(l:l))
265 if(ic.ge.97.and.ic.le.122) ic = ic-32
266 ic = ic-32 ! Limits values between 1 and 64
267 score(i) = score(i) + real(ic) * 64.**(strlen-l)
268 end do
269 end do
270
271 call sorted_index_list(dble(score), index_list)
272
273 end subroutine sort_string_array
274 !*********************************************************************************************************************************
275
276
277
278
279
280
281end module sufr_sorting
282!***********************************************************************************************************************************
Provides kinds and related constants/routines.
Definition kinds.f90:26
integer, parameter double
Double-precision float. Precision = 15, range = 307.
Definition kinds.f90:35
Procedures for sorting.
Definition sorting.f90:23
subroutine sort_string_array(narr, strlen, strarr, index_list)
Sort strings alphabetically.
Definition sorting.f90:253
subroutine sorted_index_list(array, index_list, mask, index_n)
Return a list of indices index_list that sorts the members of array to ascending value,...
Definition sorting.f90:45
subroutine sort_array_sp(array)
Sort an array to ascending value - single-precision wrapper for sort_array.
Definition sorting.f90:229
subroutine sort_array(array)
Sort an array to ascending value.
Definition sorting.f90:204
System-related procedures.
Definition system.f90:23
subroutine quit_program_error(message, status)
Print an error message to StdErr and stop the execution of the current program.
Definition system.f90:78