Split g96 I/O routines from confio.cpp
[gromacs.git] / src / gromacs / linearalgebra / gmx_lapack / dsteqr.c
blobead523ba54df4a718473e62113caf1e02d26fbf6
1 #include <math.h>
2 #include "gromacs/utility/real.h"
4 #include "../gmx_blas.h"
5 #include "../gmx_lapack.h"
6 #include "lapack_limits.h"
8 void
9 F77_FUNC(dsteqr,DSTEQR)(const char * compz,
10 int * n,
11 double * d__,
12 double * e,
13 double * z__,
14 int * ldz,
15 double * work,
16 int * info)
18 double c_b9 = 0.;
19 double c_b10 = 1.;
20 int c__0 = 0;
21 int c__1 = 1;
22 int c__2 = 2;
23 int z_dim1, z_offset, i__1, i__2;
24 double d__1, d__2;
26 double b, c__, f, g;
27 int i__, j, k, l, m;
28 double p, r__, s;
29 int l1, ii, mm, lm1, mm1, nm1;
30 double rt1, rt2, eps;
31 int lsv;
32 double tst, eps2;
33 int lend, jtot;
34 double anorm;
35 int lendm1, lendp1;
36 int iscale;
37 double safmin,minval;
38 double safmax;
39 int lendsv;
40 double ssfmin;
41 int nmaxit, icompz;
42 double ssfmax;
45 --d__;
46 --e;
47 z_dim1 = *ldz;
48 z_offset = 1 + z_dim1;
49 z__ -= z_offset;
50 --work;
52 *info = 0;
54 if (*compz=='N' || *compz=='n') {
55 icompz = 0;
56 } else if (*compz=='V' || *compz=='v') {
57 icompz = 1;
58 } else if (*compz=='I' || *compz=='i') {
59 icompz = 2;
60 } else {
61 icompz = -1;
63 if (icompz < 0) {
64 *info = -1;
65 } else if (*n < 0) {
66 *info = -2;
67 } else if (*ldz < 1 || (icompz > 0 && *ldz < ((*n>1) ? *n : 1))) {
68 *info = -6;
70 if (*info != 0) {
71 return;
75 if (*n == 0) {
76 return;
79 if (*n == 1) {
80 if (icompz == 2) {
81 z__[z_dim1 + 1] = 1.;
83 return;
86 eps = GMX_DOUBLE_EPS;
87 d__1 = eps;
88 eps2 = d__1 * d__1;
89 minval = GMX_DOUBLE_MIN;
90 safmin = minval*(1.0+GMX_DOUBLE_EPS);
92 safmax = 1. / safmin;
93 ssfmax = sqrt(safmax) / 3.;
94 ssfmin = sqrt(safmin) / eps2;
96 if (icompz == 2) {
97 F77_FUNC(dlaset,DLASET)("Full", n, n, &c_b9, &c_b10, &z__[z_offset], ldz);
100 nmaxit = *n * 30;
101 jtot = 0;
103 l1 = 1;
104 nm1 = *n - 1;
106 L10:
107 if (l1 > *n) {
108 goto L160;
110 if (l1 > 1) {
111 e[l1 - 1] = 0.;
113 if (l1 <= nm1) {
114 i__1 = nm1;
115 for (m = l1; m <= i__1; ++m) {
116 tst = fabs(e[m]);
117 if (fabs(tst)<GMX_DOUBLE_MIN) {
118 goto L30;
120 if (tst <= sqrt(fabs(d__[m])) * sqrt(fabs(d__[m + 1])) * eps) {
121 e[m] = 0.;
122 goto L30;
126 m = *n;
128 L30:
129 l = l1;
130 lsv = l;
131 lend = m;
132 lendsv = lend;
133 l1 = m + 1;
134 if (lend == l) {
135 goto L10;
138 i__1 = lend - l + 1;
139 anorm = F77_FUNC(dlanst,DLANST)("I", &i__1, &d__[l], &e[l]);
140 iscale = 0;
141 if (fabs(anorm)<GMX_DOUBLE_MIN) {
142 goto L10;
144 if (anorm > ssfmax) {
145 iscale = 1;
146 i__1 = lend - l + 1;
147 F77_FUNC(dlascl,DLASCL)("G", &c__0, &c__0, &anorm, &ssfmax, &i__1, &c__1, &d__[l], n,
148 info);
149 i__1 = lend - l;
150 F77_FUNC(dlascl,DLASCL)("G", &c__0, &c__0, &anorm, &ssfmax, &i__1, &c__1, &e[l], n,
151 info);
152 } else if (anorm < ssfmin) {
153 iscale = 2;
154 i__1 = lend - l + 1;
155 F77_FUNC(dlascl,DLASCL)("G", &c__0, &c__0, &anorm, &ssfmin, &i__1, &c__1, &d__[l], n,
156 info);
157 i__1 = lend - l;
158 F77_FUNC(dlascl,DLASCL)("G", &c__0, &c__0, &anorm, &ssfmin, &i__1, &c__1, &e[l], n,
159 info);
162 if (fabs(d__[lend]) < fabs(d__[l])) {
163 lend = lsv;
164 l = lendsv;
167 if (lend > l) {
169 L40:
170 if (l != lend) {
171 lendm1 = lend - 1;
172 i__1 = lendm1;
173 for (m = l; m <= i__1; ++m) {
174 d__2 = fabs(e[m]);
175 tst = d__2 * d__2;
176 if (tst <= eps2 * fabs(d__[m]) * fabs(d__[m+ 1]) + safmin) {
177 goto L60;
182 m = lend;
184 L60:
185 if (m < lend) {
186 e[m] = 0.;
188 p = d__[l];
189 if (m == l) {
190 goto L80;
193 if (m == l + 1) {
194 if (icompz > 0) {
195 F77_FUNC(dlaev2,DLAEV2)(&d__[l], &e[l], &d__[l + 1], &rt1, &rt2, &c__, &s);
196 work[l] = c__;
197 work[*n - 1 + l] = s;
198 F77_FUNC(dlasr,DLASR)("R", "V", "B", n, &c__2, &work[l], &work[*n - 1 + l], &
199 z__[l * z_dim1 + 1], ldz);
200 } else {
201 F77_FUNC(dlae2,DLAE2)(&d__[l], &e[l], &d__[l + 1], &rt1, &rt2);
203 d__[l] = rt1;
204 d__[l + 1] = rt2;
205 e[l] = 0.;
206 l += 2;
207 if (l <= lend) {
208 goto L40;
210 goto L140;
213 if (jtot == nmaxit) {
214 goto L140;
216 ++jtot;
218 g = (d__[l + 1] - p) / (e[l] * 2.);
219 r__ = F77_FUNC(dlapy2,DLAPY2)(&g, &c_b10);
220 g = d__[m] - p + e[l] / (g + ( (g>0) ? r__ : -r__ ) );
222 s = 1.;
223 c__ = 1.;
224 p = 0.;
226 mm1 = m - 1;
227 i__1 = l;
228 for (i__ = mm1; i__ >= i__1; --i__) {
229 f = s * e[i__];
230 b = c__ * e[i__];
231 F77_FUNC(dlartg,DLARTG)(&g, &f, &c__, &s, &r__);
232 if (i__ != m - 1) {
233 e[i__ + 1] = r__;
235 g = d__[i__ + 1] - p;
236 r__ = (d__[i__] - g) * s + c__ * 2. * b;
237 p = s * r__;
238 d__[i__ + 1] = g + p;
239 g = c__ * r__ - b;
241 if (icompz > 0) {
242 work[i__] = c__;
243 work[*n - 1 + i__] = -s;
247 if (icompz > 0) {
248 mm = m - l + 1;
249 F77_FUNC(dlasr,DLASR)("R", "V", "B", n, &mm, &work[l], &work[*n - 1 + l], &z__[l
250 * z_dim1 + 1], ldz);
253 d__[l] -= p;
254 e[l] = g;
255 goto L40;
257 L80:
258 d__[l] = p;
260 ++l;
261 if (l <= lend) {
262 goto L40;
264 goto L140;
266 } else {
268 L90:
269 if (l != lend) {
270 lendp1 = lend + 1;
271 i__1 = lendp1;
272 for (m = l; m >= i__1; --m) {
273 d__2 = fabs(e[m - 1]);
274 tst = d__2 * d__2;
275 if (tst <= eps2 * fabs(d__[m]) * fabs(d__[m- 1]) + safmin) {
276 goto L110;
281 m = lend;
283 L110:
284 if (m > lend) {
285 e[m - 1] = 0.;
287 p = d__[l];
288 if (m == l) {
289 goto L130;
291 if (m == l - 1) {
292 if (icompz > 0) {
293 F77_FUNC(dlaev2,DLAEV2)(&d__[l - 1], &e[l - 1], &d__[l], &rt1, &rt2, &c__, &s)
295 work[m] = c__;
296 work[*n - 1 + m] = s;
297 F77_FUNC(dlasr,DLASR)("R", "V", "F", n, &c__2, &work[m], &work[*n - 1 + m], &
298 z__[(l - 1) * z_dim1 + 1], ldz);
299 } else {
300 F77_FUNC(dlae2,DLAE2)(&d__[l - 1], &e[l - 1], &d__[l], &rt1, &rt2);
302 d__[l - 1] = rt1;
303 d__[l] = rt2;
304 e[l - 1] = 0.;
305 l += -2;
306 if (l >= lend) {
307 goto L90;
309 goto L140;
312 if (jtot == nmaxit) {
313 goto L140;
315 ++jtot;
317 g = (d__[l - 1] - p) / (e[l - 1] * 2.);
318 r__ = F77_FUNC(dlapy2,DLAPY2)(&g, &c_b10);
319 g = d__[m] - p + e[l - 1] / (g + ( (g>0) ? r__ : -r__ ));
321 s = 1.;
322 c__ = 1.;
323 p = 0.;
325 lm1 = l - 1;
326 i__1 = lm1;
327 for (i__ = m; i__ <= i__1; ++i__) {
328 f = s * e[i__];
329 b = c__ * e[i__];
330 F77_FUNC(dlartg,DLARTG)(&g, &f, &c__, &s, &r__);
331 if (i__ != m) {
332 e[i__ - 1] = r__;
334 g = d__[i__] - p;
335 r__ = (d__[i__ + 1] - g) * s + c__ * 2. * b;
336 p = s * r__;
337 d__[i__] = g + p;
338 g = c__ * r__ - b;
340 if (icompz > 0) {
341 work[i__] = c__;
342 work[*n - 1 + i__] = s;
346 if (icompz > 0) {
347 mm = l - m + 1;
348 F77_FUNC(dlasr,DLASR)("R", "V", "F", n, &mm, &work[m], &work[*n - 1 + m], &z__[m
349 * z_dim1 + 1], ldz);
352 d__[l] -= p;
353 e[lm1] = g;
354 goto L90;
356 L130:
357 d__[l] = p;
359 --l;
360 if (l >= lend) {
361 goto L90;
363 goto L140;
367 L140:
368 if (iscale == 1) {
369 i__1 = lendsv - lsv + 1;
370 F77_FUNC(dlascl,DLASCL)("G", &c__0, &c__0, &ssfmax, &anorm, &i__1, &c__1, &d__[lsv],
371 n, info);
372 i__1 = lendsv - lsv;
373 F77_FUNC(dlascl,DLASCL)("G", &c__0, &c__0, &ssfmax, &anorm, &i__1, &c__1, &e[lsv], n,
374 info);
375 } else if (iscale == 2) {
376 i__1 = lendsv - lsv + 1;
377 F77_FUNC(dlascl,DLASCL)("G", &c__0, &c__0, &ssfmin, &anorm, &i__1, &c__1, &d__[lsv],
378 n, info);
379 i__1 = lendsv - lsv;
380 F77_FUNC(dlascl,DLASCL)("G", &c__0, &c__0, &ssfmin, &anorm, &i__1, &c__1, &e[lsv], n,
381 info);
384 if (jtot < nmaxit) {
385 goto L10;
387 i__1 = *n - 1;
388 for (i__ = 1; i__ <= i__1; ++i__) {
389 if (fabs(e[i__])>GMX_DOUBLE_MIN) {
390 ++(*info);
393 goto L190;
395 L160:
396 if (icompz == 0) {
398 F77_FUNC(dlasrt,DLASRT)("I", n, &d__[1], info);
400 } else {
402 i__1 = *n;
403 for (ii = 2; ii <= i__1; ++ii) {
404 i__ = ii - 1;
405 k = i__;
406 p = d__[i__];
407 i__2 = *n;
408 for (j = ii; j <= i__2; ++j) {
409 if (d__[j] < p) {
410 k = j;
411 p = d__[j];
414 if (k != i__) {
415 d__[k] = d__[i__];
416 d__[i__] = p;
417 F77_FUNC(dswap,DSWAP)(n, &z__[i__ * z_dim1 + 1], &c__1, &z__[k * z_dim1 + 1],
418 &c__1);
423 L190:
424 return;