amino  1.0-beta2
Lightweight Robot Utility Library
lapack_impl.h
1 /* -*- mode: C; c-basic-offset: 4 -*- */
2 /* ex: set shiftwidth=4 tabstop=4 expandtab: */
3 /*
4  * Copyright (c) 2008-2012, Georgia Tech Research Corporation
5  * All rights reserved.
6  *
7  * Author(s): Neil T. Dantam <ntd@gatech.edu>
8  * Georgia Tech Humanoid Robotics Lab
9  * Under Direction of Prof. Mike Stilman <mstilman@cc.gatech.edu>
10  *
11  *
12  * This file is provided under the following "BSD-style" License:
13  *
14  *
15  * Redistribution and use in source and binary forms, with or
16  * without modification, are permitted provided that the following
17  * conditions are met:
18  *
19  * * Redistributions of source code must retain the above copyright
20  * notice, this list of conditions and the following disclaimer.
21  *
22  * * Redistributions in binary form must reproduce the above
23  * copyright notice, this list of conditions and the following
24  * disclaimer in the documentation and/or other materials provided
25  * with the distribution.
26  *
27  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND
28  * CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES,
29  * INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF
30  * MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
31  * DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR
32  * CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
33  * SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
34  * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF
35  * USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED
36  * AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
37  * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
38  * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
39  * POSSIBILITY OF SUCH DAMAGE.
40  *
41  */
42 
43 #include "amino/def.h"
44 
45 /* FILE lapack_impl.h
46  * BRIEF: C prototypes to various fortran lapack routines.
47  *
48  * Since there is no official c binding to lapack as there as with
49  * BLAS, the only reasonable way to interface with lapack from C is to
50  * call the fortran methods directly.
51  *
52  * Authors:
53  * Neil T. Dantam
54  */
55 
71 AA_API void AA_LAPACK_NAME(getri)
72 ( const int *N, AA_TYPE *A, const int *LDA,
73  const int *IPIV, AA_TYPE *WORK, const int *LWORK, int *INFO );
74 
113 AA_API void AA_LAPACK_NAME(getrf)
114 ( const int *M, const int *N, AA_TYPE *A, const int *LDA,
115  int *IPIV, int *INFO );
116 
184 AA_API void AA_LAPACK_NAME(geqrf)
185 ( const int *M, const int *N, AA_TYPE *A, const int *LDA,
186  AA_TYPE *TAU, AA_TYPE *WORK, int *LWORK, int *INFO );
187 
242 AA_API void AA_LAPACK_NAME(orgqr)
243 ( const int *M, const int *N, const int *K,
244  AA_TYPE *A, const int *LDA, const AA_TYPE *TAU,
245  AA_TYPE *WORK, const int *LWORK, int *INFO );
246 
247 
248 AA_API void AA_LAPACK_NAME(posv)
249 ( const char UPLO[1], const int *N, const int *NRHS,
250  AA_TYPE *A, const int *LDA,
251  AA_TYPE *B, const int *LDB,
252  int *INFO );
253 
254 
255 AA_API void AA_LAPACK_NAME(sysv)
256 ( const char uplo[1], const int *n, const int *nrhs,
257  AA_TYPE *A, const int *lda,
258  int *ipiv,
259  AA_TYPE *B, const int *ldb,
260  AA_TYPE *work, int *lwork, int *info );
261 
262 
360 AA_API void AA_LAPACK_NAME(gesvd)
361 ( const char jobu[1], const char jobvt[1],
362  const int *m, const int *n,
363  AA_TYPE *A, const int *lda,
364  AA_TYPE *S, AA_TYPE *U,
365  const int *ldu, AA_TYPE *Vt, const int *ldvt,
366  AA_TYPE *work, const int *lwork, int *info );
367 
368 AA_API void AA_LAPACK_NAME(gesdd)
369 (const char *JOBZ, const int *M, const int *N,
370  AA_TYPE *A, const int *LDA,
371  AA_TYPE *S,
372  AA_TYPE *U, const int *LDU,
373  AA_TYPE *VT, const int * LDVT,
374  AA_TYPE *WORK, const int *LWORK,
375  int *IWORK, int *INFO);
376 
377 AA_API int AA_LAPACK_NAME(geev)
378 (const char *jobvl, const char *jobvr,
379  int *n, AA_TYPE *a, int *lda,
380  AA_TYPE *wr, AA_TYPE *wi,
381  AA_TYPE *vl, int *ldvl,
382  AA_TYPE *vr, int *ldvr,
383  AA_TYPE *work, int *lwork, int *info);
384 
385 
386 
496 AA_API void AA_LAPACK_NAME(gelsd)
497 ( const int *M, const int *N, const int *NRHS,
498  AA_TYPE *A, const int *LDA, AA_TYPE *B, const int *LDB,
499  AA_TYPE *S, const AA_TYPE *RCOND, int *RANK,
500  AA_TYPE *WORK, int *LWORK, int *IWORK, int *INFO );
501 
502 
559 AA_API void AA_LAPACK_NAME(gebal)
560 ( const char JOB[1], int *N, AA_TYPE *A, const int *LDA,
561  int *ILO, int *IHI, AA_TYPE *SCALE, int *INFO );
562 
563 
680 AA_API void AA_LAPACK_NAME(gees)
681 ( const char JOBVS[1], const char SORT[1],
682  int (*SELECT)(const AA_TYPE*,const AA_TYPE*),
683  const int *N, AA_TYPE *A, const int *LDA, int *SDIM,
684  AA_TYPE *WR, AA_TYPE *WI,
685  AA_TYPE *VS, const int *LDVS,
686  AA_TYPE *WORK, const int *LWORK, int *BWORK, int *INFO );
687 
786 AA_API void AA_LAPACK_NAME(gels)
787 ( const char TRANS[1], const int *M, const int *N, const int *NRHS,
788  AA_TYPE *A, const int *LDA, AA_TYPE *B, const int *LDB, AA_TYPE *WORK,
789  const int *LWORK, int *INFO );
790 
791 /*
792 * DPOTRF computes the Cholesky factorization of a real symmetric
793 * positive definite matrix A.
794 *
795 * The factorization has the form
796 * \f[ A = U^T * U \f] if UPLO = 'U'
797 * or
798 * \f[ A = L * L^T \f] if UPLO = 'L'
799 *
800 * where U is an upper triangular matrix and L is lower triangular.
801 *
802 * This is the block version of the algorithm, calling Level 3 BLAS.
803 *
804 *
805 * \param[in] UPLO
806 * = 'U': Upper triangle of A is stored;
807 * = 'L': Lower triangle of A is stored.
808 *
809 * \param[in] N
810 * The order of the matrix A. N >= 0.
811 *
812 * \param[inout] A
813 * On entry, the symmetric matrix A. If UPLO = 'U', the leading
814 * N-by-N upper triangular part of A contains the upper
815 * triangular part of the matrix A, and the strictly lower
816 * triangular part of A is not referenced. If UPLO = 'L', the
817 * leading N-by-N lower triangular part of A contains the lower
818 * triangular part of the matrix A, and the strictly upper
819 * triangular part of A is not referenced.
820 * On exit, if INFO = 0, the factor U or L from the Cholesky
821 * factorization A = U**T*U or A = L*L**T.
822 *
823 * \param[in] LDA
824 * The leading dimension of the array A. LDA >= max(1,N).
825 *
826 * \param[out] INFO
827 * - = 0: successful exit
828 * - < 0: if INFO = -i, the i-th argument had an illegal value
829 * - > 0: if INFO = i, the leading minor of order i is not
830 * positive definite, and the factorization could not be
831 * completed.
832 */
833 
834 AA_API void AA_LAPACK_NAME(potrf)
835 ( const char UPLO[1], const int *N,
836  AA_TYPE *A, const int *LDA,
837  int *INFO );
838 
839 AA_API void AA_LAPACK_NAME(potrs)
840 ( const char UPLO[1], const int *N, const int *nrhs,
841  AA_TYPE *A, const int *LDA,
842  AA_TYPE *B, const int *LDB,
843  int *INFO );
844 
877 AA_API void AA_LAPACK_NAME(lacpy)
878 ( const char UPLO[1], const int *M, const int *N,
879  const AA_TYPE *A, const int *LDA, AA_TYPE *B, const int *LDB );
880 
881 
885 AA_API AA_TYPE AA_LAPACK_NAME(lapy2)
886 ( const AA_TYPE *x, const AA_TYPE *y );
887 
891 AA_API AA_TYPE AA_LAPACK_NAME(lapy3)
892 ( const AA_TYPE *x, const AA_TYPE *y, const AA_TYPE *z );
893 
894 
935 AA_API void AA_LAPACK_NAME(laruv)
936 ( int ISEED[4], const int *N, AA_TYPE *X );
937 
965 AA_API void AA_LAPACK_NAME(larnv)
966 ( const int *IDIST, int ISEED[4],
967  const int *N, AA_TYPE *X );
968 
1024 AA_API void AA_LAPACK_NAME(lascl)
1025 ( const char TYPE[1], const int *KL, const int *KU,
1026  const AA_TYPE *CFROM, const AA_TYPE *CTO,
1027  const int *M, const int *N, AA_TYPE *A, const int *LDA,
1028  int *INFO );
1029 
1030 
1055 AA_API double AA_LAPACK_NAME(lange)
1056 ( const char TYPE[1],
1057  const int *M, const int *N, AA_TYPE *A, const int *LDA,
1058  AA_TYPE *work );
1059 
1110 AA_API void AA_LAPACK_NAME(laset)
1111 ( const char UPLO[1], const int *M, const int *N,
1112  const AA_TYPE *ALPHA,
1113  const AA_TYPE *BETA,
1114  AA_TYPE *A, const int *LDA );
1115 
1116 
1117 
1118 
1119 
1120 #include "amino/undef.h"
1121 
1122 #if AA_TYPE == double
1153 AA_API void dlag2s_ ( const int *M, const int *N,
1154  double *A, const int *LDA,
1155  float *SA, const int *LDSA,
1156  const int *INFO );
1157 
1158 
1159 #endif // AA_TYPE == double
1160 
1161 
1162 #if AA_TYPE == float
1198 AA_API void slag2d_ ( const int *M, const int *N,
1199  float *SA, const int *LDSA,
1200  double *A, const int *LDA,
1201  const int *INFO );
1202 
1203 #endif // AA_TYPE == float
#define AA_API
calling and name mangling convention for functions
Definition: amino.h:95