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.
LU.h
Go to the documentation of this file.
1 // -*- c++ -*-
2 
3 // Copyright (C) 2005,2009 Tom Drummond (twd20@cam.ac.uk),
4 // Ed Rosten (er258@cam.ac.uk)
5 //
6 // This file is part of the TooN Library. This library is free
7 // software; you can redistribute it and/or modify it under the
8 // terms of the GNU General Public License as published by the
9 // Free Software Foundation; either version 2, or (at your option)
10 // any later version.
11 
12 // This library is distributed in the hope that it will be useful,
13 // but WITHOUT ANY WARRANTY; without even the implied warranty of
14 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 // GNU General Public License for more details.
16 
17 // You should have received a copy of the GNU General Public License along
18 // with this library; see the file COPYING. If not, write to the Free
19 // Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307,
20 // USA.
21 
22 // As a special exception, you may use this file as part of a free software
23 // library without restriction. Specifically, if other files instantiate
24 // templates or use macros or inline functions from this file, or you compile
25 // this file and link it with other files to produce an executable, this
26 // file does not by itself cause the resulting executable to be covered by
27 // the GNU General Public License. This exception does not however
28 // invalidate any other reasons why the executable file might be covered by
29 // the GNU General Public License.
30 
31 #ifndef TOON_INCLUDE_LU_H
32 #define TOON_INCLUDE_LU_H
33 
34 #include <iostream>
35 
36 #include "lapack.h"
37 
38 #include "TooN.h"
39 
40 namespace TooN {
68 template <int Size=-1, class Precision=double>
69 class LU {
70  public:
71 
74  template<int S1, int S2, class Base>
76  :my_lu(m.num_rows(),m.num_cols()),my_IPIV(m.num_rows()){
77  compute(m);
78  }
79 
81  template<int S1, int S2, class Base>
83  //check for consistency with Size
84  SizeMismatch<Size, S1>::test(my_lu.num_rows(),m.num_rows());
85  SizeMismatch<Size, S2>::test(my_lu.num_rows(),m.num_cols());
86 
87  //Make a local copy. This is guaranteed contiguous
88  my_lu=m;
89  int lda = m.num_rows();
90  int M = m.num_rows();
91  int N = m.num_rows();
92 
93  getrf_(&M,&N,&my_lu[0][0],&lda,&my_IPIV[0],&my_info);
94 
95  if(my_info < 0){
96  std::cerr << "error in LU, INFO was " << my_info << std::endl;
97  }
98  }
99 
102  template <int Rows, int NRHS, class Base>
104  //Check the number of rows is OK.
105  SizeMismatch<Size, Rows>::test(my_lu.num_rows(), rhs.num_rows());
106 
107  Matrix<Size, NRHS, Precision> result(rhs);
108 
109  int M=rhs.num_cols();
110  int N=my_lu.num_rows();
111  double alpha=1;
112  int lda=my_lu.num_rows();
113  int ldb=rhs.num_cols();
114  trsm_("R","U","N","N",&M,&N,&alpha,&my_lu[0][0],&lda,&result[0][0],&ldb);
115  trsm_("R","L","N","U",&M,&N,&alpha,&my_lu[0][0],&lda,&result[0][0],&ldb);
116 
117  // now do the row swapping (lapack dlaswp.f only shuffles fortran rows = Rowmajor cols)
118  for(int i=N-1; i>=0; i--){
119  const int swaprow = my_IPIV[i]-1; // fortran arrays start at 1
120  for(int j=0; j<NRHS; j++){
121  Precision temp = result[i][j];
122  result[i][j] = result[swaprow][j];
123  result[swaprow][j] = temp;
124  }
125  }
126  return result;
127  }
128 
131  template <int Rows, class Base>
133  //Check the number of rows is OK.
134  SizeMismatch<Size, Rows>::test(my_lu.num_rows(), rhs.size());
135 
136  Vector<Size, Precision> result(rhs);
137 
138  int M=1;
139  int N=my_lu.num_rows();
140  double alpha=1;
141  int lda=my_lu.num_rows();
142  int ldb=1;
143  trsm_("R","U","N","N",&M,&N,&alpha,&my_lu[0][0],&lda,&result[0],&ldb);
144  trsm_("R","L","N","U",&M,&N,&alpha,&my_lu[0][0],&lda,&result[0],&ldb);
145 
146  // now do the row swapping (lapack dlaswp.f only shuffles fortran rows = Rowmajor cols)
147  for(int i=N-1; i>=0; i--){
148  const int swaprow = my_IPIV[i]-1; // fortran arrays start at 1
149  Precision temp = result[i];
150  result[i] = result[swaprow];
151  result[swaprow] = temp;
152  }
153  return result;
154  }
155 
160  int N = my_lu.num_rows();
161  int lda=my_lu.num_rows();
162  int lwork=-1;
163  Precision size;
164  getri_(&N, &Inverse[0][0], &lda, &my_IPIV[0], &size, &lwork, &my_info);
165  lwork=int(size);
166  Precision* WORK = new Precision[lwork];
167  getri_(&N, &Inverse[0][0], &lda, &my_IPIV[0], WORK, &lwork, &my_info);
168  delete [] WORK;
169  return Inverse;
170  }
171 
177  const Matrix<Size,Size,Precision>& get_lu()const {return my_lu;}
178 
179  private:
180  inline int get_sign() const {
181  int result=1;
182  for(int i=0; i<my_lu.num_rows()-1; i++){
183  if(my_IPIV[i] > i+1){
184  result=-result;
185  }
186  }
187  return result;
188  }
189  public:
190 
192  inline Precision determinant() const {
193  Precision result = get_sign();
194  for (int i=0; i<my_lu.num_rows(); i++){
195  result*=my_lu(i,i);
196  }
197  return result;
198  }
199 
201  int get_info() const { return my_info; }
202 
203  private:
204 
206  int my_info;
207  Vector<Size, int> my_IPIV; //Convenient static-or-dynamic array of ints :-)
208 
209 };
210 }
211 
212 
213 #endif
LU(const Matrix< S1, S2, Precision, Base > &m)
Definition: LU.h:75
const Matrix< Size, Size, Precision > & get_lu() const
Definition: LU.h:177
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 compute(const Matrix< S1, S2, Precision, Base > &m)
Perform the LU decompsition of another matrix.
Definition: LU.h:82
Everything lives inside this namespace.
Definition: allocator.hh:48
Precision determinant() const
Calculate the determinant of the matrix.
Definition: LU.h:192
int get_sign() const
Definition: LU.h:180
Matrix< Size, NRHS, Precision > backsub(const Matrix< Rows, NRHS, Precision, Base > &rhs)
Definition: LU.h:103
static void test(int s1, int s2)
void getri_(int *N, double *A, int *lda, int *IPIV, double *WORK, int *lwork, int *INFO)
Definition: lapack.h:109
Matrix< Size, Size, Precision > my_lu
Definition: LU.h:205
int my_info
Definition: LU.h:206
int get_info() const
Get the LAPACK info.
Definition: LU.h:201
Vector< Size, Precision > backsub(const Vector< Rows, Precision, Base > &rhs)
Definition: LU.h:132
Vector< Size, int > my_IPIV
Definition: LU.h:207
Matrix< Size, Size, Precision > get_inverse()
Definition: LU.h:158
void getrf_(int *M, int *N, float *A, int *lda, int *IPIV, int *INFO)
Definition: lapack.h:93
Definition: LU.h:69