TBCI Numerical high perf. C++ Library 2.8.0
own_lapack_routines.cpp
Go to the documentation of this file.
1
5/* ev.f -- translated by f2c (version 19940714).
6 You must link the resulting object file with the libraries:
7 -lf2c -lm (in that order)
8*/
9
10#ifdef __cplusplus
11extern "C" {
12#endif
13#include "tbci/lapack/f2c.h"
14#include <stdio.h>
15
16#ifdef _DOXYGEN
17# define integer long int
18# define doublereal double
19# define doublecomplex __complex__ double
20# define real float
21# define complex __complex__ float
22#endif
23
24/* Table of constant values */
25
26static integer c__9 = 9;
27static integer c__1 = 1;
28static doublereal c_b14 = 1.;
29static doublecomplex c_b26 = {0.,0.};
30static doublecomplex c_b27 = {1.,0.};
31
169
171 ab, integer *ldab, doublecomplex *bb, integer *ldbb, doublereal *vl,
172 doublereal *vu, doublereal *w, doublecomplex *q, integer *ldq,
173 doublecomplex *work, doublereal *rwork, integer *iwork)
174{
175 /* System generated locals */
176 integer ab_dim1, ab_offset, bb_dim1, bb_offset, q_dim1, q_offset, ret_val,
177 i__1;
178 doublereal d__1, d__2;
179
180 /* Builtin functions */
181 double sqrt(doublereal);
182 integer s_wsle(cilist *), do_lio(integer *, integer *, const char *, ftnlen),
183 e_wsle();
184
185 /* Local variables */
186 static integer indd, info;
187 static doublereal anrm;
188 static char vect[1], jobz[1];
189 static doublereal rmin, rmax;
190 static char uplo[1];
191 static integer indnd;
192 static doublereal sigma;
193 static char order[1];
194 static logical upper, wantz;
195 static integer ka, kb, il;
196 extern doublereal dlamch_(const char *, ftnlen);
197 static integer iu, iscale, indibl;
198 static logical valeig;
199 static doublereal safmin;
200 extern doublereal zlanhb_(const char *, const char *, integer *, integer *,
202 extern /* Subroutine */ int xerbla_(const char *, integer *, ftnlen);
203 static doublereal abstll, bignum, abstol;
204 static integer indiwk, indisp;
205 extern /* Subroutine */ int zlascl_(const char *, integer *, integer *,
208 dstebz_(const char *, const char *, integer *, doublereal *,
212 zhbtrd_(const char *, const char *, integer *, integer *,
215 ftnlen);
216 static integer indrwk;
217 extern /* Subroutine */ int zhbgst_(const char *, const char *, integer *,
221 zpbstf_(const char *, integer *, integer *,
223 static integer nsplit;
224 static doublereal smlnum, eps, vll, vuu;
225
226 /* Fortran I/O blocks */
227 static cilist io___22 = { 0, 6, 0, 0, 0 };
228
229
230
231/* .. Executable Statements .. */
232
233/* Test the input parameters. */
234
235 /* Parameter adjustments */
236 ab_dim1 = *ldab;
237 ab_offset = ab_dim1 + 1;
238 ab -= ab_offset;
239 bb_dim1 = *ldbb;
240 bb_offset = bb_dim1 + 1;
241 bb -= bb_offset;
242 --w;
243 q_dim1 = *ldq;
244 q_offset = q_dim1 + 1;
245 q -= q_offset;
246 --work;
247 --rwork;
248 --iwork;
249
250 /* Function Body */
251 wantz = *job == 1;
252 if (wantz) {
253 *jobz = 'V';
254 } else {
255 *jobz = 'N';
256 }
257 upper = *up == 1;
258 if (upper) {
259 *uplo = 'U';
260 } else {
261 *uplo = 'L';
262 }
263 valeig = TRUE_;
264/* Werte setzen */
265 ka = *ldab - 1;
266 kb = *ldbb - 1;
267
268 info = 0;
269 if (*n < 0) {
270 info = -3;
271 } else if (ka < 0) {
272 info = -4;
273 } else if (kb < 0 || kb > ka) {
274 info = -5;
275 } else if (*ldab < ka + 1) {
276 info = -7;
277 } else if (*ldbb < kb + 1) {
278 info = -9;
279 } else if (*n <= 0) {
280 info = -110;
281 }
282 if (info != 0) {
283 i__1 = -info;
284 xerbla_("ZHBGVX", &i__1, 6L);
285 ret_val = -info;
286 return ret_val;
287 }
288
289/* Get machine constants. */
290
291 safmin = dlamch_("Safe minimum", 12L);
292 eps = dlamch_("Precision", 9L);
293 smlnum = safmin / eps;
294 bignum = 1. / smlnum;
295 rmin = sqrt(smlnum);
296/* Computing MIN */
297 d__1 = sqrt(bignum), d__2 = 1. / sqrt(sqrt(safmin));
298 rmax = min(d__1,d__2);
299
300/* Fuer maximale Genauigkeit der Eigenwerte wird gesetzt: */
301 abstol = safmin * 2.;
302
303/* Scale both matrices to allowable range, if necessary. */
304
305 iscale = 0;
306 abstll = abstol;
307 vll = *vl;
308 vuu = *vu;
309 anrm = zlanhb_("M", uplo, n, &ka, &ab[ab_offset], ldab, &rwork[1], 1L, 1L)
310 ;
311 if (anrm > 0. && anrm < rmin) {
312 iscale = 1;
313 sigma = rmin / anrm;
314 } else if (anrm > rmax) {
315 iscale = 1;
316 sigma = rmax / anrm;
317 }
318 if (iscale == 1) {
319 s_wsle(&io___22);
320 do_lio(&c__9, &c__1, "Matrix in bad condition, scale matrix ...", 49L);
321 e_wsle();
322 if (! upper) {
323 zlascl_("B", &ka, &ka, &c_b14, &sigma, n, n, &ab[ab_offset], ldab,
324 &info, 1L);
325 zlascl_("B", &kb, &kb, &c_b14, &sigma, n, n, &bb[bb_offset], ldbb,
326 &info, 1L);
327 } else {
328 zlascl_("Q", &ka, &ka, &c_b14, &sigma, n, n, &ab[ab_offset], ldab,
329 &info, 1L);
330 zlascl_("Q", &kb, &kb, &c_b14, &sigma, n, n, &bb[bb_offset], ldbb,
331 &info, 1L);
332 }
333 if (abstol > 0.) {
334 abstll = abstol * sigma;
335 }
336 vll = *vl * sigma;
337 vuu = *vu * sigma;
338 }
339
340/* Form a split Cholesky factorization of B. */
341
342 zpbstf_(uplo, n, &kb, &bb[bb_offset], ldbb, &info, 1L);
343
344 //printf("Cholesky factorization performed..\n");
345 if (info != 0) {
346 info = *n + info;
347 return info;
348 }
349
350
351/* Transform problem to standard eigenvalue problem. */
352
353 indd = 1;
354 indnd = indd + *n;
355 indrwk = *n << 1;
356 zhbgst_(jobz, uplo, n, &ka, &kb, &ab[ab_offset], ldab, &bb[bb_offset],
357 ldbb, &q[q_offset], ldq, &work[1], &rwork[indrwk], &info, 1L, 1L);
358 //printf("Transformed to standard eigenvalueform..\n");
359
360 if (info != 0) return info;
361
362
363
364
365/* Reduce to tridiagonal form. */
366
367 if (wantz) {
368 *vect = 'U';
369 } else {
370 *vect = 'N';
371 }
372 zhbtrd_(vect, uplo, n, &ka, &ab[ab_offset], ldab, &rwork[indd], &rwork[
373 indnd], &q[q_offset], ldq, &work[1], &info, 1L, 1L);
374
375
376 //printf("Reduced to tridiagonal form..\n");
377
378
379/* Bei der Rueckgabe enthält nach dieser Routine RWORK(IndD) */
380/* die Diagonalelemente,RWORK(IndND) die Nebendiagonale der */
381/* Tridiagonalmatrix */
382
383/* Call of the Eigenvaluesolver DSTEBZ */
384
385 if (wantz) {
386 *order = 'B';
387 } else {
388 *order = 'E';
389 }
390
391 indibl = 1;
392 indisp = indibl + *n;
393 indiwk = indisp + *n;
394 dstebz_("V", order, n, &vll, &vuu, &il, &iu, &abstll, &rwork[indd], &
395 rwork[indnd], m, &nsplit, &w[1], &iwork[indibl], &iwork[indisp], &
396 rwork[indrwk], &iwork[indiwk], &info, 1L, 1L);
397
398 //printf("Eigenvaluesolver DSTEBZ called..\n");
399
400
401/* End of INTEGER FUNCTION ew */
402 ret_val = info;
403 return ret_val;
404} /* ew_ */
405
406/* ********************************************************************** */
453 ldz, doublecomplex *q, integer *ldq, doublecomplex *work, doublereal *
454 rwork, integer *iwork, integer *ifail)
455{
456 /* System generated locals */
457 integer z_dim1, z_offset, q_dim1, q_offset, ret_val, i__1, i__2;
458
459 /* Local variables */
460 static integer indd, info, itmp1, i, j, indnd;
461 extern /* Subroutine */ int zgemv_(const char *, integer *, integer *,
464 zcopy_(integer *, doublecomplex *, integer *, doublecomplex *,
465 integer *),
466 zswap_(integer *, doublecomplex *, integer *, doublecomplex *,
467 integer *);
468 static integer jj, indibl, indiwk, indisp, indrwk;
469 extern /* Subroutine */ int zstein_(integer *, doublereal *, doublereal *,
472 static doublereal tmp1;
473
474
475
476/* =====================================================================
477*/
478/* .. Local Scalars .. */
479
480
481/* =====================================================================
482*/
483/* Executable statements */
484
485 /* Parameter adjustments */
486 --w;
487 z_dim1 = *ldz;
488 z_offset = z_dim1 + 1;
489 z -= z_offset;
490 q_dim1 = *ldq;
491 q_offset = q_dim1 + 1;
492 q -= q_offset;
493 --work;
494 --rwork;
495 --iwork;
496 --ifail;
497
498 /* Function Body */
499 indd = 1;
500 indnd = indd + *n;
501 indrwk = *n << 1;
502 indibl = 1;
503 indisp = indibl + *n;
504 indiwk = indisp + *n;
505
506 zstein_(n, &rwork[indd], &rwork[indnd], m, &w[1], &iwork[indibl], &iwork[
507 indisp], &z[z_offset], ldz, &rwork[indrwk], &iwork[indiwk], &
508 ifail[1], &info);
509
510/* Apply unitary matrix used in reduction to tridiagonal */
511/* form to eigenvectors returned by ZSTEIN. */
512
513 i__1 = *m;
514 for (j = 1; j <= i__1; ++j) {
515 zcopy_(n, &z[j * z_dim1 + 1], &c__1, &work[1], &c__1);
516 zgemv_("N", n, n, &c_b27, &q[q_offset], ldq, &work[1], &c__1, &c_b26,
517 &z[j * z_dim1 + 1], &c__1, 1L);
518 }
519
520/* If eigenvalues are not in order, then sort them, along with */
521/* eigenvectors. */
522
523 i__1 = *m - 1;
524 for (j = 1; j <= i__1; ++j) {
525 i = 0;
526 tmp1 = w[j];
527 i__2 = *m;
528 for (jj = j + 1; jj <= i__2; ++jj) {
529 if (w[jj] < tmp1) {
530 i = jj;
531 tmp1 = w[jj];
532 }
533 }
534
535 if (i != 0) {
536 itmp1 = iwork[indibl + i - 1];
537 w[i] = w[j];
538 iwork[indibl + i - 1] = iwork[indibl + j - 1];
539 w[j] = tmp1;
540 iwork[indibl + j - 1] = itmp1;
541 zswap_(n, &z[i * z_dim1 + 1], &c__1, &z[j * z_dim1 + 1], &c__1);
542 if (info != 0) {
543 itmp1 = ifail[i];
544 ifail[i] = ifail[j];
545 ifail[j] = itmp1;
546 }
547 }
548 }
549
550/* End of EV */
551
552 ret_val = info;
553 return ret_val;
554} /* ev_ */
555
556#ifdef __cplusplus
557 }
558#endif
int i
Definition LM_fit.h:71
NAMESPACE_END NAMESPACE_CSTD TBCI__ cplx< T > sqrt(const TBCI__ cplx< T > &z)
Definition cplx.h:751
int xerbla_(char *, int *)
long int ftnlen
Definition f2c.h:67
int logical
Definition f2c.h:35
#define min(a, b)
Definition f2c.h:180
#define TRUE_
Definition f2c.h:49
double dlamch_(const char *)
static __complex__ double c_b27
#define doublereal
static long int c__1
static __complex__ double c_b26
#define integer
long int own_ev_(long int *n, long int *m, double *w, __complex__ double *z, long int *ldz, __complex__ double *q, long int *ldq, __complex__ double *work, double *rwork, long int *iwork, long int *ifail)
– AHLAND driver routine (version 2.0) –
static double c_b14
long int own_ew_(long int *job, long int *up, long int *n, long int *m, __complex__ double *ab, long int *ldab, __complex__ double *bb, long int *ldbb, double *vl, double *vu, double *w, __complex__ double *q, long int *ldq, __complex__ double *work, double *rwork, long int *iwork)
– AHLAND driver routine (version 2.0) – Modified LAPACK-ROUTINE for the calculation of selected Eigen...
#define doublecomplex
static long int c__9
integer e_wsle()
integer do_lio(integer *, integer *, char *, ftnlen)
integer s_wsle(cilist *)
Definition f2c.h:73