-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathutils_sort.f90
260 lines (235 loc) · 6.01 KB
/
utils_sort.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
!--------------------------------------------------------------------------!
! The Phantom Smoothed Particle Hydrodynamics code, by Daniel Price et al. !
! Copyright (c) 2007-2015 The Authors (see AUTHORS) !
! See LICENCE file for usage and distribution conditions !
! http://users.monash.edu.au/~dprice/phantom !
!--------------------------------------------------------------------------!
!+
! MODULE: sortutils
!
! DESCRIPTION:
! contains low level sorting utilities
!
! REFERENCES: None
!
! OWNER: Daniel Price
!
! $Id: 2b98dc6b068598c483b2f11b028f29aa49ca867a $
!
! RUNTIME PARAMETERS: None
!
! DEPENDENCIES: None
!+
!--------------------------------------------------------------------------
module sortutils
implicit none
public :: indexx,indexxfunc,r2func,r2func_origin,set_r2func_origin
private
real, private :: x0,y0,z0
contains
!----------------------------------------------------------------
!+
! function returning radius squared given xyzh
!+
!----------------------------------------------------------------
real function r2func(xyzh)
real, intent(in) :: xyzh(4)
r2func = xyzh(1)**2 + xyzh(2)**2 + xyzh(3)**2
end function r2func
!----------------------------------------------------------------
!+
! function returning radius squared given xyzh and origin pos
!+
!----------------------------------------------------------------
subroutine set_r2func_origin(xorigin,yorigin,zorigin)
real, intent(in) :: xorigin,yorigin,zorigin
x0 = xorigin
y0 = yorigin
z0 = zorigin
end subroutine set_r2func_origin
real function r2func_origin(xyzh)
real, intent(in) :: xyzh(4)
r2func_origin = (xyzh(1)-x0)**2 + (xyzh(2)-y0)**2 + (xyzh(3)-z0)**2
end function r2func_origin
!----------------------------------------------------------------
!+
! standard sorting routine using Quicksort
!+
!----------------------------------------------------------------
subroutine indexx(n, arr, indx)
integer, parameter :: m=7, nstack=500
integer, intent(in) :: n
real, intent(in) :: arr(n)
integer, intent(out) :: indx(n)
integer :: i,j,k,l,ir,jstack,indxt,itemp
integer :: istack(nstack)
real :: a
do j = 1, n
indx(j) = j
enddo
jstack = 0
l = 1
ir = n
1 if (ir - l < m) then
do j = l + 1, ir
indxt = indx(j)
a = arr(indxt)
do i = j - 1, 1, -1
if (arr(indx(i)) <= a) goto 2
indx(i + 1) = indx(i)
end do
i = 0
2 indx(i + 1) = indxt
end do
if (jstack==0) return
ir = istack(jstack)
l = istack(jstack - 1)
jstack = jstack - 2
else
k = (l + ir)/2
itemp = indx(k)
indx(k) = indx(l + 1)
indx(l + 1) = itemp
if (arr(indx(l + 1)) > arr(indx(ir))) then
itemp = indx(l + 1)
indx(l + 1) = indx(ir)
indx(ir) = itemp
endif
if (arr(indx(l)) > arr(indx(ir))) then
itemp = indx(l)
indx(l) = indx(ir)
indx(ir) = itemp
endif
if (arr(indx(l + 1)) > arr(indx(l))) then
itemp = indx(l + 1)
indx(l + 1) = indx(l)
indx(l) = itemp
endif
i = l + 1
j = ir
indxt = indx(l)
a = arr(indxt)
3 continue
i = i + 1
if (arr(indx(i)) < a) goto 3
4 continue
j = j - 1
if (arr(indx(j)) > a) goto 4
if (j < i) goto 5
itemp = indx(i)
indx(i) = indx(j)
indx(j) = itemp
goto 3
5 indx(l) = indx(j)
indx(j) = indxt
jstack = jstack + 2
if (jstack > nstack) then
print*,'fatal error!!! stacksize exceeded in sort'
print*,'need to set parameter nstack higher in subroutine indexx '
stop
endif
if (ir - i + 1 >= j - l) then
istack(jstack) = ir
istack(jstack - 1) = i
ir = j - 1
else
istack(jstack) = j - 1
istack(jstack - 1) = l
l = i
endif
endif
goto 1
end subroutine indexx
!----------------------------------------------------------------
!+
! customised low-memory sorting routine using Quicksort
! sort key value on-the-fly by calling the function func
! which can be any function of the particle positions
!+
!----------------------------------------------------------------
subroutine indexxfunc(n, func, xyzh, indx)
integer, parameter :: m=7, nstack=500
integer, intent(in) :: n
real, external :: func
real, intent(in) :: xyzh(:,:)
integer, intent(out) :: indx(n)
integer :: i,j,k,l,ir,jstack,indxt,itemp
integer :: istack(nstack)
real :: a
do j = 1, n
indx(j) = j
enddo
jstack = 0
l = 1
ir = n
1 if (ir - l < m) then
do j = l + 1, ir
indxt = indx(j)
a = func(xyzh(:,indxt))
do i = j - 1, 1, -1
if (func(xyzh(:,indx(i))) <= a) goto 2
indx(i + 1) = indx(i)
end do
i = 0
2 indx(i + 1) = indxt
end do
if (jstack==0) return
ir = istack(jstack)
l = istack(jstack - 1)
jstack = jstack - 2
else
k = (l + ir)/2
itemp = indx(k)
indx(k) = indx(l + 1)
indx(l + 1) = itemp
if (func(xyzh(:,indx(l+1))) > func(xyzh(:,indx(ir)))) then
itemp = indx(l + 1)
indx(l + 1) = indx(ir)
indx(ir) = itemp
endif
if (func(xyzh(:,indx(l))) > func(xyzh(:,indx(ir)))) then
itemp = indx(l)
indx(l) = indx(ir)
indx(ir) = itemp
endif
if (func(xyzh(:,indx(l+1))) > func(xyzh(:,indx(l)))) then
itemp = indx(l + 1)
indx(l + 1) = indx(l)
indx(l) = itemp
endif
i = l + 1
j = ir
indxt = indx(l)
a = func(xyzh(:,indxt))
3 continue
i = i + 1
if (func(xyzh(:,indx(i))) < a) goto 3
4 continue
j = j - 1
if (func(xyzh(:,indx(j))) > a) goto 4
if (j < i) goto 5
itemp = indx(i)
indx(i) = indx(j)
indx(j) = itemp
goto 3
5 indx(l) = indx(j)
indx(j) = indxt
jstack = jstack + 2
if (jstack > nstack) then
print*,'fatal error!!! stacksize exceeded in sort'
print*,'need to set parameter nstack higher in subroutine indexxfunc '
stop
endif
if (ir - i + 1 >= j - l) then
istack(jstack) = ir
istack(jstack - 1) = i
ir = j - 1
else
istack(jstack) = j - 1
istack(jstack - 1) = l
l = i
endif
endif
goto 1
end subroutine indexxfunc
end module sortutils