You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

spbtrf.c 25 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903
  1. /* f2c.h -- Standard Fortran to C header file */
  2. /** barf [ba:rf] 2. "He suggested using FORTRAN, and everybody barfed."
  3. - From The Shogakukan DICTIONARY OF NEW ENGLISH (Second edition) */
  4. #ifndef F2C_INCLUDE
  5. #define F2C_INCLUDE
  6. #include <math.h>
  7. #include <stdlib.h>
  8. #include <string.h>
  9. #include <stdio.h>
  10. #include <complex.h>
  11. #ifdef complex
  12. #undef complex
  13. #endif
  14. #ifdef I
  15. #undef I
  16. #endif
  17. typedef int integer;
  18. typedef unsigned int uinteger;
  19. typedef char *address;
  20. typedef short int shortint;
  21. typedef float real;
  22. typedef double doublereal;
  23. typedef struct { real r, i; } complex;
  24. typedef struct { doublereal r, i; } doublecomplex;
  25. static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
  26. static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
  27. static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
  28. static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
  29. #define pCf(z) (*_pCf(z))
  30. #define pCd(z) (*_pCd(z))
  31. typedef int logical;
  32. typedef short int shortlogical;
  33. typedef char logical1;
  34. typedef char integer1;
  35. #define TRUE_ (1)
  36. #define FALSE_ (0)
  37. /* Extern is for use with -E */
  38. #ifndef Extern
  39. #define Extern extern
  40. #endif
  41. /* I/O stuff */
  42. typedef int flag;
  43. typedef int ftnlen;
  44. typedef int ftnint;
  45. /*external read, write*/
  46. typedef struct
  47. { flag cierr;
  48. ftnint ciunit;
  49. flag ciend;
  50. char *cifmt;
  51. ftnint cirec;
  52. } cilist;
  53. /*internal read, write*/
  54. typedef struct
  55. { flag icierr;
  56. char *iciunit;
  57. flag iciend;
  58. char *icifmt;
  59. ftnint icirlen;
  60. ftnint icirnum;
  61. } icilist;
  62. /*open*/
  63. typedef struct
  64. { flag oerr;
  65. ftnint ounit;
  66. char *ofnm;
  67. ftnlen ofnmlen;
  68. char *osta;
  69. char *oacc;
  70. char *ofm;
  71. ftnint orl;
  72. char *oblnk;
  73. } olist;
  74. /*close*/
  75. typedef struct
  76. { flag cerr;
  77. ftnint cunit;
  78. char *csta;
  79. } cllist;
  80. /*rewind, backspace, endfile*/
  81. typedef struct
  82. { flag aerr;
  83. ftnint aunit;
  84. } alist;
  85. /* inquire */
  86. typedef struct
  87. { flag inerr;
  88. ftnint inunit;
  89. char *infile;
  90. ftnlen infilen;
  91. ftnint *inex; /*parameters in standard's order*/
  92. ftnint *inopen;
  93. ftnint *innum;
  94. ftnint *innamed;
  95. char *inname;
  96. ftnlen innamlen;
  97. char *inacc;
  98. ftnlen inacclen;
  99. char *inseq;
  100. ftnlen inseqlen;
  101. char *indir;
  102. ftnlen indirlen;
  103. char *infmt;
  104. ftnlen infmtlen;
  105. char *inform;
  106. ftnint informlen;
  107. char *inunf;
  108. ftnlen inunflen;
  109. ftnint *inrecl;
  110. ftnint *innrec;
  111. char *inblank;
  112. ftnlen inblanklen;
  113. } inlist;
  114. #define VOID void
  115. union Multitype { /* for multiple entry points */
  116. integer1 g;
  117. shortint h;
  118. integer i;
  119. /* longint j; */
  120. real r;
  121. doublereal d;
  122. complex c;
  123. doublecomplex z;
  124. };
  125. typedef union Multitype Multitype;
  126. struct Vardesc { /* for Namelist */
  127. char *name;
  128. char *addr;
  129. ftnlen *dims;
  130. int type;
  131. };
  132. typedef struct Vardesc Vardesc;
  133. struct Namelist {
  134. char *name;
  135. Vardesc **vars;
  136. int nvars;
  137. };
  138. typedef struct Namelist Namelist;
  139. #define abs(x) ((x) >= 0 ? (x) : -(x))
  140. #define dabs(x) (fabs(x))
  141. #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
  142. #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
  143. #define dmin(a,b) (f2cmin(a,b))
  144. #define dmax(a,b) (f2cmax(a,b))
  145. #define bit_test(a,b) ((a) >> (b) & 1)
  146. #define bit_clear(a,b) ((a) & ~((uinteger)1 << (b)))
  147. #define bit_set(a,b) ((a) | ((uinteger)1 << (b)))
  148. #define abort_() { sig_die("Fortran abort routine called", 1); }
  149. #define c_abs(z) (cabsf(Cf(z)))
  150. #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
  151. #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
  152. #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
  153. #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
  154. #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
  155. #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
  156. //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
  157. #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
  158. #define d_abs(x) (fabs(*(x)))
  159. #define d_acos(x) (acos(*(x)))
  160. #define d_asin(x) (asin(*(x)))
  161. #define d_atan(x) (atan(*(x)))
  162. #define d_atn2(x, y) (atan2(*(x),*(y)))
  163. #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
  164. #define r_cnjg(R, Z) { pCf(R) = conj(Cf(Z)); }
  165. #define d_cos(x) (cos(*(x)))
  166. #define d_cosh(x) (cosh(*(x)))
  167. #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
  168. #define d_exp(x) (exp(*(x)))
  169. #define d_imag(z) (cimag(Cd(z)))
  170. #define r_imag(z) (cimag(Cf(z)))
  171. #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  172. #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  173. #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  174. #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  175. #define d_log(x) (log(*(x)))
  176. #define d_mod(x, y) (fmod(*(x), *(y)))
  177. #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
  178. #define d_nint(x) u_nint(*(x))
  179. #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
  180. #define d_sign(a,b) u_sign(*(a),*(b))
  181. #define r_sign(a,b) u_sign(*(a),*(b))
  182. #define d_sin(x) (sin(*(x)))
  183. #define d_sinh(x) (sinh(*(x)))
  184. #define d_sqrt(x) (sqrt(*(x)))
  185. #define d_tan(x) (tan(*(x)))
  186. #define d_tanh(x) (tanh(*(x)))
  187. #define i_abs(x) abs(*(x))
  188. #define i_dnnt(x) ((integer)u_nint(*(x)))
  189. #define i_len(s, n) (n)
  190. #define i_nint(x) ((integer)u_nint(*(x)))
  191. #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
  192. #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
  193. #define pow_si(B,E) spow_ui(*(B),*(E))
  194. #define pow_ri(B,E) spow_ui(*(B),*(E))
  195. #define pow_di(B,E) dpow_ui(*(B),*(E))
  196. #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
  197. #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
  198. #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
  199. #define s_cat(lpp, rpp, rnp, np, llp) { ftnlen i, nc, ll; char *f__rp, *lp; ll = (llp); lp = (lpp); for(i=0; i < (int)*(np); ++i) { nc = ll; if((rnp)[i] < nc) nc = (rnp)[i]; ll -= nc; f__rp = (rpp)[i]; while(--nc >= 0) *lp++ = *(f__rp)++; } while(--ll >= 0) *lp++ = ' '; }
  200. #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
  201. #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
  202. #define sig_die(s, kill) { exit(1); }
  203. #define s_stop(s, n) {exit(0);}
  204. static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
  205. #define z_abs(z) (cabs(Cd(z)))
  206. #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
  207. #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
  208. #define myexit_() break;
  209. #define mycycle() continue;
  210. #define myceiling(w) {ceil(w)}
  211. #define myhuge(w) {HUGE_VAL}
  212. //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
  213. #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
  214. /* procedure parameter types for -A and -C++ */
  215. #define F2C_proc_par_types 1
  216. #ifdef __cplusplus
  217. typedef logical (*L_fp)(...);
  218. #else
  219. typedef logical (*L_fp)();
  220. #endif
  221. static float spow_ui(float x, integer n) {
  222. float pow=1.0; unsigned long int u;
  223. if(n != 0) {
  224. if(n < 0) n = -n, x = 1/x;
  225. for(u = n; ; ) {
  226. if(u & 01) pow *= x;
  227. if(u >>= 1) x *= x;
  228. else break;
  229. }
  230. }
  231. return pow;
  232. }
  233. static double dpow_ui(double x, integer n) {
  234. double pow=1.0; unsigned long int u;
  235. if(n != 0) {
  236. if(n < 0) n = -n, x = 1/x;
  237. for(u = n; ; ) {
  238. if(u & 01) pow *= x;
  239. if(u >>= 1) x *= x;
  240. else break;
  241. }
  242. }
  243. return pow;
  244. }
  245. static _Complex float cpow_ui(_Complex float x, integer n) {
  246. _Complex float pow=1.0; unsigned long int u;
  247. if(n != 0) {
  248. if(n < 0) n = -n, x = 1/x;
  249. for(u = n; ; ) {
  250. if(u & 01) pow *= x;
  251. if(u >>= 1) x *= x;
  252. else break;
  253. }
  254. }
  255. return pow;
  256. }
  257. static _Complex double zpow_ui(_Complex double x, integer n) {
  258. _Complex double pow=1.0; unsigned long int u;
  259. if(n != 0) {
  260. if(n < 0) n = -n, x = 1/x;
  261. for(u = n; ; ) {
  262. if(u & 01) pow *= x;
  263. if(u >>= 1) x *= x;
  264. else break;
  265. }
  266. }
  267. return pow;
  268. }
  269. static integer pow_ii(integer x, integer n) {
  270. integer pow; unsigned long int u;
  271. if (n <= 0) {
  272. if (n == 0 || x == 1) pow = 1;
  273. else if (x != -1) pow = x == 0 ? 1/x : 0;
  274. else n = -n;
  275. }
  276. if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
  277. u = n;
  278. for(pow = 1; ; ) {
  279. if(u & 01) pow *= x;
  280. if(u >>= 1) x *= x;
  281. else break;
  282. }
  283. }
  284. return pow;
  285. }
  286. static integer dmaxloc_(double *w, integer s, integer e, integer *n)
  287. {
  288. double m; integer i, mi;
  289. for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
  290. if (w[i-1]>m) mi=i ,m=w[i-1];
  291. return mi-s+1;
  292. }
  293. static integer smaxloc_(float *w, integer s, integer e, integer *n)
  294. {
  295. float m; integer i, mi;
  296. for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
  297. if (w[i-1]>m) mi=i ,m=w[i-1];
  298. return mi-s+1;
  299. }
  300. static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
  301. integer n = *n_, incx = *incx_, incy = *incy_, i;
  302. _Complex float zdotc = 0.0;
  303. if (incx == 1 && incy == 1) {
  304. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  305. zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
  306. }
  307. } else {
  308. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  309. zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
  310. }
  311. }
  312. pCf(z) = zdotc;
  313. }
  314. static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
  315. integer n = *n_, incx = *incx_, incy = *incy_, i;
  316. _Complex double zdotc = 0.0;
  317. if (incx == 1 && incy == 1) {
  318. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  319. zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
  320. }
  321. } else {
  322. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  323. zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
  324. }
  325. }
  326. pCd(z) = zdotc;
  327. }
  328. static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
  329. integer n = *n_, incx = *incx_, incy = *incy_, i;
  330. _Complex float zdotc = 0.0;
  331. if (incx == 1 && incy == 1) {
  332. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  333. zdotc += Cf(&x[i]) * Cf(&y[i]);
  334. }
  335. } else {
  336. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  337. zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
  338. }
  339. }
  340. pCf(z) = zdotc;
  341. }
  342. static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
  343. integer n = *n_, incx = *incx_, incy = *incy_, i;
  344. _Complex double zdotc = 0.0;
  345. if (incx == 1 && incy == 1) {
  346. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  347. zdotc += Cd(&x[i]) * Cd(&y[i]);
  348. }
  349. } else {
  350. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  351. zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
  352. }
  353. }
  354. pCd(z) = zdotc;
  355. }
  356. #endif
  357. /* -- translated by f2c (version 20000121).
  358. You must link the resulting object file with the libraries:
  359. -lf2c -lm (in that order)
  360. */
  361. /* Table of constant values */
  362. static integer c__1 = 1;
  363. static integer c_n1 = -1;
  364. static real c_b18 = 1.f;
  365. static real c_b21 = -1.f;
  366. static integer c__33 = 33;
  367. /* > \brief \b SPBTRF */
  368. /* =========== DOCUMENTATION =========== */
  369. /* Online html documentation available at */
  370. /* http://www.netlib.org/lapack/explore-html/ */
  371. /* > \htmlonly */
  372. /* > Download SPBTRF + dependencies */
  373. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/spbtrf.
  374. f"> */
  375. /* > [TGZ]</a> */
  376. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/spbtrf.
  377. f"> */
  378. /* > [ZIP]</a> */
  379. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/spbtrf.
  380. f"> */
  381. /* > [TXT]</a> */
  382. /* > \endhtmlonly */
  383. /* Definition: */
  384. /* =========== */
  385. /* SUBROUTINE SPBTRF( UPLO, N, KD, AB, LDAB, INFO ) */
  386. /* CHARACTER UPLO */
  387. /* INTEGER INFO, KD, LDAB, N */
  388. /* REAL AB( LDAB, * ) */
  389. /* > \par Purpose: */
  390. /* ============= */
  391. /* > */
  392. /* > \verbatim */
  393. /* > */
  394. /* > SPBTRF computes the Cholesky factorization of a real symmetric */
  395. /* > positive definite band matrix A. */
  396. /* > */
  397. /* > The factorization has the form */
  398. /* > A = U**T * U, if UPLO = 'U', or */
  399. /* > A = L * L**T, if UPLO = 'L', */
  400. /* > where U is an upper triangular matrix and L is lower triangular. */
  401. /* > \endverbatim */
  402. /* Arguments: */
  403. /* ========== */
  404. /* > \param[in] UPLO */
  405. /* > \verbatim */
  406. /* > UPLO is CHARACTER*1 */
  407. /* > = 'U': Upper triangle of A is stored; */
  408. /* > = 'L': Lower triangle of A is stored. */
  409. /* > \endverbatim */
  410. /* > */
  411. /* > \param[in] N */
  412. /* > \verbatim */
  413. /* > N is INTEGER */
  414. /* > The order of the matrix A. N >= 0. */
  415. /* > \endverbatim */
  416. /* > */
  417. /* > \param[in] KD */
  418. /* > \verbatim */
  419. /* > KD is INTEGER */
  420. /* > The number of superdiagonals of the matrix A if UPLO = 'U', */
  421. /* > or the number of subdiagonals if UPLO = 'L'. KD >= 0. */
  422. /* > \endverbatim */
  423. /* > */
  424. /* > \param[in,out] AB */
  425. /* > \verbatim */
  426. /* > AB is REAL array, dimension (LDAB,N) */
  427. /* > On entry, the upper or lower triangle of the symmetric band */
  428. /* > matrix A, stored in the first KD+1 rows of the array. The */
  429. /* > j-th column of A is stored in the j-th column of the array AB */
  430. /* > as follows: */
  431. /* > if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for f2cmax(1,j-kd)<=i<=j; */
  432. /* > if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=f2cmin(n,j+kd). */
  433. /* > */
  434. /* > On exit, if INFO = 0, the triangular factor U or L from the */
  435. /* > Cholesky factorization A = U**T*U or A = L*L**T of the band */
  436. /* > matrix A, in the same storage format as A. */
  437. /* > \endverbatim */
  438. /* > */
  439. /* > \param[in] LDAB */
  440. /* > \verbatim */
  441. /* > LDAB is INTEGER */
  442. /* > The leading dimension of the array AB. LDAB >= KD+1. */
  443. /* > \endverbatim */
  444. /* > */
  445. /* > \param[out] INFO */
  446. /* > \verbatim */
  447. /* > INFO is INTEGER */
  448. /* > = 0: successful exit */
  449. /* > < 0: if INFO = -i, the i-th argument had an illegal value */
  450. /* > > 0: if INFO = i, the leading minor of order i is not */
  451. /* > positive definite, and the factorization could not be */
  452. /* > completed. */
  453. /* > \endverbatim */
  454. /* Authors: */
  455. /* ======== */
  456. /* > \author Univ. of Tennessee */
  457. /* > \author Univ. of California Berkeley */
  458. /* > \author Univ. of Colorado Denver */
  459. /* > \author NAG Ltd. */
  460. /* > \date December 2016 */
  461. /* > \ingroup realOTHERcomputational */
  462. /* > \par Further Details: */
  463. /* ===================== */
  464. /* > */
  465. /* > \verbatim */
  466. /* > */
  467. /* > The band storage scheme is illustrated by the following example, when */
  468. /* > N = 6, KD = 2, and UPLO = 'U': */
  469. /* > */
  470. /* > On entry: On exit: */
  471. /* > */
  472. /* > * * a13 a24 a35 a46 * * u13 u24 u35 u46 */
  473. /* > * a12 a23 a34 a45 a56 * u12 u23 u34 u45 u56 */
  474. /* > a11 a22 a33 a44 a55 a66 u11 u22 u33 u44 u55 u66 */
  475. /* > */
  476. /* > Similarly, if UPLO = 'L' the format of A is as follows: */
  477. /* > */
  478. /* > On entry: On exit: */
  479. /* > */
  480. /* > a11 a22 a33 a44 a55 a66 l11 l22 l33 l44 l55 l66 */
  481. /* > a21 a32 a43 a54 a65 * l21 l32 l43 l54 l65 * */
  482. /* > a31 a42 a53 a64 * * l31 l42 l53 l64 * * */
  483. /* > */
  484. /* > Array elements marked * are not used by the routine. */
  485. /* > \endverbatim */
  486. /* > \par Contributors: */
  487. /* ================== */
  488. /* > */
  489. /* > Peter Mayes and Giuseppe Radicati, IBM ECSEC, Rome, March 23, 1989 */
  490. /* ===================================================================== */
  491. /* Subroutine */ int spbtrf_(char *uplo, integer *n, integer *kd, real *ab,
  492. integer *ldab, integer *info)
  493. {
  494. /* System generated locals */
  495. integer ab_dim1, ab_offset, i__1, i__2, i__3, i__4;
  496. /* Local variables */
  497. real work[1056] /* was [33][32] */;
  498. integer i__, j;
  499. extern logical lsame_(char *, char *);
  500. extern /* Subroutine */ int sgemm_(char *, char *, integer *, integer *,
  501. integer *, real *, real *, integer *, real *, integer *, real *,
  502. real *, integer *);
  503. integer i2, i3;
  504. extern /* Subroutine */ int strsm_(char *, char *, char *, char *,
  505. integer *, integer *, real *, real *, integer *, real *, integer *
  506. ), ssyrk_(char *, char *, integer
  507. *, integer *, real *, real *, integer *, real *, real *, integer *
  508. ), spbtf2_(char *, integer *, integer *, real *,
  509. integer *, integer *);
  510. integer ib;
  511. extern /* Subroutine */ int spotf2_(char *, integer *, real *, integer *,
  512. integer *);
  513. integer nb, ii, jj;
  514. extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
  515. extern integer ilaenv_(integer *, char *, char *, integer *, integer *,
  516. integer *, integer *, ftnlen, ftnlen);
  517. /* -- LAPACK computational routine (version 3.7.0) -- */
  518. /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
  519. /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
  520. /* December 2016 */
  521. /* ===================================================================== */
  522. /* Test the input parameters. */
  523. /* Parameter adjustments */
  524. ab_dim1 = *ldab;
  525. ab_offset = 1 + ab_dim1 * 1;
  526. ab -= ab_offset;
  527. /* Function Body */
  528. *info = 0;
  529. if (! lsame_(uplo, "U") && ! lsame_(uplo, "L")) {
  530. *info = -1;
  531. } else if (*n < 0) {
  532. *info = -2;
  533. } else if (*kd < 0) {
  534. *info = -3;
  535. } else if (*ldab < *kd + 1) {
  536. *info = -5;
  537. }
  538. if (*info != 0) {
  539. i__1 = -(*info);
  540. xerbla_("SPBTRF", &i__1, (ftnlen)6);
  541. return 0;
  542. }
  543. /* Quick return if possible */
  544. if (*n == 0) {
  545. return 0;
  546. }
  547. /* Determine the block size for this environment */
  548. nb = ilaenv_(&c__1, "SPBTRF", uplo, n, kd, &c_n1, &c_n1, (ftnlen)6, (
  549. ftnlen)1);
  550. /* The block size must not exceed the semi-bandwidth KD, and must not */
  551. /* exceed the limit set by the size of the local array WORK. */
  552. nb = f2cmin(nb,32);
  553. if (nb <= 1 || nb > *kd) {
  554. /* Use unblocked code */
  555. spbtf2_(uplo, n, kd, &ab[ab_offset], ldab, info);
  556. } else {
  557. /* Use blocked code */
  558. if (lsame_(uplo, "U")) {
  559. /* Compute the Cholesky factorization of a symmetric band */
  560. /* matrix, given the upper triangle of the matrix in band */
  561. /* storage. */
  562. /* Zero the upper triangle of the work array. */
  563. i__1 = nb;
  564. for (j = 1; j <= i__1; ++j) {
  565. i__2 = j - 1;
  566. for (i__ = 1; i__ <= i__2; ++i__) {
  567. work[i__ + j * 33 - 34] = 0.f;
  568. /* L10: */
  569. }
  570. /* L20: */
  571. }
  572. /* Process the band matrix one diagonal block at a time. */
  573. i__1 = *n;
  574. i__2 = nb;
  575. for (i__ = 1; i__2 < 0 ? i__ >= i__1 : i__ <= i__1; i__ += i__2) {
  576. /* Computing MIN */
  577. i__3 = nb, i__4 = *n - i__ + 1;
  578. ib = f2cmin(i__3,i__4);
  579. /* Factorize the diagonal block */
  580. i__3 = *ldab - 1;
  581. spotf2_(uplo, &ib, &ab[*kd + 1 + i__ * ab_dim1], &i__3, &ii);
  582. if (ii != 0) {
  583. *info = i__ + ii - 1;
  584. goto L150;
  585. }
  586. if (i__ + ib <= *n) {
  587. /* Update the relevant part of the trailing submatrix. */
  588. /* If A11 denotes the diagonal block which has just been */
  589. /* factorized, then we need to update the remaining */
  590. /* blocks in the diagram: */
  591. /* A11 A12 A13 */
  592. /* A22 A23 */
  593. /* A33 */
  594. /* The numbers of rows and columns in the partitioning */
  595. /* are IB, I2, I3 respectively. The blocks A12, A22 and */
  596. /* A23 are empty if IB = KD. The upper triangle of A13 */
  597. /* lies outside the band. */
  598. /* Computing MIN */
  599. i__3 = *kd - ib, i__4 = *n - i__ - ib + 1;
  600. i2 = f2cmin(i__3,i__4);
  601. /* Computing MIN */
  602. i__3 = ib, i__4 = *n - i__ - *kd + 1;
  603. i3 = f2cmin(i__3,i__4);
  604. if (i2 > 0) {
  605. /* Update A12 */
  606. i__3 = *ldab - 1;
  607. i__4 = *ldab - 1;
  608. strsm_("Left", "Upper", "Transpose", "Non-unit", &ib,
  609. &i2, &c_b18, &ab[*kd + 1 + i__ * ab_dim1], &
  610. i__3, &ab[*kd + 1 - ib + (i__ + ib) * ab_dim1]
  611. , &i__4);
  612. /* Update A22 */
  613. i__3 = *ldab - 1;
  614. i__4 = *ldab - 1;
  615. ssyrk_("Upper", "Transpose", &i2, &ib, &c_b21, &ab[*
  616. kd + 1 - ib + (i__ + ib) * ab_dim1], &i__3, &
  617. c_b18, &ab[*kd + 1 + (i__ + ib) * ab_dim1], &
  618. i__4);
  619. }
  620. if (i3 > 0) {
  621. /* Copy the lower triangle of A13 into the work array. */
  622. i__3 = i3;
  623. for (jj = 1; jj <= i__3; ++jj) {
  624. i__4 = ib;
  625. for (ii = jj; ii <= i__4; ++ii) {
  626. work[ii + jj * 33 - 34] = ab[ii - jj + 1 + (
  627. jj + i__ + *kd - 1) * ab_dim1];
  628. /* L30: */
  629. }
  630. /* L40: */
  631. }
  632. /* Update A13 (in the work array). */
  633. i__3 = *ldab - 1;
  634. strsm_("Left", "Upper", "Transpose", "Non-unit", &ib,
  635. &i3, &c_b18, &ab[*kd + 1 + i__ * ab_dim1], &
  636. i__3, work, &c__33);
  637. /* Update A23 */
  638. if (i2 > 0) {
  639. i__3 = *ldab - 1;
  640. i__4 = *ldab - 1;
  641. sgemm_("Transpose", "No Transpose", &i2, &i3, &ib,
  642. &c_b21, &ab[*kd + 1 - ib + (i__ + ib) *
  643. ab_dim1], &i__3, work, &c__33, &c_b18, &
  644. ab[ib + 1 + (i__ + *kd) * ab_dim1], &i__4);
  645. }
  646. /* Update A33 */
  647. i__3 = *ldab - 1;
  648. ssyrk_("Upper", "Transpose", &i3, &ib, &c_b21, work, &
  649. c__33, &c_b18, &ab[*kd + 1 + (i__ + *kd) *
  650. ab_dim1], &i__3);
  651. /* Copy the lower triangle of A13 back into place. */
  652. i__3 = i3;
  653. for (jj = 1; jj <= i__3; ++jj) {
  654. i__4 = ib;
  655. for (ii = jj; ii <= i__4; ++ii) {
  656. ab[ii - jj + 1 + (jj + i__ + *kd - 1) *
  657. ab_dim1] = work[ii + jj * 33 - 34];
  658. /* L50: */
  659. }
  660. /* L60: */
  661. }
  662. }
  663. }
  664. /* L70: */
  665. }
  666. } else {
  667. /* Compute the Cholesky factorization of a symmetric band */
  668. /* matrix, given the lower triangle of the matrix in band */
  669. /* storage. */
  670. /* Zero the lower triangle of the work array. */
  671. i__2 = nb;
  672. for (j = 1; j <= i__2; ++j) {
  673. i__1 = nb;
  674. for (i__ = j + 1; i__ <= i__1; ++i__) {
  675. work[i__ + j * 33 - 34] = 0.f;
  676. /* L80: */
  677. }
  678. /* L90: */
  679. }
  680. /* Process the band matrix one diagonal block at a time. */
  681. i__2 = *n;
  682. i__1 = nb;
  683. for (i__ = 1; i__1 < 0 ? i__ >= i__2 : i__ <= i__2; i__ += i__1) {
  684. /* Computing MIN */
  685. i__3 = nb, i__4 = *n - i__ + 1;
  686. ib = f2cmin(i__3,i__4);
  687. /* Factorize the diagonal block */
  688. i__3 = *ldab - 1;
  689. spotf2_(uplo, &ib, &ab[i__ * ab_dim1 + 1], &i__3, &ii);
  690. if (ii != 0) {
  691. *info = i__ + ii - 1;
  692. goto L150;
  693. }
  694. if (i__ + ib <= *n) {
  695. /* Update the relevant part of the trailing submatrix. */
  696. /* If A11 denotes the diagonal block which has just been */
  697. /* factorized, then we need to update the remaining */
  698. /* blocks in the diagram: */
  699. /* A11 */
  700. /* A21 A22 */
  701. /* A31 A32 A33 */
  702. /* The numbers of rows and columns in the partitioning */
  703. /* are IB, I2, I3 respectively. The blocks A21, A22 and */
  704. /* A32 are empty if IB = KD. The lower triangle of A31 */
  705. /* lies outside the band. */
  706. /* Computing MIN */
  707. i__3 = *kd - ib, i__4 = *n - i__ - ib + 1;
  708. i2 = f2cmin(i__3,i__4);
  709. /* Computing MIN */
  710. i__3 = ib, i__4 = *n - i__ - *kd + 1;
  711. i3 = f2cmin(i__3,i__4);
  712. if (i2 > 0) {
  713. /* Update A21 */
  714. i__3 = *ldab - 1;
  715. i__4 = *ldab - 1;
  716. strsm_("Right", "Lower", "Transpose", "Non-unit", &i2,
  717. &ib, &c_b18, &ab[i__ * ab_dim1 + 1], &i__3, &
  718. ab[ib + 1 + i__ * ab_dim1], &i__4);
  719. /* Update A22 */
  720. i__3 = *ldab - 1;
  721. i__4 = *ldab - 1;
  722. ssyrk_("Lower", "No Transpose", &i2, &ib, &c_b21, &ab[
  723. ib + 1 + i__ * ab_dim1], &i__3, &c_b18, &ab[(
  724. i__ + ib) * ab_dim1 + 1], &i__4);
  725. }
  726. if (i3 > 0) {
  727. /* Copy the upper triangle of A31 into the work array. */
  728. i__3 = ib;
  729. for (jj = 1; jj <= i__3; ++jj) {
  730. i__4 = f2cmin(jj,i3);
  731. for (ii = 1; ii <= i__4; ++ii) {
  732. work[ii + jj * 33 - 34] = ab[*kd + 1 - jj +
  733. ii + (jj + i__ - 1) * ab_dim1];
  734. /* L100: */
  735. }
  736. /* L110: */
  737. }
  738. /* Update A31 (in the work array). */
  739. i__3 = *ldab - 1;
  740. strsm_("Right", "Lower", "Transpose", "Non-unit", &i3,
  741. &ib, &c_b18, &ab[i__ * ab_dim1 + 1], &i__3,
  742. work, &c__33);
  743. /* Update A32 */
  744. if (i2 > 0) {
  745. i__3 = *ldab - 1;
  746. i__4 = *ldab - 1;
  747. sgemm_("No transpose", "Transpose", &i3, &i2, &ib,
  748. &c_b21, work, &c__33, &ab[ib + 1 + i__ *
  749. ab_dim1], &i__3, &c_b18, &ab[*kd + 1 - ib
  750. + (i__ + ib) * ab_dim1], &i__4);
  751. }
  752. /* Update A33 */
  753. i__3 = *ldab - 1;
  754. ssyrk_("Lower", "No Transpose", &i3, &ib, &c_b21,
  755. work, &c__33, &c_b18, &ab[(i__ + *kd) *
  756. ab_dim1 + 1], &i__3);
  757. /* Copy the upper triangle of A31 back into place. */
  758. i__3 = ib;
  759. for (jj = 1; jj <= i__3; ++jj) {
  760. i__4 = f2cmin(jj,i3);
  761. for (ii = 1; ii <= i__4; ++ii) {
  762. ab[*kd + 1 - jj + ii + (jj + i__ - 1) *
  763. ab_dim1] = work[ii + jj * 33 - 34];
  764. /* L120: */
  765. }
  766. /* L130: */
  767. }
  768. }
  769. }
  770. /* L140: */
  771. }
  772. }
  773. }
  774. return 0;
  775. L150:
  776. return 0;
  777. /* End of SPBTRF */
  778. } /* spbtrf_ */