1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: m; c-basic-offset: 4 -*-
4 * This file is part of Gromacs Copyright (c) 1991-2004
5 * David van der Spoel, Erik Lindahl, University of Groningen.
7 * This file contains a subset of ARPACK functions to perform
8 * diagonalization and SVD for sparse matrices in Gromacs.
10 * The code has been translated to C to avoid being dependent on
11 * a Fotran compiler, and it has been made threadsafe by using
12 * additional workspace arrays to store data during reverse communication.
14 * You might prefer the original ARPACK library for general use, but
15 * in case you want to this version can be redistributed freely, just
16 * as the original library. However, please make clear that it is the
17 * hacked version from Gromacs so any bugs are blamed on us and not
18 * the original authors. You should also be aware that the double
19 * precision work array workd needs to be of size (3*N+4) here
20 * (4 more than the general library), and there is an extra argument
21 * iwork, which should be an integer work array of length 80.
23 * ARPACK was written by
25 * Danny Sorensen Phuong Vu
26 * Rihard Lehoucq CRPC / Rice University
27 * Dept. of Computational & Houston, Texas
39 #include "types/simple.h"
40 #include "gmx_lapack.h"
44 /* Default Fortran name mangling */
46 #define F77_FUNC(name,NAME) name ## _
52 F77_FUNC(dstqrb
,DSTQRB
)(int * n
,
68 int l1
, ii
, mm
, lm1
, mm1
, nm1
;
72 int lend
, jtot
, lendm1
, lendp1
, iscale
;
74 int lendsv
, nmaxit
, icompz
;
75 double ssfmax
, ssfmin
,safmin
,minval
,safmax
,anorm
;
102 minval
= GMX_DOUBLE_MIN
;
103 safmin
= minval
/ GMX_DOUBLE_EPS
;
104 safmax
= 1. / safmin
;
105 ssfmax
= sqrt(safmax
) / 3.;
106 ssfmin
= sqrt(safmin
) / eps2
;
110 for (j
= 1; j
<= i__1
; ++j
) {
132 for (m
= l1
; m
<= i__1
; ++m
) {
137 if (tst
<= sqrt(fabs(d__
[m
])) * sqrt(fabs(d__
[m
+1])) * eps
) {
156 anorm
=F77_FUNC(dlanst
,DLANST
)("i", &i__1
, &d__
[l
], &e
[l
]);
161 if (anorm
> ssfmax
) {
164 F77_FUNC(dlascl
,DLASCL
)("g", &c__0
, &c__0
, &anorm
, &ssfmax
, &i__1
, &c__1
, &d__
[l
], n
,
167 F77_FUNC(dlascl
,DLASCL
)("g", &c__0
, &c__0
, &anorm
, &ssfmax
, &i__1
, &c__1
, &e
[l
], n
,
169 } else if (anorm
< ssfmin
) {
172 F77_FUNC(dlascl
,DLASCL
)("g", &c__0
, &c__0
, &anorm
, &ssfmin
, &i__1
, &c__1
, &d__
[l
], n
,
175 F77_FUNC(dlascl
,DLASCL
)("g", &c__0
, &c__0
, &anorm
, &ssfmin
, &i__1
, &c__1
, &e
[l
], n
,
179 if (fabs(d__
[lend
]) < fabs(d__
[l
])) {
190 for (m
= l
; m
<= i__1
; ++m
) {
193 if (tst
<= eps2
* fabs(d__
[m
]) * fabs(d__
[m
+ 1]) + safmin
) {
212 F77_FUNC(dlaev2
,DLAEV2
)(&d__
[l
], &e
[l
], &d__
[l
+ 1], &rt1
, &rt2
, &c__
, &s
);
214 work
[*n
- 1 + l
] = s
;
217 z__
[l
+ 1] = c__
* tst
- s
* z__
[l
];
218 z__
[l
] = s
* tst
+ c__
* z__
[l
];
220 F77_FUNC(dlae2
,DLAE2
)(&d__
[l
], &e
[l
], &d__
[l
+ 1], &rt1
, &rt2
);
232 if (jtot
== nmaxit
) {
237 g
= (d__
[l
+ 1] - p
) / (e
[l
] * 2.);
238 r__
=F77_FUNC(dlapy2
,DLAPY2
)(&g
, &c_b31
);
239 g
= d__
[m
] - p
+ e
[l
] / (g
+ ((g
>0) ? r__
: -r__
));
247 for (i__
= mm1
; i__
>= i__1
; --i__
) {
250 F77_FUNC(dlartg
,DLARTG
)(&g
, &f
, &c__
, &s
, &r__
);
254 g
= d__
[i__
+ 1] - p
;
255 r__
= (d__
[i__
] - g
) * s
+ c__
* 2. * b
;
257 d__
[i__
+ 1] = g
+ p
;
262 work
[*n
- 1 + i__
] = -s
;
270 F77_FUNC(dlasr
,DLASR
)("r", "v", "b", &c__1
, &mm
, &work
[l
], &work
[*n
- 1 + l
], &
293 for (m
= l
; m
>= i__1
; --m
) {
294 d__2
= fabs(e
[m
- 1]);
296 if (tst
<= eps2
* fabs(d__
[m
]) * fabs(d__
[m
- 1]) + safmin
) {
315 F77_FUNC(dlaev2
,DLAEV2
)(&d__
[l
- 1], &e
[l
- 1], &d__
[l
], &rt1
, &rt2
, &c__
, &s
)
319 z__
[l
] = c__
* tst
- s
* z__
[l
- 1];
320 z__
[l
- 1] = s
* tst
+ c__
* z__
[l
- 1];
322 F77_FUNC(dlae2
,DLAE2
)(&d__
[l
- 1], &e
[l
- 1], &d__
[l
], &rt1
, &rt2
);
334 if (jtot
== nmaxit
) {
340 g
= (d__
[l
- 1] - p
) / (e
[l
- 1] * 2.);
341 r__
=F77_FUNC(dlapy2
,DLAPY2
)(&g
, &c_b31
);
342 g
= d__
[m
] - p
+ e
[l
- 1] / (g
+ ((g
>0) ? r__
: -r__
));
350 for (i__
= m
; i__
<= i__1
; ++i__
) {
353 F77_FUNC(dlartg
,DLARTG
)(&g
, &f
, &c__
, &s
, &r__
);
358 r__
= (d__
[i__
+ 1] - g
) * s
+ c__
* 2. * b
;
365 work
[*n
- 1 + i__
] = s
;
373 F77_FUNC(dlasr
,DLASR
)("r", "v", "f", &c__1
, &mm
, &work
[m
], &work
[*n
- 1 + m
], &
394 i__1
= lendsv
- lsv
+ 1;
395 F77_FUNC(dlascl
,DLASCL
)("g", &c__0
, &c__0
, &ssfmax
, &anorm
, &i__1
, &c__1
, &d__
[lsv
],
398 F77_FUNC(dlascl
,DLASCL
)("g", &c__0
, &c__0
, &ssfmax
, &anorm
, &i__1
, &c__1
, &e
[lsv
], n
,
400 } else if (iscale
== 2) {
401 i__1
= lendsv
- lsv
+ 1;
402 F77_FUNC(dlascl
,DLASCL
)("g", &c__0
, &c__0
, &ssfmin
, &anorm
, &i__1
, &c__1
, &d__
[lsv
],
405 F77_FUNC(dlascl
,DLASCL
)("g", &c__0
, &c__0
, &ssfmin
, &anorm
, &i__1
, &c__1
, &e
[lsv
], n
,
413 for (i__
= 1; i__
<= i__1
; ++i__
) {
423 F77_FUNC(dlasrt
,DLASRT
)("i", n
, &d__
[1], info
);
428 for (ii
= 2; ii
<= i__1
; ++ii
) {
433 for (j
= ii
; j
<= i__2
; ++j
) {
456 F77_FUNC(dgetv0
,DGETV0
)(int * ido
,
475 int v_dim1
, v_offset
, i__1
;
483 v_offset
= 1 + v_dim1
;
497 F77_FUNC(dlarnv
,DLARNV
)(&idist
, &iwork
[1], n
, &resid
[1]);
503 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
519 F77_FUNC(dcopy
,DCOPY
)(n
, &workd
[*n
+ 1], &c__1
, &resid
[1], &c__1
);
524 } else if (*bmat
== 'I') {
525 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
533 workd
[*n
* 3 + 4] =F77_FUNC(ddot
,DDOT
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
534 workd
[*n
* 3 + 4] = sqrt(fabs(workd
[*n
* 3 + 4]));
535 } else if (*bmat
== 'I') {
536 workd
[*n
* 3 + 4] =F77_FUNC(dnrm2
,DNRM2
)(n
, &resid
[1], &c__1
);
538 *rnorm
= workd
[*n
* 3 + 4];
547 F77_FUNC(dgemv
,DGEMV
)("T", n
, &i__1
, &c_b22
, &v
[v_offset
], ldv
, &workd
[1], &c__1
, &c_b24
,
548 &workd
[*n
+ 1], &c__1
);
550 F77_FUNC(dgemv
,DGEMV
)("N", n
, &i__1
, &c_b27
, &v
[v_offset
], ldv
, &workd
[*n
+ 1], &c__1
, &
551 c_b22
, &resid
[1], &c__1
);
554 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[*n
+ 1], &c__1
);
559 } else if (*bmat
== 'I') {
560 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
566 *rnorm
=F77_FUNC(ddot
,DDOT
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
567 *rnorm
= sqrt(fabs(*rnorm
));
568 } else if (*bmat
== 'I') {
569 *rnorm
=F77_FUNC(dnrm2
,DNRM2
)(n
, &resid
[1], &c__1
);
572 if (*rnorm
> workd
[*n
* 3 + 4] * .717f
) {
579 workd
[*n
* 3 + 4] = *rnorm
;
584 for (jj
= 1; jj
<= i__1
; ++jj
) {
604 F77_FUNC(dsapps
,DSAPPS
)(int * n
,
621 int h_dim1
, h_offset
, q_dim1
, q_offset
, v_dim1
, v_offset
, i__1
, i__2
,
625 double r__
, s
, a1
, a2
, a3
, a4
;
636 v_offset
= 1 + v_dim1
;
639 h_offset
= 1 + h_dim1
;
642 q_offset
= 1 + q_dim1
;
645 epsmch
= GMX_DOUBLE_EPS
;
651 F77_FUNC(dlaset
,DLASET
)("All", &kplusp
, &kplusp
, &c_b4
, &c_b5
, &q
[q_offset
], ldq
);
658 for (jj
= 1; jj
<= i__1
; ++jj
) {
665 for (i__
= istart
; i__
<= i__2
; ++i__
) {
666 big
= fabs(h__
[i__
+ (h_dim1
*2)]) + fabs(h__
[i__
+ 1 + (h_dim1
*2)]);
667 if (h__
[i__
+ 1 + h_dim1
] <= epsmch
* big
) {
668 h__
[i__
+ 1 + h_dim1
] = 0.;
678 f
= h__
[istart
+ (h_dim1
<< 1)] - shift
[jj
];
679 g
= h__
[istart
+ 1 + h_dim1
];
680 F77_FUNC(dlartg
,DLARTG
)(&f
, &g
, &c__
, &s
, &r__
);
682 a1
= c__
* h__
[istart
+ (h_dim1
<< 1)] + s
* h__
[istart
+ 1 +
684 a2
= c__
* h__
[istart
+ 1 + h_dim1
] + s
* h__
[istart
+ 1 + (
686 a4
= c__
* h__
[istart
+ 1 + (h_dim1
<< 1)] - s
* h__
[istart
+ 1 +
688 a3
= c__
* h__
[istart
+ 1 + h_dim1
] - s
* h__
[istart
+ (h_dim1
<<
690 h__
[istart
+ (h_dim1
<< 1)] = c__
* a1
+ s
* a2
;
691 h__
[istart
+ 1 + (h_dim1
<< 1)] = c__
* a4
- s
* a3
;
692 h__
[istart
+ 1 + h_dim1
] = c__
* a3
+ s
* a4
;
695 i__2
= (i__3
<kplusp
) ? i__3
: kplusp
;
696 for (j
= 1; j
<= i__2
; ++j
) {
697 a1
= c__
* q
[j
+ istart
* q_dim1
] + s
* q
[j
+ (istart
+ 1) *
699 q
[j
+ (istart
+ 1) * q_dim1
] = -s
* q
[j
+ istart
* q_dim1
] +
700 c__
* q
[j
+ (istart
+ 1) * q_dim1
];
701 q
[j
+ istart
* q_dim1
] = a1
;
706 for (i__
= istart
+ 1; i__
<= i__2
; ++i__
) {
708 f
= h__
[i__
+ h_dim1
];
709 g
= s
* h__
[i__
+ 1 + h_dim1
];
711 h__
[i__
+ 1 + h_dim1
] = c__
* h__
[i__
+ 1 + h_dim1
];
712 F77_FUNC(dlartg
,DLARTG
)(&f
, &g
, &c__
, &s
, &r__
);
720 h__
[i__
+ h_dim1
] = r__
;
722 a1
= c__
* h__
[i__
+ (h_dim1
<< 1)] + s
* h__
[i__
+ 1 +
724 a2
= c__
* h__
[i__
+ 1 + h_dim1
] + s
* h__
[i__
+ 1 + (h_dim1
726 a3
= c__
* h__
[i__
+ 1 + h_dim1
] - s
* h__
[i__
+ (h_dim1
<< 1)
728 a4
= c__
* h__
[i__
+ 1 + (h_dim1
<< 1)] - s
* h__
[i__
+ 1 +
731 h__
[i__
+ (h_dim1
<< 1)] = c__
* a1
+ s
* a2
;
732 h__
[i__
+ 1 + (h_dim1
<< 1)] = c__
* a4
- s
* a3
;
733 h__
[i__
+ 1 + h_dim1
] = c__
* a3
+ s
* a4
;
736 i__3
= (i__4
<kplusp
) ? i__4
: kplusp
;
737 for (j
= 1; j
<= i__3
; ++j
) {
738 a1
= c__
* q
[j
+ i__
* q_dim1
] + s
* q
[j
+ (i__
+ 1) *
740 q
[j
+ (i__
+ 1) * q_dim1
] = -s
* q
[j
+ i__
* q_dim1
] +
741 c__
* q
[j
+ (i__
+ 1) * q_dim1
];
742 q
[j
+ i__
* q_dim1
] = a1
;
751 if (h__
[iend
+ h_dim1
] < 0.) {
752 h__
[iend
+ h_dim1
] = -h__
[iend
+ h_dim1
];
753 F77_FUNC(dscal
,DSCAL
)(&kplusp
, &c_b14
, &q
[iend
* q_dim1
+ 1], &c__1
);
761 for (i__
= itop
; i__
<= i__2
; ++i__
) {
762 if (h__
[i__
+ 1 + h_dim1
] > 0.) {
773 for (i__
= itop
; i__
<= i__1
; ++i__
) {
774 big
= fabs(h__
[i__
+ (h_dim1
*2)]) + fabs(h__
[i__
+ 1 + (h_dim1
*2)]);
775 if (h__
[i__
+ 1 + h_dim1
] <= epsmch
* big
) {
776 h__
[i__
+ 1 + h_dim1
] = 0.;
781 if (h__
[*kev
+ 1 + h_dim1
] > 0.) {
782 F77_FUNC(dgemv
,DGEMV
)("N", n
, &kplusp
, &c_b5
, &v
[v_offset
], ldv
, &q
[(*kev
+ 1) *
783 q_dim1
+ 1], &c__1
, &c_b4
, &workd
[*n
+ 1], &c__1
);
787 for (i__
= 1; i__
<= i__1
; ++i__
) {
788 i__2
= kplusp
- i__
+ 1;
789 F77_FUNC(dgemv
,DGEMV
)("N", n
, &i__2
, &c_b5
, &v
[v_offset
], ldv
, &q
[(*kev
- i__
+ 1) *
790 q_dim1
+ 1], &c__1
, &c_b4
, &workd
[1], &c__1
);
791 F77_FUNC(dcopy
,DCOPY
)(n
, &workd
[1], &c__1
, &v
[(kplusp
- i__
+ 1) * v_dim1
+ 1], &
796 F77_FUNC(dlacpy
,DLACPY
)("All", n
, kev
, &v
[(*np
+ 1) * v_dim1
+ 1], ldv
, &v
[v_offset
], ldv
);
798 if (h__
[*kev
+ 1 + h_dim1
] > 0.) {
799 F77_FUNC(dcopy
,DCOPY
)(n
, &workd
[*n
+ 1], &c__1
, &v
[(*kev
+ 1) * v_dim1
+ 1], &c__1
);
802 F77_FUNC(dscal
,DSCAL
)(n
, &q
[kplusp
+ *kev
* q_dim1
], &resid
[1], &c__1
);
803 if (h__
[*kev
+ 1 + h_dim1
] > 0.) {
804 F77_FUNC(daxpy
,DAXPY
)(n
, &h__
[*kev
+ 1 + h_dim1
], &v
[(*kev
+ 1) * v_dim1
+ 1], &c__1
,
818 F77_FUNC(dsortr
,DSORTR
)(const char * which
,
833 if (!strncmp(which
, "SA", 2)) {
840 for (i__
= igap
; i__
<= i__1
; ++i__
) {
848 if (x1
[j
] < x1
[j
+ igap
]) {
850 x1
[j
] = x1
[j
+ igap
];
854 x2
[j
] = x2
[j
+ igap
];
868 } else if (!strncmp(which
, "SM", 2)) {
875 for (i__
= igap
; i__
<= i__1
; ++i__
) {
883 if (fabs(x1
[j
]) < fabs(x1
[j
+ igap
])) {
885 x1
[j
] = x1
[j
+ igap
];
889 x2
[j
] = x2
[j
+ igap
];
903 } else if (!strncmp(which
, "LA", 2)) {
910 for (i__
= igap
; i__
<= i__1
; ++i__
) {
918 if (x1
[j
] > x1
[j
+ igap
]) {
920 x1
[j
] = x1
[j
+ igap
];
924 x2
[j
] = x2
[j
+ igap
];
938 } else if (!strncmp(which
, "LM", 2)) {
946 for (i__
= igap
; i__
<= i__1
; ++i__
) {
954 if (fabs(x1
[j
]) > fabs(x1
[j
+ igap
])) {
956 x1
[j
] = x1
[j
+ igap
];
960 x2
[j
] = x2
[j
+ igap
];
984 F77_FUNC(dsesrt
,DSESRT
)(const char * which
,
992 int a_dim1
, a_offset
, i__1
;
999 a_offset
= 1 + a_dim1
* 0;
1004 if (!strncmp(which
, "SA", 2)) {
1011 for (i__
= igap
; i__
<= i__1
; ++i__
) {
1019 if (x
[j
] < x
[j
+ igap
]) {
1024 F77_FUNC(dswap
,DSWAP
)(na
, &a
[j
* a_dim1
+ 1], &c__1
, &a
[(j
+ igap
) *
1025 a_dim1
+ 1], &c__1
);
1038 } else if (!strncmp(which
, "SM", 2)) {
1045 for (i__
= igap
; i__
<= i__1
; ++i__
) {
1053 if (fabs(x
[j
]) < fabs(x
[j
+ igap
])) {
1058 F77_FUNC(dswap
,DSWAP
)(na
, &a
[j
* a_dim1
+ 1], &c__1
, &a
[(j
+ igap
) *
1059 a_dim1
+ 1], &c__1
);
1072 } else if (!strncmp(which
, "LA", 2)) {
1079 for (i__
= igap
; i__
<= i__1
; ++i__
) {
1087 if (x
[j
] > x
[j
+ igap
]) {
1092 F77_FUNC(dswap
,DSWAP
)(na
, &a
[j
* a_dim1
+ 1], &c__1
, &a
[(j
+ igap
) *
1093 a_dim1
+ 1], &c__1
);
1106 } else if (!strncmp(which
, "LM", 2)) {
1113 for (i__
= igap
; i__
<= i__1
; ++i__
) {
1121 if (fabs(x
[j
]) > fabs(x
[j
+ igap
])) {
1126 F77_FUNC(dswap
,DSWAP
)(na
, &a
[j
* a_dim1
+ 1], &c__1
, &a
[(j
+ igap
) *
1127 a_dim1
+ 1], &c__1
);
1150 F77_FUNC(dsgets
,DSGETS
)(int * ishift
,
1166 if (!strncmp(which
, "BE", 2)) {
1168 F77_FUNC(dsortr
,DSORTR
)("LA", &c__1
, &i__1
, &ritz
[1], &bounds
[1]);
1171 i__1
= (kevd2
<*np
) ? kevd2
: *np
;
1172 i__2
= (kevd2
>*np
) ? kevd2
: *np
;
1173 F77_FUNC(dswap
,DSWAP
)(&i__1
, &ritz
[1], &c__1
,
1174 &ritz
[i__2
+ 1], &c__1
);
1175 i__1
= (kevd2
<*np
) ? kevd2
: *np
;
1176 i__2
= (kevd2
>*np
) ? kevd2
: *np
;
1177 F77_FUNC(dswap
,DSWAP
)(&i__1
, &bounds
[1], &c__1
,
1178 &bounds
[i__2
+ 1], &c__1
);
1183 F77_FUNC(dsortr
,DSORTR
)(which
, &c__1
, &i__1
, &ritz
[1], &bounds
[1]);
1186 if (*ishift
== 1 && *np
> 0) {
1188 F77_FUNC(dsortr
,DSORTR
)("SM", &c__1
, np
, &bounds
[1], &ritz
[1]);
1189 F77_FUNC(dcopy
,DCOPY
)(np
, &ritz
[1], &c__1
, &shifts
[1], &c__1
);
1199 F77_FUNC(dsconv
,DSCONV
)(int * n
,
1215 eps23
= GMX_DOUBLE_EPS
;
1216 eps23
= pow(eps23
, c_b3
);
1220 for (i__
= 1; i__
<= i__1
; ++i__
) {
1223 d__3
= fabs(ritz
[i__
]);
1224 temp
= (d__2
> d__3
) ? d__2
: d__3
;
1225 if (bounds
[i__
] <= *tol
* temp
) {
1235 F77_FUNC(dseigt
,DSEIGT
)(double * rnorm
,
1245 int h_dim1
, h_offset
, i__1
;
1254 h_offset
= 1 + h_dim1
;
1257 F77_FUNC(dcopy
,DCOPY
)(n
, &h__
[(h_dim1
<< 1) + 1], &c__1
, &eig
[1], &c__1
);
1259 F77_FUNC(dcopy
,DCOPY
)(&i__1
, &h__
[h_dim1
+ 2], &c__1
, &workl
[1], &c__1
);
1260 F77_FUNC(dstqrb
,DSTQRB
)(n
, &eig
[1], &workl
[1], &bounds
[1], &workl
[*n
+ 1], ierr
);
1266 for (k
= 1; k
<= i__1
; ++k
) {
1267 bounds
[k
] = *rnorm
* fabs(bounds
[k
]);
1280 F77_FUNC(dsaitr
,DSAITR
)(int * ido
,
1304 int h_dim1
, h_offset
, v_dim1
, v_offset
, i__1
;
1308 double safmin
,minval
;
1314 v_offset
= 1 + v_dim1
;
1317 h_offset
= 1 + h_dim1
;
1321 minval
= GMX_DOUBLE_MIN
;
1322 safmin
= minval
/ GMX_DOUBLE_EPS
;
1335 iwork
[9] = iwork
[8] + *n
;
1336 iwork
[10] = iwork
[9] + *n
;
1339 if (iwork
[5] == 1) {
1342 if (iwork
[6] == 1) {
1345 if (iwork
[2] == 1) {
1348 if (iwork
[3] == 1) {
1351 if (iwork
[4] == 1) {
1368 F77_FUNC(dgetv0
,DGETV0
)(ido
, bmat
, &iwork
[11], &c__0
, n
, &iwork
[12], &v
[v_offset
], ldv
,
1369 &resid
[1], rnorm
, &ipntr
[1], &workd
[1], &iwork
[21], &iwork
[7]);
1375 if (iwork
[11] <= 3) {
1379 *info
= iwork
[12] - 1;
1386 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &v
[iwork
[12] * v_dim1
+ 1], &c__1
);
1387 if (*rnorm
>= safmin
) {
1388 temp1
= 1. / *rnorm
;
1389 F77_FUNC(dscal
,DSCAL
)(n
, &temp1
, &v
[iwork
[12] * v_dim1
+ 1], &c__1
);
1390 F77_FUNC(dscal
,DSCAL
)(n
, &temp1
, &workd
[iwork
[8]], &c__1
);
1393 F77_FUNC(dlascl
,DLASCL
)("General", &i__
, &i__
, rnorm
, &c_b18
, n
, &c__1
, &v
[iwork
[12] *
1394 v_dim1
+ 1], n
, &infol
);
1395 F77_FUNC(dlascl
,DLASCL
)("General", &i__
, &i__
, rnorm
, &c_b18
, n
, &c__1
, &workd
[iwork
[
1400 F77_FUNC(dcopy
,DCOPY
)(n
, &v
[iwork
[12] * v_dim1
+ 1], &c__1
, &workd
[iwork
[10]], &c__1
);
1401 ipntr
[1] = iwork
[10];
1402 ipntr
[2] = iwork
[9];
1403 ipntr
[3] = iwork
[8];
1412 F77_FUNC(dcopy
,DCOPY
)(n
, &workd
[iwork
[9]], &c__1
, &resid
[1], &c__1
);
1419 ipntr
[1] = iwork
[9];
1420 ipntr
[2] = iwork
[8];
1424 } else if (*bmat
== 'I') {
1425 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &c__1
);
1434 workd
[*n
* 3 + 3] =F77_FUNC(ddot
,DDOT
)(n
, &resid
[1], &c__1
, &workd
[iwork
[10]], &
1436 workd
[*n
* 3 + 3] = sqrt(fabs(workd
[*n
* 3 + 3]));
1437 } else if (*bmat
== 'G') {
1438 workd
[*n
* 3 + 3] =F77_FUNC(ddot
,DDOT
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &
1440 workd
[*n
* 3 + 3] = sqrt(fabs(workd
[*n
* 3 + 3]));
1441 } else if (*bmat
== 'I') {
1442 workd
[*n
* 3 + 3] =F77_FUNC(dnrm2
,DNRM2
)(n
, &resid
[1], &c__1
);
1446 F77_FUNC(dgemv
,DGEMV
)("T", n
, &iwork
[12], &c_b18
, &v
[v_offset
], ldv
, &workd
[iwork
[8]]
1447 , &c__1
, &c_b42
, &workd
[iwork
[9]], &c__1
);
1448 } else if (*mode
== 2) {
1449 F77_FUNC(dgemv
,DGEMV
)("T", n
, &iwork
[12], &c_b18
, &v
[v_offset
], ldv
, &workd
[iwork
[10]
1450 ], &c__1
, &c_b42
, &workd
[iwork
[9]], &c__1
);
1453 F77_FUNC(dgemv
,DGEMV
)("N", n
, &iwork
[12], &c_b50
, &v
[v_offset
], ldv
, &workd
[iwork
[9]], &
1454 c__1
, &c_b18
, &resid
[1], &c__1
);
1456 h__
[iwork
[12] + (h_dim1
<< 1)] = workd
[iwork
[9] + iwork
[12] - 1];
1457 if (iwork
[12] == 1 || iwork
[4] == 1) {
1458 h__
[iwork
[12] + h_dim1
] = 0.;
1460 h__
[iwork
[12] + h_dim1
] = *rnorm
;
1467 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[9]], &c__1
);
1468 ipntr
[1] = iwork
[9];
1469 ipntr
[2] = iwork
[8];
1473 } else if (*bmat
== 'I') {
1474 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &c__1
);
1481 *rnorm
=F77_FUNC(ddot
,DDOT
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &c__1
);
1482 *rnorm
= sqrt(fabs(*rnorm
));
1483 } else if (*bmat
== 'I') {
1484 *rnorm
=F77_FUNC(dnrm2
,DNRM2
)(n
, &resid
[1], &c__1
);
1487 if (*rnorm
> workd
[*n
* 3 + 3] * .717f
) {
1493 F77_FUNC(dgemv
,DGEMV
)("T", n
, &iwork
[12], &c_b18
, &v
[v_offset
], ldv
, &workd
[iwork
[8]], &
1494 c__1
, &c_b42
, &workd
[iwork
[9]], &c__1
);
1496 F77_FUNC(dgemv
,DGEMV
)("N", n
, &iwork
[12], &c_b50
, &v
[v_offset
], ldv
, &workd
[iwork
[9]], &
1497 c__1
, &c_b18
, &resid
[1], &c__1
);
1499 if (iwork
[12] == 1 || iwork
[4] == 1) {
1500 h__
[iwork
[12] + h_dim1
] = 0.;
1502 h__
[iwork
[12] + (h_dim1
<< 1)] += workd
[iwork
[9] + iwork
[12] - 1];
1506 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[9]], &c__1
);
1507 ipntr
[1] = iwork
[9];
1508 ipntr
[2] = iwork
[8];
1512 } else if (*bmat
== 'I') {
1513 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &c__1
);
1519 workd
[*n
* 3 + 2] =F77_FUNC(ddot
,DDOT
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &
1521 workd
[*n
* 3 + 2] = sqrt(fabs(workd
[*n
* 3 + 2]));
1522 } else if (*bmat
== 'I') {
1523 workd
[*n
* 3 + 2] =F77_FUNC(dnrm2
,DNRM2
)(n
, &resid
[1], &c__1
);
1527 if (workd
[*n
* 3 + 2] > *rnorm
* .717f
) {
1529 *rnorm
= workd
[*n
* 3 + 2];
1533 *rnorm
= workd
[*n
* 3 + 2];
1535 if (iwork
[1] <= 1) {
1540 for (jj
= 1; jj
<= i__1
; ++jj
) {
1551 if (h__
[iwork
[12] + h_dim1
] < 0.) {
1552 h__
[iwork
[12] + h_dim1
] = -h__
[iwork
[12] + h_dim1
];
1553 if (iwork
[12] < *k
+ *np
) {
1554 F77_FUNC(dscal
,DSCAL
)(n
, &c_b50
, &v
[(iwork
[12] + 1) * v_dim1
+ 1], &c__1
);
1556 F77_FUNC(dscal
,DSCAL
)(n
, &c_b50
, &resid
[1], &c__1
);
1561 if (iwork
[12] > *k
+ *np
) {
1580 F77_FUNC(dsaup2
,DSAUP2
)(int * ido
,
1610 int h_dim1
, h_offset
, q_dim1
, q_offset
, v_dim1
, v_offset
, i__1
, i__2
,
1629 v_offset
= 1 + v_dim1
;
1632 h_offset
= 1 + h_dim1
;
1635 q_offset
= 1 + q_dim1
;
1639 eps23
= GMX_DOUBLE_EPS
;
1640 eps23
= pow(eps23
, c_b3
);
1652 iwork
[7] = iwork
[9] + iwork
[10];
1670 if (iwork
[2] == 1) {
1671 F77_FUNC(dgetv0
,DGETV0
)(ido
, bmat
, &c__1
, &iwork
[3], n
, &c__1
, &v
[v_offset
], ldv
, &
1672 resid
[1], &workd
[*n
* 3 + 1], &ipntr
[1], &workd
[1], &iwork
[41]
1679 if (workd
[*n
* 3 + 1] == 0.) {
1688 if (iwork
[4] == 1) {
1692 if (iwork
[5] == 1) {
1696 if (iwork
[1] == 1) {
1700 F77_FUNC(dsaitr
,DSAITR
)(ido
, bmat
, n
, &c__0
, &iwork
[9], mode
, &resid
[1], &workd
[*n
* 3 +
1701 1], &v
[v_offset
], ldv
, &h__
[h_offset
], ldh
, &ipntr
[1], &workd
[1],
1725 F77_FUNC(dsaitr
,DSAITR
)(ido
, bmat
, n
, nev
, np
, mode
, &resid
[1], &workd
[*n
* 3 + 1],
1726 &v
[v_offset
], ldv
, &h__
[h_offset
], ldh
, &ipntr
[1], &workd
[1],
1742 F77_FUNC(dseigt
,DSEIGT
)(&workd
[*n
* 3 + 1], &iwork
[7], &h__
[h_offset
], ldh
, &ritz
[1], &
1743 bounds
[1], &workl
[1], &ierr
);
1750 F77_FUNC(dcopy
,DCOPY
)(&iwork
[7], &ritz
[1], &c__1
, &workl
[iwork
[7] + 1], &c__1
);
1751 F77_FUNC(dcopy
,DCOPY
)(&iwork
[7], &bounds
[1], &c__1
, &workl
[(iwork
[7] << 1) + 1], &c__1
);
1755 F77_FUNC(dsgets
,DSGETS
)(ishift
, which
, nev
, np
, &ritz
[1], &bounds
[1], &workl
[1]);
1757 F77_FUNC(dcopy
,DCOPY
)(nev
, &bounds
[*np
+ 1], &c__1
, &workl
[*np
+ 1], &c__1
);
1758 F77_FUNC(dsconv
,DSCONV
)(nev
, &ritz
[*np
+ 1], &workl
[*np
+ 1], tol
, &iwork
[8]);
1762 for (j
= 1; j
<= i__1
; ++j
) {
1763 if (bounds
[j
] == 0.) {
1769 if (iwork
[8] >= iwork
[9] || iwork
[6] > *mxiter
|| *np
== 0) {
1771 if (!strncmp(which
, "BE", 2)) {
1773 strncpy(wprime
, "SA",2);
1774 F77_FUNC(dsortr
,DSORTR
)(wprime
, &c__1
, &iwork
[7], &ritz
[1], &bounds
[1]);
1776 nevm2
= *nev
- nevd2
;
1778 i__1
= (nevd2
< *np
) ? nevd2
: *np
;
1779 i__2
= iwork
[7] - nevd2
+ 1, i__3
= iwork
[7] - *np
+ 1;
1780 F77_FUNC(dswap
,DSWAP
)(&i__1
, &ritz
[nevm2
+ 1], &c__1
,
1781 &ritz
[((i__2
>i__3
) ? i__2
: i__3
)],
1783 i__1
= (nevd2
< *np
) ? nevd2
: *np
;
1784 i__2
= iwork
[7] - nevd2
+ 1, i__3
= iwork
[7] - *np
;
1785 F77_FUNC(dswap
,DSWAP
)(&i__1
, &bounds
[nevm2
+ 1], &c__1
,
1786 &bounds
[((i__2
>i__3
) ? i__2
: i__3
) + 1],
1792 if (!strncmp(which
, "LM", 2)) {
1793 strncpy(wprime
, "SM", 2);
1795 if (!strncmp(which
, "SM", 2)) {
1796 strncpy(wprime
, "LM", 2);
1798 if (!strncmp(which
, "LA", 2)) {
1799 strncpy(wprime
, "SA", 2);
1801 if (!strncmp(which
, "SA", 2)) {
1802 strncpy(wprime
, "LA", 2);
1805 F77_FUNC(dsortr
,DSORTR
)(wprime
, &c__1
, &iwork
[7], &ritz
[1], &bounds
[1]);
1810 for (j
= 1; j
<= i__1
; ++j
) {
1812 d__3
= fabs(ritz
[j
]);
1813 temp
= (d__2
>d__3
) ? d__2
: d__3
;
1817 strncpy(wprime
, "LA",2);
1818 F77_FUNC(dsortr
,DSORTR
)(wprime
, &c__1
, &iwork
[9], &bounds
[1], &ritz
[1]);
1821 for (j
= 1; j
<= i__1
; ++j
) {
1823 d__3
= fabs(ritz
[j
]);
1824 temp
= (d__2
>d__3
) ? d__2
: d__3
;
1828 if (!strncmp(which
, "BE", 2)) {
1830 strncpy(wprime
, "LA", 2);
1831 F77_FUNC(dsortr
,DSORTR
)(wprime
, &c__1
, &iwork
[8], &ritz
[1], &bounds
[1]);
1834 F77_FUNC(dsortr
,DSORTR
)(which
, &c__1
, &iwork
[8], &ritz
[1], &bounds
[1]);
1838 h__
[h_dim1
+ 1] = workd
[*n
* 3 + 1];
1841 if (iwork
[6] > *mxiter
&& iwork
[8] < *nev
) {
1844 if (*np
== 0 && iwork
[8] < iwork
[9]) {
1851 } else if (iwork
[8] < *nev
&& *ishift
== 1) {
1853 i__1
= iwork
[8], i__2
= *np
/ 2;
1854 *nev
+= (i__1
< i__2
) ? i__1
: i__2
;
1855 if (*nev
== 1 && iwork
[7] >= 6) {
1856 *nev
= iwork
[7] / 2;
1857 } else if (*nev
== 1 && iwork
[7] > 2) {
1860 *np
= iwork
[7] - *nev
;
1863 if (nevbef
< *nev
) {
1864 F77_FUNC(dsgets
,DSGETS
)(ishift
, which
, nev
, np
, &ritz
[1], &bounds
[1], &workl
[1]);
1882 F77_FUNC(dcopy
,DCOPY
)(np
, &workl
[1], &c__1
, &ritz
[1], &c__1
);
1885 F77_FUNC(dsapps
,DSAPPS
)(n
, nev
, np
, &ritz
[1], &v
[v_offset
], ldv
, &h__
[h_offset
], ldh
, &
1886 resid
[1], &q
[q_offset
], ldq
, &workd
[1]);
1890 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[*n
+ 1], &c__1
);
1896 } else if (*bmat
== 'I') {
1897 F77_FUNC(dcopy
,DCOPY
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
1903 workd
[*n
* 3 + 1] =F77_FUNC(ddot
,DDOT
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
1904 workd
[*n
* 3 + 1] = sqrt(fabs(workd
[*n
* 3 + 1]));
1905 } else if (*bmat
== 'I') {
1906 workd
[*n
* 3 + 1] =F77_FUNC(dnrm2
,DNRM2
)(n
, &resid
[1], &c__1
);
1928 F77_FUNC(dsaupd
,DSAUPD
)(int * ido
,
1946 int v_dim1
, v_offset
, i__1
, i__2
;
1952 v_offset
= 1 + v_dim1
;
1963 iwork
[5] = iparam
[1];
1964 iwork
[10] = iparam
[3];
1965 iwork
[12] = iparam
[4];
1968 iwork
[11] = iparam
[7];
1973 } else if (*nev
<= 0) {
1975 } else if (*ncv
<= *nev
|| *ncv
> *n
) {
1980 iwork
[15] = *ncv
- *nev
;
1982 if (iwork
[10] <= 0) {
1985 if (strncmp(which
,"LM",2) && strncmp(which
,"SM",2) &&
1986 strncmp(which
,"LA",2) && strncmp(which
,"SA",2) &&
1987 strncmp(which
,"BE",2)) {
1990 if (*bmat
!= 'I' && *bmat
!= 'G') {
1995 if (*lworkl
< i__1
* i__1
+ (*ncv
<< 3)) {
1998 if (iwork
[11] < 1 || iwork
[11] > 5) {
2000 } else if (iwork
[11] == 1 && *bmat
== 'G') {
2002 } else if (iwork
[5] < 0 || iwork
[5] > 1) {
2004 } else if (*nev
== 1 && !strncmp(which
, "BE", 2)) {
2008 if (iwork
[2] != 0) {
2014 if (iwork
[12] <= 0) {
2018 *tol
= GMX_DOUBLE_EPS
;
2021 iwork
[15] = *ncv
- *nev
;
2024 i__1
= i__2
* i__2
+ (*ncv
<< 3);
2025 for (j
= 1; j
<= i__1
; ++j
) {
2032 iwork
[16] = iwork
[3] + (iwork
[8] << 1);
2033 iwork
[1] = iwork
[16] + *ncv
;
2034 iwork
[4] = iwork
[1] + *ncv
;
2036 iwork
[7] = iwork
[4] + i__1
* i__1
;
2037 iwork
[14] = iwork
[7] + *ncv
* 3;
2039 ipntr
[4] = iwork
[14];
2040 ipntr
[5] = iwork
[3];
2041 ipntr
[6] = iwork
[16];
2042 ipntr
[7] = iwork
[1];
2043 ipntr
[11] = iwork
[7];
2046 F77_FUNC(dsaup2
,DSAUP2
)(ido
, bmat
, n
, which
, &iwork
[13], &iwork
[15], tol
, &resid
[1], &
2047 iwork
[11], &iwork
[6], &iwork
[5], &iwork
[10], &v
[v_offset
], ldv
, &
2048 workl
[iwork
[3]], &iwork
[8], &workl
[iwork
[16]], &workl
[iwork
[1]], &
2049 workl
[iwork
[4]], &iwork
[9], &workl
[iwork
[7]], &ipntr
[1], &workd
[1]
2050 , &iwork
[21], info
);
2053 iparam
[8] = iwork
[15];
2059 iparam
[3] = iwork
[10];
2060 iparam
[5] = iwork
[15];
2078 F77_FUNC(dseupd
,DSEUPD
)(int * rvec
,
2079 const char * howmny
,
2104 int v_dim1
, v_offset
, z_dim1
, z_offset
, i__1
;
2105 double d__1
, d__2
, d__3
;
2107 int j
, k
, ih
, iq
, iw
, ibd
, ihb
, ihd
, ldh
, ilg
, ldq
, ism
, irz
;
2119 double thres1
=0, thres2
=0;
2123 int leftptr
, rghtptr
;
2129 z_offset
= 1 + z_dim1
;
2134 v_offset
= 1 + v_dim1
;
2158 if (*ncv
<= *nev
|| *ncv
> *n
) {
2161 if (strncmp(which
,"LM",2) && strncmp(which
,"SM",2) &&
2162 strncmp(which
,"LA",2) && strncmp(which
,"SA",2) &&
2163 strncmp(which
,"BE",2)) {
2166 if (*bmat
!= 'I' && *bmat
!= 'G') {
2169 if (*howmny
!= 'A' && *howmny
!= 'P' &&
2170 *howmny
!= 'S' && *rvec
) {
2173 if (*rvec
&& *howmny
== 'S') {
2177 if (*rvec
&& *lworkl
< i__1
* i__1
+ (*ncv
<< 3)) {
2181 if (mode
== 1 || mode
== 2) {
2182 strncpy(type__
, "REGULR",6);
2183 } else if (mode
== 3) {
2184 strncpy(type__
, "SHIFTI",6);
2185 } else if (mode
== 4) {
2186 strncpy(type__
, "BUCKLE",6);
2187 } else if (mode
== 5) {
2188 strncpy(type__
, "CAYLEY",6);
2192 if (mode
== 1 && *bmat
== 'G') {
2195 if (*nev
== 1 && !strncmp(which
, "BE",2)) {
2212 iw
= iq
+ ldh
* *ncv
;
2213 next
= iw
+ (*ncv
<< 1);
2219 irz
= ipntr
[11] + *ncv
;
2223 eps23
= GMX_DOUBLE_EPS
;
2224 eps23
= pow(eps23
, c_b21
);
2229 } else if (*bmat
== 'G') {
2230 bnorm2
=F77_FUNC(dnrm2
,DNRM2
)(n
, &workd
[1], &c__1
);
2235 if (!strncmp(which
,"LM",2) || !strncmp(which
,"SM",2) ||
2236 !strncmp(which
,"LA",2) || !strncmp(which
,"SA",2)) {
2238 } else if (!strncmp(which
,"BE",2)) {
2241 ism
= (*nev
>nconv
) ? *nev
: nconv
;
2244 thres1
= workl
[ism
];
2245 thres2
= workl
[ilg
];
2253 for (j
= 0; j
<= i__1
; ++j
) {
2255 if (!strncmp(which
,"LM",2)) {
2256 if (fabs(workl
[irz
+ j
]) >= fabs(thres1
)) {
2258 d__3
= fabs(workl
[irz
+ j
]);
2259 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
2260 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
2264 } else if (!strncmp(which
,"SM",2)) {
2265 if (fabs(workl
[irz
+ j
]) <= fabs(thres1
)) {
2267 d__3
= fabs(workl
[irz
+ j
]);
2268 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
2269 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
2273 } else if (!strncmp(which
,"LA",2)) {
2274 if (workl
[irz
+ j
] >= thres1
) {
2276 d__3
= fabs(workl
[irz
+ j
]);
2277 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
2278 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
2282 } else if (!strncmp(which
,"SA",2)) {
2283 if (workl
[irz
+ j
] <= thres1
) {
2285 d__3
= fabs(workl
[irz
+ j
]);
2286 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
2287 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
2291 } else if (!strncmp(which
,"BE",2)) {
2292 if (workl
[irz
+ j
] <= thres1
|| workl
[irz
+ j
] >= thres2
) {
2294 d__3
= fabs(workl
[irz
+ j
]);
2295 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
2296 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
2301 if (j
+ 1 > nconv
) {
2302 reord
= select
[j
+ 1] || reord
;
2304 if (select
[j
+ 1]) {
2310 F77_FUNC(dcopy
,DCOPY
)(&i__1
, &workl
[ih
+ 1], &c__1
, &workl
[ihb
], &c__1
);
2311 F77_FUNC(dcopy
,DCOPY
)(ncv
, &workl
[ih
+ ldh
], &c__1
, &workl
[ihd
], &c__1
);
2313 F77_FUNC(dsteqr
,DSTEQR
)("Identity", ncv
, &workl
[ihd
], &workl
[ihb
], &workl
[iq
], &ldq
, &
2332 if (select
[leftptr
]) {
2336 } else if (! select
[rghtptr
]) {
2342 temp
= workl
[ihd
+ leftptr
- 1];
2343 workl
[ihd
+ leftptr
- 1] = workl
[ihd
+ rghtptr
- 1];
2344 workl
[ihd
+ rghtptr
- 1] = temp
;
2345 F77_FUNC(dcopy
,DCOPY
)(ncv
, &workl
[iq
+ *ncv
* (leftptr
- 1)], &c__1
, &workl
[
2347 F77_FUNC(dcopy
,DCOPY
)(ncv
, &workl
[iq
+ *ncv
* (rghtptr
- 1)], &c__1
, &workl
[
2348 iq
+ *ncv
* (leftptr
- 1)], &c__1
);
2349 F77_FUNC(dcopy
,DCOPY
)(ncv
, &workl
[iw
], &c__1
, &workl
[iq
+ *ncv
* (rghtptr
-
2356 if (leftptr
< rghtptr
) {
2364 F77_FUNC(dcopy
,DCOPY
)(&nconv
, &workl
[ihd
], &c__1
, &d__
[1], &c__1
);
2368 F77_FUNC(dcopy
,DCOPY
)(&nconv
, &workl
[ritz
], &c__1
, &d__
[1], &c__1
);
2369 F77_FUNC(dcopy
,DCOPY
)(ncv
, &workl
[ritz
], &c__1
, &workl
[ihd
], &c__1
);
2372 if (!strncmp(type__
, "REGULR",6)) {
2375 F77_FUNC(dsesrt
,DSESRT
)("LA", rvec
, &nconv
, &d__
[1], ncv
, &workl
[iq
], &ldq
);
2377 F77_FUNC(dcopy
,DCOPY
)(ncv
, &workl
[bounds
], &c__1
, &workl
[ihb
], &c__1
);
2382 F77_FUNC(dcopy
,DCOPY
)(ncv
, &workl
[ihd
], &c__1
, &workl
[iw
], &c__1
);
2383 if (!strncmp(type__
, "SHIFTI", 6)) {
2385 for (k
= 1; k
<= i__1
; ++k
) {
2386 workl
[ihd
+ k
- 1] = 1. / workl
[ihd
+ k
- 1] + *sigma
;
2388 } else if (!strncmp(type__
, "BUCKLE",6)) {
2390 for (k
= 1; k
<= i__1
; ++k
) {
2391 workl
[ihd
+ k
- 1] = *sigma
* workl
[ihd
+ k
- 1] / (workl
[ihd
2394 } else if (!strncmp(type__
, "CAYLEY",6)) {
2396 for (k
= 1; k
<= i__1
; ++k
) {
2397 workl
[ihd
+ k
- 1] = *sigma
* (workl
[ihd
+ k
- 1] + 1.) / (
2398 workl
[ihd
+ k
- 1] - 1.);
2402 F77_FUNC(dcopy
,DCOPY
)(&nconv
, &workl
[ihd
], &c__1
, &d__
[1], &c__1
);
2403 F77_FUNC(dsortr
,DSORTR
)("LA", &c__1
, &nconv
, &workl
[ihd
], &workl
[iw
]);
2405 F77_FUNC(dsesrt
,DSESRT
)("LA", rvec
, &nconv
, &d__
[1], ncv
, &workl
[iq
], &ldq
);
2407 F77_FUNC(dcopy
,DCOPY
)(ncv
, &workl
[bounds
], &c__1
, &workl
[ihb
], &c__1
);
2408 d__1
= bnorm2
/ rnorm
;
2409 F77_FUNC(dscal
,DSCAL
)(ncv
, &d__1
, &workl
[ihb
], &c__1
);
2410 F77_FUNC(dsortr
,DSORTR
)("LA", &c__1
, &nconv
, &d__
[1], &workl
[ihb
]);
2415 if (*rvec
&& *howmny
== 'A') {
2417 F77_FUNC(dgeqr2
,DGEQR2
)(ncv
, &nconv
, &workl
[iq
], &ldq
, &workl
[iw
+ *ncv
], &workl
[ihb
],
2420 F77_FUNC(dorm2r
,DORM2R
)("Right", "Notranspose", n
, ncv
, &nconv
, &workl
[iq
], &ldq
, &
2421 workl
[iw
+ *ncv
], &v
[v_offset
], ldv
, &workd
[*n
+ 1], &ierr
);
2422 F77_FUNC(dlacpy
,DLACPY
)("All", n
, &nconv
, &v
[v_offset
], ldv
, &z__
[z_offset
], ldz
);
2425 for (j
= 1; j
<= i__1
; ++j
) {
2426 workl
[ihb
+ j
- 1] = 0.;
2428 workl
[ihb
+ *ncv
- 1] = 1.;
2429 F77_FUNC(dorm2r
,DORM2R
)("Left", "Transpose", ncv
, &c__1
, &nconv
, &workl
[iq
], &ldq
, &
2430 workl
[iw
+ *ncv
], &workl
[ihb
], ncv
, &temp
, &ierr
);
2432 } else if (*rvec
&& *howmny
== 'S') {
2436 if (!strncmp(type__
, "REGULR",6) && *rvec
) {
2439 for (j
= 1; j
<= i__1
; ++j
) {
2440 workl
[ihb
+ j
- 1] = rnorm
* fabs(workl
[ihb
+ j
- 1]);
2443 } else if (strncmp(type__
, "REGULR",6) && *rvec
) {
2445 F77_FUNC(dscal
,DSCAL
)(ncv
, &bnorm2
, &workl
[ihb
], &c__1
);
2446 if (!strncmp(type__
, "SHIFTI",6)) {
2449 for (k
= 1; k
<= i__1
; ++k
) {
2450 d__2
= workl
[iw
+ k
- 1];
2451 workl
[ihb
+ k
- 1] = fabs(workl
[ihb
+ k
- 1])/(d__2
* d__2
);
2454 } else if (!strncmp(type__
, "BUCKLE",6)) {
2457 for (k
= 1; k
<= i__1
; ++k
) {
2458 d__2
= workl
[iw
+ k
- 1] - 1.;
2459 workl
[ihb
+ k
- 1] = *sigma
* fabs(workl
[ihb
+ k
- 1])/(d__2
* d__2
);
2462 } else if (!strncmp(type__
, "CAYLEY",6)) {
2465 for (k
= 1; k
<= i__1
; ++k
) {
2466 workl
[ihb
+ k
- 1] = fabs(workl
[ihb
+ k
- 1] / workl
[iw
+ k
- 1] * (workl
[iw
+ k
- 1] - 1.));
2474 if (*rvec
&& (!strncmp(type__
, "SHIFTI",6) || !strncmp(type__
, "CAYLEY",6))) {
2477 for (k
= 0; k
<= i__1
; ++k
) {
2478 workl
[iw
+ k
] = workl
[iq
+ k
* ldq
+ *ncv
- 1] / workl
[iw
+ k
];
2481 } else if (*rvec
&& !strncmp(type__
, "BUCKLE", 6)) {
2484 for (k
= 0; k
<= i__1
; ++k
) {
2485 workl
[iw
+ k
] = workl
[iq
+ k
* ldq
+ *ncv
- 1] / (workl
[iw
+ k
] -
2491 if (strncmp(type__
, "REGULR",6)) {
2492 F77_FUNC(dger
,DGER
)(n
, &nconv
, &c_b102
, &resid
[1], &c__1
, &workl
[iw
], &c__1
, &z__
[
2506 /* Selected single precision arpack routines */
2510 F77_FUNC(sstqrb
,SSTQRB
)(int * n
,
2524 int i__
, j
, k
, l
, m
;
2526 int l1
, ii
, mm
, lm1
, mm1
, nm1
;
2527 float rt1
, rt2
, eps
;
2530 int lend
, jtot
, lendm1
, lendp1
, iscale
;
2532 int lendsv
, nmaxit
, icompz
;
2533 float ssfmax
, ssfmin
,safmin
,minval
,safmax
,anorm
;
2556 eps
= GMX_FLOAT_EPS
;
2560 minval
= GMX_FLOAT_MIN
;
2561 safmin
= minval
/ GMX_FLOAT_EPS
;
2562 safmax
= 1. / safmin
;
2563 ssfmax
= sqrt(safmax
) / 3.;
2564 ssfmin
= sqrt(safmin
) / eps2
;
2568 for (j
= 1; j
<= i__1
; ++j
) {
2590 for (m
= l1
; m
<= i__1
; ++m
) {
2595 if (tst
<= sqrt(fabs(d__
[m
])) * sqrt(fabs(d__
[m
+1])) * eps
) {
2613 i__1
= lend
- l
+ 1;
2614 anorm
=F77_FUNC(slanst
,SLANST
)("i", &i__1
, &d__
[l
], &e
[l
]);
2619 if (anorm
> ssfmax
) {
2621 i__1
= lend
- l
+ 1;
2622 F77_FUNC(slascl
,SLASCL
)("g", &c__0
, &c__0
, &anorm
, &ssfmax
, &i__1
, &c__1
, &d__
[l
], n
,
2625 F77_FUNC(slascl
,SLASCL
)("g", &c__0
, &c__0
, &anorm
, &ssfmax
, &i__1
, &c__1
, &e
[l
], n
,
2627 } else if (anorm
< ssfmin
) {
2629 i__1
= lend
- l
+ 1;
2630 F77_FUNC(slascl
,SLASCL
)("g", &c__0
, &c__0
, &anorm
, &ssfmin
, &i__1
, &c__1
, &d__
[l
], n
,
2633 F77_FUNC(slascl
,SLASCL
)("g", &c__0
, &c__0
, &anorm
, &ssfmin
, &i__1
, &c__1
, &e
[l
], n
,
2637 if (fabs(d__
[lend
]) < fabs(d__
[l
])) {
2648 for (m
= l
; m
<= i__1
; ++m
) {
2651 if (tst
<= eps2
* fabs(d__
[m
]) * fabs(d__
[m
+ 1]) + safmin
) {
2670 F77_FUNC(slaev2
,SLAEV2
)(&d__
[l
], &e
[l
], &d__
[l
+ 1], &rt1
, &rt2
, &c__
, &s
);
2672 work
[*n
- 1 + l
] = s
;
2675 z__
[l
+ 1] = c__
* tst
- s
* z__
[l
];
2676 z__
[l
] = s
* tst
+ c__
* z__
[l
];
2678 F77_FUNC(slae2
,SLAE2
)(&d__
[l
], &e
[l
], &d__
[l
+ 1], &rt1
, &rt2
);
2690 if (jtot
== nmaxit
) {
2695 g
= (d__
[l
+ 1] - p
) / (e
[l
] * 2.);
2696 r__
=F77_FUNC(slapy2
,SLAPY2
)(&g
, &c_b31
);
2697 g
= d__
[m
] - p
+ e
[l
] / (g
+ ((g
>0) ? r__
: -r__
));
2705 for (i__
= mm1
; i__
>= i__1
; --i__
) {
2708 F77_FUNC(slartg
,SLARTG
)(&g
, &f
, &c__
, &s
, &r__
);
2712 g
= d__
[i__
+ 1] - p
;
2713 r__
= (d__
[i__
] - g
) * s
+ c__
* 2. * b
;
2715 d__
[i__
+ 1] = g
+ p
;
2720 work
[*n
- 1 + i__
] = -s
;
2728 F77_FUNC(slasr
,SLASR
)("r", "v", "b", &c__1
, &mm
, &work
[l
], &work
[*n
- 1 + l
], &
2751 for (m
= l
; m
>= i__1
; --m
) {
2752 d__2
= fabs(e
[m
- 1]);
2754 if (tst
<= eps2
* fabs(d__
[m
]) * fabs(d__
[m
- 1]) + safmin
) {
2773 F77_FUNC(slaev2
,SLAEV2
)(&d__
[l
- 1], &e
[l
- 1], &d__
[l
], &rt1
, &rt2
, &c__
, &s
)
2777 z__
[l
] = c__
* tst
- s
* z__
[l
- 1];
2778 z__
[l
- 1] = s
* tst
+ c__
* z__
[l
- 1];
2780 F77_FUNC(slae2
,SLAE2
)(&d__
[l
- 1], &e
[l
- 1], &d__
[l
], &rt1
, &rt2
);
2792 if (jtot
== nmaxit
) {
2798 g
= (d__
[l
- 1] - p
) / (e
[l
- 1] * 2.);
2799 r__
=F77_FUNC(slapy2
,SLAPY2
)(&g
, &c_b31
);
2800 g
= d__
[m
] - p
+ e
[l
- 1] / (g
+ ((g
>0) ? r__
: -r__
));
2808 for (i__
= m
; i__
<= i__1
; ++i__
) {
2811 F77_FUNC(slartg
,SLARTG
)(&g
, &f
, &c__
, &s
, &r__
);
2816 r__
= (d__
[i__
+ 1] - g
) * s
+ c__
* 2. * b
;
2823 work
[*n
- 1 + i__
] = s
;
2831 F77_FUNC(slasr
,SLASR
)("r", "v", "f", &c__1
, &mm
, &work
[m
], &work
[*n
- 1 + m
], &
2852 i__1
= lendsv
- lsv
+ 1;
2853 F77_FUNC(slascl
,SLASCL
)("g", &c__0
, &c__0
, &ssfmax
, &anorm
, &i__1
, &c__1
, &d__
[lsv
],
2855 i__1
= lendsv
- lsv
;
2856 F77_FUNC(slascl
,SLASCL
)("g", &c__0
, &c__0
, &ssfmax
, &anorm
, &i__1
, &c__1
, &e
[lsv
], n
,
2858 } else if (iscale
== 2) {
2859 i__1
= lendsv
- lsv
+ 1;
2860 F77_FUNC(slascl
,SLASCL
)("g", &c__0
, &c__0
, &ssfmin
, &anorm
, &i__1
, &c__1
, &d__
[lsv
],
2862 i__1
= lendsv
- lsv
;
2863 F77_FUNC(slascl
,SLASCL
)("g", &c__0
, &c__0
, &ssfmin
, &anorm
, &i__1
, &c__1
, &e
[lsv
], n
,
2867 if (jtot
< nmaxit
) {
2871 for (i__
= 1; i__
<= i__1
; ++i__
) {
2881 F77_FUNC(slasrt
,SLASRT
)("i", n
, &d__
[1], info
);
2886 for (ii
= 2; ii
<= i__1
; ++ii
) {
2891 for (j
= ii
; j
<= i__2
; ++j
) {
2914 F77_FUNC(sgetv0
,SGETV0
)(int * ido
,
2933 int v_dim1
, v_offset
, i__1
;
2941 v_offset
= 1 + v_dim1
;
2955 F77_FUNC(slarnv
,SLARNV
)(&idist
, &iwork
[1], n
, &resid
[1]);
2961 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
2967 if (iwork
[5] == 1) {
2971 if (iwork
[6] == 1) {
2977 F77_FUNC(scopy
,SCOPY
)(n
, &workd
[*n
+ 1], &c__1
, &resid
[1], &c__1
);
2982 } else if (*bmat
== 'I') {
2983 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
2991 workd
[*n
* 3 + 4] =F77_FUNC(sdot
,SDOT
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
2992 workd
[*n
* 3 + 4] = sqrt(fabs(workd
[*n
* 3 + 4]));
2993 } else if (*bmat
== 'I') {
2994 workd
[*n
* 3 + 4] =F77_FUNC(snrm2
,SNRM2
)(n
, &resid
[1], &c__1
);
2996 *rnorm
= workd
[*n
* 3 + 4];
3005 F77_FUNC(sgemv
,SGEMV
)("T", n
, &i__1
, &c_b22
, &v
[v_offset
], ldv
, &workd
[1], &c__1
, &c_b24
,
3006 &workd
[*n
+ 1], &c__1
);
3008 F77_FUNC(sgemv
,SGEMV
)("N", n
, &i__1
, &c_b27
, &v
[v_offset
], ldv
, &workd
[*n
+ 1], &c__1
, &
3009 c_b22
, &resid
[1], &c__1
);
3012 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[*n
+ 1], &c__1
);
3017 } else if (*bmat
== 'I') {
3018 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
3024 *rnorm
=F77_FUNC(sdot
,SDOT
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
3025 *rnorm
= sqrt(fabs(*rnorm
));
3026 } else if (*bmat
== 'I') {
3027 *rnorm
=F77_FUNC(snrm2
,SNRM2
)(n
, &resid
[1], &c__1
);
3030 if (*rnorm
> workd
[*n
* 3 + 4] * .717f
) {
3035 if (iwork
[7] <= 1) {
3037 workd
[*n
* 3 + 4] = *rnorm
;
3042 for (jj
= 1; jj
<= i__1
; ++jj
) {
3062 F77_FUNC(ssapps
,SSAPPS
)(int * n
,
3079 int h_dim1
, h_offset
, q_dim1
, q_offset
, v_dim1
, v_offset
, i__1
, i__2
,
3083 float r__
, s
, a1
, a2
, a3
, a4
;
3094 v_offset
= 1 + v_dim1
;
3097 h_offset
= 1 + h_dim1
;
3100 q_offset
= 1 + q_dim1
;
3103 epsmch
= GMX_FLOAT_EPS
;
3107 kplusp
= *kev
+ *np
;
3109 F77_FUNC(slaset
,SLASET
)("All", &kplusp
, &kplusp
, &c_b4
, &c_b5
, &q
[q_offset
], ldq
);
3116 for (jj
= 1; jj
<= i__1
; ++jj
) {
3123 for (i__
= istart
; i__
<= i__2
; ++i__
) {
3124 big
= fabs(h__
[i__
+ (h_dim1
*2)]) + fabs(h__
[i__
+ 1 + (h_dim1
*2)]);
3125 if (h__
[i__
+ 1 + h_dim1
] <= epsmch
* big
) {
3126 h__
[i__
+ 1 + h_dim1
] = 0.;
3134 if (istart
< iend
) {
3136 f
= h__
[istart
+ (h_dim1
<< 1)] - shift
[jj
];
3137 g
= h__
[istart
+ 1 + h_dim1
];
3138 F77_FUNC(slartg
,SLARTG
)(&f
, &g
, &c__
, &s
, &r__
);
3140 a1
= c__
* h__
[istart
+ (h_dim1
<< 1)] + s
* h__
[istart
+ 1 +
3142 a2
= c__
* h__
[istart
+ 1 + h_dim1
] + s
* h__
[istart
+ 1 + (
3144 a4
= c__
* h__
[istart
+ 1 + (h_dim1
<< 1)] - s
* h__
[istart
+ 1 +
3146 a3
= c__
* h__
[istart
+ 1 + h_dim1
] - s
* h__
[istart
+ (h_dim1
<<
3148 h__
[istart
+ (h_dim1
<< 1)] = c__
* a1
+ s
* a2
;
3149 h__
[istart
+ 1 + (h_dim1
<< 1)] = c__
* a4
- s
* a3
;
3150 h__
[istart
+ 1 + h_dim1
] = c__
* a3
+ s
* a4
;
3153 i__2
= (i__3
<kplusp
) ? i__3
: kplusp
;
3154 for (j
= 1; j
<= i__2
; ++j
) {
3155 a1
= c__
* q
[j
+ istart
* q_dim1
] + s
* q
[j
+ (istart
+ 1) *
3157 q
[j
+ (istart
+ 1) * q_dim1
] = -s
* q
[j
+ istart
* q_dim1
] +
3158 c__
* q
[j
+ (istart
+ 1) * q_dim1
];
3159 q
[j
+ istart
* q_dim1
] = a1
;
3164 for (i__
= istart
+ 1; i__
<= i__2
; ++i__
) {
3166 f
= h__
[i__
+ h_dim1
];
3167 g
= s
* h__
[i__
+ 1 + h_dim1
];
3169 h__
[i__
+ 1 + h_dim1
] = c__
* h__
[i__
+ 1 + h_dim1
];
3170 F77_FUNC(slartg
,SLARTG
)(&f
, &g
, &c__
, &s
, &r__
);
3178 h__
[i__
+ h_dim1
] = r__
;
3180 a1
= c__
* h__
[i__
+ (h_dim1
<< 1)] + s
* h__
[i__
+ 1 +
3182 a2
= c__
* h__
[i__
+ 1 + h_dim1
] + s
* h__
[i__
+ 1 + (h_dim1
3184 a3
= c__
* h__
[i__
+ 1 + h_dim1
] - s
* h__
[i__
+ (h_dim1
<< 1)
3186 a4
= c__
* h__
[i__
+ 1 + (h_dim1
<< 1)] - s
* h__
[i__
+ 1 +
3189 h__
[i__
+ (h_dim1
<< 1)] = c__
* a1
+ s
* a2
;
3190 h__
[i__
+ 1 + (h_dim1
<< 1)] = c__
* a4
- s
* a3
;
3191 h__
[i__
+ 1 + h_dim1
] = c__
* a3
+ s
* a4
;
3194 i__3
= (i__4
<kplusp
) ? i__4
: kplusp
;
3195 for (j
= 1; j
<= i__3
; ++j
) {
3196 a1
= c__
* q
[j
+ i__
* q_dim1
] + s
* q
[j
+ (i__
+ 1) *
3198 q
[j
+ (i__
+ 1) * q_dim1
] = -s
* q
[j
+ i__
* q_dim1
] +
3199 c__
* q
[j
+ (i__
+ 1) * q_dim1
];
3200 q
[j
+ i__
* q_dim1
] = a1
;
3209 if (h__
[iend
+ h_dim1
] < 0.) {
3210 h__
[iend
+ h_dim1
] = -h__
[iend
+ h_dim1
];
3211 F77_FUNC(sscal
,SSCAL
)(&kplusp
, &c_b14
, &q
[iend
* q_dim1
+ 1], &c__1
);
3214 if (iend
< kplusp
) {
3219 for (i__
= itop
; i__
<= i__2
; ++i__
) {
3220 if (h__
[i__
+ 1 + h_dim1
] > 0.) {
3231 for (i__
= itop
; i__
<= i__1
; ++i__
) {
3232 big
= fabs(h__
[i__
+ (h_dim1
*2)]) + fabs(h__
[i__
+ 1 + (h_dim1
*2)]);
3233 if (h__
[i__
+ 1 + h_dim1
] <= epsmch
* big
) {
3234 h__
[i__
+ 1 + h_dim1
] = 0.;
3239 if (h__
[*kev
+ 1 + h_dim1
] > 0.) {
3240 F77_FUNC(sgemv
,SGEMV
)("N", n
, &kplusp
, &c_b5
, &v
[v_offset
], ldv
, &q
[(*kev
+ 1) *
3241 q_dim1
+ 1], &c__1
, &c_b4
, &workd
[*n
+ 1], &c__1
);
3245 for (i__
= 1; i__
<= i__1
; ++i__
) {
3246 i__2
= kplusp
- i__
+ 1;
3247 F77_FUNC(sgemv
,SGEMV
)("N", n
, &i__2
, &c_b5
, &v
[v_offset
], ldv
, &q
[(*kev
- i__
+ 1) *
3248 q_dim1
+ 1], &c__1
, &c_b4
, &workd
[1], &c__1
);
3249 F77_FUNC(scopy
,SCOPY
)(n
, &workd
[1], &c__1
, &v
[(kplusp
- i__
+ 1) * v_dim1
+ 1], &
3254 F77_FUNC(slacpy
,SLACPY
)("All", n
, kev
, &v
[(*np
+ 1) * v_dim1
+ 1], ldv
, &v
[v_offset
], ldv
);
3256 if (h__
[*kev
+ 1 + h_dim1
] > 0.) {
3257 F77_FUNC(scopy
,SCOPY
)(n
, &workd
[*n
+ 1], &c__1
, &v
[(*kev
+ 1) * v_dim1
+ 1], &c__1
);
3260 F77_FUNC(sscal
,SSCAL
)(n
, &q
[kplusp
+ *kev
* q_dim1
], &resid
[1], &c__1
);
3261 if (h__
[*kev
+ 1 + h_dim1
] > 0.) {
3262 F77_FUNC(saxpy
,SAXPY
)(n
, &h__
[*kev
+ 1 + h_dim1
], &v
[(*kev
+ 1) * v_dim1
+ 1], &c__1
,
3276 F77_FUNC(ssortr
,SSORTR
)(const char * which
,
3291 if (!strncmp(which
, "SA", 2)) {
3298 for (i__
= igap
; i__
<= i__1
; ++i__
) {
3306 if (x1
[j
] < x1
[j
+ igap
]) {
3308 x1
[j
] = x1
[j
+ igap
];
3309 x1
[j
+ igap
] = temp
;
3312 x2
[j
] = x2
[j
+ igap
];
3313 x2
[j
+ igap
] = temp
;
3326 } else if (!strncmp(which
, "SM", 2)) {
3333 for (i__
= igap
; i__
<= i__1
; ++i__
) {
3341 if (fabs(x1
[j
]) < fabs(x1
[j
+ igap
])) {
3343 x1
[j
] = x1
[j
+ igap
];
3344 x1
[j
+ igap
] = temp
;
3347 x2
[j
] = x2
[j
+ igap
];
3348 x2
[j
+ igap
] = temp
;
3361 } else if (!strncmp(which
, "LA", 2)) {
3368 for (i__
= igap
; i__
<= i__1
; ++i__
) {
3376 if (x1
[j
] > x1
[j
+ igap
]) {
3378 x1
[j
] = x1
[j
+ igap
];
3379 x1
[j
+ igap
] = temp
;
3382 x2
[j
] = x2
[j
+ igap
];
3383 x2
[j
+ igap
] = temp
;
3396 } else if (!strncmp(which
, "LM", 2)) {
3404 for (i__
= igap
; i__
<= i__1
; ++i__
) {
3412 if (fabs(x1
[j
]) > fabs(x1
[j
+ igap
])) {
3414 x1
[j
] = x1
[j
+ igap
];
3415 x1
[j
+ igap
] = temp
;
3418 x2
[j
] = x2
[j
+ igap
];
3419 x2
[j
+ igap
] = temp
;
3442 F77_FUNC(ssesrt
,SSESRT
)(const char * which
,
3450 int a_dim1
, a_offset
, i__1
;
3457 a_offset
= 1 + a_dim1
* 0;
3462 if (!strncmp(which
, "SA", 2)) {
3469 for (i__
= igap
; i__
<= i__1
; ++i__
) {
3477 if (x
[j
] < x
[j
+ igap
]) {
3482 F77_FUNC(sswap
,SSWAP
)(na
, &a
[j
* a_dim1
+ 1], &c__1
, &a
[(j
+ igap
) *
3483 a_dim1
+ 1], &c__1
);
3496 } else if (!strncmp(which
, "SM", 2)) {
3503 for (i__
= igap
; i__
<= i__1
; ++i__
) {
3511 if (fabs(x
[j
]) < fabs(x
[j
+ igap
])) {
3516 F77_FUNC(sswap
,SSWAP
)(na
, &a
[j
* a_dim1
+ 1], &c__1
, &a
[(j
+ igap
) *
3517 a_dim1
+ 1], &c__1
);
3530 } else if (!strncmp(which
, "LA", 2)) {
3537 for (i__
= igap
; i__
<= i__1
; ++i__
) {
3545 if (x
[j
] > x
[j
+ igap
]) {
3550 F77_FUNC(sswap
,SSWAP
)(na
, &a
[j
* a_dim1
+ 1], &c__1
, &a
[(j
+ igap
) *
3551 a_dim1
+ 1], &c__1
);
3564 } else if (!strncmp(which
, "LM", 2)) {
3571 for (i__
= igap
; i__
<= i__1
; ++i__
) {
3579 if (fabs(x
[j
]) > fabs(x
[j
+ igap
])) {
3584 F77_FUNC(sswap
,SSWAP
)(na
, &a
[j
* a_dim1
+ 1], &c__1
, &a
[(j
+ igap
) *
3585 a_dim1
+ 1], &c__1
);
3608 F77_FUNC(ssgets
,SSGETS
)(int * ishift
,
3624 if (!strncmp(which
, "BE", 2)) {
3626 F77_FUNC(ssortr
,SSORTR
)("LA", &c__1
, &i__1
, &ritz
[1], &bounds
[1]);
3629 i__1
= (kevd2
<*np
) ? kevd2
: *np
;
3630 i__2
= (kevd2
>*np
) ? kevd2
: *np
;
3631 F77_FUNC(sswap
,SSWAP
)(&i__1
, &ritz
[1], &c__1
,
3632 &ritz
[i__2
+ 1], &c__1
);
3633 i__1
= (kevd2
<*np
) ? kevd2
: *np
;
3634 i__2
= (kevd2
>*np
) ? kevd2
: *np
;
3635 F77_FUNC(sswap
,SSWAP
)(&i__1
, &bounds
[1], &c__1
,
3636 &bounds
[i__2
+ 1], &c__1
);
3641 F77_FUNC(ssortr
,SSORTR
)(which
, &c__1
, &i__1
, &ritz
[1], &bounds
[1]);
3644 if (*ishift
== 1 && *np
> 0) {
3646 F77_FUNC(ssortr
,SSORTR
)("SM", &c__1
, np
, &bounds
[1], &ritz
[1]);
3647 F77_FUNC(scopy
,SCOPY
)(np
, &ritz
[1], &c__1
, &shifts
[1], &c__1
);
3657 F77_FUNC(ssconv
,SSCONV
)(int * n
,
3673 eps23
= GMX_FLOAT_EPS
;
3674 eps23
= pow(eps23
, c_b3
);
3678 for (i__
= 1; i__
<= i__1
; ++i__
) {
3681 d__3
= fabs(ritz
[i__
]);
3682 temp
= (d__2
> d__3
) ? d__2
: d__3
;
3683 if (bounds
[i__
] <= *tol
* temp
) {
3693 F77_FUNC(sseigt
,SSEIGT
)(float * rnorm
,
3703 int h_dim1
, h_offset
, i__1
;
3712 h_offset
= 1 + h_dim1
;
3715 F77_FUNC(scopy
,SCOPY
)(n
, &h__
[(h_dim1
<< 1) + 1], &c__1
, &eig
[1], &c__1
);
3717 F77_FUNC(scopy
,SCOPY
)(&i__1
, &h__
[h_dim1
+ 2], &c__1
, &workl
[1], &c__1
);
3718 F77_FUNC(sstqrb
,SSTQRB
)(n
, &eig
[1], &workl
[1], &bounds
[1], &workl
[*n
+ 1], ierr
);
3724 for (k
= 1; k
<= i__1
; ++k
) {
3725 bounds
[k
] = *rnorm
* fabs(bounds
[k
]);
3738 F77_FUNC(ssaitr
,SSAITR
)(int * ido
,
3762 int h_dim1
, h_offset
, v_dim1
, v_offset
, i__1
;
3766 float safmin
,minval
;
3772 v_offset
= 1 + v_dim1
;
3775 h_offset
= 1 + h_dim1
;
3779 minval
= GMX_FLOAT_MIN
;
3780 safmin
= minval
/ GMX_FLOAT_EPS
;
3793 iwork
[9] = iwork
[8] + *n
;
3794 iwork
[10] = iwork
[9] + *n
;
3797 if (iwork
[5] == 1) {
3800 if (iwork
[6] == 1) {
3803 if (iwork
[2] == 1) {
3806 if (iwork
[3] == 1) {
3809 if (iwork
[4] == 1) {
3826 F77_FUNC(sgetv0
,sgetv0
)(ido
, bmat
, &iwork
[11], &c__0
, n
, &iwork
[12], &v
[v_offset
], ldv
,
3827 &resid
[1], rnorm
, &ipntr
[1], &workd
[1], &iwork
[21], &iwork
[7]);
3833 if (iwork
[11] <= 3) {
3837 *info
= iwork
[12] - 1;
3844 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &v
[iwork
[12] * v_dim1
+ 1], &c__1
);
3845 if (*rnorm
>= safmin
) {
3846 temp1
= 1. / *rnorm
;
3847 F77_FUNC(sscal
,SSCAL
)(n
, &temp1
, &v
[iwork
[12] * v_dim1
+ 1], &c__1
);
3848 F77_FUNC(sscal
,SSCAL
)(n
, &temp1
, &workd
[iwork
[8]], &c__1
);
3851 F77_FUNC(slascl
,SLASCL
)("General", &i__
, &i__
, rnorm
, &c_b18
, n
, &c__1
, &v
[iwork
[12] *
3852 v_dim1
+ 1], n
, &infol
);
3853 F77_FUNC(slascl
,SLASCL
)("General", &i__
, &i__
, rnorm
, &c_b18
, n
, &c__1
, &workd
[iwork
[
3858 F77_FUNC(scopy
,SCOPY
)(n
, &v
[iwork
[12] * v_dim1
+ 1], &c__1
, &workd
[iwork
[10]], &c__1
);
3859 ipntr
[1] = iwork
[10];
3860 ipntr
[2] = iwork
[9];
3861 ipntr
[3] = iwork
[8];
3870 F77_FUNC(scopy
,SCOPY
)(n
, &workd
[iwork
[9]], &c__1
, &resid
[1], &c__1
);
3877 ipntr
[1] = iwork
[9];
3878 ipntr
[2] = iwork
[8];
3882 } else if (*bmat
== 'I') {
3883 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &c__1
);
3892 workd
[*n
* 3 + 3] =F77_FUNC(sdot
,SDOT
)(n
, &resid
[1], &c__1
, &workd
[iwork
[10]], &
3894 workd
[*n
* 3 + 3] = sqrt(fabs(workd
[*n
* 3 + 3]));
3895 } else if (*bmat
== 'G') {
3896 workd
[*n
* 3 + 3] =F77_FUNC(sdot
,SDOT
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &
3898 workd
[*n
* 3 + 3] = sqrt(fabs(workd
[*n
* 3 + 3]));
3899 } else if (*bmat
== 'I') {
3900 workd
[*n
* 3 + 3] =F77_FUNC(snrm2
,SNRM2
)(n
, &resid
[1], &c__1
);
3904 F77_FUNC(sgemv
,SGEMV
)("T", n
, &iwork
[12], &c_b18
, &v
[v_offset
], ldv
, &workd
[iwork
[8]]
3905 , &c__1
, &c_b42
, &workd
[iwork
[9]], &c__1
);
3906 } else if (*mode
== 2) {
3907 F77_FUNC(sgemv
,SGEMV
)("T", n
, &iwork
[12], &c_b18
, &v
[v_offset
], ldv
, &workd
[iwork
[10]
3908 ], &c__1
, &c_b42
, &workd
[iwork
[9]], &c__1
);
3911 F77_FUNC(sgemv
,SGEMV
)("N", n
, &iwork
[12], &c_b50
, &v
[v_offset
], ldv
, &workd
[iwork
[9]], &
3912 c__1
, &c_b18
, &resid
[1], &c__1
);
3914 h__
[iwork
[12] + (h_dim1
<< 1)] = workd
[iwork
[9] + iwork
[12] - 1];
3915 if (iwork
[12] == 1 || iwork
[4] == 1) {
3916 h__
[iwork
[12] + h_dim1
] = 0.;
3918 h__
[iwork
[12] + h_dim1
] = *rnorm
;
3925 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[9]], &c__1
);
3926 ipntr
[1] = iwork
[9];
3927 ipntr
[2] = iwork
[8];
3931 } else if (*bmat
== 'I') {
3932 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &c__1
);
3939 *rnorm
=F77_FUNC(sdot
,SDOT
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &c__1
);
3940 *rnorm
= sqrt(fabs(*rnorm
));
3941 } else if (*bmat
== 'I') {
3942 *rnorm
=F77_FUNC(snrm2
,SNRM2
)(n
, &resid
[1], &c__1
);
3945 if (*rnorm
> workd
[*n
* 3 + 3] * .717f
) {
3951 F77_FUNC(sgemv
,SGEMV
)("T", n
, &iwork
[12], &c_b18
, &v
[v_offset
], ldv
, &workd
[iwork
[8]], &
3952 c__1
, &c_b42
, &workd
[iwork
[9]], &c__1
);
3954 F77_FUNC(sgemv
,SGEMV
)("N", n
, &iwork
[12], &c_b50
, &v
[v_offset
], ldv
, &workd
[iwork
[9]], &
3955 c__1
, &c_b18
, &resid
[1], &c__1
);
3957 if (iwork
[12] == 1 || iwork
[4] == 1) {
3958 h__
[iwork
[12] + h_dim1
] = 0.;
3960 h__
[iwork
[12] + (h_dim1
<< 1)] += workd
[iwork
[9] + iwork
[12] - 1];
3964 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[9]], &c__1
);
3965 ipntr
[1] = iwork
[9];
3966 ipntr
[2] = iwork
[8];
3970 } else if (*bmat
== 'I') {
3971 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &c__1
);
3977 workd
[*n
* 3 + 2] =F77_FUNC(sdot
,SDOT
)(n
, &resid
[1], &c__1
, &workd
[iwork
[8]], &
3979 workd
[*n
* 3 + 2] = sqrt(fabs(workd
[*n
* 3 + 2]));
3980 } else if (*bmat
== 'I') {
3981 workd
[*n
* 3 + 2] =F77_FUNC(snrm2
,SNRM2
)(n
, &resid
[1], &c__1
);
3985 if (workd
[*n
* 3 + 2] > *rnorm
* .717f
) {
3987 *rnorm
= workd
[*n
* 3 + 2];
3991 *rnorm
= workd
[*n
* 3 + 2];
3993 if (iwork
[1] <= 1) {
3998 for (jj
= 1; jj
<= i__1
; ++jj
) {
4009 if (h__
[iwork
[12] + h_dim1
] < 0.) {
4010 h__
[iwork
[12] + h_dim1
] = -h__
[iwork
[12] + h_dim1
];
4011 if (iwork
[12] < *k
+ *np
) {
4012 F77_FUNC(sscal
,SSCAL
)(n
, &c_b50
, &v
[(iwork
[12] + 1) * v_dim1
+ 1], &c__1
);
4014 F77_FUNC(sscal
,SSCAL
)(n
, &c_b50
, &resid
[1], &c__1
);
4019 if (iwork
[12] > *k
+ *np
) {
4038 F77_FUNC(ssaup2
,SSAUP2
)(int * ido
,
4068 int h_dim1
, h_offset
, q_dim1
, q_offset
, v_dim1
, v_offset
, i__1
, i__2
,
4087 v_offset
= 1 + v_dim1
;
4090 h_offset
= 1 + h_dim1
;
4093 q_offset
= 1 + q_dim1
;
4097 eps23
= GMX_FLOAT_EPS
;
4098 eps23
= pow(eps23
, c_b3
);
4110 iwork
[7] = iwork
[9] + iwork
[10];
4128 if (iwork
[2] == 1) {
4129 F77_FUNC(sgetv0
,SGETV0
)(ido
, bmat
, &c__1
, &iwork
[3], n
, &c__1
, &v
[v_offset
], ldv
, &
4130 resid
[1], &workd
[*n
* 3 + 1], &ipntr
[1], &workd
[1], &iwork
[41]
4137 if (workd
[*n
* 3 + 1] == 0.) {
4146 if (iwork
[4] == 1) {
4150 if (iwork
[5] == 1) {
4154 if (iwork
[1] == 1) {
4158 F77_FUNC(ssaitr
,SSAITR
)(ido
, bmat
, n
, &c__0
, &iwork
[9], mode
, &resid
[1], &workd
[*n
* 3 +
4159 1], &v
[v_offset
], ldv
, &h__
[h_offset
], ldh
, &ipntr
[1], &workd
[1],
4183 F77_FUNC(ssaitr
,SSAITR
)(ido
, bmat
, n
, nev
, np
, mode
, &resid
[1], &workd
[*n
* 3 + 1], &v
[
4184 v_offset
], ldv
, &h__
[h_offset
], ldh
, &ipntr
[1], &workd
[1], &iwork
[
4200 F77_FUNC(sseigt
,SSEIGT
)(&workd
[*n
* 3 + 1], &iwork
[7], &h__
[h_offset
], ldh
, &ritz
[1], &
4201 bounds
[1], &workl
[1], &ierr
);
4208 F77_FUNC(scopy
,SCOPY
)(&iwork
[7], &ritz
[1], &c__1
, &workl
[iwork
[7] + 1], &c__1
);
4209 F77_FUNC(scopy
,SCOPY
)(&iwork
[7], &bounds
[1], &c__1
, &workl
[(iwork
[7] << 1) + 1], &c__1
);
4213 F77_FUNC(ssgets
,SSGETS
)(ishift
, which
, nev
, np
, &ritz
[1], &bounds
[1], &workl
[1]);
4215 F77_FUNC(scopy
,SCOPY
)(nev
, &bounds
[*np
+ 1], &c__1
, &workl
[*np
+ 1], &c__1
);
4216 F77_FUNC(ssconv
,SSCONV
)(nev
, &ritz
[*np
+ 1], &workl
[*np
+ 1], tol
, &iwork
[8]);
4221 for (j
= 1; j
<= i__1
; ++j
) {
4222 if (bounds
[j
] == 0.) {
4228 if (iwork
[8] >= iwork
[9] || iwork
[6] > *mxiter
|| *np
== 0) {
4230 if (!strncmp(which
, "BE", 2)) {
4232 strncpy(wprime
, "SA",2);
4233 F77_FUNC(ssortr
,SSORTR
)(wprime
, &c__1
, &iwork
[7], &ritz
[1], &bounds
[1]);
4235 nevm2
= *nev
- nevd2
;
4237 i__1
= (nevd2
< *np
) ? nevd2
: *np
;
4238 i__2
= iwork
[7] - nevd2
+ 1, i__3
= iwork
[7] - *np
+ 1;
4239 F77_FUNC(sswap
,SSWAP
)(&i__1
, &ritz
[nevm2
+ 1], &c__1
,
4240 &ritz
[((i__2
>i__3
) ? i__2
: i__3
)],
4242 i__1
= (nevd2
< *np
) ? nevd2
: *np
;
4243 i__2
= iwork
[7] - nevd2
+ 1, i__3
= iwork
[7] - *np
;
4244 F77_FUNC(sswap
,SSWAP
)(&i__1
, &bounds
[nevm2
+ 1], &c__1
,
4245 &bounds
[((i__2
>i__3
) ? i__2
: i__3
) + 1],
4251 if (!strncmp(which
, "LM", 2)) {
4252 strncpy(wprime
, "SM", 2);
4254 if (!strncmp(which
, "SM", 2)) {
4255 strncpy(wprime
, "LM", 2);
4257 if (!strncmp(which
, "LA", 2)) {
4258 strncpy(wprime
, "SA", 2);
4260 if (!strncmp(which
, "SA", 2)) {
4261 strncpy(wprime
, "LA", 2);
4264 F77_FUNC(ssortr
,SSORTR
)(wprime
, &c__1
, &iwork
[7], &ritz
[1], &bounds
[1]);
4269 for (j
= 1; j
<= i__1
; ++j
) {
4271 d__3
= fabs(ritz
[j
]);
4272 temp
= (d__2
>d__3
) ? d__2
: d__3
;
4276 strncpy(wprime
, "LA",2);
4277 F77_FUNC(ssortr
,SSORTR
)(wprime
, &c__1
, &iwork
[9], &bounds
[1], &ritz
[1]);
4280 for (j
= 1; j
<= i__1
; ++j
) {
4282 d__3
= fabs(ritz
[j
]);
4283 temp
= (d__2
>d__3
) ? d__2
: d__3
;
4287 if (!strncmp(which
, "BE", 2)) {
4289 strncpy(wprime
, "LA", 2);
4290 F77_FUNC(ssortr
,SSORTR
)(wprime
, &c__1
, &iwork
[8], &ritz
[1], &bounds
[1]);
4293 F77_FUNC(ssortr
,SSORTR
)(which
, &c__1
, &iwork
[8], &ritz
[1], &bounds
[1]);
4297 h__
[h_dim1
+ 1] = workd
[*n
* 3 + 1];
4300 if (iwork
[6] > *mxiter
&& iwork
[8] < *nev
) {
4303 if (*np
== 0 && iwork
[8] < iwork
[9]) {
4310 } else if (iwork
[8] < *nev
&& *ishift
== 1) {
4312 i__1
= iwork
[8], i__2
= *np
/ 2;
4313 *nev
+= (i__1
< i__2
) ? i__1
: i__2
;
4314 if (*nev
== 1 && iwork
[7] >= 6) {
4315 *nev
= iwork
[7] / 2;
4316 } else if (*nev
== 1 && iwork
[7] > 2) {
4319 *np
= iwork
[7] - *nev
;
4322 if (nevbef
< *nev
) {
4323 F77_FUNC(ssgets
,SSGETS
)(ishift
, which
, nev
, np
, &ritz
[1], &bounds
[1], &workl
[1]);
4341 F77_FUNC(scopy
,SCOPY
)(np
, &workl
[1], &c__1
, &ritz
[1], &c__1
);
4344 F77_FUNC(ssapps
,SSAPPS
)(n
, nev
, np
, &ritz
[1], &v
[v_offset
], ldv
, &h__
[h_offset
], ldh
, &
4345 resid
[1], &q
[q_offset
], ldq
, &workd
[1]);
4349 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[*n
+ 1], &c__1
);
4355 } else if (*bmat
== 'I') {
4356 F77_FUNC(scopy
,SCOPY
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
4362 workd
[*n
* 3 + 1] =F77_FUNC(sdot
,SDOT
)(n
, &resid
[1], &c__1
, &workd
[1], &c__1
);
4363 workd
[*n
* 3 + 1] = sqrt(fabs(workd
[*n
* 3 + 1]));
4364 } else if (*bmat
== 'I') {
4365 workd
[*n
* 3 + 1] =F77_FUNC(snrm2
,SNRM2
)(n
, &resid
[1], &c__1
);
4387 F77_FUNC(ssaupd
,SSAUPD
)(int * ido
,
4405 int v_dim1
, v_offset
, i__1
, i__2
;
4411 v_offset
= 1 + v_dim1
;
4422 iwork
[5] = iparam
[1];
4423 iwork
[10] = iparam
[3];
4424 iwork
[12] = iparam
[4];
4427 iwork
[11] = iparam
[7];
4432 } else if (*nev
<= 0) {
4434 } else if (*ncv
<= *nev
|| *ncv
> *n
) {
4439 iwork
[15] = *ncv
- *nev
;
4441 if (iwork
[10] <= 0) {
4444 if (strncmp(which
,"LM",2) && strncmp(which
,"SM",2) &&
4445 strncmp(which
,"LA",2) && strncmp(which
,"SA",2) &&
4446 strncmp(which
,"BE",2)) {
4449 if (*bmat
!= 'I' && *bmat
!= 'G') {
4454 if (*lworkl
< i__1
* i__1
+ (*ncv
<< 3)) {
4457 if (iwork
[11] < 1 || iwork
[11] > 5) {
4459 } else if (iwork
[11] == 1 && *bmat
== 'G') {
4461 } else if (iwork
[5] < 0 || iwork
[5] > 1) {
4463 } else if (*nev
== 1 && !strncmp(which
, "BE", 2)) {
4467 if (iwork
[2] != 0) {
4473 if (iwork
[12] <= 0) {
4477 *tol
= GMX_FLOAT_EPS
;
4480 iwork
[15] = *ncv
- *nev
;
4483 i__1
= i__2
* i__2
+ (*ncv
<< 3);
4484 for (j
= 1; j
<= i__1
; ++j
) {
4491 iwork
[16] = iwork
[3] + (iwork
[8] << 1);
4492 iwork
[1] = iwork
[16] + *ncv
;
4493 iwork
[4] = iwork
[1] + *ncv
;
4495 iwork
[7] = iwork
[4] + i__1
* i__1
;
4496 iwork
[14] = iwork
[7] + *ncv
* 3;
4498 ipntr
[4] = iwork
[14];
4499 ipntr
[5] = iwork
[3];
4500 ipntr
[6] = iwork
[16];
4501 ipntr
[7] = iwork
[1];
4502 ipntr
[11] = iwork
[7];
4505 F77_FUNC(ssaup2
,SSAUP2
)(ido
, bmat
, n
, which
, &iwork
[13], &iwork
[15], tol
, &resid
[1], &
4506 iwork
[11], &iwork
[6], &iwork
[5], &iwork
[10], &v
[v_offset
], ldv
, &
4507 workl
[iwork
[3]], &iwork
[8], &workl
[iwork
[16]], &workl
[iwork
[1]], &
4508 workl
[iwork
[4]], &iwork
[9], &workl
[iwork
[7]], &ipntr
[1], &workd
[1]
4509 , &iwork
[21], info
);
4512 iparam
[8] = iwork
[15];
4518 iparam
[3] = iwork
[10];
4519 iparam
[5] = iwork
[15];
4537 F77_FUNC(sseupd
,SSEUPD
)(int * rvec
,
4538 const char * howmny
,
4563 int v_dim1
, v_offset
, z_dim1
, z_offset
, i__1
;
4564 float d__1
, d__2
, d__3
;
4566 int j
, k
, ih
, iq
, iw
, ibd
, ihb
, ihd
, ldh
, ilg
, ldq
, ism
, irz
;
4578 float thres1
=0, thres2
=0;
4582 int leftptr
, rghtptr
;
4588 z_offset
= 1 + z_dim1
;
4593 v_offset
= 1 + v_dim1
;
4617 if (*ncv
<= *nev
|| *ncv
> *n
) {
4620 if (strncmp(which
,"LM",2) && strncmp(which
,"SM",2) &&
4621 strncmp(which
,"LA",2) && strncmp(which
,"SA",2) &&
4622 strncmp(which
,"BE",2)) {
4625 if (*bmat
!= 'I' && *bmat
!= 'G') {
4628 if (*howmny
!= 'A' && *howmny
!= 'P' &&
4629 *howmny
!= 'S' && *rvec
) {
4632 if (*rvec
&& *howmny
== 'S') {
4636 if (*rvec
&& *lworkl
< i__1
* i__1
+ (*ncv
<< 3)) {
4640 if (mode
== 1 || mode
== 2) {
4641 strncpy(type__
, "REGULR",6);
4642 } else if (mode
== 3) {
4643 strncpy(type__
, "SHIFTI",6);
4644 } else if (mode
== 4) {
4645 strncpy(type__
, "BUCKLE",6);
4646 } else if (mode
== 5) {
4647 strncpy(type__
, "CAYLEY",6);
4651 if (mode
== 1 && *bmat
== 'G') {
4654 if (*nev
== 1 && !strncmp(which
, "BE",2)) {
4671 iw
= iq
+ ldh
* *ncv
;
4672 next
= iw
+ (*ncv
<< 1);
4678 irz
= ipntr
[11] + *ncv
;
4682 eps23
= GMX_FLOAT_EPS
;
4683 eps23
= pow(eps23
, c_b21
);
4688 } else if (*bmat
== 'G') {
4689 bnorm2
=F77_FUNC(snrm2
,SNRM2
)(n
, &workd
[1], &c__1
);
4694 if (!strncmp(which
,"LM",2) || !strncmp(which
,"SM",2) ||
4695 !strncmp(which
,"LA",2) || !strncmp(which
,"SA",2)) {
4697 } else if (!strncmp(which
,"BE",2)) {
4700 ism
= (*nev
>nconv
) ? *nev
: nconv
;
4703 thres1
= workl
[ism
];
4704 thres2
= workl
[ilg
];
4712 for (j
= 0; j
<= i__1
; ++j
) {
4714 if (!strncmp(which
,"LM",2)) {
4715 if (fabs(workl
[irz
+ j
]) >= fabs(thres1
)) {
4717 d__3
= fabs(workl
[irz
+ j
]);
4718 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
4719 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
4723 } else if (!strncmp(which
,"SM",2)) {
4724 if (fabs(workl
[irz
+ j
]) <= fabs(thres1
)) {
4726 d__3
= fabs(workl
[irz
+ j
]);
4727 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
4728 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
4732 } else if (!strncmp(which
,"LA",2)) {
4733 if (workl
[irz
+ j
] >= thres1
) {
4735 d__3
= fabs(workl
[irz
+ j
]);
4736 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
4737 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
4741 } else if (!strncmp(which
,"SA",2)) {
4742 if (workl
[irz
+ j
] <= thres1
) {
4744 d__3
= fabs(workl
[irz
+ j
]);
4745 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
4746 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
4750 } else if (!strncmp(which
,"BE",2)) {
4751 if (workl
[irz
+ j
] <= thres1
|| workl
[irz
+ j
] >= thres2
) {
4753 d__3
= fabs(workl
[irz
+ j
]);
4754 tempbnd
= (d__2
>d__3
) ? d__2
: d__3
;
4755 if (workl
[ibd
+ j
] <= *tol
* tempbnd
) {
4760 if (j
+ 1 > nconv
) {
4761 reord
= select
[j
+ 1] || reord
;
4763 if (select
[j
+ 1]) {
4769 F77_FUNC(scopy
,SCOPY
)(&i__1
, &workl
[ih
+ 1], &c__1
, &workl
[ihb
], &c__1
);
4770 F77_FUNC(scopy
,SCOPY
)(ncv
, &workl
[ih
+ ldh
], &c__1
, &workl
[ihd
], &c__1
);
4772 F77_FUNC(ssteqr
,SSTEQR
)("Identity", ncv
, &workl
[ihd
], &workl
[ihb
], &workl
[iq
], &ldq
, &
4791 if (select
[leftptr
]) {
4795 } else if (! select
[rghtptr
]) {
4801 temp
= workl
[ihd
+ leftptr
- 1];
4802 workl
[ihd
+ leftptr
- 1] = workl
[ihd
+ rghtptr
- 1];
4803 workl
[ihd
+ rghtptr
- 1] = temp
;
4804 F77_FUNC(scopy
,SCOPY
)(ncv
, &workl
[iq
+ *ncv
* (leftptr
- 1)], &c__1
, &workl
[
4806 F77_FUNC(scopy
,SCOPY
)(ncv
, &workl
[iq
+ *ncv
* (rghtptr
- 1)], &c__1
, &workl
[
4807 iq
+ *ncv
* (leftptr
- 1)], &c__1
);
4808 F77_FUNC(scopy
,SCOPY
)(ncv
, &workl
[iw
], &c__1
, &workl
[iq
+ *ncv
* (rghtptr
-
4815 if (leftptr
< rghtptr
) {
4823 F77_FUNC(scopy
,SCOPY
)(&nconv
, &workl
[ihd
], &c__1
, &d__
[1], &c__1
);
4827 F77_FUNC(scopy
,SCOPY
)(&nconv
, &workl
[ritz
], &c__1
, &d__
[1], &c__1
);
4828 F77_FUNC(scopy
,SCOPY
)(ncv
, &workl
[ritz
], &c__1
, &workl
[ihd
], &c__1
);
4831 if (!strncmp(type__
, "REGULR",6)) {
4834 F77_FUNC(ssesrt
,SSESRT
)("LA", rvec
, &nconv
, &d__
[1], ncv
, &workl
[iq
], &ldq
);
4836 F77_FUNC(scopy
,SCOPY
)(ncv
, &workl
[bounds
], &c__1
, &workl
[ihb
], &c__1
);
4841 F77_FUNC(scopy
,SCOPY
)(ncv
, &workl
[ihd
], &c__1
, &workl
[iw
], &c__1
);
4842 if (!strncmp(type__
, "SHIFTI", 6)) {
4844 for (k
= 1; k
<= i__1
; ++k
) {
4845 workl
[ihd
+ k
- 1] = 1. / workl
[ihd
+ k
- 1] + *sigma
;
4847 } else if (!strncmp(type__
, "BUCKLE",6)) {
4849 for (k
= 1; k
<= i__1
; ++k
) {
4850 workl
[ihd
+ k
- 1] = *sigma
* workl
[ihd
+ k
- 1] / (workl
[ihd
4853 } else if (!strncmp(type__
, "CAYLEY",6)) {
4855 for (k
= 1; k
<= i__1
; ++k
) {
4856 workl
[ihd
+ k
- 1] = *sigma
* (workl
[ihd
+ k
- 1] + 1.) / (
4857 workl
[ihd
+ k
- 1] - 1.);
4861 F77_FUNC(scopy
,SCOPY
)(&nconv
, &workl
[ihd
], &c__1
, &d__
[1], &c__1
);
4862 F77_FUNC(ssortr
,SSORTR
)("LA", &c__1
, &nconv
, &workl
[ihd
], &workl
[iw
]);
4864 F77_FUNC(ssesrt
,SSESRT
)("LA", rvec
, &nconv
, &d__
[1], ncv
, &workl
[iq
], &ldq
);
4866 F77_FUNC(scopy
,SCOPY
)(ncv
, &workl
[bounds
], &c__1
, &workl
[ihb
], &c__1
);
4867 d__1
= bnorm2
/ rnorm
;
4868 F77_FUNC(sscal
,SSCAL
)(ncv
, &d__1
, &workl
[ihb
], &c__1
);
4869 F77_FUNC(ssortr
,SSORTR
)("LA", &c__1
, &nconv
, &d__
[1], &workl
[ihb
]);
4874 if (*rvec
&& *howmny
== 'A') {
4876 F77_FUNC(sgeqr2
,SGEQR2
)(ncv
, &nconv
, &workl
[iq
], &ldq
, &workl
[iw
+ *ncv
], &workl
[ihb
],
4879 F77_FUNC(sorm2r
,SORM2R
)("Right", "Notranspose", n
, ncv
, &nconv
, &workl
[iq
], &ldq
, &
4880 workl
[iw
+ *ncv
], &v
[v_offset
], ldv
, &workd
[*n
+ 1], &ierr
);
4881 F77_FUNC(slacpy
,SLACPY
)("All", n
, &nconv
, &v
[v_offset
], ldv
, &z__
[z_offset
], ldz
);
4884 for (j
= 1; j
<= i__1
; ++j
) {
4885 workl
[ihb
+ j
- 1] = 0.;
4887 workl
[ihb
+ *ncv
- 1] = 1.;
4888 F77_FUNC(sorm2r
,SORM2R
)("Left", "Transpose", ncv
, &c__1
, &nconv
, &workl
[iq
], &ldq
, &
4889 workl
[iw
+ *ncv
], &workl
[ihb
], ncv
, &temp
, &ierr
);
4891 } else if (*rvec
&& *howmny
== 'S') {
4895 if (!strncmp(type__
, "REGULR",6) && *rvec
) {
4898 for (j
= 1; j
<= i__1
; ++j
) {
4899 workl
[ihb
+ j
- 1] = rnorm
* fabs(workl
[ihb
+ j
- 1]);
4902 } else if (strncmp(type__
, "REGULR",6) && *rvec
) {
4904 F77_FUNC(sscal
,SSCAL
)(ncv
, &bnorm2
, &workl
[ihb
], &c__1
);
4905 if (!strncmp(type__
, "SHIFTI",6)) {
4908 for (k
= 1; k
<= i__1
; ++k
) {
4909 d__2
= workl
[iw
+ k
- 1];
4910 workl
[ihb
+ k
- 1] = fabs(workl
[ihb
+ k
- 1])/(d__2
* d__2
);
4913 } else if (!strncmp(type__
, "BUCKLE",6)) {
4916 for (k
= 1; k
<= i__1
; ++k
) {
4917 d__2
= workl
[iw
+ k
- 1] - 1.;
4918 workl
[ihb
+ k
- 1] = *sigma
* fabs(workl
[ihb
+ k
- 1])/(d__2
* d__2
);
4921 } else if (!strncmp(type__
, "CAYLEY",6)) {
4924 for (k
= 1; k
<= i__1
; ++k
) {
4925 workl
[ihb
+ k
- 1] = fabs(workl
[ihb
+ k
- 1] / workl
[iw
+ k
- 1] * (workl
[iw
+ k
- 1] - 1.));
4933 if (*rvec
&& (!strncmp(type__
, "SHIFTI",6) || !strncmp(type__
, "CAYLEY",6))) {
4936 for (k
= 0; k
<= i__1
; ++k
) {
4937 workl
[iw
+ k
] = workl
[iq
+ k
* ldq
+ *ncv
- 1] / workl
[iw
+ k
];
4940 } else if (*rvec
&& !strncmp(type__
, "BUCKLE", 6)) {
4943 for (k
= 0; k
<= i__1
; ++k
) {
4944 workl
[iw
+ k
] = workl
[iq
+ k
* ldq
+ *ncv
- 1] / (workl
[iw
+ k
] -
4950 if (strncmp(type__
, "REGULR",6)) {
4951 F77_FUNC(sger
,SGER
)(n
, &nconv
, &c_b102
, &resid
[1], &c__1
, &workl
[iw
], &c__1
, &z__
[