SLAMflex SE  0.1.0
SLAMflex provides detection and tracking of dominant planes for smartphone devices. This plane can then be used to show AR content relative to the plane orientation. The detection of plane is performed in the field of view of the smartphone camera. In subsequent frames it is tracked. The interface returns the plane position and orientation.
lapack.h
Go to the documentation of this file.
1 // -*- c++ -*-
2 
3 // Copyright (C) 2005,2009,2010 Tom Drummond (twd20@cam.ac.uk), E. Rosten
4 //
5 // This file is part of the TooN Library. This library is free
6 // software; you can redistribute it and/or modify it under the
7 // terms of the GNU General Public License as published by the
8 // Free Software Foundation; either version 2, or (at your option)
9 // any later version.
10 
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 // GNU General Public License for more details.
15 
16 // You should have received a copy of the GNU General Public License along
17 // with this library; see the file COPYING. If not, write to the Free
18 // Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307,
19 // USA.
20 
21 // As a special exception, you may use this file as part of a free software
22 // library without restriction. Specifically, if other files instantiate
23 // templates or use macros or inline functions from this file, or you compile
24 // this file and link it with other files to produce an executable, this
25 // file does not by itself cause the resulting executable to be covered by
26 // the GNU General Public License. This exception does not however
27 // invalidate any other reasons why the executable file might be covered by
28 // the GNU General Public License.
29 
30 #ifndef TOON_INCLUDE_LAPCK_H
31 #define TOON_INCLUDE_LAPCK_H
32 
33 // LAPACK and BLAS routines
34 namespace TooN {
35 
36  extern "C" {
37  // LU decomoposition of a general matrix
38  void dgetrf_(int* M, int *N, double* A, int* lda, int* IPIV, int* INFO);
39  void sgetrf_(int* M, int *N, float* A, int* lda, int* IPIV, int* INFO);
40 
41  // generate inverse of a matrix given its LU decomposition
42  void dgetri_(int* N, double* A, int* lda, int* IPIV, double* WORK, int* lwork, int* INFO);
43  void sgetri_(int* N, float* A, int* lda, int* IPIV, float* WORK, int* lwork, int* INFO);
44 
45  // inverse of a triangular matrix * a vector (BLAS level 2)
46  void dtrsm_(char* SIDE, char* UPLO, char* TRANSA, char* DIAG, int* M, int* N, double* alpha, double* A, int* lda, double* B, int* ldb);
47  void strsm_(char* SIDE, char* UPLO, char* TRANSA, char* DIAG, int* M, int* N, float* alpha, float* A, int* lda, float* B, int* ldb);
48 
49 
50  // SVD of a general matrix
51  void dgesvd_(const char* JOBU, const char* JOBVT, int* M, int *N, double* A, int* lda,
52  double* S, double *U, int* ldu, double* VT, int* ldvt,
53  double* WORK, int* lwork, int* INFO);
54 
55  void sgesvd_(const char* JOBU, const char* JOBVT, int* M, int *N, float* A, int* lda,
56  float* S, float *U, int* ldu, float* VT, int* ldvt,
57  float* WORK, int* lwork, int* INFO);
58 
59  // Eigen decomposition of a symmetric matrix
60  void dsyev_(const char* JOBZ, const char* UPLO, int* N, double* A, int* lda, double* W, double* WORK, int* LWORK, int* INFO);
61  void ssyev_(const char* JOBZ, const char* UPLO, int* N, float* A, int* lda, float* W, float* WORK, int* LWORK, int* INFO);
62 
63  // Eigen decomposition of a non-symmetric matrix
64  void dgeev_(const char* JOBVL, const char* JOBVR, int* N, double* A, int* lda, double* WR, double* WI, double* VL, int* LDVL, double* VR, int* LDVR , double* WORK, int* LWORK, int* INFO);
65  void sgeev_(const char* JOBVL, const char* JOBVR, int* N, float* A, int* lda, float* WR, float* WI, float* VL, int* LDVL, float* VR, int* LDVR , float* WORK, int* LWORK, int* INFO);
66 
67  // Cholesky decomposition
68  void dpotrf_(const char* UPLO, const int* N, double* A, const int* LDA, int* INFO);
69  void spotrf_(const char* UPLO, const int* N, float* A, const int* LDA, int* INFO);
70 
71  // Cholesky solve AX=B given decomposition
72  void dpotrs_(const char* UPLO, const int* N, const int* NRHS, const double* A, const int* LDA, double* B, const int* LDB, int* INFO);
73  void spotrs_(const char* UPLO, const int* N, const int* NRHS, const float* A, const int* LDA, float* B, const int* LDB, int* INFO);
74 
75  // Cholesky inverse given decomposition
76  void dpotri_(const char* UPLO, const int* N, double* A, const int* LDA, int* INFO);
77  void spotri_(const char* UPLO, const int* N, float* A, const int* LDA, int* INFO);
78 
79  // Computes a QR decomposition of a general rectangular matrix with column pivoting
80  void sgeqp3_(int* M, int* N, float* A, int* LDA, int* JPVT, float* TAU, float* WORK, int* LWORK, int* INFO );
81  void dgeqp3_(int* M, int* N, double* A, int* LDA, int* JPVT, double* TAU, double* WORK, int* LWORK, int* INFO );
82 
83  //Reconstruct Q from a QR decomposition
84  void sorgqr_(int* M,int* N,int* K, float* A, int* LDA, float* TAU, float* WORK, int* LWORK, int* INFO );
85  void dorgqr_(int* M,int* N,int* K, double* A, int* LDA, double* TAU, double* WORK, int* LWORK, int* INFO );
86  }
87 
88 
90  // C++ overloaded functions to access single and double precision automatically //
92 
93  inline void getrf_(int* M, int *N, float* A, int* lda, int* IPIV, int* INFO){
94  sgetrf_(M, N, A, lda, IPIV, INFO);
95  }
96 
97  inline void getrf_(int* M, int *N, double* A, int* lda, int* IPIV, int* INFO){
98  dgetrf_(M, N, A, lda, IPIV, INFO);
99  }
100 
101  inline void trsm_(const char* SIDE, const char* UPLO, const char* TRANSA, const char* DIAG, int* M, int* N, float* alpha, float* A, int* lda, float* B, int* ldb) {
102  strsm_(const_cast<char*>(SIDE), const_cast<char*>(UPLO), const_cast<char*>(TRANSA), const_cast<char*>(DIAG), M, N, alpha, A, lda, B, ldb);
103  }
104 
105  inline void trsm_(const char* SIDE, const char* UPLO, const char* TRANSA, const char* DIAG, int* M, int* N, double* alpha, double* A, int* lda, double* B, int* ldb) {
106  dtrsm_(const_cast<char*>(SIDE), const_cast<char*>(UPLO), const_cast<char*>(TRANSA), const_cast<char*>(DIAG), M, N, alpha, A, lda, B, ldb);
107  }
108 
109  inline void getri_(int* N, double* A, int* lda, int* IPIV, double* WORK, int* lwork, int* INFO){
110  dgetri_(N, A, lda, IPIV, WORK, lwork, INFO);
111  }
112 
113  inline void getri_(int* N, float* A, int* lda, int* IPIV, float* WORK, int* lwork, int* INFO){
114  sgetri_(N, A, lda, IPIV, WORK, lwork, INFO);
115  }
116 
117  inline void potrf_(const char * UPLO, const int* N, double* A, const int* LDA, int* INFO){
118  dpotrf_(UPLO, N, A, LDA, INFO);
119  }
120 
121  inline void potrf_(const char * UPLO, const int* N, float* A, const int* LDA, int* INFO){
122  spotrf_(UPLO, N, A, LDA, INFO);
123  }
124 
125  // SVD
126  inline void gesvd_(const char* JOBU, const char* JOBVT, int* M, int *N, double* A, int* lda,
127  double* S, double *U, int* ldu, double* VT, int* ldvt,
128  double* WORK, int* lwork, int* INFO){
129  dgesvd_(JOBU, JOBVT, M, N, A, lda, S, U, ldu, VT, ldvt, WORK, lwork, INFO);
130  }
131 
132  inline void gesvd_(const char* JOBU, const char* JOBVT, int* M, int *N, float* A, int* lda,
133  float* S, float *U, int* ldu, float* VT, int* ldvt,
134  float* WORK, int* lwork, int* INFO){
135  sgesvd_(JOBU, JOBVT, M, N, A, lda, S, U, ldu, VT, ldvt, WORK, lwork, INFO);
136  }
137 
138  // Cholesky solve AX=B given decomposition
139  inline void potrs_(const char* UPLO, const int* N, const int* NRHS, const double* A, const int* LDA, double* B, const int* LDB, int* INFO){
140  dpotrs_(UPLO, N, NRHS, A, LDA, B, LDB, INFO);
141  }
142 
143  inline void potrs_(const char* UPLO, const int* N, const int* NRHS, const float* A, const int* LDA, float* B, const int* LDB, int* INFO){
144  spotrs_(UPLO, N, NRHS, A, LDA, B, LDB, INFO);
145  }
146 
147  // Cholesky inverse given decomposition
148  inline void potri_(const char* UPLO, const int* N, double* A, const int* LDA, int* INFO){
149  dpotri_(UPLO, N, A, LDA, INFO);
150  }
151 
152  inline void potri_(const char* UPLO, const int* N, float* A, const int* LDA, int* INFO){
153  spotri_(UPLO, N, A, LDA, INFO);
154  }
155 
156  inline void syev_(const char* JOBZ, const char* UPLO, int* N, double* A, int* lda, double* W, double* WORK, int* LWORK, int* INFO){
157  dsyev_(JOBZ, UPLO, N, A, lda, W, WORK, LWORK, INFO);
158  }
159  inline void syev_(const char* JOBZ, const char* UPLO, int* N, float* A, int* lda, float* W, float* WORK, int* LWORK, int* INFO){
160  ssyev_(JOBZ, UPLO, N, A, lda, W, WORK, LWORK, INFO);
161  }
162 
163  //QR decomposition
164  inline void geqp3_(int* M, int* N, float* A, int* LDA, int* JPVT, float* TAU, float* WORK, int* LWORK, int* INFO )
165  {
166  sgeqp3_(M, N, A, LDA, JPVT, TAU, WORK, LWORK, INFO);
167  }
168 
169  inline void geqp3_(int* M, int* N, double* A, int* LDA, int* JPVT, double* TAU, double* WORK, int* LWORK, int* INFO )
170  {
171  dgeqp3_(M, N, A, LDA, JPVT, TAU, WORK, LWORK, INFO);
172  }
173 
174  inline void orgqr_(int* M,int* N,int* K, float* A, int* LDA, float* TAU, float* WORK, int* LWORK, int* INFO )
175  {
176  sorgqr_(M, N, K, A, LDA, TAU, WORK, LWORK, INFO);
177  }
178 
179  inline void orgqr_(int* M,int* N,int* K, double* A, int* LDA, double* TAU, double* WORK, int* LWORK, int* INFO )
180  {
181  dorgqr_(M, N, K, A, LDA, TAU, WORK, LWORK, INFO);
182  }
183 
184  //Non symmetric (general) eigen decomposition
185  inline void geev_(const char* JOBVL, const char* JOBVR, int* N, double* A, int* lda, double* WR, double* WI, double* VL, int* LDVL, double* VR, int* LDVR , double* WORK, int* LWORK, int* INFO){
186  dgeev_(JOBVL, JOBVR, N, A, lda, WR, WI, VL, LDVL, VR, LDVR , WORK, LWORK, INFO);
187  }
188 
189  inline void geev_(const char* JOBVL, const char* JOBVR, int* N, float* A, int* lda, float* WR, float* WI, float* VL, int* LDVL, float* VR, int* LDVR , float* WORK, int* LWORK, int* INFO){
190  sgeev_(JOBVL, JOBVR, N, A, lda, WR, WI, VL, LDVL, VR, LDVR , WORK, LWORK, INFO);
191  }
192 }
193 #endif
void dpotri_(const char *UPLO, const int *N, double *A, const int *LDA, int *INFO)
void dpotrf_(const char *UPLO, const int *N, double *A, const int *LDA, int *INFO)
void potrf_(const char *UPLO, const int *N, double *A, const int *LDA, int *INFO)
Definition: lapack.h:117
void dgeev_(const char *JOBVL, const char *JOBVR, int *N, double *A, int *lda, double *WR, double *WI, double *VL, int *LDVL, double *VR, int *LDVR, double *WORK, int *LWORK, int *INFO)
void sgeqp3_(int *M, int *N, float *A, int *LDA, int *JPVT, float *TAU, float *WORK, int *LWORK, int *INFO)
void trsm_(const char *SIDE, const char *UPLO, const char *TRANSA, const char *DIAG, int *M, int *N, float *alpha, float *A, int *lda, float *B, int *ldb)
Definition: lapack.h:101
void dgetrf_(int *M, int *N, double *A, int *lda, int *IPIV, int *INFO)
void spotrs_(const char *UPLO, const int *N, const int *NRHS, const float *A, const int *LDA, float *B, const int *LDB, int *INFO)
void sorgqr_(int *M, int *N, int *K, float *A, int *LDA, float *TAU, float *WORK, int *LWORK, int *INFO)
Everything lives inside this namespace.
Definition: allocator.hh:48
void orgqr_(int *M, int *N, int *K, float *A, int *LDA, float *TAU, float *WORK, int *LWORK, int *INFO)
Definition: lapack.h:174
void sgeev_(const char *JOBVL, const char *JOBVR, int *N, float *A, int *lda, float *WR, float *WI, float *VL, int *LDVL, float *VR, int *LDVR, float *WORK, int *LWORK, int *INFO)
void dorgqr_(int *M, int *N, int *K, double *A, int *LDA, double *TAU, double *WORK, int *LWORK, int *INFO)
void potri_(const char *UPLO, const int *N, double *A, const int *LDA, int *INFO)
Definition: lapack.h:148
void spotrf_(const char *UPLO, const int *N, float *A, const int *LDA, int *INFO)
void dtrsm_(char *SIDE, char *UPLO, char *TRANSA, char *DIAG, int *M, int *N, double *alpha, double *A, int *lda, double *B, int *ldb)
void dpotrs_(const char *UPLO, const int *N, const int *NRHS, const double *A, const int *LDA, double *B, const int *LDB, int *INFO)
void potrs_(const char *UPLO, const int *N, const int *NRHS, const double *A, const int *LDA, double *B, const int *LDB, int *INFO)
Definition: lapack.h:139
void dsyev_(const char *JOBZ, const char *UPLO, int *N, double *A, int *lda, double *W, double *WORK, int *LWORK, int *INFO)
void syev_(const char *JOBZ, const char *UPLO, int *N, double *A, int *lda, double *W, double *WORK, int *LWORK, int *INFO)
Definition: lapack.h:156
void strsm_(char *SIDE, char *UPLO, char *TRANSA, char *DIAG, int *M, int *N, float *alpha, float *A, int *lda, float *B, int *ldb)
void getri_(int *N, double *A, int *lda, int *IPIV, double *WORK, int *lwork, int *INFO)
Definition: lapack.h:109
void sgesvd_(const char *JOBU, const char *JOBVT, int *M, int *N, float *A, int *lda, float *S, float *U, int *ldu, float *VT, int *ldvt, float *WORK, int *lwork, int *INFO)
void geev_(const char *JOBVL, const char *JOBVR, int *N, double *A, int *lda, double *WR, double *WI, double *VL, int *LDVL, double *VR, int *LDVR, double *WORK, int *LWORK, int *INFO)
Definition: lapack.h:185
void sgetri_(int *N, float *A, int *lda, int *IPIV, float *WORK, int *lwork, int *INFO)
void dgeqp3_(int *M, int *N, double *A, int *LDA, int *JPVT, double *TAU, double *WORK, int *LWORK, int *INFO)
void dgetri_(int *N, double *A, int *lda, int *IPIV, double *WORK, int *lwork, int *INFO)
void spotri_(const char *UPLO, const int *N, float *A, const int *LDA, int *INFO)
void sgetrf_(int *M, int *N, float *A, int *lda, int *IPIV, int *INFO)
void getrf_(int *M, int *N, float *A, int *lda, int *IPIV, int *INFO)
Definition: lapack.h:93
void dgesvd_(const char *JOBU, const char *JOBVT, int *M, int *N, double *A, int *lda, double *S, double *U, int *ldu, double *VT, int *ldvt, double *WORK, int *lwork, int *INFO)
void ssyev_(const char *JOBZ, const char *UPLO, int *N, float *A, int *lda, float *W, float *WORK, int *LWORK, int *INFO)
void gesvd_(const char *JOBU, const char *JOBVT, int *M, int *N, double *A, int *lda, double *S, double *U, int *ldu, double *VT, int *ldvt, double *WORK, int *lwork, int *INFO)
Definition: lapack.h:126
void geqp3_(int *M, int *N, float *A, int *LDA, int *JPVT, float *TAU, float *WORK, int *LWORK, int *INFO)
Definition: lapack.h:164