1 /* 2 * This file is part of the GROMACS molecular simulation package. 3 * 4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands. 5 * Copyright (c) 2001-2008, The GROMACS development team. 6 * Copyright (c) 2012,2013,2014,2019, by the GROMACS development team, led by 7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl, 8 * and including many others, as listed in the AUTHORS file in the 9 * top-level source directory and at http://www.gromacs.org. 10 * 11 * GROMACS is free software; you can redistribute it and/or 12 * modify it under the terms of the GNU Lesser General Public License 13 * as published by the Free Software Foundation; either version 2.1 14 * of the License, or (at your option) any later version. 15 * 16 * GROMACS is distributed in the hope that it will be useful, 17 * but WITHOUT ANY WARRANTY; without even the implied warranty of 18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 19 * Lesser General Public License for more details. 20 * 21 * You should have received a copy of the GNU Lesser General Public 22 * License along with GROMACS; if not, see 23 * http://www.gnu.org/licenses, or write to the Free Software Foundation, 24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. 25 * 26 * If you want to redistribute modifications to GROMACS, please 27 * consider that scientific software is very special. Version 28 * control is crucial - bugs must be traceable. We will be happy to 29 * consider code for inclusion in the official distribution, but 30 * derived work must not be called official GROMACS. Details are found 31 * in the README & COPYING files - if they are missing, get the 32 * official version at http://www.gromacs.org. 33 * 34 * To help us fund GROMACS development, we humbly ask that you cite 35 * the research papers on the package. Check out http://www.gromacs.org. 36 */ 37 /*! \internal \file 38 * \brief 39 * Header definitions for the standard LAPACK library. 40 * 41 * This is the subset of LAPACK routines used for the 42 * linear algebra operations in Gromacs. Most of the execution time 43 * will be spent in the BLAS routines, which you hopefully have an 44 * optimized version of. Gromacs includes reference implementations 45 * of both BLAS and LAPACK so it compiles everywhere, but you should 46 * really try to find a vendor or otherwise optimized version at least 47 * of BLAS for better performance. 48 * 49 * Do NOT use this code for other purposes - we only provide this as a 50 * simple fallback/reference implementation when no optimized BLAS 51 * is present. If you need an implementation for your own code 52 * there are several much faster versions out there. 53 * 54 * All routines are compatible with the LAPACK/BLAS reference implementations, 55 * meaning they assume fortran-style matrix row/column organization. 56 * 57 * There is plenty of documentation for these routines available 58 * at http://www.netlib.org/lapack , so there is no point in repeating 59 * it here. 60 */ 61 #ifndef GMX_LAPACK_H 62 #define GMX_LAPACK_H 63 64 /*! \cond */ 65 66 #include "config.h" 67 68 /* These are not required by this file, but by the internal LAPACK 69 * implementation. In principle, they could be included in each file 70 * that requires them, but this is simpler. Since the header is internal 71 * to the linearyalgebra/ module, the added complexity may not be worth it. */ 72 #include "gromacs/utility/basedefinitions.h" 73 #include "gromacs/utility/real.h" 74 75 #ifdef __cplusplus 76 extern "C" 77 { 78 #endif 79 #if 0 80 } 81 #endif 82 /* Double precision */ 83 84 void F77_FUNC(dbdsdc, DBDSDC)(const char* uplo, 85 const char* compq, 86 int* n, 87 double* d, 88 double* e, 89 double* u, 90 int* ldu, 91 double* vt, 92 int* ldvt, 93 double* q, 94 int* iq, 95 double* work, 96 int* iwork, 97 int* info); 98 99 void F77_FUNC(dgetf2, DGETF2)(int* m, int* n, double* a, int* lda, int* ipiv, int* info); 100 101 void F77_FUNC(dlamrg, DLAMRG)(int* n1, int* n2, double* a, int* dtrd1, int* dtrd2, int* index); 102 103 void F77_FUNC(dlarnv, DLARNV)(int* idist, int* iseed, int* n, double* x); 104 105 void F77_FUNC(dlasd0, DLASD0)(int* n, 106 int* sqre, 107 double* d, 108 double* e, 109 double* u, 110 int* ldu, 111 double* vt, 112 int* ldvt, 113 int* smlsiz, 114 int* iwork, 115 double* work, 116 int* info); 117 118 void F77_FUNC(dlasda, DLASDA)(int* icompq, 119 int* smlsiz, 120 int* n, 121 int* sqre, 122 double* d, 123 double* e, 124 double* u, 125 int* ldu, 126 double* vt, 127 int* k, 128 double* difl, 129 double* difr, 130 double* z, 131 double* poles, 132 int* givptr, 133 int* givcol, 134 int* ldgcol, 135 int* perm, 136 double* givnum, 137 double* c, 138 double* s, 139 double* work, 140 int* iwork, 141 int* info); 142 143 void F77_FUNC(dlasq6, DLASQ6)(int* i0, 144 int* n0, 145 double* z, 146 int* pp, 147 double* dmin, 148 double* dmin1, 149 double* dmin2, 150 double* dn, 151 double* dnm1, 152 double* dnm2); 153 154 void F77_FUNC(dorgl2, 155 DORGL2)(int* m, int* n, int* k, double* a, int* lda, double* tau, double* work, int* info); 156 157 void F77_FUNC(dbdsqr, DBDSQR)(const char* uplo, 158 int* n, 159 int* ncvt, 160 int* nru, 161 int* ncc, 162 double* d, 163 double* e, 164 double* vt, 165 int* ldvt, 166 double* u, 167 int* ldu, 168 double* c, 169 int* ldc, 170 double* work, 171 int* info); 172 173 void F77_FUNC(dgetrf, DGETRF)(int* m, int* n, double* a, int* lda, int* ipiv, int* info); 174 175 void F77_FUNC(dgetri, 176 DGETRI)(int* n, double* a, int* lda, int* ipiv, double* work, int* lwork, int* info); 177 178 void F77_FUNC(dgetrs, DGETRS)(const char* trans, 179 int* n, 180 int* nrhs, 181 double* a, 182 int* lda, 183 int* ipiv, 184 double* b, 185 int* ldb, 186 int* info); 187 188 void F77_FUNC(dtrtri, 189 DTRTRI)(const char* uplo, const char* diag, int* n, double* a, int* lda, int* info); 190 191 void F77_FUNC(dtrti2, 192 DTRTI2)(const char* uplo, const char* diag, int* n, double* a, int* lda, int* info); 193 194 double F77_FUNC(dlange, DLANGE)(const char* norm, int* m, int* n, double* a, int* lda, double* work); 195 196 void F77_FUNC(dlarrbx, DLARRBX)(int* n, 197 double* d, 198 double* l, 199 double* ld, 200 double* lld, 201 int* ifirst, 202 int* ilast, 203 double* rtol1, 204 double* rtol2, 205 int* offset, 206 double* w, 207 double* wgap, 208 double* werr, 209 double* work, 210 int* iwork, 211 int* info); 212 213 void F77_FUNC(dlasd1, DLASD1)(int* nl, 214 int* nr, 215 int* sqre, 216 double* d, 217 double* alpha, 218 double* beta, 219 double* u, 220 int* ldu, 221 double* vt, 222 int* ldvt, 223 int* idxq, 224 int* iwork, 225 double* work, 226 int* info); 227 228 void F77_FUNC(dlasdq, DLASDQ)(const char* uplo, 229 int* sqre, 230 int* n, 231 int* ncvt, 232 int* nru, 233 int* ncc, 234 double* d, 235 double* e, 236 double* vt, 237 int* ldvt, 238 double* u, 239 int* ldu, 240 double* c, 241 int* ldc, 242 double* work, 243 int* info); 244 245 void F77_FUNC(dlasr, DLASR)(const char* side, 246 const char* pivot, 247 const char* direct, 248 int* m, 249 int* n, 250 double* c, 251 double* s, 252 double* a, 253 int* lda); 254 255 void F77_FUNC(dorglq, 256 DORGLQ)(int* m, int* n, int* k, double* a, int* lda, double* tau, double* work, int* lwork, int* info); 257 258 void F77_FUNC(dormtr, DORMTR)(const char* side, 259 const char* uplo, 260 const char* trans, 261 int* m, 262 int* n, 263 double* a, 264 int* lda, 265 double* tau, 266 double* c, 267 int* ldc, 268 double* work, 269 int* lwork, 270 int* info); 271 272 void F77_FUNC(dgebd2, DGEBD2)(int* m, 273 int* n, 274 double* a, 275 int* lda, 276 double* d, 277 double* e, 278 double* tauq, 279 double* taup, 280 double* work, 281 int* info); 282 283 void F77_FUNC(dlabrd, DLABRD)(int* m, 284 int* n, 285 int* nb, 286 double* a, 287 int* lda, 288 double* d, 289 double* e, 290 double* tauq, 291 double* taup, 292 double* x, 293 int* ldx, 294 double* y, 295 int* ldy); 296 297 double F77_FUNC(dlanst, DLANST)(const char* norm, int* n, double* d, double* e); 298 299 double F77_FUNC(dlansy, 300 DLANSY)(const char* norm, const char* uplo, int* n, double* a, int* lda, double* work); 301 302 void F77_FUNC(dlarrex, DLARREX)(const char* range, 303 int* n, 304 double* vl, 305 double* vu, 306 int* il, 307 int* iu, 308 double* d, 309 double* e, 310 double* tol, 311 int* nsplit, 312 int* isplit, 313 int* m, 314 double* w, 315 int* iblock, 316 int* indexw, 317 double* gersch, 318 double* work, 319 int* iwork, 320 int* info); 321 322 void F77_FUNC(dlasd2, DLASD2)(int* nl, 323 int* nr, 324 int* sqre, 325 int* k, 326 double* d, 327 double* z, 328 double* alpha, 329 double* beta, 330 double* u, 331 int* ldu, 332 double* vt, 333 int* ldvt, 334 double* dsigma, 335 double* u2, 336 int* ldu2, 337 double* vt2, 338 int* ldvt2, 339 int* idxp, 340 int* idx, 341 int* idxc, 342 int* idxq, 343 int* coltyp, 344 int* info); 345 346 void F77_FUNC(dlasdt, DLASDT)(int* n, int* lvl, int* nd, int* inode, int* ndiml, int* ndimr, int* msub); 347 348 void F77_FUNC(dlasrt, DLASRT)(const char* id, int* n, double* d, int* info); 349 350 void F77_FUNC(dlasrt2, DLASRT2)(const char* id, int* n, double* d, int* key, int* info); 351 352 void F77_FUNC(ilasrt2, ILASRT2)(const char* id, int* n, int* d, int* key, int* info); 353 354 void F77_FUNC(dorgqr, 355 DORGQR)(int* m, int* n, int* k, double* a, int* lda, double* tau, double* work, int* lwork, int* info); 356 357 void F77_FUNC(dstebz, DSTEBZ)(const char* range, 358 const char* order, 359 int* n, 360 double* vl, 361 double* vu, 362 int* il, 363 int* iu, 364 double* abstol, 365 double* d, 366 double* e, 367 int* m, 368 int* nsplit, 369 double* w, 370 int* iblock, 371 int* isplit, 372 double* work, 373 int* iwork, 374 int* info); 375 376 void F77_FUNC(dsteqr, DSTEQR)(const char* compz, 377 int* n, 378 double* d__, 379 double* e, 380 double* z__, 381 int* ldz, 382 double* work, 383 int* info); 384 385 void F77_FUNC(dgebrd, DGEBRD)(int* m, 386 int* n, 387 double* a, 388 int* lda, 389 double* d, 390 double* e, 391 double* tauq, 392 double* taup, 393 double* work, 394 int* lwork, 395 int* info); 396 397 void F77_FUNC(dlacpy, 398 DLACPY)(const char* uplo, int* m, int* n, double* a, int* lda, double* b, int* ldb); 399 400 double F77_FUNC(dlapy2, DLAPY2)(double* x, double* y); 401 402 403 void F77_FUNC(dlarrfx, DLARRFX)(int* n, 404 double* d, 405 double* l, 406 double* ld, 407 double* lld, 408 int* ifirst, 409 int* ilast, 410 double* w, 411 double* sigma, 412 double* dplus, 413 double* lplus, 414 double* work, 415 int* info); 416 417 void F77_FUNC(dlasd3, DLASD3)(int* nl, 418 int* nr, 419 int* sqre, 420 int* k, 421 double* d, 422 double* q, 423 int* ldq, 424 double* dsigma, 425 double* u, 426 int* ldu, 427 double* u2, 428 int* ldu2, 429 double* vt, 430 int* ldvt, 431 double* vt2, 432 int* ldvt2, 433 int* idxc, 434 int* ctot, 435 double* z, 436 int* info); 437 438 void F77_FUNC(dlaset, 439 DLASET)(const char* uplo, int* m, int* n, double* alpha, double* beta, double* a, int* lda); 440 441 void F77_FUNC(dlassq, DLASSQ)(int* n, double* x, int* incx, double* scale, double* sumsq); 442 443 void F77_FUNC(dorm2l, DORM2L)(const char* side, 444 const char* trans, 445 int* m, 446 int* n, 447 int* k, 448 double* a, 449 int* lda, 450 double* tau, 451 double* c, 452 int* ldc, 453 double* work, 454 int* info); 455 456 void F77_FUNC(dstegr, DSTEGR)(const char* jobz, 457 const char* range, 458 int* n, 459 double* d, 460 double* e, 461 double* vl, 462 double* vu, 463 int* il, 464 int* iu, 465 double* abstol, 466 int* m, 467 double* w, 468 double* z, 469 int* ldz, 470 int* isuppz, 471 double* work, 472 int* lwork, 473 int* iwork, 474 int* liwork, 475 int* info); 476 477 void F77_FUNC(ssteqr, 478 SSTEQR)(const char* compz, int* n, float* d__, float* e, float* z__, int* ldz, float* work, int* info); 479 480 void F77_FUNC(dgelq2, DGELQ2)(int* m, int* n, double* a, int* lda, double* tau, double* work, int* info); 481 482 void F77_FUNC(dlae2, DLAE2)(double* a, double* b, double* c, double* rt1, double* rt2); 483 484 void F77_FUNC(dlaev2, 485 DLAEV2)(double* a, double* b, double* c, double* rt1, double* rt2, double* cs1, double* cs2); 486 487 void F77_FUNC(dlar1vx, DLAR1VX)(int* n, 488 int* b1, 489 int* bn, 490 double* sigma, 491 double* d, 492 double* l, 493 double* ld, 494 double* lld, 495 double* eval, 496 double* gersch, 497 double* z, 498 double* ztz, 499 double* mingma, 500 int* r, 501 int* isuppz, 502 double* work); 503 504 void F77_FUNC(dlarrvx, DLARRVX)(int* n, 505 double* d, 506 double* l, 507 int* isplit, 508 int* m, 509 double* w, 510 int* iblock, 511 int* indexw, 512 double* gersch, 513 double* tol, 514 double* z, 515 int* ldz, 516 int* isuppz, 517 double* work, 518 int* iwork, 519 int* info); 520 521 void F77_FUNC(dlasd4, DLASD4)(int* n, 522 int* i, 523 double* d, 524 double* z, 525 double* delta, 526 double* rho, 527 double* sigma, 528 double* work, 529 int* info); 530 531 void F77_FUNC(dlasq1, DLASQ1)(int* n, double* d, double* e, double* work, int* info); 532 533 534 void F77_FUNC(dlasv2, DLASV2)(double* f, 535 double* g, 536 double* h, 537 double* ssmin, 538 double* ssmax, 539 double* snr, 540 double* csr, 541 double* snl, 542 double* csl); 543 544 void F77_FUNC(dorm2r, DORM2R)(const char* side, 545 const char* trans, 546 int* m, 547 int* n, 548 int* k, 549 double* a, 550 int* lda, 551 double* tau, 552 double* c, 553 int* ldc, 554 double* work, 555 int* info); 556 557 void F77_FUNC(dstein, DSTEIN)(int* n, 558 double* d, 559 double* e, 560 int* m, 561 double* w, 562 int* iblock, 563 int* isplit, 564 double* z, 565 int* ldz, 566 double* work, 567 int* iwork, 568 int* ifail, 569 int* info); 570 571 void F77_FUNC(dgelqf, 572 DGELQF)(int* m, int* n, double* a, int* lda, double* tau, double* work, int* lwork, int* info); 573 574 void F77_FUNC(dlaebz, DLAEBZ)(int* ijob, 575 int* nitmax, 576 int* n, 577 int* mmax, 578 int* minp, 579 int* nbmin, 580 double* abstol, 581 double* reltol, 582 double* pivmin, 583 double* d, 584 double* e, 585 double* e2, 586 int* nval, 587 double* ab, 588 double* c, 589 int* mout, 590 int* nab, 591 double* work, 592 int* iwork, 593 int* info); 594 595 void F77_FUNC(dlarf, DLARF)(const char* side, 596 int* m, 597 int* n, 598 double* v, 599 int* incv, 600 double* tau, 601 double* c, 602 int* ldc, 603 double* work); 604 605 void F77_FUNC(dlartg, DLARTG)(double* f, double* g, double* cs, double* sn, double* r); 606 607 void F77_FUNC(dlasd5, 608 DLASD5)(int* i, double* d, double* z, double* delta, double* rho, double* dsigma, double* work); 609 610 void F77_FUNC(dlasq2, DLASQ2)(int* n, double* z, int* info); 611 612 void F77_FUNC(dlasq3, DLASQ3)(int* i0, 613 int* n0, 614 double* z, 615 int* pp, 616 double* dmin, 617 double* sigma, 618 double* desig, 619 double* qmax, 620 int* nfail, 621 int* iter, 622 int* ndiv, 623 int* ieee); 624 625 void F77_FUNC(dlaswp, DLASWP)(int* n, double* a, int* lda, int* k1, int* k2, int* ipiv, int* incx); 626 627 void F77_FUNC(dormbr, DORMBR)(const char* vect, 628 const char* side, 629 const char* trans, 630 int* m, 631 int* n, 632 int* k, 633 double* a, 634 int* lda, 635 double* tau, 636 double* c, 637 int* ldc, 638 double* work, 639 int* lwork, 640 int* info); 641 642 void F77_FUNC(dsterf, DSTERF)(int* n, double* d, double* e, int* info); 643 644 void F77_FUNC(dgeqr2, DGEQR2)(int* m, int* n, double* a, int* lda, double* tau, double* work, int* info); 645 646 void F77_FUNC(dlaed6, DLAED6)(int* kniter, 647 int* orgati, 648 double* rho, 649 double* d, 650 double* z, 651 double* finit, 652 double* tau, 653 int* info); 654 655 void F77_FUNC(dlarfb, DLARFB)(const char* side, 656 const char* trans, 657 const char* direct, 658 const char* storev, 659 int* m, 660 int* n, 661 int* k, 662 double* v, 663 int* ldv, 664 double* t, 665 int* ldt, 666 double* c, 667 int* ldc, 668 double* work, 669 int* ldwork); 670 671 void F77_FUNC(dlaruv, DLARUV)(int* iseed, int* n, double* x); 672 673 void F77_FUNC(dlasd6, DLASD6)(int* icompq, 674 int* nl, 675 int* nr, 676 int* sqre, 677 double* d, 678 double* vf, 679 double* vl, 680 double* alpha, 681 double* beta, 682 int* idxq, 683 int* perm, 684 int* givptr, 685 int* givcol, 686 int* ldgcol, 687 double* givnum, 688 int* ldgnum, 689 double* poles, 690 double* difl, 691 double* difr, 692 double* z, 693 int* k, 694 double* c, 695 double* s, 696 double* work, 697 int* iwork, 698 int* info); 699 700 void F77_FUNC(dlatrd, DLATRD)(const char* uplo, 701 int* n, 702 int* nb, 703 double* a, 704 int* lda, 705 double* e, 706 double* tau, 707 double* w, 708 int* ldw); 709 710 void F77_FUNC(dorml2, DORML2)(const char* side, 711 const char* trans, 712 int* m, 713 int* n, 714 int* k, 715 double* a, 716 int* lda, 717 double* tau, 718 double* c, 719 int* ldc, 720 double* work, 721 int* info); 722 723 void F77_FUNC(dstevr, DSTEVR)(const char* jobz, 724 const char* range, 725 int* n, 726 double* d, 727 double* e, 728 double* vl, 729 double* vu, 730 int* il, 731 int* iu, 732 double* abstol, 733 int* m, 734 double* w, 735 double* z, 736 int* ldz, 737 int* isuppz, 738 double* work, 739 int* lwork, 740 int* iwork, 741 int* liwork, 742 int* info); 743 744 void F77_FUNC(dsytrd, DSYTRD)(const char* uplo, 745 int* n, 746 double* a, 747 int* lda, 748 double* d, 749 double* e, 750 double* tau, 751 double* work, 752 int* lwork, 753 int* info); 754 755 void F77_FUNC(dsyevr, DSYEVR)(const char* jobz, 756 const char* range, 757 const char* uplo, 758 int* n, 759 double* a, 760 int* lda, 761 double* vl, 762 double* vu, 763 int* il, 764 int* iu, 765 double* abstol, 766 int* m, 767 double* w, 768 double* z__, 769 int* ldz, 770 int* isuppz, 771 double* work, 772 int* lwork, 773 int* iwork, 774 int* liwork, 775 int* info); 776 777 void F77_FUNC(dormql, DORMQL)(const char* side, 778 const char* trans, 779 int* m, 780 int* n, 781 int* k, 782 double* a, 783 int* lda, 784 double* tau, 785 double* c, 786 int* ldc, 787 double* work, 788 int* lwork, 789 int* info); 790 791 void F77_FUNC(dormqr, DORMQR)(const char* side, 792 const char* trans, 793 int* m, 794 int* n, 795 int* k, 796 double* a, 797 int* lda, 798 double* tau, 799 double* c, 800 int* ldc, 801 double* work, 802 int* lwork, 803 int* info); 804 805 void F77_FUNC(dorgbr, DORGBR)(const char* vect, 806 int* m, 807 int* n, 808 int* k, 809 double* a, 810 int* lda, 811 double* tau, 812 double* work, 813 int* lwork, 814 int* info); 815 816 void F77_FUNC(dlasq5, DLASQ5)(int* i0, 817 int* n0, 818 double* z, 819 int* pp, 820 double* tau, 821 double* dmin, 822 double* dmin1, 823 double* dmin2, 824 double* dn, 825 double* dnm1, 826 double* dnm2, 827 int* ieee); 828 829 void F77_FUNC(dlasd8, DLASD8)(int* icompq, 830 int* k, 831 double* d, 832 double* z, 833 double* vf, 834 double* vl, 835 double* difl, 836 double* difr, 837 int* lddifr, 838 double* dsigma, 839 double* work, 840 int* info); 841 842 void F77_FUNC(dlascl, DLASCL)(const char* type, 843 int* kl, 844 int* ku, 845 double* cfrom, 846 double* cto, 847 int* m, 848 int* n, 849 double* a, 850 int* lda, 851 int* info); 852 853 void F77_FUNC(dlarft, DLARFT)(const char* direct, 854 const char* storev, 855 int* n, 856 int* k, 857 double* v, 858 int* ldv, 859 double* tau, 860 double* t, 861 int* ldt); 862 863 void F77_FUNC(dlagts, DLAGTS)(int* job, 864 int* n, 865 double* a, 866 double* b, 867 double* c, 868 double* d, 869 int* in, 870 double* y, 871 double* tol, 872 int* info); 873 874 void F77_FUNC(dgesdd, DGESDD)(const char* jobz, 875 int* m, 876 int* n, 877 double* a, 878 int* lda, 879 double* s, 880 double* u, 881 int* ldu, 882 double* vt, 883 int* ldvt, 884 double* work, 885 int* lwork, 886 int* iwork, 887 int* info); 888 889 void F77_FUNC(dsytd2, 890 DSYTD2)(const char* uplo, int* n, double* a, int* lda, double* d, double* e, double* tau, int* info); 891 892 void F77_FUNC(dormlq, DORMLQ)(const char* side, 893 const char* trans, 894 int* m, 895 int* n, 896 int* k, 897 double* a, 898 int* lda, 899 double* tau, 900 double* c, 901 int* ldc, 902 double* work, 903 int* lwork, 904 int* info); 905 906 void F77_FUNC(dorg2r, 907 DORG2R)(int* m, int* n, int* k, double* a, int* lda, double* tau, double* work, int* info); 908 909 void F77_FUNC(dlasq4, DLASQ4)(int* i0, 910 int* n0, 911 double* z, 912 int* pp, 913 int* n0in, 914 double* dmin, 915 double* dmin1, 916 double* dmin2, 917 double* dn, 918 double* dn1, 919 double* dn2, 920 double* tau, 921 int* ttype); 922 923 void F77_FUNC(dlasd7, DLASD7)(int* icompq, 924 int* nl, 925 int* nr, 926 int* sqre, 927 int* k, 928 double* d, 929 double* z, 930 double* zw, 931 double* vf, 932 double* vfw, 933 double* vl, 934 double* vlw, 935 double* alpha, 936 double* beta, 937 double* dsigma, 938 int* idx, 939 int* idxp, 940 int* idxq, 941 int* perm, 942 int* givptr, 943 int* givcol, 944 int* ldgcol, 945 double* givnum, 946 int* ldgnum, 947 double* c, 948 double* s, 949 int* info); 950 951 void F77_FUNC(dlas2, DLAS2)(double* f, double* g, double* h, double* ssmin, double* ssmax); 952 953 void F77_FUNC(dlarfg, DLARFG)(int* n, double* alpha, double* x, int* incx, double* tau); 954 955 void F77_FUNC(dlagtf, DLAGTF)(int* n, 956 double* a, 957 double* lambda, 958 double* b, 959 double* c, 960 double* tol, 961 double* d, 962 int* in, 963 int* info); 964 965 void F77_FUNC(dgeqrf, 966 DGEQRF)(int* m, int* n, double* a, int* lda, double* tau, double* work, int* lwork, int* info); 967 968 969 /* Single precision */ 970 971 void F77_FUNC(sbdsdc, SBDSDC)(const char* uplo, 972 const char* compq, 973 int* n, 974 float* d, 975 float* e, 976 float* u, 977 int* ldu, 978 float* vt, 979 int* ldvt, 980 float* q, 981 int* iq, 982 float* work, 983 int* iwork, 984 int* info); 985 986 void F77_FUNC(sgetf2, SGETF2)(int* m, int* n, float* a, int* lda, int* ipiv, int* info); 987 988 void F77_FUNC(slamrg, SLAMRG)(int* n1, int* n2, float* a, int* dtrd1, int* dtrd2, int* index); 989 990 void F77_FUNC(slarnv, SLARNV)(int* idist, int* iseed, int* n, float* x); 991 992 void F77_FUNC(slasd0, SLASD0)(int* n, 993 int* sqre, 994 float* d, 995 float* e, 996 float* u, 997 int* ldu, 998 float* vt, 999 int* ldvt, 1000 int* smlsiz, 1001 int* iwork, 1002 float* work, 1003 int* info); 1004 1005 void F77_FUNC(slasda, SLASDA)(int* icompq, 1006 int* smlsiz, 1007 int* n, 1008 int* sqre, 1009 float* d, 1010 float* e, 1011 float* u, 1012 int* ldu, 1013 float* vt, 1014 int* k, 1015 float* difl, 1016 float* difr, 1017 float* z, 1018 float* poles, 1019 int* givptr, 1020 int* givcol, 1021 int* ldgcol, 1022 int* perm, 1023 float* givnum, 1024 float* c, 1025 float* s, 1026 float* work, 1027 int* iwork, 1028 int* info); 1029 1030 void F77_FUNC(slasq6, SLASQ6)(int* i0, 1031 int* n0, 1032 float* z, 1033 int* pp, 1034 float* dmin, 1035 float* dmin1, 1036 float* dmin2, 1037 float* dn, 1038 float* dnm1, 1039 float* dnm2); 1040 1041 void F77_FUNC(sorgl2, 1042 SORGL2)(int* m, int* n, int* k, float* a, int* lda, float* tau, float* work, int* info); 1043 1044 void F77_FUNC(sbdsqr, SBDSQR)(const char* uplo, 1045 int* n, 1046 int* ncvt, 1047 int* nru, 1048 int* ncc, 1049 float* d, 1050 float* e, 1051 float* vt, 1052 int* ldvt, 1053 float* u, 1054 int* ldu, 1055 float* c, 1056 int* ldc, 1057 float* work, 1058 int* info); 1059 1060 void F77_FUNC(sgetrf, SGETRF)(int* m, int* n, float* a, int* lda, int* ipiv, int* info); 1061 1062 void F77_FUNC(sgetri, SGETRI)(int* n, float* a, int* lda, int* ipiv, float* work, int* lwork, int* info); 1063 1064 void F77_FUNC(sgetrs, SGETRS)(const char* trans, 1065 int* n, 1066 int* nrhs, 1067 float* a, 1068 int* lda, 1069 int* ipiv, 1070 float* b, 1071 int* ldb, 1072 int* info); 1073 1074 void F77_FUNC(strtri, STRTRI)(const char* uplo, const char* diag, int* n, float* a, int* lda, int* info); 1075 1076 void F77_FUNC(strti2, STRTI2)(const char* uplo, const char* diag, int* n, float* a, int* lda, int* info); 1077 1078 float F77_FUNC(slange, SLANGE)(const char* norm, int* m, int* n, float* a, int* lda, float* work); 1079 1080 void F77_FUNC(slarrbx, SLARRBX)(int* n, 1081 float* d, 1082 float* l, 1083 float* ld, 1084 float* lld, 1085 int* ifirst, 1086 int* ilast, 1087 float* rtol1, 1088 float* rtol2, 1089 int* offset, 1090 float* w, 1091 float* wgap, 1092 float* werr, 1093 float* work, 1094 int* iwork, 1095 int* info); 1096 1097 void F77_FUNC(slasd1, SLASD1)(int* nl, 1098 int* nr, 1099 int* sqre, 1100 float* d, 1101 float* alpha, 1102 float* beta, 1103 float* u, 1104 int* ldu, 1105 float* vt, 1106 int* ldvt, 1107 int* idxq, 1108 int* iwork, 1109 float* work, 1110 int* info); 1111 1112 void F77_FUNC(slasdq, SLASDQ)(const char* uplo, 1113 int* sqre, 1114 int* n, 1115 int* ncvt, 1116 int* nru, 1117 int* ncc, 1118 float* d, 1119 float* e, 1120 float* vt, 1121 int* ldvt, 1122 float* u, 1123 int* ldu, 1124 float* c, 1125 int* ldc, 1126 float* work, 1127 int* info); 1128 1129 void F77_FUNC(slasr, SLASR)(const char* side, 1130 const char* pivot, 1131 const char* direct, 1132 int* m, 1133 int* n, 1134 float* c, 1135 float* s, 1136 float* a, 1137 int* lda); 1138 1139 void F77_FUNC(sorglq, 1140 SORGLQ)(int* m, int* n, int* k, float* a, int* lda, float* tau, float* work, int* lwork, int* info); 1141 1142 void F77_FUNC(sormtr, SORMTR)(const char* side, 1143 const char* uplo, 1144 const char* trans, 1145 int* m, 1146 int* n, 1147 float* a, 1148 int* lda, 1149 float* tau, 1150 float* c, 1151 int* ldc, 1152 float* work, 1153 int* lwork, 1154 int* info); 1155 1156 void F77_FUNC(sgebd2, SGEBD2)(int* m, 1157 int* n, 1158 float* a, 1159 int* lda, 1160 float* d, 1161 float* e, 1162 float* tauq, 1163 float* taup, 1164 float* work, 1165 int* info); 1166 1167 void F77_FUNC(slabrd, SLABRD)(int* m, 1168 int* n, 1169 int* nb, 1170 float* a, 1171 int* lda, 1172 float* d, 1173 float* e, 1174 float* tauq, 1175 float* taup, 1176 float* x, 1177 int* ldx, 1178 float* y, 1179 int* ldy); 1180 1181 float F77_FUNC(slanst, SLANST)(const char* norm, int* n, float* d, float* e); 1182 1183 float F77_FUNC(slansy, 1184 SLANSY)(const char* norm, const char* uplo, int* n, float* a, int* lda, float* work); 1185 1186 void F77_FUNC(slarrex, SLARREX)(const char* range, 1187 int* n, 1188 float* vl, 1189 float* vu, 1190 int* il, 1191 int* iu, 1192 float* d, 1193 float* e, 1194 float* tol, 1195 int* nsplit, 1196 int* isplit, 1197 int* m, 1198 float* w, 1199 int* iblock, 1200 int* indexw, 1201 float* gersch, 1202 float* work, 1203 int* iwork, 1204 int* info); 1205 1206 void F77_FUNC(slasd2, SLASD2)(int* nl, 1207 int* nr, 1208 int* sqre, 1209 int* k, 1210 float* d, 1211 float* z, 1212 float* alpha, 1213 float* beta, 1214 float* u, 1215 int* ldu, 1216 float* vt, 1217 int* ldvt, 1218 float* dsigma, 1219 float* u2, 1220 int* ldu2, 1221 float* vt2, 1222 int* ldvt2, 1223 int* idxp, 1224 int* idx, 1225 int* idxc, 1226 int* idxq, 1227 int* coltyp, 1228 int* info); 1229 1230 void F77_FUNC(slasdt, SLASDT)(int* n, int* lvl, int* nd, int* inode, int* ndiml, int* ndimr, int* msub); 1231 1232 void F77_FUNC(slasrt, SLASRT)(const char* id, int* n, float* d, int* info); 1233 1234 void F77_FUNC(slasrt2, SLASRT2)(const char* id, int* n, float* d, int* key, int* info); 1235 1236 void F77_FUNC(sorgqr, 1237 SORGQR)(int* m, int* n, int* k, float* a, int* lda, float* tau, float* work, int* lwork, int* info); 1238 1239 void F77_FUNC(sstebz, SSTEBZ)(const char* range, 1240 const char* order, 1241 int* n, 1242 float* vl, 1243 float* vu, 1244 int* il, 1245 int* iu, 1246 float* abstol, 1247 float* d, 1248 float* e, 1249 int* m, 1250 int* nsplit, 1251 float* w, 1252 int* iblock, 1253 int* isplit, 1254 float* work, 1255 int* iwork, 1256 int* info); 1257 1258 void F77_FUNC(sgebrd, SGEBRD)(int* m, 1259 int* n, 1260 float* a, 1261 int* lda, 1262 float* d, 1263 float* e, 1264 float* tauq, 1265 float* taup, 1266 float* work, 1267 int* lwork, 1268 int* info); 1269 1270 void F77_FUNC(slacpy, SLACPY)(const char* uplo, int* m, int* n, float* a, int* lda, float* b, int* ldb); 1271 1272 float F77_FUNC(slapy2, SLAPY2)(float* x, float* y); 1273 1274 void F77_FUNC(slarrfx, SLARRFX)(int* n, 1275 float* d, 1276 float* l, 1277 float* ld, 1278 float* lld, 1279 int* ifirst, 1280 int* ilast, 1281 float* w, 1282 float* sigma, 1283 float* dplus, 1284 float* lplus, 1285 float* work, 1286 int* info); 1287 1288 void F77_FUNC(slasd3, SLASD3)(int* nl, 1289 int* nr, 1290 int* sqre, 1291 int* k, 1292 float* d, 1293 float* q, 1294 int* ldq, 1295 float* dsigma, 1296 float* u, 1297 int* ldu, 1298 float* u2, 1299 int* ldu2, 1300 float* vt, 1301 int* ldvt, 1302 float* vt2, 1303 int* ldvt2, 1304 int* idxc, 1305 int* ctot, 1306 float* z, 1307 int* info); 1308 1309 void F77_FUNC(slaset, 1310 SLASET)(const char* uplo, int* m, int* n, float* alpha, float* beta, float* a, int* lda); 1311 1312 void F77_FUNC(slassq, SLASSQ)(int* n, float* x, int* incx, float* scale, float* sumsq); 1313 1314 void F77_FUNC(sorm2l, SORM2L)(const char* side, 1315 const char* trans, 1316 int* m, 1317 int* n, 1318 int* k, 1319 float* a, 1320 int* lda, 1321 float* tau, 1322 float* c, 1323 int* ldc, 1324 float* work, 1325 int* info); 1326 1327 void F77_FUNC(sstegr, SSTEGR)(const char* jobz, 1328 const char* range, 1329 int* n, 1330 float* d, 1331 float* e, 1332 float* vl, 1333 float* vu, 1334 int* il, 1335 int* iu, 1336 float* abstol, 1337 int* m, 1338 float* w, 1339 float* z, 1340 int* ldz, 1341 int* isuppz, 1342 float* work, 1343 int* lwork, 1344 int* iwork, 1345 int* liwork, 1346 int* info); 1347 1348 void F77_FUNC(sgelq2, SGELQ2)(int* m, int* n, float* a, int* lda, float* tau, float* work, int* info); 1349 1350 void F77_FUNC(slae2, SLAE2)(float* a, float* b, float* c, float* rt1, float* rt2); 1351 1352 void F77_FUNC(slaev2, 1353 SLAEV2)(float* a, float* b, float* c, float* rt1, float* rt2, float* cs1, float* cs2); 1354 1355 void F77_FUNC(slar1vx, SLAR1VX)(int* n, 1356 int* b1, 1357 int* bn, 1358 float* sigma, 1359 float* d, 1360 float* l, 1361 float* ld, 1362 float* lld, 1363 float* eval, 1364 float* gersch, 1365 float* z, 1366 float* ztz, 1367 float* mingma, 1368 int* r, 1369 int* isuppz, 1370 float* work); 1371 1372 void F77_FUNC(slarrvx, SLARRVX)(int* n, 1373 float* d, 1374 float* l, 1375 int* isplit, 1376 int* m, 1377 float* w, 1378 int* iblock, 1379 int* indexw, 1380 float* gersch, 1381 float* tol, 1382 float* z, 1383 int* ldz, 1384 int* isuppz, 1385 float* work, 1386 int* iwork, 1387 int* info); 1388 1389 void F77_FUNC(slasd4, SLASD4)(int* n, 1390 int* i, 1391 float* d, 1392 float* z, 1393 float* delta, 1394 float* rho, 1395 float* sigma, 1396 float* work, 1397 int* info); 1398 1399 void F77_FUNC(slasq1, SLASQ1)(int* n, float* d, float* e, float* work, int* info); 1400 1401 1402 void F77_FUNC(slasv2, SLASV2)(float* f, 1403 float* g, 1404 float* h, 1405 float* ssmin, 1406 float* ssmax, 1407 float* snr, 1408 float* csr, 1409 float* snl, 1410 float* csl); 1411 1412 void F77_FUNC(sorm2r, SORM2R)(const char* side, 1413 const char* trans, 1414 int* m, 1415 int* n, 1416 int* k, 1417 float* a, 1418 int* lda, 1419 float* tau, 1420 float* c, 1421 int* ldc, 1422 float* work, 1423 int* info); 1424 1425 void F77_FUNC(sstein, SSTEIN)(int* n, 1426 float* d, 1427 float* e, 1428 int* m, 1429 float* w, 1430 int* iblock, 1431 int* isplit, 1432 float* z, 1433 int* ldz, 1434 float* work, 1435 int* iwork, 1436 int* ifail, 1437 int* info); 1438 1439 void F77_FUNC(sgelqf, 1440 SGELQF)(int* m, int* n, float* a, int* lda, float* tau, float* work, int* lwork, int* info); 1441 1442 void F77_FUNC(slaebz, SLAEBZ)(int* ijob, 1443 int* nitmax, 1444 int* n, 1445 int* mmax, 1446 int* minp, 1447 int* nbmin, 1448 float* abstol, 1449 float* reltol, 1450 float* pivmin, 1451 float* d, 1452 float* e, 1453 float* e2, 1454 int* nval, 1455 float* ab, 1456 float* c, 1457 int* mout, 1458 int* nab, 1459 float* work, 1460 int* iwork, 1461 int* info); 1462 1463 void F77_FUNC(slarf, 1464 SLARF)(const char* side, int* m, int* n, float* v, int* incv, float* tau, float* c, int* ldc, float* work); 1465 1466 void F77_FUNC(slartg, SLARTG)(float* f, float* g, float* cs, float* sn, float* r); 1467 1468 void F77_FUNC(slasd5, 1469 SLASD5)(int* i, float* d, float* z, float* delta, float* rho, float* dsigma, float* work); 1470 1471 void F77_FUNC(slasq2, SLASQ2)(int* n, float* z, int* info); 1472 1473 void F77_FUNC(slasq3, SLASQ3)(int* i0, 1474 int* n0, 1475 float* z, 1476 int* pp, 1477 float* dmin, 1478 float* sigma, 1479 float* desig, 1480 float* qmax, 1481 int* nfail, 1482 int* iter, 1483 int* ndiv, 1484 int* ieee); 1485 1486 void F77_FUNC(slaswp, SLASWP)(int* n, float* a, int* lda, int* k1, int* k2, int* ipiv, int* incx); 1487 1488 void F77_FUNC(sormbr, SORMBR)(const char* vect, 1489 const char* side, 1490 const char* trans, 1491 int* m, 1492 int* n, 1493 int* k, 1494 float* a, 1495 int* lda, 1496 float* tau, 1497 float* c, 1498 int* ldc, 1499 float* work, 1500 int* lwork, 1501 int* info); 1502 1503 void F77_FUNC(ssterf, SSTERF)(int* n, float* d, float* e, int* info); 1504 1505 void F77_FUNC(sgeqr2, SGEQR2)(int* m, int* n, float* a, int* lda, float* tau, float* work, int* info); 1506 1507 void F77_FUNC(slaed6, 1508 SLAED6)(int* kniter, int* orgati, float* rho, float* d, float* z, float* finit, float* tau, int* info); 1509 1510 void F77_FUNC(slarfb, SLARFB)(const char* side, 1511 const char* trans, 1512 const char* direct, 1513 const char* storev, 1514 int* m, 1515 int* n, 1516 int* k, 1517 float* v, 1518 int* ldv, 1519 float* t, 1520 int* ldt, 1521 float* c, 1522 int* ldc, 1523 float* work, 1524 int* ldwork); 1525 1526 void F77_FUNC(slaruv, SLARUV)(int* iseed, int* n, float* x); 1527 1528 void F77_FUNC(slasd6, SLASD6)(int* icompq, 1529 int* nl, 1530 int* nr, 1531 int* sqre, 1532 float* d, 1533 float* vf, 1534 float* vl, 1535 float* alpha, 1536 float* beta, 1537 int* idxq, 1538 int* perm, 1539 int* givptr, 1540 int* givcol, 1541 int* ldgcol, 1542 float* givnum, 1543 int* ldgnum, 1544 float* poles, 1545 float* difl, 1546 float* difr, 1547 float* z, 1548 int* k, 1549 float* c, 1550 float* s, 1551 float* work, 1552 int* iwork, 1553 int* info); 1554 1555 void F77_FUNC(slatrd, 1556 SLATRD)(const char* uplo, int* n, int* nb, float* a, int* lda, float* e, float* tau, float* w, int* ldw); 1557 1558 void F77_FUNC(sorml2, SORML2)(const char* side, 1559 const char* trans, 1560 int* m, 1561 int* n, 1562 int* k, 1563 float* a, 1564 int* lda, 1565 float* tau, 1566 float* c, 1567 int* ldc, 1568 float* work, 1569 int* info); 1570 1571 void F77_FUNC(sstevr, SSTEVR)(const char* jobz, 1572 const char* range, 1573 int* n, 1574 float* d, 1575 float* e, 1576 float* vl, 1577 float* vu, 1578 int* il, 1579 int* iu, 1580 float* abstol, 1581 int* m, 1582 float* w, 1583 float* z, 1584 int* ldz, 1585 int* isuppz, 1586 float* work, 1587 int* lwork, 1588 int* iwork, 1589 int* liwork, 1590 int* info); 1591 1592 void F77_FUNC(ssytrd, SSYTRD)(const char* uplo, 1593 int* n, 1594 float* a, 1595 int* lda, 1596 float* d, 1597 float* e, 1598 float* tau, 1599 float* work, 1600 int* lwork, 1601 int* info); 1602 1603 void F77_FUNC(ssyevr, SSYEVR)(const char* jobz, 1604 const char* range, 1605 const char* uplo, 1606 int* n, 1607 float* a, 1608 int* lda, 1609 float* vl, 1610 float* vu, 1611 int* il, 1612 int* iu, 1613 float* abstol, 1614 int* m, 1615 float* w, 1616 float* z__, 1617 int* ldz, 1618 int* isuppz, 1619 float* work, 1620 int* lwork, 1621 int* iwork, 1622 int* liwork, 1623 int* info); 1624 1625 void F77_FUNC(sormql, SORMQL)(const char* side, 1626 const char* trans, 1627 int* m, 1628 int* n, 1629 int* k, 1630 float* a, 1631 int* lda, 1632 float* tau, 1633 float* c, 1634 int* ldc, 1635 float* work, 1636 int* lwork, 1637 int* info); 1638 1639 void F77_FUNC(sormqr, SORMQR)(const char* side, 1640 const char* trans, 1641 int* m, 1642 int* n, 1643 int* k, 1644 float* a, 1645 int* lda, 1646 float* tau, 1647 float* c, 1648 int* ldc, 1649 float* work, 1650 int* lwork, 1651 int* info); 1652 1653 void F77_FUNC(sorgbr, SORGBR)(const char* vect, 1654 int* m, 1655 int* n, 1656 int* k, 1657 float* a, 1658 int* lda, 1659 float* tau, 1660 float* work, 1661 int* lwork, 1662 int* info); 1663 1664 void F77_FUNC(slasq5, SLASQ5)(int* i0, 1665 int* n0, 1666 float* z, 1667 int* pp, 1668 float* tau, 1669 float* dmin, 1670 float* dmin1, 1671 float* dmin2, 1672 float* dn, 1673 float* dnm1, 1674 float* dnm2, 1675 int* ieee); 1676 1677 void F77_FUNC(slasd8, SLASD8)(int* icompq, 1678 int* k, 1679 float* d, 1680 float* z, 1681 float* vf, 1682 float* vl, 1683 float* difl, 1684 float* difr, 1685 int* lddifr, 1686 float* dsigma, 1687 float* work, 1688 int* info); 1689 1690 void F77_FUNC(slascl, SLASCL)(const char* type, 1691 int* kl, 1692 int* ku, 1693 float* cfrom, 1694 float* cto, 1695 int* m, 1696 int* n, 1697 float* a, 1698 int* lda, 1699 int* info); 1700 1701 void F77_FUNC(slarft, SLARFT)(const char* direct, 1702 const char* storev, 1703 int* n, 1704 int* k, 1705 float* v, 1706 int* ldv, 1707 float* tau, 1708 float* t, 1709 int* ldt); 1710 1711 void F77_FUNC(slagts, 1712 SLAGTS)(int* job, int* n, float* a, float* b, float* c, float* d, int* in, float* y, float* tol, int* info); 1713 1714 void F77_FUNC(sgesdd, SGESDD)(const char* jobz, 1715 int* m, 1716 int* n, 1717 float* a, 1718 int* lda, 1719 float* s, 1720 float* u, 1721 int* ldu, 1722 float* vt, 1723 int* ldvt, 1724 float* work, 1725 int* lwork, 1726 int* iwork, 1727 int* info); 1728 1729 void F77_FUNC(ssytd2, 1730 SSYTD2)(const char* uplo, int* n, float* a, int* lda, float* d, float* e, float* tau, int* info); 1731 1732 void F77_FUNC(sormlq, SORMLQ)(const char* side, 1733 const char* trans, 1734 int* m, 1735 int* n, 1736 int* k, 1737 float* a, 1738 int* lda, 1739 float* tau, 1740 float* c, 1741 int* ldc, 1742 float* work, 1743 int* lwork, 1744 int* info); 1745 1746 void F77_FUNC(sorg2r, 1747 SORG2R)(int* m, int* n, int* k, float* a, int* lda, float* tau, float* work, int* info); 1748 1749 void F77_FUNC(slasq4, SLASQ4)(int* i0, 1750 int* n0, 1751 float* z, 1752 int* pp, 1753 int* n0in, 1754 float* dmin, 1755 float* dmin1, 1756 float* dmin2, 1757 float* dn, 1758 float* dn1, 1759 float* dn2, 1760 float* tau, 1761 int* ttype); 1762 1763 void F77_FUNC(slasd7, SLASD7)(int* icompq, 1764 int* nl, 1765 int* nr, 1766 int* sqre, 1767 int* k, 1768 float* d, 1769 float* z, 1770 float* zw, 1771 float* vf, 1772 float* vfw, 1773 float* vl, 1774 float* vlw, 1775 float* alpha, 1776 float* beta, 1777 float* dsigma, 1778 int* idx, 1779 int* idxp, 1780 int* idxq, 1781 int* perm, 1782 int* givptr, 1783 int* givcol, 1784 int* ldgcol, 1785 float* givnum, 1786 int* ldgnum, 1787 float* c, 1788 float* s, 1789 int* info); 1790 1791 void F77_FUNC(slas2, SLAS2)(float* f, float* g, float* h, float* ssmin, float* ssmax); 1792 1793 void F77_FUNC(slarfg, SLARFG)(int* n, float* alpha, float* x, int* incx, float* tau); 1794 1795 void F77_FUNC(slagtf, 1796 SLAGTF)(int* n, float* a, float* lambda, float* b, float* c, float* tol, float* d, int* in, int* info); 1797 1798 void F77_FUNC(sgeqrf, 1799 SGEQRF)(int* m, int* n, float* a, int* lda, float* tau, float* work, int* lwork, int* info); 1800 1801 1802 #ifdef __cplusplus 1803 } 1804 #endif 1805 1806 /*! \endcond */ 1807 1808 #endif /* GMX_LAPACK_H */ 1809