zla_gbrcond_c.c
Go to the documentation of this file.
00001 /* zla_gbrcond_c.f -- translated by f2c (version 20061008).
00002    You must link the resulting object file with libf2c:
00003         on Microsoft Windows system, link with libf2c.lib;
00004         on Linux or Unix systems, link with .../path/to/libf2c.a -lm
00005         or, if you install libf2c.a in a standard place, with -lf2c -lm
00006         -- in that order, at the end of the command line, as in
00007                 cc *.o -lf2c -lm
00008         Source for libf2c is in /netlib/f2c/libf2c.zip, e.g.,
00009 
00010                 http://www.netlib.org/f2c/libf2c.zip
00011 */
00012 
00013 #include "f2c.h"
00014 #include "blaswrap.h"
00015 
00016 /* Table of constant values */
00017 
00018 static integer c__1 = 1;
00019 
00020 doublereal zla_gbrcond_c__(char *trans, integer *n, integer *kl, integer *ku, 
00021         doublecomplex *ab, integer *ldab, doublecomplex *afb, integer *ldafb, 
00022         integer *ipiv, doublereal *c__, logical *capply, integer *info, 
00023         doublecomplex *work, doublereal *rwork, ftnlen trans_len)
00024 {
00025     /* System generated locals */
00026     integer ab_dim1, ab_offset, afb_dim1, afb_offset, i__1, i__2, i__3, i__4;
00027     doublereal ret_val, d__1, d__2;
00028     doublecomplex z__1;
00029 
00030     /* Builtin functions */
00031     double d_imag(doublecomplex *);
00032 
00033     /* Local variables */
00034     integer i__, j, kd, ke;
00035     doublereal tmp;
00036     integer kase;
00037     extern logical lsame_(char *, char *);
00038     integer isave[3];
00039     doublereal anorm;
00040     extern /* Subroutine */ int zlacn2_(integer *, doublecomplex *, 
00041             doublecomplex *, doublereal *, integer *, integer *), xerbla_(
00042             char *, integer *);
00043     doublereal ainvnm;
00044     extern /* Subroutine */ int zgbtrs_(char *, integer *, integer *, integer 
00045             *, integer *, doublecomplex *, integer *, integer *, 
00046             doublecomplex *, integer *, integer *);
00047     logical notrans;
00048 
00049 
00050 /*     -- LAPACK routine (version 3.2.1)                               -- */
00051 /*     -- Contributed by James Demmel, Deaglan Halligan, Yozo Hida and -- */
00052 /*     -- Jason Riedy of Univ. of California Berkeley.                 -- */
00053 /*     -- April 2009                                                   -- */
00054 
00055 /*     -- LAPACK is a software package provided by Univ. of Tennessee, -- */
00056 /*     -- Univ. of California Berkeley and NAG Ltd.                    -- */
00057 
00058 /*     .. */
00059 /*     .. Scalar Arguments .. */
00060 /*     .. */
00061 /*     .. Array Arguments .. */
00062 
00063 
00064 /*  Purpose */
00065 /*  ======= */
00066 
00067 /*     ZLA_GBRCOND_C Computes the infinity norm condition number of */
00068 /*     op(A) * inv(diag(C)) where C is a DOUBLE PRECISION vector. */
00069 
00070 /*  Arguments */
00071 /*  ========= */
00072 
00073 /*     TRANS   (input) CHARACTER*1 */
00074 /*     Specifies the form of the system of equations: */
00075 /*       = 'N':  A * X = B     (No transpose) */
00076 /*       = 'T':  A**T * X = B  (Transpose) */
00077 /*       = 'C':  A**H * X = B  (Conjugate Transpose = Transpose) */
00078 
00079 /*     N       (input) INTEGER */
00080 /*     The number of linear equations, i.e., the order of the */
00081 /*     matrix A.  N >= 0. */
00082 
00083 /*     KL      (input) INTEGER */
00084 /*     The number of subdiagonals within the band of A.  KL >= 0. */
00085 
00086 /*     KU      (input) INTEGER */
00087 /*     The number of superdiagonals within the band of A.  KU >= 0. */
00088 
00089 /*     AB      (input) COMPLEX*16 array, dimension (LDAB,N) */
00090 /*     On entry, the matrix A in band storage, in rows 1 to KL+KU+1. */
00091 /*     The j-th column of A is stored in the j-th column of the */
00092 /*     array AB as follows: */
00093 /*     AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl) */
00094 
00095 /*     LDAB    (input) INTEGER */
00096 /*     The leading dimension of the array AB.  LDAB >= KL+KU+1. */
00097 
00098 /*     AFB     (input) COMPLEX*16 array, dimension (LDAFB,N) */
00099 /*     Details of the LU factorization of the band matrix A, as */
00100 /*     computed by ZGBTRF.  U is stored as an upper triangular */
00101 /*     band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, */
00102 /*     and the multipliers used during the factorization are stored */
00103 /*     in rows KL+KU+2 to 2*KL+KU+1. */
00104 
00105 /*     LDAFB   (input) INTEGER */
00106 /*     The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1. */
00107 
00108 /*     IPIV    (input) INTEGER array, dimension (N) */
00109 /*     The pivot indices from the factorization A = P*L*U */
00110 /*     as computed by ZGBTRF; row i of the matrix was interchanged */
00111 /*     with row IPIV(i). */
00112 
00113 /*     C       (input) DOUBLE PRECISION array, dimension (N) */
00114 /*     The vector C in the formula op(A) * inv(diag(C)). */
00115 
00116 /*     CAPPLY  (input) LOGICAL */
00117 /*     If .TRUE. then access the vector C in the formula above. */
00118 
00119 /*     INFO    (output) INTEGER */
00120 /*       = 0:  Successful exit. */
00121 /*     i > 0:  The ith argument is invalid. */
00122 
00123 /*     WORK    (input) COMPLEX*16 array, dimension (2*N). */
00124 /*     Workspace. */
00125 
00126 /*     RWORK   (input) DOUBLE PRECISION array, dimension (N). */
00127 /*     Workspace. */
00128 
00129 /*  ===================================================================== */
00130 
00131 /*     .. Local Scalars .. */
00132 /*     .. */
00133 /*     .. Local Arrays .. */
00134 /*     .. */
00135 /*     .. External Functions .. */
00136 /*     .. */
00137 /*     .. External Subroutines .. */
00138 /*     .. */
00139 /*     .. Intrinsic Functions .. */
00140 /*     .. */
00141 /*     .. Statement Functions .. */
00142 /*     .. */
00143 /*     .. Statement Function Definitions .. */
00144 /*     .. */
00145 /*     .. Executable Statements .. */
00146     /* Parameter adjustments */
00147     ab_dim1 = *ldab;
00148     ab_offset = 1 + ab_dim1;
00149     ab -= ab_offset;
00150     afb_dim1 = *ldafb;
00151     afb_offset = 1 + afb_dim1;
00152     afb -= afb_offset;
00153     --ipiv;
00154     --c__;
00155     --work;
00156     --rwork;
00157 
00158     /* Function Body */
00159     ret_val = 0.;
00160 
00161     *info = 0;
00162     notrans = lsame_(trans, "N");
00163     if (! notrans && ! lsame_(trans, "T") && ! lsame_(
00164             trans, "C")) {
00165         *info = -1;
00166     } else if (*n < 0) {
00167         *info = -2;
00168     } else if (*kl < 0 || *kl > *n - 1) {
00169         *info = -3;
00170     } else if (*ku < 0 || *ku > *n - 1) {
00171         *info = -4;
00172     } else if (*ldab < *kl + *ku + 1) {
00173         *info = -6;
00174     } else if (*ldafb < (*kl << 1) + *ku + 1) {
00175         *info = -8;
00176     }
00177     if (*info != 0) {
00178         i__1 = -(*info);
00179         xerbla_("ZLA_GBRCOND_C", &i__1);
00180         return ret_val;
00181     }
00182 
00183 /*     Compute norm of op(A)*op2(C). */
00184 
00185     anorm = 0.;
00186     kd = *ku + 1;
00187     ke = *kl + 1;
00188     if (notrans) {
00189         i__1 = *n;
00190         for (i__ = 1; i__ <= i__1; ++i__) {
00191             tmp = 0.;
00192             if (*capply) {
00193 /* Computing MAX */
00194                 i__2 = i__ - *kl;
00195 /* Computing MIN */
00196                 i__4 = i__ + *ku;
00197                 i__3 = min(i__4,*n);
00198                 for (j = max(i__2,1); j <= i__3; ++j) {
00199                     i__2 = kd + i__ - j + j * ab_dim1;
00200                     tmp += ((d__1 = ab[i__2].r, abs(d__1)) + (d__2 = d_imag(&
00201                             ab[kd + i__ - j + j * ab_dim1]), abs(d__2))) / 
00202                             c__[j];
00203                 }
00204             } else {
00205 /* Computing MAX */
00206                 i__3 = i__ - *kl;
00207 /* Computing MIN */
00208                 i__4 = i__ + *ku;
00209                 i__2 = min(i__4,*n);
00210                 for (j = max(i__3,1); j <= i__2; ++j) {
00211                     i__3 = kd + i__ - j + j * ab_dim1;
00212                     tmp += (d__1 = ab[i__3].r, abs(d__1)) + (d__2 = d_imag(&
00213                             ab[kd + i__ - j + j * ab_dim1]), abs(d__2));
00214                 }
00215             }
00216             rwork[i__] = tmp;
00217             anorm = max(anorm,tmp);
00218         }
00219     } else {
00220         i__1 = *n;
00221         for (i__ = 1; i__ <= i__1; ++i__) {
00222             tmp = 0.;
00223             if (*capply) {
00224 /* Computing MAX */
00225                 i__2 = i__ - *kl;
00226 /* Computing MIN */
00227                 i__4 = i__ + *ku;
00228                 i__3 = min(i__4,*n);
00229                 for (j = max(i__2,1); j <= i__3; ++j) {
00230                     i__2 = ke - i__ + j + i__ * ab_dim1;
00231                     tmp += ((d__1 = ab[i__2].r, abs(d__1)) + (d__2 = d_imag(&
00232                             ab[ke - i__ + j + i__ * ab_dim1]), abs(d__2))) / 
00233                             c__[j];
00234                 }
00235             } else {
00236 /* Computing MAX */
00237                 i__3 = i__ - *kl;
00238 /* Computing MIN */
00239                 i__4 = i__ + *ku;
00240                 i__2 = min(i__4,*n);
00241                 for (j = max(i__3,1); j <= i__2; ++j) {
00242                     i__3 = ke - i__ + j + i__ * ab_dim1;
00243                     tmp += (d__1 = ab[i__3].r, abs(d__1)) + (d__2 = d_imag(&
00244                             ab[ke - i__ + j + i__ * ab_dim1]), abs(d__2));
00245                 }
00246             }
00247             rwork[i__] = tmp;
00248             anorm = max(anorm,tmp);
00249         }
00250     }
00251 
00252 /*     Quick return if possible. */
00253 
00254     if (*n == 0) {
00255         ret_val = 1.;
00256         return ret_val;
00257     } else if (anorm == 0.) {
00258         return ret_val;
00259     }
00260 
00261 /*     Estimate the norm of inv(op(A)). */
00262 
00263     ainvnm = 0.;
00264 
00265     kase = 0;
00266 L10:
00267     zlacn2_(n, &work[*n + 1], &work[1], &ainvnm, &kase, isave);
00268     if (kase != 0) {
00269         if (kase == 2) {
00270 
00271 /*           Multiply by R. */
00272 
00273             i__1 = *n;
00274             for (i__ = 1; i__ <= i__1; ++i__) {
00275                 i__2 = i__;
00276                 i__3 = i__;
00277                 i__4 = i__;
00278                 z__1.r = rwork[i__4] * work[i__3].r, z__1.i = rwork[i__4] * 
00279                         work[i__3].i;
00280                 work[i__2].r = z__1.r, work[i__2].i = z__1.i;
00281             }
00282 
00283             if (notrans) {
00284                 zgbtrs_("No transpose", n, kl, ku, &c__1, &afb[afb_offset], 
00285                         ldafb, &ipiv[1], &work[1], n, info);
00286             } else {
00287                 zgbtrs_("Conjugate transpose", n, kl, ku, &c__1, &afb[
00288                         afb_offset], ldafb, &ipiv[1], &work[1], n, info);
00289             }
00290 
00291 /*           Multiply by inv(C). */
00292 
00293             if (*capply) {
00294                 i__1 = *n;
00295                 for (i__ = 1; i__ <= i__1; ++i__) {
00296                     i__2 = i__;
00297                     i__3 = i__;
00298                     i__4 = i__;
00299                     z__1.r = c__[i__4] * work[i__3].r, z__1.i = c__[i__4] * 
00300                             work[i__3].i;
00301                     work[i__2].r = z__1.r, work[i__2].i = z__1.i;
00302                 }
00303             }
00304         } else {
00305 
00306 /*           Multiply by inv(C'). */
00307 
00308             if (*capply) {
00309                 i__1 = *n;
00310                 for (i__ = 1; i__ <= i__1; ++i__) {
00311                     i__2 = i__;
00312                     i__3 = i__;
00313                     i__4 = i__;
00314                     z__1.r = c__[i__4] * work[i__3].r, z__1.i = c__[i__4] * 
00315                             work[i__3].i;
00316                     work[i__2].r = z__1.r, work[i__2].i = z__1.i;
00317                 }
00318             }
00319 
00320             if (notrans) {
00321                 zgbtrs_("Conjugate transpose", n, kl, ku, &c__1, &afb[
00322                         afb_offset], ldafb, &ipiv[1], &work[1], n, info);
00323             } else {
00324                 zgbtrs_("No transpose", n, kl, ku, &c__1, &afb[afb_offset], 
00325                         ldafb, &ipiv[1], &work[1], n, info);
00326             }
00327 
00328 /*           Multiply by R. */
00329 
00330             i__1 = *n;
00331             for (i__ = 1; i__ <= i__1; ++i__) {
00332                 i__2 = i__;
00333                 i__3 = i__;
00334                 i__4 = i__;
00335                 z__1.r = rwork[i__4] * work[i__3].r, z__1.i = rwork[i__4] * 
00336                         work[i__3].i;
00337                 work[i__2].r = z__1.r, work[i__2].i = z__1.i;
00338             }
00339         }
00340         goto L10;
00341     }
00342 
00343 /*     Compute the estimate of the reciprocal condition number. */
00344 
00345     if (ainvnm != 0.) {
00346         ret_val = 1. / ainvnm;
00347     }
00348 
00349     return ret_val;
00350 
00351 } /* zla_gbrcond_c__ */


swiftnav
Author(s):
autogenerated on Sat Jun 8 2019 18:56:39