OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
myqsort.F
Go to the documentation of this file.
1Copyright> OpenRadioss
2Copyright> Copyright (C) 1986-2025 Altair Engineering Inc.
3Copyright>
4Copyright> This program is free software: you can redistribute it and/or modify
5Copyright> it under the terms of the GNU Affero General Public License as published by
6Copyright> the Free Software Foundation, either version 3 of the License, or
7Copyright> (at your option) any later version.
8Copyright>
9Copyright> This program is distributed in the hope that it will be useful,
10Copyright> but WITHOUT ANY WARRANTY; without even the implied warranty of
11Copyright> MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12Copyright> GNU Affero General Public License for more details.
13Copyright>
14Copyright> You should have received a copy of the GNU Affero General Public License
15Copyright> along with this program. If not, see <https://www.gnu.org/licenses/>.
16Copyright>
17Copyright>
18Copyright> Commercial Alternative: Altair Radioss Software
19Copyright>
20Copyright> As an alternative to this open-source version, Altair also offers Altair Radioss
21Copyright> software under a commercial license. Contact Altair to discuss further if the
22Copyright> commercial version may interest you: https://www.altair.com/radioss/.
23!||====================================================================
24!|| myqsort ../common_source/tools/sort/myqsort.F
25!||--- called by ------------------------------------------------------
26!|| add_mass_stat ../starter/source/tools/admas/add_mass_stat.F
27!|| find_dt_engine ../starter/source/coupling/rad2rad/r2r_speedup.F
28!|| find_dt_for_targeted_added_mass ../engine/source/time_step/find_dt_for_targeted_added_mass.F
29!|| hm_read_mat70 ../starter/source/materials/mat/mat070/hm_read_mat70.F
30!|| hm_read_slipring ../starter/source/tools/seatbelts/hm_read_slipring.F
31!|| hm_read_table2_1 ../starter/source/tools/curve/hm_read_table2_1.F
32!|| i20normnp ../engine/source/interfaces/int20/i20rcurv.F
33!|| i20normp ../engine/source/interfaces/int20/i20curv.F
34!|| i20normsp ../engine/source/interfaces/int20/i20curv.F
35!|| i7normnp ../engine/source/interfaces/int07/i7rcurv.F
36!|| i7normp ../engine/source/interfaces/int07/i7curv.F
37!|| inibcs_cy ../starter/source/constraints/general/bcs/lecbcscyc.F
38!|| myqsort3d ../starter/source/elements/ige3d/searchigeo3d.F
39!|| nloc_dmg_init ../starter/source/materials/fail/nloc_dmg_init.F
40!|| outri ../starter/source/materials/time_step/outri.F
41!|| outrin ../starter/source/materials/time_step/outri.F
42!|| phase_detection ../starter/source/initial_conditions/inivol/phase_detection.F
43!|| spbuc3 ../engine/source/elements/sph/spbuc3.F
44!|| spclasv ../engine/source/elements/sph/spclasv.F
45!|| sppro3 ../engine/source/elements/sph/sppro3.F
46!|| sppro31 ../starter/source/elements/sph/spbuc31.F
47!|| spsymp ../engine/source/elements/sph/spsym.F
48!|| unify_abscissa_2d ../starter/source/materials/tools/unify_abscissas_2d.F
49!||====================================================================
50 SUBROUTINE myqsort(N, A, PERM, ERROR)
51C-----------------------------------------------
52c q u i c k s o r t
53C Sedgewick algorithm from "Implementing Quicksort Programs"
54C A: data
55C N: len
56C PERM: permutations
57C-----------------------------------------------
58C I m p l i c i t T y p e s
59C-----------------------------------------------
60#include "implicit_f.inc"
61C-----------------------------------------------
62C-----------------------------------------------
63C D u m m y A r g u m e n t s
64C-----------------------------------------------
65 INTEGER N,ERROR,PERM(N)
66C REAL
68 . a(n)
69C-----------------------------------------------
70C L o c a l V a r i a b l e s
71C-----------------------------------------------
72 INTEGER :: STACKLEN
73 INTEGER :: TRESHOLD
74 INTEGER :: DONE
75C the max STACKLEN <= 1 + 2 x log2 (N+1)/(TRESHOLD + 2)
76 parameter( stacklen = 128 ,
77 . treshold = 9 )
78C
79 INTEGER :: I
80 INTEGER :: IPLUS1
81 INTEGER :: J
82 INTEGER :: JMINUS1
83 INTEGER :: K
84 INTEGER :: LEFT
85 INTEGER :: LLEN
86 INTEGER :: RIGHT
87 INTEGER :: RLEN
88 INTEGER :: TOP
89 INTEGER :: STACK(STACKLEN)
90C REAL
92 . rk, rv
93C
94 error = 0
95C
96 IF (n < 1) THEN
97 error = -1
98 RETURN
99 ENDIF
100
101 IF (n == 1) THEN
102 perm(1)=1
103 RETURN
104 ENDIF
105
106 DO i = 1, n
107 perm(i) = i
108 ENDDO
109C
110 top = 1
111 left = 1
112 right = n
113 IF (n <= treshold) THEN
114 done = 1
115 ELSE
116 done = 0
117 ENDIF
118
119c QUICKSORT
120c
121 DO WHILE (done /= 1)
122 rk = a((left+right)/2)
123 a((left+right)/2) = a(left)
124 a(left) = rk
125C
126 k = perm((left+right)/2)
127 perm((left+right)/2) = perm(left)
128 perm(left) = k
129
130 IF( a(left+1) > a(right) ) THEN
131 rk = a(left+1)
132 a(left+1) = a(right)
133 a(right) = rk
134 k = perm(left+1)
135 perm(left+1) = perm(right)
136 perm(right) = k
137 ENDIF
138 IF( a(left) > a(right) ) THEN
139 rk = a(left)
140 a(left) = a(right)
141 a(right) = rk
142 k = perm(left)
143 perm(left) = perm(right)
144 perm(right) = k
145 ENDIF
146 IF( a(left+1) > a(left) ) THEN
147 rk = a(left+1)
148 a(left+1) = a(left)
149 a(left) = rk
150 k = perm(left+1)
151 perm(left+1) = perm(left)
152 perm(left) = k
153 ENDIF
154
155 rv = a(left)
156 i = left+1
157 j = right
158
159 DO WHILE(j >= i)
160 i = i + 1
161 DO WHILE(a(i) < rv)
162 i = i +1
163 ENDDO
164 j = j - 1
165 DO WHILE(a(j) > rv)
166 j = j - 1
167 ENDDO
168 IF (j >= i) THEN
169 rk = a(i)
170 a(i) = a(j)
171 a(j) = rk
172 k = perm(i)
173 perm(i) = perm(j)
174 perm(j) = k
175 ENDIF
176 ENDDO
177C
178 rk = a(left)
179 a(left) = a(j)
180 a(j) = rk
181C
182 k = perm(left)
183 perm(left) = perm(j)
184 perm(j) = k
185C
186 llen = j-left
187 rlen = right - i + 1
188
189 IF(max(llen, rlen) <= treshold ) THEN
190 IF (top == 1) THEN
191 done = 1
192 ELSE
193 top = top - 2
194 left = stack(top)
195 right = stack(top+1)
196 ENDIF
197 ELSE IF(min(llen, rlen) <= treshold) THEN
198 IF( llen > rlen ) THEN
199 right = j - 1
200 ELSE
201 left = i
202 ENDIF
203 ELSE
204 IF( llen > rlen ) THEN
205 stack(top) = left
206 stack(top+1) = j-1
207 left = i
208 ELSE
209 stack(top) = i
210 stack(top+1) = right
211 right = j-1
212 ENDIF
213 top = top + 2
214 ENDIF
215 END DO
216c
217c INSERTION SORT
218c
219 i = n - 1
220 iplus1 = n
221 DO WHILE (i > 0)
222 IF( a(i) > a(iplus1) ) THEN
223 rk = a(i)
224 k = perm(i)
225 j = iplus1
226 jminus1 = i
227 DO WHILE(a(j) < rk)
228 a(jminus1) = a(j)
229 perm(jminus1) = perm(j)
230 jminus1 = j
231 j = j + 1
232 IF ( j > n ) EXIT
233 ENDDO
234 a(jminus1) = rk
235 perm(jminus1) = k
236 ENDIF
237C
238 iplus1 = i
239 i = i - 1
240 ENDDO
241c
242 RETURN
243c
244c -------------------
245c
246 end
#define my_real
Definition cppsort.cpp:32
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine myqsort(n, a, perm, error)
Definition myqsort.F:51