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.

dsysvx.c 27 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866
  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. #if defined(_WIN64)
  18. typedef long long BLASLONG;
  19. typedef unsigned long long BLASULONG;
  20. #else
  21. typedef long BLASLONG;
  22. typedef unsigned long BLASULONG;
  23. #endif
  24. #ifdef LAPACK_ILP64
  25. typedef BLASLONG blasint;
  26. #if defined(_WIN64)
  27. #define blasabs(x) llabs(x)
  28. #else
  29. #define blasabs(x) labs(x)
  30. #endif
  31. #else
  32. typedef int blasint;
  33. #define blasabs(x) abs(x)
  34. #endif
  35. typedef blasint integer;
  36. typedef unsigned int uinteger;
  37. typedef char *address;
  38. typedef short int shortint;
  39. typedef float real;
  40. typedef double doublereal;
  41. typedef struct { real r, i; } complex;
  42. typedef struct { doublereal r, i; } doublecomplex;
  43. static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
  44. static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
  45. static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
  46. static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
  47. #define pCf(z) (*_pCf(z))
  48. #define pCd(z) (*_pCd(z))
  49. typedef int logical;
  50. typedef short int shortlogical;
  51. typedef char logical1;
  52. typedef char integer1;
  53. #define TRUE_ (1)
  54. #define FALSE_ (0)
  55. /* Extern is for use with -E */
  56. #ifndef Extern
  57. #define Extern extern
  58. #endif
  59. /* I/O stuff */
  60. typedef int flag;
  61. typedef int ftnlen;
  62. typedef int ftnint;
  63. /*external read, write*/
  64. typedef struct
  65. { flag cierr;
  66. ftnint ciunit;
  67. flag ciend;
  68. char *cifmt;
  69. ftnint cirec;
  70. } cilist;
  71. /*internal read, write*/
  72. typedef struct
  73. { flag icierr;
  74. char *iciunit;
  75. flag iciend;
  76. char *icifmt;
  77. ftnint icirlen;
  78. ftnint icirnum;
  79. } icilist;
  80. /*open*/
  81. typedef struct
  82. { flag oerr;
  83. ftnint ounit;
  84. char *ofnm;
  85. ftnlen ofnmlen;
  86. char *osta;
  87. char *oacc;
  88. char *ofm;
  89. ftnint orl;
  90. char *oblnk;
  91. } olist;
  92. /*close*/
  93. typedef struct
  94. { flag cerr;
  95. ftnint cunit;
  96. char *csta;
  97. } cllist;
  98. /*rewind, backspace, endfile*/
  99. typedef struct
  100. { flag aerr;
  101. ftnint aunit;
  102. } alist;
  103. /* inquire */
  104. typedef struct
  105. { flag inerr;
  106. ftnint inunit;
  107. char *infile;
  108. ftnlen infilen;
  109. ftnint *inex; /*parameters in standard's order*/
  110. ftnint *inopen;
  111. ftnint *innum;
  112. ftnint *innamed;
  113. char *inname;
  114. ftnlen innamlen;
  115. char *inacc;
  116. ftnlen inacclen;
  117. char *inseq;
  118. ftnlen inseqlen;
  119. char *indir;
  120. ftnlen indirlen;
  121. char *infmt;
  122. ftnlen infmtlen;
  123. char *inform;
  124. ftnint informlen;
  125. char *inunf;
  126. ftnlen inunflen;
  127. ftnint *inrecl;
  128. ftnint *innrec;
  129. char *inblank;
  130. ftnlen inblanklen;
  131. } inlist;
  132. #define VOID void
  133. union Multitype { /* for multiple entry points */
  134. integer1 g;
  135. shortint h;
  136. integer i;
  137. /* longint j; */
  138. real r;
  139. doublereal d;
  140. complex c;
  141. doublecomplex z;
  142. };
  143. typedef union Multitype Multitype;
  144. struct Vardesc { /* for Namelist */
  145. char *name;
  146. char *addr;
  147. ftnlen *dims;
  148. int type;
  149. };
  150. typedef struct Vardesc Vardesc;
  151. struct Namelist {
  152. char *name;
  153. Vardesc **vars;
  154. int nvars;
  155. };
  156. typedef struct Namelist Namelist;
  157. #define abs(x) ((x) >= 0 ? (x) : -(x))
  158. #define dabs(x) (fabs(x))
  159. #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
  160. #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
  161. #define dmin(a,b) (f2cmin(a,b))
  162. #define dmax(a,b) (f2cmax(a,b))
  163. #define bit_test(a,b) ((a) >> (b) & 1)
  164. #define bit_clear(a,b) ((a) & ~((uinteger)1 << (b)))
  165. #define bit_set(a,b) ((a) | ((uinteger)1 << (b)))
  166. #define abort_() { sig_die("Fortran abort routine called", 1); }
  167. #define c_abs(z) (cabsf(Cf(z)))
  168. #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
  169. #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
  170. #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
  171. #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
  172. #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
  173. #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
  174. //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
  175. #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
  176. #define d_abs(x) (fabs(*(x)))
  177. #define d_acos(x) (acos(*(x)))
  178. #define d_asin(x) (asin(*(x)))
  179. #define d_atan(x) (atan(*(x)))
  180. #define d_atn2(x, y) (atan2(*(x),*(y)))
  181. #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
  182. #define r_cnjg(R, Z) { pCf(R) = conj(Cf(Z)); }
  183. #define d_cos(x) (cos(*(x)))
  184. #define d_cosh(x) (cosh(*(x)))
  185. #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
  186. #define d_exp(x) (exp(*(x)))
  187. #define d_imag(z) (cimag(Cd(z)))
  188. #define r_imag(z) (cimag(Cf(z)))
  189. #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  190. #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  191. #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  192. #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  193. #define d_log(x) (log(*(x)))
  194. #define d_mod(x, y) (fmod(*(x), *(y)))
  195. #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
  196. #define d_nint(x) u_nint(*(x))
  197. #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
  198. #define d_sign(a,b) u_sign(*(a),*(b))
  199. #define r_sign(a,b) u_sign(*(a),*(b))
  200. #define d_sin(x) (sin(*(x)))
  201. #define d_sinh(x) (sinh(*(x)))
  202. #define d_sqrt(x) (sqrt(*(x)))
  203. #define d_tan(x) (tan(*(x)))
  204. #define d_tanh(x) (tanh(*(x)))
  205. #define i_abs(x) abs(*(x))
  206. #define i_dnnt(x) ((integer)u_nint(*(x)))
  207. #define i_len(s, n) (n)
  208. #define i_nint(x) ((integer)u_nint(*(x)))
  209. #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
  210. #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
  211. #define pow_si(B,E) spow_ui(*(B),*(E))
  212. #define pow_ri(B,E) spow_ui(*(B),*(E))
  213. #define pow_di(B,E) dpow_ui(*(B),*(E))
  214. #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
  215. #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
  216. #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
  217. #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++ = ' '; }
  218. #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
  219. #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]; }
  220. #define sig_die(s, kill) { exit(1); }
  221. #define s_stop(s, n) {exit(0);}
  222. static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
  223. #define z_abs(z) (cabs(Cd(z)))
  224. #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
  225. #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
  226. #define myexit_() break;
  227. #define mycycle() continue;
  228. #define myceiling(w) {ceil(w)}
  229. #define myhuge(w) {HUGE_VAL}
  230. //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
  231. #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
  232. /* procedure parameter types for -A and -C++ */
  233. #define F2C_proc_par_types 1
  234. #ifdef __cplusplus
  235. typedef logical (*L_fp)(...);
  236. #else
  237. typedef logical (*L_fp)();
  238. #endif
  239. static float spow_ui(float x, integer n) {
  240. float pow=1.0; unsigned long int u;
  241. if(n != 0) {
  242. if(n < 0) n = -n, x = 1/x;
  243. for(u = n; ; ) {
  244. if(u & 01) pow *= x;
  245. if(u >>= 1) x *= x;
  246. else break;
  247. }
  248. }
  249. return pow;
  250. }
  251. static double dpow_ui(double x, integer n) {
  252. double pow=1.0; unsigned long int u;
  253. if(n != 0) {
  254. if(n < 0) n = -n, x = 1/x;
  255. for(u = n; ; ) {
  256. if(u & 01) pow *= x;
  257. if(u >>= 1) x *= x;
  258. else break;
  259. }
  260. }
  261. return pow;
  262. }
  263. static _Complex float cpow_ui(_Complex float x, integer n) {
  264. _Complex float pow=1.0; unsigned long int u;
  265. if(n != 0) {
  266. if(n < 0) n = -n, x = 1/x;
  267. for(u = n; ; ) {
  268. if(u & 01) pow *= x;
  269. if(u >>= 1) x *= x;
  270. else break;
  271. }
  272. }
  273. return pow;
  274. }
  275. static _Complex double zpow_ui(_Complex double x, integer n) {
  276. _Complex double pow=1.0; unsigned long int u;
  277. if(n != 0) {
  278. if(n < 0) n = -n, x = 1/x;
  279. for(u = n; ; ) {
  280. if(u & 01) pow *= x;
  281. if(u >>= 1) x *= x;
  282. else break;
  283. }
  284. }
  285. return pow;
  286. }
  287. static integer pow_ii(integer x, integer n) {
  288. integer pow; unsigned long int u;
  289. if (n <= 0) {
  290. if (n == 0 || x == 1) pow = 1;
  291. else if (x != -1) pow = x == 0 ? 1/x : 0;
  292. else n = -n;
  293. }
  294. if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
  295. u = n;
  296. for(pow = 1; ; ) {
  297. if(u & 01) pow *= x;
  298. if(u >>= 1) x *= x;
  299. else break;
  300. }
  301. }
  302. return pow;
  303. }
  304. static integer dmaxloc_(double *w, integer s, integer e, integer *n)
  305. {
  306. double m; integer i, mi;
  307. for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
  308. if (w[i-1]>m) mi=i ,m=w[i-1];
  309. return mi-s+1;
  310. }
  311. static integer smaxloc_(float *w, integer s, integer e, integer *n)
  312. {
  313. float m; integer i, mi;
  314. for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
  315. if (w[i-1]>m) mi=i ,m=w[i-1];
  316. return mi-s+1;
  317. }
  318. static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
  319. integer n = *n_, incx = *incx_, incy = *incy_, i;
  320. _Complex float zdotc = 0.0;
  321. if (incx == 1 && incy == 1) {
  322. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  323. zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
  324. }
  325. } else {
  326. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  327. zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
  328. }
  329. }
  330. pCf(z) = zdotc;
  331. }
  332. static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
  333. integer n = *n_, incx = *incx_, incy = *incy_, i;
  334. _Complex double zdotc = 0.0;
  335. if (incx == 1 && incy == 1) {
  336. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  337. zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
  338. }
  339. } else {
  340. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  341. zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
  342. }
  343. }
  344. pCd(z) = zdotc;
  345. }
  346. static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
  347. integer n = *n_, incx = *incx_, incy = *incy_, i;
  348. _Complex float zdotc = 0.0;
  349. if (incx == 1 && incy == 1) {
  350. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  351. zdotc += Cf(&x[i]) * Cf(&y[i]);
  352. }
  353. } else {
  354. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  355. zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
  356. }
  357. }
  358. pCf(z) = zdotc;
  359. }
  360. static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
  361. integer n = *n_, incx = *incx_, incy = *incy_, i;
  362. _Complex double zdotc = 0.0;
  363. if (incx == 1 && incy == 1) {
  364. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  365. zdotc += Cd(&x[i]) * Cd(&y[i]);
  366. }
  367. } else {
  368. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  369. zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
  370. }
  371. }
  372. pCd(z) = zdotc;
  373. }
  374. #endif
  375. /* -- translated by f2c (version 20000121).
  376. You must link the resulting object file with the libraries:
  377. -lf2c -lm (in that order)
  378. */
  379. /* Table of constant values */
  380. static integer c__1 = 1;
  381. static integer c_n1 = -1;
  382. /* > \brief <b> DSYSVX computes the solution to system of linear equations A * X = B for SY matrices</b> */
  383. /* =========== DOCUMENTATION =========== */
  384. /* Online html documentation available at */
  385. /* http://www.netlib.org/lapack/explore-html/ */
  386. /* > \htmlonly */
  387. /* > Download DSYSVX + dependencies */
  388. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsysvx.
  389. f"> */
  390. /* > [TGZ]</a> */
  391. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsysvx.
  392. f"> */
  393. /* > [ZIP]</a> */
  394. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsysvx.
  395. f"> */
  396. /* > [TXT]</a> */
  397. /* > \endhtmlonly */
  398. /* Definition: */
  399. /* =========== */
  400. /* SUBROUTINE DSYSVX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV, B, */
  401. /* LDB, X, LDX, RCOND, FERR, BERR, WORK, LWORK, */
  402. /* IWORK, INFO ) */
  403. /* CHARACTER FACT, UPLO */
  404. /* INTEGER INFO, LDA, LDAF, LDB, LDX, LWORK, N, NRHS */
  405. /* DOUBLE PRECISION RCOND */
  406. /* INTEGER IPIV( * ), IWORK( * ) */
  407. /* DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ), */
  408. /* $ BERR( * ), FERR( * ), WORK( * ), X( LDX, * ) */
  409. /* > \par Purpose: */
  410. /* ============= */
  411. /* > */
  412. /* > \verbatim */
  413. /* > */
  414. /* > DSYSVX uses the diagonal pivoting factorization to compute the */
  415. /* > solution to a real system of linear equations A * X = B, */
  416. /* > where A is an N-by-N symmetric matrix and X and B are N-by-NRHS */
  417. /* > matrices. */
  418. /* > */
  419. /* > Error bounds on the solution and a condition estimate are also */
  420. /* > provided. */
  421. /* > \endverbatim */
  422. /* > \par Description: */
  423. /* ================= */
  424. /* > */
  425. /* > \verbatim */
  426. /* > */
  427. /* > The following steps are performed: */
  428. /* > */
  429. /* > 1. If FACT = 'N', the diagonal pivoting method is used to factor A. */
  430. /* > The form of the factorization is */
  431. /* > A = U * D * U**T, if UPLO = 'U', or */
  432. /* > A = L * D * L**T, if UPLO = 'L', */
  433. /* > where U (or L) is a product of permutation and unit upper (lower) */
  434. /* > triangular matrices, and D is symmetric and block diagonal with */
  435. /* > 1-by-1 and 2-by-2 diagonal blocks. */
  436. /* > */
  437. /* > 2. If some D(i,i)=0, so that D is exactly singular, then the routine */
  438. /* > returns with INFO = i. Otherwise, the factored form of A is used */
  439. /* > to estimate the condition number of the matrix A. If the */
  440. /* > reciprocal of the condition number is less than machine precision, */
  441. /* > INFO = N+1 is returned as a warning, but the routine still goes on */
  442. /* > to solve for X and compute error bounds as described below. */
  443. /* > */
  444. /* > 3. The system of equations is solved for X using the factored form */
  445. /* > of A. */
  446. /* > */
  447. /* > 4. Iterative refinement is applied to improve the computed solution */
  448. /* > matrix and calculate error bounds and backward error estimates */
  449. /* > for it. */
  450. /* > \endverbatim */
  451. /* Arguments: */
  452. /* ========== */
  453. /* > \param[in] FACT */
  454. /* > \verbatim */
  455. /* > FACT is CHARACTER*1 */
  456. /* > Specifies whether or not the factored form of A has been */
  457. /* > supplied on entry. */
  458. /* > = 'F': On entry, AF and IPIV contain the factored form of */
  459. /* > A. AF and IPIV will not be modified. */
  460. /* > = 'N': The matrix A will be copied to AF and factored. */
  461. /* > \endverbatim */
  462. /* > */
  463. /* > \param[in] UPLO */
  464. /* > \verbatim */
  465. /* > UPLO is CHARACTER*1 */
  466. /* > = 'U': Upper triangle of A is stored; */
  467. /* > = 'L': Lower triangle of A is stored. */
  468. /* > \endverbatim */
  469. /* > */
  470. /* > \param[in] N */
  471. /* > \verbatim */
  472. /* > N is INTEGER */
  473. /* > The number of linear equations, i.e., the order of the */
  474. /* > matrix A. N >= 0. */
  475. /* > \endverbatim */
  476. /* > */
  477. /* > \param[in] NRHS */
  478. /* > \verbatim */
  479. /* > NRHS is INTEGER */
  480. /* > The number of right hand sides, i.e., the number of columns */
  481. /* > of the matrices B and X. NRHS >= 0. */
  482. /* > \endverbatim */
  483. /* > */
  484. /* > \param[in] A */
  485. /* > \verbatim */
  486. /* > A is DOUBLE PRECISION array, dimension (LDA,N) */
  487. /* > The symmetric matrix A. If UPLO = 'U', the leading N-by-N */
  488. /* > upper triangular part of A contains the upper triangular part */
  489. /* > of the matrix A, and the strictly lower triangular part of A */
  490. /* > is not referenced. If UPLO = 'L', the leading N-by-N lower */
  491. /* > triangular part of A contains the lower triangular part of */
  492. /* > the matrix A, and the strictly upper triangular part of A is */
  493. /* > not referenced. */
  494. /* > \endverbatim */
  495. /* > */
  496. /* > \param[in] LDA */
  497. /* > \verbatim */
  498. /* > LDA is INTEGER */
  499. /* > The leading dimension of the array A. LDA >= f2cmax(1,N). */
  500. /* > \endverbatim */
  501. /* > */
  502. /* > \param[in,out] AF */
  503. /* > \verbatim */
  504. /* > AF is DOUBLE PRECISION array, dimension (LDAF,N) */
  505. /* > If FACT = 'F', then AF is an input argument and on entry */
  506. /* > contains the block diagonal matrix D and the multipliers used */
  507. /* > to obtain the factor U or L from the factorization */
  508. /* > A = U*D*U**T or A = L*D*L**T as computed by DSYTRF. */
  509. /* > */
  510. /* > If FACT = 'N', then AF is an output argument and on exit */
  511. /* > returns the block diagonal matrix D and the multipliers used */
  512. /* > to obtain the factor U or L from the factorization */
  513. /* > A = U*D*U**T or A = L*D*L**T. */
  514. /* > \endverbatim */
  515. /* > */
  516. /* > \param[in] LDAF */
  517. /* > \verbatim */
  518. /* > LDAF is INTEGER */
  519. /* > The leading dimension of the array AF. LDAF >= f2cmax(1,N). */
  520. /* > \endverbatim */
  521. /* > */
  522. /* > \param[in,out] IPIV */
  523. /* > \verbatim */
  524. /* > IPIV is INTEGER array, dimension (N) */
  525. /* > If FACT = 'F', then IPIV is an input argument and on entry */
  526. /* > contains details of the interchanges and the block structure */
  527. /* > of D, as determined by DSYTRF. */
  528. /* > If IPIV(k) > 0, then rows and columns k and IPIV(k) were */
  529. /* > interchanged and D(k,k) is a 1-by-1 diagonal block. */
  530. /* > If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and */
  531. /* > columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) */
  532. /* > is a 2-by-2 diagonal block. If UPLO = 'L' and IPIV(k) = */
  533. /* > IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were */
  534. /* > interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. */
  535. /* > */
  536. /* > If FACT = 'N', then IPIV is an output argument and on exit */
  537. /* > contains details of the interchanges and the block structure */
  538. /* > of D, as determined by DSYTRF. */
  539. /* > \endverbatim */
  540. /* > */
  541. /* > \param[in] B */
  542. /* > \verbatim */
  543. /* > B is DOUBLE PRECISION array, dimension (LDB,NRHS) */
  544. /* > The N-by-NRHS right hand side matrix B. */
  545. /* > \endverbatim */
  546. /* > */
  547. /* > \param[in] LDB */
  548. /* > \verbatim */
  549. /* > LDB is INTEGER */
  550. /* > The leading dimension of the array B. LDB >= f2cmax(1,N). */
  551. /* > \endverbatim */
  552. /* > */
  553. /* > \param[out] X */
  554. /* > \verbatim */
  555. /* > X is DOUBLE PRECISION array, dimension (LDX,NRHS) */
  556. /* > If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X. */
  557. /* > \endverbatim */
  558. /* > */
  559. /* > \param[in] LDX */
  560. /* > \verbatim */
  561. /* > LDX is INTEGER */
  562. /* > The leading dimension of the array X. LDX >= f2cmax(1,N). */
  563. /* > \endverbatim */
  564. /* > */
  565. /* > \param[out] RCOND */
  566. /* > \verbatim */
  567. /* > RCOND is DOUBLE PRECISION */
  568. /* > The estimate of the reciprocal condition number of the matrix */
  569. /* > A. If RCOND is less than the machine precision (in */
  570. /* > particular, if RCOND = 0), the matrix is singular to working */
  571. /* > precision. This condition is indicated by a return code of */
  572. /* > INFO > 0. */
  573. /* > \endverbatim */
  574. /* > */
  575. /* > \param[out] FERR */
  576. /* > \verbatim */
  577. /* > FERR is DOUBLE PRECISION array, dimension (NRHS) */
  578. /* > The estimated forward error bound for each solution vector */
  579. /* > X(j) (the j-th column of the solution matrix X). */
  580. /* > If XTRUE is the true solution corresponding to X(j), FERR(j) */
  581. /* > is an estimated upper bound for the magnitude of the largest */
  582. /* > element in (X(j) - XTRUE) divided by the magnitude of the */
  583. /* > largest element in X(j). The estimate is as reliable as */
  584. /* > the estimate for RCOND, and is almost always a slight */
  585. /* > overestimate of the true error. */
  586. /* > \endverbatim */
  587. /* > */
  588. /* > \param[out] BERR */
  589. /* > \verbatim */
  590. /* > BERR is DOUBLE PRECISION array, dimension (NRHS) */
  591. /* > The componentwise relative backward error of each solution */
  592. /* > vector X(j) (i.e., the smallest relative change in */
  593. /* > any element of A or B that makes X(j) an exact solution). */
  594. /* > \endverbatim */
  595. /* > */
  596. /* > \param[out] WORK */
  597. /* > \verbatim */
  598. /* > WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)) */
  599. /* > On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
  600. /* > \endverbatim */
  601. /* > */
  602. /* > \param[in] LWORK */
  603. /* > \verbatim */
  604. /* > LWORK is INTEGER */
  605. /* > The length of WORK. LWORK >= f2cmax(1,3*N), and for best */
  606. /* > performance, when FACT = 'N', LWORK >= f2cmax(1,3*N,N*NB), where */
  607. /* > NB is the optimal blocksize for DSYTRF. */
  608. /* > */
  609. /* > If LWORK = -1, then a workspace query is assumed; the routine */
  610. /* > only calculates the optimal size of the WORK array, returns */
  611. /* > this value as the first entry of the WORK array, and no error */
  612. /* > message related to LWORK is issued by XERBLA. */
  613. /* > \endverbatim */
  614. /* > */
  615. /* > \param[out] IWORK */
  616. /* > \verbatim */
  617. /* > IWORK is INTEGER array, dimension (N) */
  618. /* > \endverbatim */
  619. /* > */
  620. /* > \param[out] INFO */
  621. /* > \verbatim */
  622. /* > INFO is INTEGER */
  623. /* > = 0: successful exit */
  624. /* > < 0: if INFO = -i, the i-th argument had an illegal value */
  625. /* > > 0: if INFO = i, and i is */
  626. /* > <= N: D(i,i) is exactly zero. The factorization */
  627. /* > has been completed but the factor D is exactly */
  628. /* > singular, so the solution and error bounds could */
  629. /* > not be computed. RCOND = 0 is returned. */
  630. /* > = N+1: D is nonsingular, but RCOND is less than machine */
  631. /* > precision, meaning that the matrix is singular */
  632. /* > to working precision. Nevertheless, the */
  633. /* > solution and error bounds are computed because */
  634. /* > there are a number of situations where the */
  635. /* > computed solution can be more accurate than the */
  636. /* > value of RCOND would suggest. */
  637. /* > \endverbatim */
  638. /* Authors: */
  639. /* ======== */
  640. /* > \author Univ. of Tennessee */
  641. /* > \author Univ. of California Berkeley */
  642. /* > \author Univ. of Colorado Denver */
  643. /* > \author NAG Ltd. */
  644. /* > \date April 2012 */
  645. /* > \ingroup doubleSYsolve */
  646. /* ===================================================================== */
  647. /* Subroutine */ int dsysvx_(char *fact, char *uplo, integer *n, integer *
  648. nrhs, doublereal *a, integer *lda, doublereal *af, integer *ldaf,
  649. integer *ipiv, doublereal *b, integer *ldb, doublereal *x, integer *
  650. ldx, doublereal *rcond, doublereal *ferr, doublereal *berr,
  651. doublereal *work, integer *lwork, integer *iwork, integer *info)
  652. {
  653. /* System generated locals */
  654. integer a_dim1, a_offset, af_dim1, af_offset, b_dim1, b_offset, x_dim1,
  655. x_offset, i__1, i__2;
  656. /* Local variables */
  657. extern logical lsame_(char *, char *);
  658. doublereal anorm;
  659. integer nb;
  660. extern doublereal dlamch_(char *);
  661. logical nofact;
  662. extern /* Subroutine */ int dlacpy_(char *, integer *, integer *,
  663. doublereal *, integer *, doublereal *, integer *),
  664. xerbla_(char *, integer *, ftnlen);
  665. extern integer ilaenv_(integer *, char *, char *, integer *, integer *,
  666. integer *, integer *, ftnlen, ftnlen);
  667. extern doublereal dlansy_(char *, char *, integer *, doublereal *,
  668. integer *, doublereal *);
  669. extern /* Subroutine */ int dsycon_(char *, integer *, doublereal *,
  670. integer *, integer *, doublereal *, doublereal *, doublereal *,
  671. integer *, integer *), dsyrfs_(char *, integer *, integer
  672. *, doublereal *, integer *, doublereal *, integer *, integer *,
  673. doublereal *, integer *, doublereal *, integer *, doublereal *,
  674. doublereal *, doublereal *, integer *, integer *),
  675. dsytrf_(char *, integer *, doublereal *, integer *, integer *,
  676. doublereal *, integer *, integer *);
  677. integer lwkopt;
  678. logical lquery;
  679. extern /* Subroutine */ int dsytrs_(char *, integer *, integer *,
  680. doublereal *, integer *, integer *, doublereal *, integer *,
  681. integer *);
  682. /* -- LAPACK driver routine (version 3.7.0) -- */
  683. /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
  684. /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
  685. /* April 2012 */
  686. /* ===================================================================== */
  687. /* Test the input parameters. */
  688. /* Parameter adjustments */
  689. a_dim1 = *lda;
  690. a_offset = 1 + a_dim1 * 1;
  691. a -= a_offset;
  692. af_dim1 = *ldaf;
  693. af_offset = 1 + af_dim1 * 1;
  694. af -= af_offset;
  695. --ipiv;
  696. b_dim1 = *ldb;
  697. b_offset = 1 + b_dim1 * 1;
  698. b -= b_offset;
  699. x_dim1 = *ldx;
  700. x_offset = 1 + x_dim1 * 1;
  701. x -= x_offset;
  702. --ferr;
  703. --berr;
  704. --work;
  705. --iwork;
  706. /* Function Body */
  707. *info = 0;
  708. nofact = lsame_(fact, "N");
  709. lquery = *lwork == -1;
  710. if (! nofact && ! lsame_(fact, "F")) {
  711. *info = -1;
  712. } else if (! lsame_(uplo, "U") && ! lsame_(uplo,
  713. "L")) {
  714. *info = -2;
  715. } else if (*n < 0) {
  716. *info = -3;
  717. } else if (*nrhs < 0) {
  718. *info = -4;
  719. } else if (*lda < f2cmax(1,*n)) {
  720. *info = -6;
  721. } else if (*ldaf < f2cmax(1,*n)) {
  722. *info = -8;
  723. } else if (*ldb < f2cmax(1,*n)) {
  724. *info = -11;
  725. } else if (*ldx < f2cmax(1,*n)) {
  726. *info = -13;
  727. } else /* if(complicated condition) */ {
  728. /* Computing MAX */
  729. i__1 = 1, i__2 = *n * 3;
  730. if (*lwork < f2cmax(i__1,i__2) && ! lquery) {
  731. *info = -18;
  732. }
  733. }
  734. if (*info == 0) {
  735. /* Computing MAX */
  736. i__1 = 1, i__2 = *n * 3;
  737. lwkopt = f2cmax(i__1,i__2);
  738. if (nofact) {
  739. nb = ilaenv_(&c__1, "DSYTRF", uplo, n, &c_n1, &c_n1, &c_n1, (
  740. ftnlen)6, (ftnlen)1);
  741. /* Computing MAX */
  742. i__1 = lwkopt, i__2 = *n * nb;
  743. lwkopt = f2cmax(i__1,i__2);
  744. }
  745. work[1] = (doublereal) lwkopt;
  746. }
  747. if (*info != 0) {
  748. i__1 = -(*info);
  749. xerbla_("DSYSVX", &i__1, (ftnlen)6);
  750. return 0;
  751. } else if (lquery) {
  752. return 0;
  753. }
  754. if (nofact) {
  755. /* Compute the factorization A = U*D*U**T or A = L*D*L**T. */
  756. dlacpy_(uplo, n, n, &a[a_offset], lda, &af[af_offset], ldaf);
  757. dsytrf_(uplo, n, &af[af_offset], ldaf, &ipiv[1], &work[1], lwork,
  758. info);
  759. /* Return if INFO is non-zero. */
  760. if (*info > 0) {
  761. *rcond = 0.;
  762. return 0;
  763. }
  764. }
  765. /* Compute the norm of the matrix A. */
  766. anorm = dlansy_("I", uplo, n, &a[a_offset], lda, &work[1]);
  767. /* Compute the reciprocal of the condition number of A. */
  768. dsycon_(uplo, n, &af[af_offset], ldaf, &ipiv[1], &anorm, rcond, &work[1],
  769. &iwork[1], info);
  770. /* Compute the solution vectors X. */
  771. dlacpy_("Full", n, nrhs, &b[b_offset], ldb, &x[x_offset], ldx);
  772. dsytrs_(uplo, n, nrhs, &af[af_offset], ldaf, &ipiv[1], &x[x_offset], ldx,
  773. info);
  774. /* Use iterative refinement to improve the computed solutions and */
  775. /* compute error bounds and backward error estimates for them. */
  776. dsyrfs_(uplo, n, nrhs, &a[a_offset], lda, &af[af_offset], ldaf, &ipiv[1],
  777. &b[b_offset], ldb, &x[x_offset], ldx, &ferr[1], &berr[1], &work[1]
  778. , &iwork[1], info);
  779. /* Set INFO = N+1 if the matrix is singular to working precision. */
  780. if (*rcond < dlamch_("Epsilon")) {
  781. *info = *n + 1;
  782. }
  783. work[1] = (doublereal) lwkopt;
  784. return 0;
  785. /* End of DSYSVX */
  786. } /* dsysvx_ */