00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 #include "f2c.h"
00014 #include "blaswrap.h"
00015
00016
00017
00018 static doublereal c_b7 = -1.;
00019 static integer c__1 = 1;
00020 static doublereal c_b19 = 1.;
00021
00022 int dsptrs_(char *uplo, integer *n, integer *nrhs,
00023 doublereal *ap, integer *ipiv, doublereal *b, integer *ldb, integer *
00024 info)
00025 {
00026
00027 integer b_dim1, b_offset, i__1;
00028 doublereal d__1;
00029
00030
00031 integer j, k;
00032 doublereal ak, bk;
00033 integer kc, kp;
00034 doublereal akm1, bkm1;
00035 extern int dger_(integer *, integer *, doublereal *,
00036 doublereal *, integer *, doublereal *, integer *, doublereal *,
00037 integer *);
00038 doublereal akm1k;
00039 extern int dscal_(integer *, doublereal *, doublereal *,
00040 integer *);
00041 extern logical lsame_(char *, char *);
00042 doublereal denom;
00043 extern int dgemv_(char *, integer *, integer *,
00044 doublereal *, doublereal *, integer *, doublereal *, integer *,
00045 doublereal *, doublereal *, integer *), dswap_(integer *,
00046 doublereal *, integer *, doublereal *, integer *);
00047 logical upper;
00048 extern int xerbla_(char *, integer *);
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118 --ap;
00119 --ipiv;
00120 b_dim1 = *ldb;
00121 b_offset = 1 + b_dim1;
00122 b -= b_offset;
00123
00124
00125 *info = 0;
00126 upper = lsame_(uplo, "U");
00127 if (! upper && ! lsame_(uplo, "L")) {
00128 *info = -1;
00129 } else if (*n < 0) {
00130 *info = -2;
00131 } else if (*nrhs < 0) {
00132 *info = -3;
00133 } else if (*ldb < max(1,*n)) {
00134 *info = -7;
00135 }
00136 if (*info != 0) {
00137 i__1 = -(*info);
00138 xerbla_("DSPTRS", &i__1);
00139 return 0;
00140 }
00141
00142
00143
00144 if (*n == 0 || *nrhs == 0) {
00145 return 0;
00146 }
00147
00148 if (upper) {
00149
00150
00151
00152
00153
00154
00155
00156
00157 k = *n;
00158 kc = *n * (*n + 1) / 2 + 1;
00159 L10:
00160
00161
00162
00163 if (k < 1) {
00164 goto L30;
00165 }
00166
00167 kc -= k;
00168 if (ipiv[k] > 0) {
00169
00170
00171
00172
00173
00174 kp = ipiv[k];
00175 if (kp != k) {
00176 dswap_(nrhs, &b[k + b_dim1], ldb, &b[kp + b_dim1], ldb);
00177 }
00178
00179
00180
00181
00182 i__1 = k - 1;
00183 dger_(&i__1, nrhs, &c_b7, &ap[kc], &c__1, &b[k + b_dim1], ldb, &b[
00184 b_dim1 + 1], ldb);
00185
00186
00187
00188 d__1 = 1. / ap[kc + k - 1];
00189 dscal_(nrhs, &d__1, &b[k + b_dim1], ldb);
00190 --k;
00191 } else {
00192
00193
00194
00195
00196
00197 kp = -ipiv[k];
00198 if (kp != k - 1) {
00199 dswap_(nrhs, &b[k - 1 + b_dim1], ldb, &b[kp + b_dim1], ldb);
00200 }
00201
00202
00203
00204
00205 i__1 = k - 2;
00206 dger_(&i__1, nrhs, &c_b7, &ap[kc], &c__1, &b[k + b_dim1], ldb, &b[
00207 b_dim1 + 1], ldb);
00208 i__1 = k - 2;
00209 dger_(&i__1, nrhs, &c_b7, &ap[kc - (k - 1)], &c__1, &b[k - 1 +
00210 b_dim1], ldb, &b[b_dim1 + 1], ldb);
00211
00212
00213
00214 akm1k = ap[kc + k - 2];
00215 akm1 = ap[kc - 1] / akm1k;
00216 ak = ap[kc + k - 1] / akm1k;
00217 denom = akm1 * ak - 1.;
00218 i__1 = *nrhs;
00219 for (j = 1; j <= i__1; ++j) {
00220 bkm1 = b[k - 1 + j * b_dim1] / akm1k;
00221 bk = b[k + j * b_dim1] / akm1k;
00222 b[k - 1 + j * b_dim1] = (ak * bkm1 - bk) / denom;
00223 b[k + j * b_dim1] = (akm1 * bk - bkm1) / denom;
00224
00225 }
00226 kc = kc - k + 1;
00227 k += -2;
00228 }
00229
00230 goto L10;
00231 L30:
00232
00233
00234
00235
00236
00237
00238 k = 1;
00239 kc = 1;
00240 L40:
00241
00242
00243
00244 if (k > *n) {
00245 goto L50;
00246 }
00247
00248 if (ipiv[k] > 0) {
00249
00250
00251
00252
00253
00254
00255 i__1 = k - 1;
00256 dgemv_("Transpose", &i__1, nrhs, &c_b7, &b[b_offset], ldb, &ap[kc]
00257 , &c__1, &c_b19, &b[k + b_dim1], ldb);
00258
00259
00260
00261 kp = ipiv[k];
00262 if (kp != k) {
00263 dswap_(nrhs, &b[k + b_dim1], ldb, &b[kp + b_dim1], ldb);
00264 }
00265 kc += k;
00266 ++k;
00267 } else {
00268
00269
00270
00271
00272
00273
00274 i__1 = k - 1;
00275 dgemv_("Transpose", &i__1, nrhs, &c_b7, &b[b_offset], ldb, &ap[kc]
00276 , &c__1, &c_b19, &b[k + b_dim1], ldb);
00277 i__1 = k - 1;
00278 dgemv_("Transpose", &i__1, nrhs, &c_b7, &b[b_offset], ldb, &ap[kc
00279 + k], &c__1, &c_b19, &b[k + 1 + b_dim1], ldb);
00280
00281
00282
00283 kp = -ipiv[k];
00284 if (kp != k) {
00285 dswap_(nrhs, &b[k + b_dim1], ldb, &b[kp + b_dim1], ldb);
00286 }
00287 kc = kc + (k << 1) + 1;
00288 k += 2;
00289 }
00290
00291 goto L40;
00292 L50:
00293
00294 ;
00295 } else {
00296
00297
00298
00299
00300
00301
00302
00303
00304 k = 1;
00305 kc = 1;
00306 L60:
00307
00308
00309
00310 if (k > *n) {
00311 goto L80;
00312 }
00313
00314 if (ipiv[k] > 0) {
00315
00316
00317
00318
00319
00320 kp = ipiv[k];
00321 if (kp != k) {
00322 dswap_(nrhs, &b[k + b_dim1], ldb, &b[kp + b_dim1], ldb);
00323 }
00324
00325
00326
00327
00328 if (k < *n) {
00329 i__1 = *n - k;
00330 dger_(&i__1, nrhs, &c_b7, &ap[kc + 1], &c__1, &b[k + b_dim1],
00331 ldb, &b[k + 1 + b_dim1], ldb);
00332 }
00333
00334
00335
00336 d__1 = 1. / ap[kc];
00337 dscal_(nrhs, &d__1, &b[k + b_dim1], ldb);
00338 kc = kc + *n - k + 1;
00339 ++k;
00340 } else {
00341
00342
00343
00344
00345
00346 kp = -ipiv[k];
00347 if (kp != k + 1) {
00348 dswap_(nrhs, &b[k + 1 + b_dim1], ldb, &b[kp + b_dim1], ldb);
00349 }
00350
00351
00352
00353
00354 if (k < *n - 1) {
00355 i__1 = *n - k - 1;
00356 dger_(&i__1, nrhs, &c_b7, &ap[kc + 2], &c__1, &b[k + b_dim1],
00357 ldb, &b[k + 2 + b_dim1], ldb);
00358 i__1 = *n - k - 1;
00359 dger_(&i__1, nrhs, &c_b7, &ap[kc + *n - k + 2], &c__1, &b[k +
00360 1 + b_dim1], ldb, &b[k + 2 + b_dim1], ldb);
00361 }
00362
00363
00364
00365 akm1k = ap[kc + 1];
00366 akm1 = ap[kc] / akm1k;
00367 ak = ap[kc + *n - k + 1] / akm1k;
00368 denom = akm1 * ak - 1.;
00369 i__1 = *nrhs;
00370 for (j = 1; j <= i__1; ++j) {
00371 bkm1 = b[k + j * b_dim1] / akm1k;
00372 bk = b[k + 1 + j * b_dim1] / akm1k;
00373 b[k + j * b_dim1] = (ak * bkm1 - bk) / denom;
00374 b[k + 1 + j * b_dim1] = (akm1 * bk - bkm1) / denom;
00375
00376 }
00377 kc = kc + (*n - k << 1) + 1;
00378 k += 2;
00379 }
00380
00381 goto L60;
00382 L80:
00383
00384
00385
00386
00387
00388
00389 k = *n;
00390 kc = *n * (*n + 1) / 2 + 1;
00391 L90:
00392
00393
00394
00395 if (k < 1) {
00396 goto L100;
00397 }
00398
00399 kc -= *n - k + 1;
00400 if (ipiv[k] > 0) {
00401
00402
00403
00404
00405
00406
00407 if (k < *n) {
00408 i__1 = *n - k;
00409 dgemv_("Transpose", &i__1, nrhs, &c_b7, &b[k + 1 + b_dim1],
00410 ldb, &ap[kc + 1], &c__1, &c_b19, &b[k + b_dim1], ldb);
00411 }
00412
00413
00414
00415 kp = ipiv[k];
00416 if (kp != k) {
00417 dswap_(nrhs, &b[k + b_dim1], ldb, &b[kp + b_dim1], ldb);
00418 }
00419 --k;
00420 } else {
00421
00422
00423
00424
00425
00426
00427 if (k < *n) {
00428 i__1 = *n - k;
00429 dgemv_("Transpose", &i__1, nrhs, &c_b7, &b[k + 1 + b_dim1],
00430 ldb, &ap[kc + 1], &c__1, &c_b19, &b[k + b_dim1], ldb);
00431 i__1 = *n - k;
00432 dgemv_("Transpose", &i__1, nrhs, &c_b7, &b[k + 1 + b_dim1],
00433 ldb, &ap[kc - (*n - k)], &c__1, &c_b19, &b[k - 1 +
00434 b_dim1], ldb);
00435 }
00436
00437
00438
00439 kp = -ipiv[k];
00440 if (kp != k) {
00441 dswap_(nrhs, &b[k + b_dim1], ldb, &b[kp + b_dim1], ldb);
00442 }
00443 kc -= *n - k + 2;
00444 k += -2;
00445 }
00446
00447 goto L90;
00448 L100:
00449 ;
00450 }
00451
00452 return 0;
00453
00454
00455
00456 }