NUMcblas.h 35 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776
  1. #ifndef _NUMcblas_h_
  2. #define _NUMcblas_h_
  3. /* NUMcblas.h
  4. *
  5. * Copyright (C) 1994-2011 David Weenink
  6. *
  7. * This code is free software; you can redistribute it and/or modify
  8. * it under the terms of the GNU General Public License as published by
  9. * the Free Software Foundation; either version 2 of the License, or (at
  10. * your option) any later version.
  11. *
  12. * This code is distributed in the hope that it will be useful, but
  13. * WITHOUT ANY WARRANTY; without even the implied warranty of
  14. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  15. * General Public License for more details.
  16. *
  17. * You should have received a copy of the GNU General Public License
  18. * along with this work. If not, see <http://www.gnu.org/licenses/>.
  19. */
  20. /*
  21. djmw 20020923 GPL header
  22. djmw 20110308 Latest modification
  23. */
  24. #define xerbla_(src,info) Melder_throw (Melder_peek8to32 (src), U": parameter ", *info, U" not correct!")
  25. int NUMblas_daxpy (integer *n, double *da, double *dx, integer *incx, double *dy, integer *incy);
  26. /* constant times a vector plus a vector.
  27. uses unrolled loops for increments equal to one. */
  28. int NUMblas_dcopy (integer *n, double *dx, integer *incx, double *dy, integer *incy);
  29. /* copies a vector, x, to a vector, y.
  30. uses unrolled loops for increments equal to one.
  31. */
  32. double NUMblas_ddot (integer *n, double *dx, integer *incx, double *dy, integer *incy);
  33. /* forms the dot product of two vectors.
  34. uses unrolled loops for increments equal to one. */
  35. int NUMblas_dgemm (const char *transa, const char *transb, integer *m, integer *n, integer *k, double *alpha,
  36. double *a, integer *lda, double *b, integer *ldb, double *beta, double *c, integer *ldc);
  37. /* Purpose
  38. =======
  39. NUMblas_dgemm performs one of the matrix-matrix operations
  40. C := alpha*op( A )*op( B ) + beta*C,
  41. where op( X ) is one of
  42. op( X ) = X or op( X ) = X',
  43. alpha and beta are scalars, and A, B and C are matrices, with op( A )
  44. an m by k matrix, op( B ) a k by n matrix and C an m by n matrix.
  45. Parameters
  46. ==========
  47. TRANSA - char*.
  48. On entry, TRANSA specifies the form of op( A ) to be used in
  49. the matrix multiplication as follows:
  50. TRANSA = 'N' or 'n', op( A ) = A.
  51. TRANSA = 'T' or 't', op( A ) = A'.
  52. TRANSA = 'C' or 'c', op( A ) = A'.
  53. Unchanged on exit.
  54. TRANSB - char*.
  55. On entry, TRANSB specifies the form of op( B ) to be used in
  56. the matrix multiplication as follows:
  57. TRANSB = 'N' or 'n', op( B ) = B.
  58. TRANSB = 'T' or 't', op( B ) = B'.
  59. TRANSB = 'C' or 'c', op( B ) = B'.
  60. Unchanged on exit.
  61. M - long.
  62. On entry, M specifies the number of rows of the matrix
  63. op( A ) and of the matrix C. M must be at least zero.
  64. Unchanged on exit.
  65. N - long.
  66. On entry, N specifies the number of columns of the matrix
  67. op( B ) and the number of columns of the matrix C. N must be
  68. at least zero.
  69. Unchanged on exit.
  70. K - long.
  71. On entry, K specifies the number of columns of the matrix
  72. op( A ) and the number of rows of the matrix op( B ). K must
  73. be at least zero.
  74. Unchanged on exit.
  75. ALPHA - double.
  76. On entry, ALPHA specifies the scalar alpha.
  77. Unchanged on exit.
  78. A - double array of DIMENSION ( LDA, ka ), where ka is
  79. k when TRANSA = 'N' or 'n', and is m otherwise.
  80. Before entry with TRANSA = 'N' or 'n', the leading m by k
  81. part of the array A must contain the matrix A, otherwise
  82. the leading k by m part of the array A must contain the
  83. matrix A.
  84. Unchanged on exit.
  85. LDA - long.
  86. On entry, LDA specifies the first dimension of A as declared
  87. in the calling (sub) program. When TRANSA = 'N' or 'n' then
  88. LDA must be at least max( 1, m ), otherwise LDA must be at
  89. least max( 1, k ).
  90. Unchanged on exit.
  91. B - double array of DIMENSION ( LDB, kb ), where kb is
  92. n when TRANSB = 'N' or 'n', and is k otherwise.
  93. Before entry with TRANSB = 'N' or 'n', the leading k by n
  94. part of the array B must contain the matrix B, otherwise
  95. the leading n by k part of the array B must contain the
  96. matrix B.
  97. Unchanged on exit.
  98. LDB - long.
  99. On entry, LDB specifies the first dimension of B as declared
  100. in the calling (sub) program. When TRANSB = 'N' or 'n' then
  101. LDB must be at least max( 1, k ), otherwise LDB must be at
  102. least max( 1, n ).
  103. Unchanged on exit.
  104. BETA - double.
  105. On entry, BETA specifies the scalar beta. When BETA is
  106. supplied as zero then C need not be set on input.
  107. Unchanged on exit.
  108. C - double array of DIMENSION ( LDC, n ).
  109. Before entry, the leading m by n part of the array C must
  110. contain the matrix C, except when beta is zero, in which
  111. case C need not be set on entry.
  112. On exit, the array C is overwritten by the m by n matrix
  113. ( alpha*op( A )*op( B ) + beta*C ).
  114. LDC - long.
  115. On entry, LDC specifies the first dimension of C as declared
  116. in the calling (sub) program. LDC must be at least
  117. max( 1, m ).
  118. Unchanged on exit.
  119. Level 3 Blas routine.
  120. -- Written on 8-February-1989.
  121. Jack Dongarra, Argonne National Laboratory.
  122. Iain Duff, AERE Harwell.
  123. Jeremy Du Croz, Numerical Algorithms Group Ltd.
  124. Sven Hammarling, Numerical Algorithms Group Ltd.
  125. */
  126. int NUMblas_dger (integer *m, integer *n, double *alpha, double *x, integer *incx, double *y,
  127. integer *incy, double *a, integer *lda);
  128. /* Purpose
  129. =======
  130. NUMblas_dger performs the rank 1 operation
  131. A := alpha*x*y' + A,
  132. where alpha is a scalar, x is an m element vector, y is an n element
  133. vector and A is an m by n matrix.
  134. Parameters
  135. ==========
  136. M - long.
  137. On entry, M specifies the number of rows of the matrix A.
  138. M must be at least zero.
  139. Unchanged on exit.
  140. N - long.
  141. On entry, N specifies the number of columns of the matrix A.
  142. N must be at least zero.
  143. Unchanged on exit.
  144. ALPHA - double.
  145. On entry, ALPHA specifies the scalar alpha.
  146. Unchanged on exit.
  147. X - double array of dimension at least
  148. ( 1 + ( m - 1 )*abs( INCX ) ).
  149. Before entry, the incremented array X must contain the m
  150. element vector x.
  151. Unchanged on exit.
  152. INCX - long.
  153. On entry, INCX specifies the increment for the elements of
  154. X. INCX must not be zero.
  155. Unchanged on exit.
  156. Y - double array of dimension at least
  157. ( 1 + ( n - 1 )*abs( INCY ) ).
  158. Before entry, the incremented array Y must contain the n
  159. element vector y.
  160. Unchanged on exit.
  161. INCY - long.
  162. On entry, INCY specifies the increment for the elements of
  163. Y. INCY must not be zero.
  164. Unchanged on exit.
  165. A - double array of DIMENSION ( LDA, n ).
  166. Before entry, the leading m by n part of the array A must
  167. contain the matrix of coefficients. On exit, A is
  168. overwritten by the updated matrix.
  169. LDA - long.
  170. On entry, LDA specifies the first dimension of A as declared
  171. in the calling (sub) program. LDA must be at least
  172. max( 1, m ).
  173. Unchanged on exit.
  174. Level 2 Blas routine.
  175. -- Written on 22-October-1986.
  176. Jack Dongarra, Argonne National Lab.
  177. Jeremy Du Croz, Nag Central Office.
  178. Sven Hammarling, Nag Central Office.
  179. Richard Hanson, Sandia National Labs.
  180. */
  181. int NUMblas_dgemv (const char *trans, integer *m, integer *n, double *alpha, double *a, integer *lda,
  182. double *x, integer *incx, double *beta, double *y, integer *incy);
  183. /* Purpose
  184. =======
  185. NUMblas_dgemv performs one of the matrix-vector operations
  186. y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y,
  187. where alpha and beta are scalars, x and y are vectors and A is an
  188. m by n matrix.
  189. Parameters
  190. ==========
  191. TRANS - char*.
  192. On entry, TRANS specifies the operation to be performed as
  193. follows:
  194. TRANS = 'N' or 'n' y := alpha*A*x + beta*y.
  195. TRANS = 'T' or 't' y := alpha*A'*x + beta*y.
  196. TRANS = 'C' or 'c' y := alpha*A'*x + beta*y.
  197. Unchanged on exit.
  198. M - long.
  199. On entry, M specifies the number of rows of the matrix A.
  200. M must be at least zero.
  201. Unchanged on exit.
  202. N - long.
  203. On entry, N specifies the number of columns of the matrix A.
  204. N must be at least zero.
  205. Unchanged on exit.
  206. ALPHA - double.
  207. On entry, ALPHA specifies the scalar alpha.
  208. Unchanged on exit.
  209. A - double array of DIMENSION ( LDA, n ).
  210. Before entry, the leading m by n part of the array A must
  211. contain the matrix of coefficients.
  212. Unchanged on exit.
  213. LDA - long.
  214. On entry, LDA specifies the first dimension of A as declared
  215. in the calling (sub) program. LDA must be at least
  216. max( 1, m ).
  217. Unchanged on exit.
  218. X - double array of DIMENSION at least
  219. ( 1 + ( n - 1 )*abs( INCX ) ) when TRANS = 'N' or 'n'
  220. and at least
  221. ( 1 + ( m - 1 )*abs( INCX ) ) otherwise.
  222. Before entry, the incremented array X must contain the
  223. vector x.
  224. Unchanged on exit.
  225. INCX - long.
  226. On entry, INCX specifies the increment for the elements of
  227. X. INCX must not be zero.
  228. Unchanged on exit.
  229. BETA - double.
  230. On entry, BETA specifies the scalar beta. When BETA is
  231. supplied as zero then Y need not be set on input.
  232. Unchanged on exit.
  233. Y - double array of DIMENSION at least
  234. ( 1 + ( m - 1 )*abs( INCY ) ) when TRANS = 'N' or 'n'
  235. and at least
  236. ( 1 + ( n - 1 )*abs( INCY ) ) otherwise.
  237. Before entry with BETA non-zero, the incremented array Y
  238. must contain the vector y. On exit, Y is overwritten by the
  239. updated vector y.
  240. INCY - long.
  241. On entry, INCY specifies the increment for the elements of
  242. Y. INCY must not be zero.
  243. Unchanged on exit.
  244. Level 2 Blas routine.
  245. -- Written on 22-October-1986.
  246. Jack Dongarra, Argonne National Lab.
  247. Jeremy Du Croz, Nag Central Office.
  248. Sven Hammarling, Nag Central Office.
  249. Richard Hanson, Sandia National Labs.
  250. */
  251. double NUMblas_dlamch (const char *cmach);
  252. /* Purpose
  253. =======
  254. NUMblas_dlamch determines double machine parameters.
  255. Arguments
  256. =========
  257. CMACH (input) char*
  258. Specifies the value to be returned by DLAMCH:
  259. = 'E' or 'e', DLAMCH := eps
  260. = 'S' or 's , DLAMCH := sfmin
  261. = 'B' or 'b', DLAMCH := base
  262. = 'P' or 'p', DLAMCH := eps*base
  263. = 'N' or 'n', DLAMCH := t
  264. = 'R' or 'r', DLAMCH := rnd
  265. = 'M' or 'm', DLAMCH := emin
  266. = 'U' or 'u', DLAMCH := rmin
  267. = 'L' or 'l', DLAMCH := emax
  268. = 'O' or 'o', DLAMCH := rmax
  269. where
  270. eps = relative machine precision
  271. sfmin = safe minimum, such that 1/sfmin does not overflow
  272. base = base of the machine
  273. prec = eps*base
  274. t = number of (base) digits in the mantissa
  275. rnd = 1.0 when rounding occurs in addition, 0.0 otherwise
  276. emin = minimum exponent before (gradual) underflow
  277. rmin = underflow threshold - base**(emin-1)
  278. emax = largest exponent before overflow
  279. rmax = overflow threshold - (base**emax)*(1-eps)
  280. =====================================================================
  281. */
  282. double NUMblas_dnrm2 (integer *n, double *x, integer *incx);
  283. /*
  284. NUMblas_dnrm2 returns the euclidean norm of a vector via the function
  285. name, so that
  286. DNRM2 := sqrt( x'*x )
  287. -- This version written on 25-October-1982.
  288. Modified on 14-October-1993 to inline the call to DLASSQ.
  289. Sven Hammarling, Nag Ltd.
  290. Parameter adjustments
  291. */
  292. int NUMblas_drot (integer *n, double *dx, integer *incx, double *dy, integer *incy,
  293. double *c__, double *s);
  294. /* applies a plane rotation. */
  295. int NUMblas_dscal (integer *n, double *da, double *dx, integer *incx);
  296. /* scales a vector by a constant.
  297. uses unrolled loops for increment equal to one.
  298. */
  299. int NUMblas_dswap (integer *n, double *dx, integer *incx, double *dy, integer *incy);
  300. /* interchanges two vectors.
  301. uses unrolled loops for increments equal one. */
  302. int NUMblas_dsymv (const char *uplo, integer *n, double *alpha, double *a, integer *lda,
  303. double *x, integer *incx, double *beta, double *y, integer *incy);
  304. /* Purpose
  305. =======
  306. NUMblas_dsymv performs the matrix-vector operation
  307. y := alpha*A*x + beta*y,
  308. where alpha and beta are scalars, x and y are n element vectors and
  309. A is an n by n symmetric matrix.
  310. Parameters
  311. ==========
  312. UPLO - char*.
  313. On entry, UPLO specifies whether the upper or lower
  314. triangular part of the array A is to be referenced as
  315. follows:
  316. UPLO = 'U' or 'u' Only the upper triangular part of A
  317. is to be referenced.
  318. UPLO = 'L' or 'l' Only the lower triangular part of A
  319. is to be referenced.
  320. Unchanged on exit.
  321. N - long.
  322. On entry, N specifies the order of the matrix A.
  323. N must be at least zero.
  324. Unchanged on exit.
  325. ALPHA - double.
  326. On entry, ALPHA specifies the scalar alpha.
  327. Unchanged on exit.
  328. A - double array of DIMENSION ( LDA, n ).
  329. Before entry with UPLO = 'U' or 'u', the leading n by n
  330. upper triangular part of the array A must contain the upper
  331. triangular part of the symmetric matrix and the strictly
  332. lower triangular part of A is not referenced.
  333. Before entry with UPLO = 'L' or 'l', the leading n by n
  334. lower triangular part of the array A must contain the lower
  335. triangular part of the symmetric matrix and the strictly
  336. upper triangular part of A is not referenced.
  337. Unchanged on exit.
  338. LDA - long.
  339. On entry, LDA specifies the first dimension of A as declared
  340. in the calling (sub) program. LDA must be at least
  341. max( 1, n ).
  342. Unchanged on exit.
  343. X - double array of dimension at least
  344. ( 1 + ( n - 1 )*abs( INCX ) ).
  345. Before entry, the incremented array X must contain the n
  346. element vector x.
  347. Unchanged on exit.
  348. INCX - long.
  349. On entry, INCX specifies the increment for the elements of
  350. X. INCX must not be zero.
  351. Unchanged on exit.
  352. BETA - double.
  353. On entry, BETA specifies the scalar beta. When BETA is
  354. supplied as zero then Y need not be set on input.
  355. Unchanged on exit.
  356. Y - double array of dimension at least
  357. ( 1 + ( n - 1 )*abs( INCY ) ).
  358. Before entry, the incremented array Y must contain the n
  359. element vector y. On exit, Y is overwritten by the updated
  360. vector y.
  361. INCY - long.
  362. On entry, INCY specifies the increment for the elements of
  363. Y. INCY must not be zero.
  364. Unchanged on exit.
  365. */
  366. int NUMblas_dsyr2 (const char *uplo, integer *n, double *alpha, double *x, integer *incx,
  367. double *y, integer *incy, double *a, integer *lda);
  368. /* Purpose
  369. =======
  370. NUMblas_dsyr2 performs the symmetric rank 2 operation
  371. A := alpha*x*y' + alpha*y*x' + A,
  372. where alpha is a scalar, x and y are n element vectors and A is an n
  373. by n symmetric matrix.
  374. Parameters
  375. ==========
  376. UPLO - char*.
  377. On entry, UPLO specifies whether the upper or lower
  378. triangular part of the array A is to be referenced as
  379. follows:
  380. UPLO = 'U' or 'u' Only the upper triangular part of A
  381. is to be referenced.
  382. UPLO = 'L' or 'l' Only the lower triangular part of A
  383. is to be referenced.
  384. Unchanged on exit.
  385. N - long.
  386. On entry, N specifies the order of the matrix A.
  387. N must be at least zero.
  388. Unchanged on exit.
  389. ALPHA - double.
  390. On entry, ALPHA specifies the scalar alpha.
  391. Unchanged on exit.
  392. X - double array of dimension at least
  393. ( 1 + ( n - 1 )*abs( INCX ) ).
  394. Before entry, the incremented array X must contain the n
  395. element vector x.
  396. Unchanged on exit.
  397. INCX - long.
  398. On entry, INCX specifies the increment for the elements of
  399. X. INCX must not be zero.
  400. Unchanged on exit.
  401. Y - double array of dimension at least
  402. ( 1 + ( n - 1 )*abs( INCY ) ).
  403. Before entry, the incremented array Y must contain the n
  404. element vector y.
  405. Unchanged on exit.
  406. INCY - long.
  407. On entry, INCY specifies the increment for the elements of
  408. Y. INCY must not be zero.
  409. Unchanged on exit.
  410. A - double array of DIMENSION ( LDA, n ).
  411. Before entry with UPLO = 'U' or 'u', the leading n by n
  412. upper triangular part of the array A must contain the upper
  413. triangular part of the symmetric matrix and the strictly
  414. lower triangular part of A is not referenced. On exit, the
  415. upper triangular part of the array A is overwritten by the
  416. upper triangular part of the updated matrix.
  417. Before entry with UPLO = 'L' or 'l', the leading n by n
  418. lower triangular part of the array A must contain the lower
  419. triangular part of the symmetric matrix and the strictly
  420. upper triangular part of A is not referenced. On exit, the
  421. lower triangular part of the array A is overwritten by the
  422. lower triangular part of the updated matrix.
  423. LDA - long.
  424. On entry, LDA specifies the first dimension of A as declared
  425. in the calling (sub) program. LDA must be at least
  426. max( 1, n ).
  427. Unchanged on exit.
  428. */
  429. int NUMblas_dsyr2k (const char *uplo, const char *trans, integer *n, integer *k, double *alpha, double *a,
  430. integer *lda, double *b, integer *ldb, double *beta, double *c, integer *ldc);
  431. /* Purpose
  432. =======
  433. NUMblas_dsyr2k performs one of the symmetric rank 2k operations
  434. C := alpha*A*B' + alpha*B*A' + beta*C,
  435. or
  436. C := alpha*A'*B + alpha*B'*A + beta*C,
  437. where alpha and beta are scalars, C is an n by n symmetric matrix
  438. and A and B are n by k matrices in the first case and k by n
  439. matrices in the second case.
  440. Parameters
  441. ==========
  442. UPLO - char*.
  443. On entry, UPLO specifies whether the upper or lower
  444. triangular part of the array C is to be referenced as
  445. follows:
  446. UPLO = 'U' or 'u' Only the upper triangular part of C
  447. is to be referenced.
  448. UPLO = 'L' or 'l' Only the lower triangular part of C
  449. is to be referenced.
  450. Unchanged on exit.
  451. TRANS - char*.
  452. On entry, TRANS specifies the operation to be performed as
  453. follows:
  454. TRANS = 'N' or 'n' C := alpha*A*B' + alpha*B*A' +
  455. beta*C.
  456. TRANS = 'T' or 't' C := alpha*A'*B + alpha*B'*A +
  457. beta*C.
  458. TRANS = 'C' or 'c' C := alpha*A'*B + alpha*B'*A +
  459. beta*C.
  460. Unchanged on exit.
  461. N - long.
  462. On entry, N specifies the order of the matrix C. N must be
  463. at least zero.
  464. Unchanged on exit.
  465. K - long.
  466. On entry with TRANS = 'N' or 'n', K specifies the number
  467. of columns of the matrices A and B, and on entry with
  468. TRANS = 'T' or 't' or 'C' or 'c', K specifies the number
  469. of rows of the matrices A and B. K must be at least zero.
  470. Unchanged on exit.
  471. ALPHA - double.
  472. On entry, ALPHA specifies the scalar alpha.
  473. Unchanged on exit.
  474. A - double array of DIMENSION ( LDA, ka ), where ka is
  475. k when TRANS = 'N' or 'n', and is n otherwise.
  476. Before entry with TRANS = 'N' or 'n', the leading n by k
  477. part of the array A must contain the matrix A, otherwise
  478. the leading k by n part of the array A must contain the
  479. matrix A.
  480. Unchanged on exit.
  481. LDA - long.
  482. On entry, LDA specifies the first dimension of A as declared
  483. in the calling (sub) program. When TRANS = 'N' or 'n'
  484. then LDA must be at least max( 1, n ), otherwise LDA must
  485. be at least max( 1, k ).
  486. Unchanged on exit.
  487. B - double array of DIMENSION ( LDB, kb ), where kb is
  488. k when TRANS = 'N' or 'n', and is n otherwise.
  489. Before entry with TRANS = 'N' or 'n', the leading n by k
  490. part of the array B must contain the matrix B, otherwise
  491. the leading k by n part of the array B must contain the
  492. matrix B.
  493. Unchanged on exit.
  494. LDB - long.
  495. On entry, LDB specifies the first dimension of B as declared
  496. in the calling (sub) program. When TRANS = 'N' or 'n'
  497. then LDB must be at least max( 1, n ), otherwise LDB must
  498. be at least max( 1, k ).
  499. Unchanged on exit.
  500. BETA - double.
  501. On entry, BETA specifies the scalar beta.
  502. Unchanged on exit.
  503. C - double array of DIMENSION ( LDC, n ).
  504. Before entry with UPLO = 'U' or 'u', the leading n by n
  505. upper triangular part of the array C must contain the upper
  506. triangular part of the symmetric matrix and the strictly
  507. lower triangular part of C is not referenced. On exit, the
  508. upper triangular part of the array C is overwritten by the
  509. upper triangular part of the updated matrix.
  510. Before entry with UPLO = 'L' or 'l', the leading n by n
  511. lower triangular part of the array C must contain the lower
  512. triangular part of the symmetric matrix and the strictly
  513. upper triangular part of C is not referenced. On exit, the
  514. lower triangular part of the array C is overwritten by the
  515. lower triangular part of the updated matrix.
  516. LDC - long.
  517. On entry, LDC specifies the first dimension of C as declared
  518. in the calling (sub) program. LDC must be at least
  519. max( 1, n ).
  520. Unchanged on exit.
  521. Level 3 Blas routine.
  522. */
  523. int NUMblas_dtrmm (const char *side, const char *uplo, const char *transa, const char *diag,
  524. integer *m, integer *n, double *alpha, double *a, integer *lda,
  525. double *b, integer *ldb);
  526. /* Purpose
  527. =======
  528. NUMblas_dtrmm performs one of the matrix-matrix operations
  529. B := alpha*op( A )*B, or B := alpha*B*op( A ),
  530. where alpha is a scalar, B is an m by n matrix, A is a unit, or
  531. non-unit, upper or lower triangular matrix and op( A ) is one of
  532. op( A ) = A or op( A ) = A'.
  533. Parameters
  534. ==========
  535. SIDE - char*.
  536. On entry, SIDE specifies whether op( A ) multiplies B from
  537. the left or right as follows:
  538. SIDE = 'L' or 'l' B := alpha*op( A )*B.
  539. SIDE = 'R' or 'r' B := alpha*B*op( A ).
  540. Unchanged on exit.
  541. UPLO - char*.
  542. On entry, UPLO specifies whether the matrix A is an upper or
  543. lower triangular matrix as follows:
  544. UPLO = 'U' or 'u' A is an upper triangular matrix.
  545. UPLO = 'L' or 'l' A is a lower triangular matrix.
  546. Unchanged on exit.
  547. TRANSA - char*.
  548. On entry, TRANSA specifies the form of op( A ) to be used in
  549. the matrix multiplication as follows:
  550. TRANSA = 'N' or 'n' op( A ) = A.
  551. TRANSA = 'T' or 't' op( A ) = A'.
  552. TRANSA = 'C' or 'c' op( A ) = A'.
  553. Unchanged on exit.
  554. DIAG - char*.
  555. On entry, DIAG specifies whether or not A is unit triangular
  556. as follows:
  557. DIAG = 'U' or 'u' A is assumed to be unit triangular.
  558. DIAG = 'N' or 'n' A is not assumed to be unit
  559. triangular.
  560. Unchanged on exit.
  561. M - long.
  562. On entry, M specifies the number of rows of B. M must be at
  563. least zero.
  564. Unchanged on exit.
  565. N - long.
  566. On entry, N specifies the number of columns of B. N must be
  567. at least zero.
  568. Unchanged on exit.
  569. ALPHA - double.
  570. On entry, ALPHA specifies the scalar alpha. When alpha is
  571. zero then A is not referenced and B need not be set before
  572. entry.
  573. Unchanged on exit.
  574. A - double array of DIMENSION ( LDA, k ), where k is m
  575. when SIDE = 'L' or 'l' and is n when SIDE = 'R' or 'r'.
  576. Before entry with UPLO = 'U' or 'u', the leading k by k
  577. upper triangular part of the array A must contain the upper
  578. triangular matrix and the strictly lower triangular part of
  579. A is not referenced.
  580. Before entry with UPLO = 'L' or 'l', the leading k by k
  581. lower triangular part of the array A must contain the lower
  582. triangular matrix and the strictly upper triangular part of
  583. A is not referenced.
  584. Note that when DIAG = 'U' or 'u', the diagonal elements of
  585. A are not referenced either, but are assumed to be unity.
  586. Unchanged on exit.
  587. LDA - long.
  588. On entry, LDA specifies the first dimension of A as declared
  589. in the calling (sub) program. When SIDE = 'L' or 'l' then
  590. LDA must be at least max( 1, m ), when SIDE = 'R' or 'r'
  591. then LDA must be at least max( 1, n ).
  592. Unchanged on exit.
  593. B - double array of DIMENSION ( LDB, n ).
  594. Before entry, the leading m by n part of the array B must
  595. contain the matrix B, and on exit is overwritten by the
  596. transformed matrix.
  597. LDB - long.
  598. On entry, LDB specifies the first dimension of B as declared
  599. in the calling (sub) program. LDB must be at least
  600. max( 1, m ).
  601. Unchanged on exit.
  602. Level 3 Blas routine.
  603. */
  604. int NUMblas_dtrmv (const char *uplo, const char *trans, const char *diag, integer *n,
  605. double *a, integer *lda, double *x, integer *incx);
  606. /* Purpose
  607. =======
  608. NUMblas_dtrmv performs one of the matrix-vector operations
  609. x := A*x, or x := A'*x,
  610. where x is an n element vector and A is an n by n unit, or non-unit,
  611. upper or lower triangular matrix.
  612. Parameters
  613. ==========
  614. UPLO - char*.
  615. On entry, UPLO specifies whether the matrix is an upper or
  616. lower triangular matrix as follows:
  617. UPLO = 'U' or 'u' A is an upper triangular matrix.
  618. UPLO = 'L' or 'l' A is a lower triangular matrix.
  619. Unchanged on exit.
  620. TRANS - char*.
  621. On entry, TRANS specifies the operation to be performed as
  622. follows:
  623. TRANS = 'N' or 'n' x := A*x.
  624. TRANS = 'T' or 't' x := A'*x.
  625. TRANS = 'C' or 'c' x := A'*x.
  626. Unchanged on exit.
  627. DIAG - char*.
  628. On entry, DIAG specifies whether or not A is unit
  629. triangular as follows:
  630. DIAG = 'U' or 'u' A is assumed to be unit triangular.
  631. DIAG = 'N' or 'n' A is not assumed to be unit
  632. triangular.
  633. Unchanged on exit.
  634. N - long.
  635. On entry, N specifies the order of the matrix A.
  636. N must be at least zero.
  637. Unchanged on exit.
  638. A - double array of DIMENSION ( LDA, n ).
  639. Before entry with UPLO = 'U' or 'u', the leading n by n
  640. upper triangular part of the array A must contain the upper
  641. triangular matrix and the strictly lower triangular part of
  642. A is not referenced.
  643. Before entry with UPLO = 'L' or 'l', the leading n by n
  644. lower triangular part of the array A must contain the lower
  645. triangular matrix and the strictly upper triangular part of
  646. A is not referenced.
  647. Note that when DIAG = 'U' or 'u', the diagonal elements of
  648. A are not referenced either, but are assumed to be unity.
  649. Unchanged on exit.
  650. LDA - long.
  651. On entry, LDA specifies the first dimension of A as declared
  652. in the calling (sub) program. LDA must be at least
  653. max( 1, n ).
  654. Unchanged on exit.
  655. X - double array of dimension at least
  656. ( 1 + ( n - 1 )*abs( INCX ) ).
  657. Before entry, the incremented array X must contain the n
  658. element vector x. On exit, X is overwritten with the
  659. tranformed vector x.
  660. INCX - long.
  661. On entry, INCX specifies the increment for the elements of
  662. X. INCX must not be zero.
  663. Unchanged on exit.
  664. Level 2 Blas routine.
  665. */
  666. int NUMblas_dtrsm (const char *side, const char *uplo, const char *transa, const char *diag, integer *m, integer *n,
  667. double *alpha, double *a, integer *lda, double *b, integer *ldb);
  668. /* Purpose
  669. =======
  670. NUMblas_dtrsm solves one of the matrix equations
  671. op( A )*X = alpha*B, or X*op( A ) = alpha*B,
  672. where alpha is a scalar, X and B are m by n matrices, A is a unit, or
  673. non-unit, upper or lower triangular matrix and op( A ) is one of
  674. op( A ) = A or op( A ) = A'.
  675. The matrix X is overwritten on B.
  676. Parameters
  677. ==========
  678. SIDE - char*.
  679. On entry, SIDE specifies whether op( A ) appears on the left
  680. or right of X as follows:
  681. SIDE = 'L' or 'l' op( A )*X = alpha*B.
  682. SIDE = 'R' or 'r' X*op( A ) = alpha*B.
  683. Unchanged on exit.
  684. UPLO - char*.
  685. On entry, UPLO specifies whether the matrix A is an upper or
  686. lower triangular matrix as follows:
  687. UPLO = 'U' or 'u' A is an upper triangular matrix.
  688. UPLO = 'L' or 'l' A is a lower triangular matrix.
  689. Unchanged on exit.
  690. TRANSA - char*.
  691. On entry, TRANSA specifies the form of op( A ) to be used in
  692. the matrix multiplication as follows:
  693. TRANSA = 'N' or 'n' op( A ) = A.
  694. TRANSA = 'T' or 't' op( A ) = A'.
  695. TRANSA = 'C' or 'c' op( A ) = A'.
  696. Unchanged on exit.
  697. DIAG - char*.
  698. On entry, DIAG specifies whether or not A is unit triangular
  699. as follows:
  700. DIAG = 'U' or 'u' A is assumed to be unit triangular.
  701. DIAG = 'N' or 'n' A is not assumed to be unit
  702. triangular.
  703. Unchanged on exit.
  704. M - long.
  705. On entry, M specifies the number of rows of B. M must be at
  706. least zero.
  707. Unchanged on exit.
  708. N - long.
  709. On entry, N specifies the number of columns of B. N must be
  710. at least zero.
  711. Unchanged on exit.
  712. ALPHA - double.
  713. On entry, ALPHA specifies the scalar alpha. When alpha is
  714. zero then A is not referenced and B need not be set before
  715. entry.
  716. Unchanged on exit.
  717. A - double array of DIMENSION ( LDA, k ), where k is m
  718. when SIDE = 'L' or 'l' and is n when SIDE = 'R' or 'r'.
  719. Before entry with UPLO = 'U' or 'u', the leading k by k
  720. upper triangular part of the array A must contain the upper
  721. triangular matrix and the strictly lower triangular part of
  722. A is not referenced.
  723. Before entry with UPLO = 'L' or 'l', the leading k by k
  724. lower triangular part of the array A must contain the lower
  725. triangular matrix and the strictly upper triangular part of
  726. A is not referenced.
  727. Note that when DIAG = 'U' or 'u', the diagonal elements of
  728. A are not referenced either, but are assumed to be unity.
  729. Unchanged on exit.
  730. LDA - long.
  731. On entry, LDA specifies the first dimension of A as declared
  732. in the calling (sub) program. When SIDE = 'L' or 'l' then
  733. LDA must be at least max( 1, m ), when SIDE = 'R' or 'r'
  734. then LDA must be at least max( 1, n ).
  735. Unchanged on exit.
  736. B - double array of DIMENSION ( LDB, n ).
  737. Before entry, the leading m by n part of the array B must
  738. contain the right-hand side matrix B, and on exit is
  739. overwritten by the solution matrix X.
  740. LDB - long.
  741. On entry, LDB specifies the first dimension of B as declared
  742. in the calling (sub) program. LDB must be at least
  743. max( 1, m ).
  744. Unchanged on exit.
  745. */
  746. integer NUMblas_idamax (integer *n, double *dx, integer *incx);
  747. /* finds the index of element having max. absolute value.*/
  748. #endif /* _NUMcblas_h_ */