14 # include <vcl_msvc_warnings.h> 23 #ifndef DOXYGEN_SHOULD_SKIP_THIS 25 inline void vnl_linpack_qrdc(vnl_netlib_qrdc_proto(T)) \ 26 { v3p_netlib_##p##qrdc_(vnl_netlib_qrdc_params); } \ 27 inline void vnl_linpack_qrsl(vnl_netlib_qrsl_proto(T)) \ 28 { v3p_netlib_##p##qrsl_(vnl_netlib_qrsl_params); } 31 macro(c, std::complex<float>);
32 macro(z, std::complex<double>);
38 qrdc_out_(M.columns(), M.rows()),
49 for (
int i = 0; i < r; ++i)
50 for (
int j = 0; j < c; ++j)
80 int m =
std::min((
int)qrdc_out_.columns(), (int)qrdc_out_.rows());
81 T det = qrdc_out_(0,0);
83 for (
int i = 1; i <
m; ++i)
84 det *= -qrdc_out_(i,i);
93 int m = qrdc_out_.columns();
94 int n = qrdc_out_.rows();
102 std::cerr << __FILE__
": vnl_qr<T>::Q()\n" 103 <<
" m,n = " <<
m <<
", " << n <<
'\n' 104 <<
" qr0 = [" << qrdc_out_ <<
"];\n" 105 <<
" aux = [" << qraux_ <<
"];\n";
117 for (
int k = n-1; k >= 0; --k) {
118 if (k >=
m)
continue;
122 for (
int j = k+1; j <
m; ++j) {
123 v[j] = qrdc_out_(k,j);
127 #ifndef DOXYGEN_SHOULD_SKIP_THIS 128 # define c vnl_complex_traits<T>::conjugate 134 abs_t scale = abs_t(2)/sq;
136 for (
int i = k; i <
m; ++i) {
138 for (
int j = k; j <
m; ++j)
139 w[i] += scale * c(
v[j]) * matrQ(j, i);
144 for (
int i = k; i <
m; ++i)
145 for (
int j = k; j <
m; ++j)
146 matrQ(i,j) -= (
v[i]) * (w[j]);
159 int m = qrdc_out_.columns();
160 int n = qrdc_out_.rows();
164 for (
int i = 0; i <
m; ++i)
165 for (
int j = 0; j < n; ++j)
169 Rmatr(i, j) = qrdc_out_(j,i);
191 long n = qrdc_out_.columns();
192 long p = qrdc_out_.rows();
201 vnl_linpack_qrsl(qrdc_out_.data_block(),
212 std::cerr << __FILE__
": vnl_qr<T>::solve() : matrix is rank-deficient by " 222 long n = qrdc_out_.columns();
223 long p = qrdc_out_.rows();
231 vnl_linpack_qrsl(qrdc_out_.data_block(),
244 std::cerr << __FILE__
": vnl_qr<T>::QtB() -- matrix is rank-deficient by " 253 unsigned int r = qrdc_out_.columns();
254 assert(r > 0 && r == qrdc_out_.rows());
259 for (
unsigned int i=0; i<r; ++i)
263 inv.set_column(i,col);
272 unsigned int r = qrdc_out_.columns();
273 assert(r > 0 && r == qrdc_out_.rows());
278 for (
unsigned int i=0; i<r; ++i)
291 assert(rhs.
rows() == qrdc_out_.columns());
292 int c = qrdc_out_.rows();
296 for (
int i=0; i<n; ++i)
307 #define VNL_QR_INSTANTIATE(T) \ 308 template class VNL_ALGO_EXPORT vnl_qr<T >; \ vnl_matrix< T > inverse() const
return the inverse matrix of M.
Print matrices and vectors in nice MATLAB format.
vnl_matrix< T > solve(const vnl_matrix< T > &rhs) const
Solve equation M x = rhs for x using the computed decomposition.
Complex additions to vnl_math.
vnl_matrix< T > tinverse() const
return the transpose of the inverse matrix of M.
VNL_EXPORT std::ostream & vnl_matlab_print(std::ostream &, vnl_diag_matrix< T > const &, char const *variable_name=nullptr, vnl_matlab_print_format=vnl_matlab_print_format_default)
print a vnl_diagonal_matrix<T>.
To allow templated algorithms to determine appropriate actions of conjugation, complexification etc.
T const * data_block() const
Access the contiguous block storing the elements in the vector. O(1).
vnl_matrix< T > qrdc_out_
Namespace with standard math functions.
vnl_matrix< T > const & Q() const
Unpack and return unitary part Q.
vnl_matrix< T > const & R() const
Unpack and return R.
Calculate inverse of a matrix using QR.
vnl_bignum squared_magnitude(vnl_bignum const &x)
vnl_matrix & set_column(unsigned i, T const *v)
Set the elements of the i'th column to v[i] (No bounds checking).
T determinant() const
Return the determinant of M. This is computed from M = Q R as follows:.
An ordinary mathematical matrix.
Declare in a central place the list of symbols from netlib.
vnl_vector< T > QtB(const vnl_vector< T > &b) const
Return residual vector d of M x = b -> d = Q'b.
vnl_vector & fill(T const &v)
Set all values to v.
bool empty() const
Return true iff the size is zero.
Mathematical vector class, templated by type of element.
vnl_vector< T > get_column(unsigned c) const
Get a vector equal to the given column.
unsigned int rows() const
Return the number of rows.
vnl_qr(vnl_matrix< T > const &M)
vnl_decnum min(vnl_decnum const &x, vnl_decnum const &y)
unsigned int columns() const
Return the number of columns.
vnl_matrix< T > recompose() const
return the original matrix M.