00001
00002 #ifndef vnl_matrix_h_
00003 #define vnl_matrix_h_
00004 #ifdef VCL_NEEDS_PRAGMA_INTERFACE
00005 #pragma interface
00006 #endif
00007
00008
00009
00010
00011 #include <vcl_iosfwd.h>
00012 #include <vnl/vnl_tag.h>
00013 #include <vnl/vnl_c_vector.h>
00014 #include <vnl/vnl_config.h>
00015 #ifndef NDEBUG
00016 # if VNL_CONFIG_CHECK_BOUNDS
00017 # include <vnl/vnl_error.h>
00018 # include <vcl_cassert.h>
00019 # endif
00020 #else
00021 # undef VNL_CONFIG_CHECK_BOUNDS
00022 # define VNL_CONFIG_CHECK_BOUNDS 0
00023 # undef ERROR_CHECKING
00024 #endif
00025
00026 export template <class T> class vnl_vector;
00027 export template <class T> class vnl_matrix;
00028
00029
00030
00031 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00032 #define v vnl_vector<T>
00033 #define m vnl_matrix<T>
00034 #endif // DOXYGEN_SHOULD_SKIP_THIS
00035 template <class T> m operator+(T const&, m const&);
00036 template <class T> m operator-(T const&, m const&);
00037 template <class T> m operator*(T const&, m const&);
00038 template <class T> m element_product(m const&, m const&);
00039 template <class T> m element_quotient(m const&, m const&);
00040 template <class T> T dot_product(m const&, m const&);
00041 template <class T> T inner_product(m const&, m const&);
00042 template <class T> T cos_angle(m const&, m const& );
00043 template <class T> vcl_ostream& operator<<(vcl_ostream&, m const&);
00044 template <class T> vcl_istream& operator>>(vcl_istream&, m&);
00045 #undef v
00046 #undef m
00047
00048
00049
00050 enum vnl_matrix_type
00051 {
00052 vnl_matrix_null,
00053 vnl_matrix_identity
00054 };
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082 template<class T>
00083 class vnl_matrix
00084 {
00085 public:
00086
00087 vnl_matrix() :
00088 num_rows(0),
00089 num_cols(0),
00090 data(0)
00091 {
00092 }
00093
00094
00095
00096
00097 vnl_matrix(unsigned r, unsigned c);
00098
00099
00100
00101 vnl_matrix(unsigned r, unsigned c, T const& v0);
00102
00103
00104
00105
00106 vnl_matrix(unsigned r, unsigned c, vnl_matrix_type t);
00107
00108
00109
00110
00111 vnl_matrix(unsigned r, unsigned c, unsigned n, T const values[]);
00112
00113
00114
00115
00116 vnl_matrix(T const* data_block, unsigned r, unsigned c);
00117
00118
00119
00120 vnl_matrix(vnl_matrix<T> const&);
00121
00122 #ifndef VXL_DOXYGEN_SHOULD_SKIP_THIS
00123
00124
00125
00126 vnl_matrix(vnl_matrix<T> const &, vnl_matrix<T> const &, vnl_tag_add);
00127 vnl_matrix(vnl_matrix<T> const &, vnl_matrix<T> const &, vnl_tag_sub);
00128 vnl_matrix(vnl_matrix<T> const &, T, vnl_tag_mul);
00129 vnl_matrix(vnl_matrix<T> const &, T, vnl_tag_div);
00130 vnl_matrix(vnl_matrix<T> const &, T, vnl_tag_add);
00131 vnl_matrix(vnl_matrix<T> const &, T, vnl_tag_sub);
00132 vnl_matrix(vnl_matrix<T> const &, vnl_matrix<T> const &, vnl_tag_mul);
00133 vnl_matrix(vnl_matrix<T> &that, vnl_tag_grab)
00134 : num_rows(that.num_rows), num_cols(that.num_cols), data(that.data)
00135 { that.num_cols=that.num_rows=0; that.data=0; }
00136
00137 #endif
00138
00139
00140 ~vnl_matrix();
00141
00142
00143
00144
00145 unsigned rows() const { return num_rows; }
00146
00147
00148
00149 unsigned columns() const { return num_cols; }
00150
00151
00152
00153 unsigned cols() const { return num_cols; }
00154
00155
00156
00157 unsigned size() const { return rows()*cols(); }
00158
00159
00160 void put(unsigned r, unsigned c, T const&);
00161
00162
00163 T get(unsigned r, unsigned c) const;
00164
00165
00166
00167 T * operator[](unsigned r) { return data[r]; }
00168
00169
00170
00171 T const * operator[](unsigned r) const { return data[r]; }
00172
00173
00174
00175 T & operator()(unsigned r, unsigned c)
00176 {
00177 #if VNL_CONFIG_CHECK_BOUNDS
00178 assert(r<rows());
00179 assert(c<cols());
00180 #endif
00181 return this->data[r][c];
00182 }
00183
00184
00185
00186 T const & operator()(unsigned r, unsigned c) const
00187 {
00188 #if VNL_CONFIG_CHECK_BOUNDS
00189 assert(r<rows());
00190 assert(c<cols());
00191 #endif
00192 return this->data[r][c];
00193 }
00194
00195
00196
00197
00198
00199
00200 void fill(T const&);
00201
00202
00203
00204 void fill_diagonal(T const&);
00205
00206
00207
00208 void copy_in(T const *);
00209
00210
00211
00212 void set(T const *d) { copy_in(d); }
00213
00214
00215
00216
00217 void copy_out(T *) const;
00218
00219
00220
00221
00222 vnl_matrix<T>& operator=(T const&v) { fill(v); return *this; }
00223
00224
00225
00226 vnl_matrix<T>& operator=(vnl_matrix<T> const&);
00227
00228
00229
00230
00231
00232
00233 vnl_matrix<T>& operator+=(T value);
00234
00235
00236 vnl_matrix<T>& operator-=(T value);
00237
00238
00239 vnl_matrix<T>& operator*=(T value);
00240
00241
00242 vnl_matrix<T>& operator/=(T value);
00243
00244
00245 vnl_matrix<T>& operator+=(vnl_matrix<T> const&);
00246
00247 vnl_matrix<T>& operator-=(vnl_matrix<T> const&);
00248
00249 vnl_matrix<T>& operator*=(vnl_matrix<T> const&rhs) { return *this = (*this) * rhs; }
00250
00251
00252 vnl_matrix<T> operator-() const;
00253
00254
00255
00256 vnl_matrix<T> operator+(T const& v) const { return vnl_matrix<T>(*this, v, vnl_tag_add()); }
00257
00258
00259 vnl_matrix<T> operator-(T const& v) const { return vnl_matrix<T>(*this, v, vnl_tag_sub()); }
00260
00261
00262 vnl_matrix<T> operator*(T const& v) const { return vnl_matrix<T>(*this, v, vnl_tag_mul()); }
00263
00264
00265 vnl_matrix<T> operator/(T const& v) const { return vnl_matrix<T>(*this, v, vnl_tag_div()); }
00266
00267
00268 vnl_matrix<T> operator+(vnl_matrix<T> const& rhs) const { return vnl_matrix<T>(*this, rhs, vnl_tag_add()); }
00269
00270 vnl_matrix<T> operator-(vnl_matrix<T> const& rhs) const { return vnl_matrix<T>(*this, rhs, vnl_tag_sub()); }
00271
00272 vnl_matrix<T> operator*(vnl_matrix<T> const& rhs) const { return vnl_matrix<T>(*this, rhs, vnl_tag_mul()); }
00273
00274
00275
00276
00277 vnl_matrix<T> apply(T (*f)(T)) const;
00278
00279
00280 vnl_matrix<T> apply(T (*f)(T const&)) const;
00281
00282
00283 vnl_matrix<T> transpose() const;
00284
00285
00286 vnl_matrix<T> conjugate_transpose() const;
00287
00288
00289 vnl_matrix<T>& update(vnl_matrix<T> const&, unsigned top=0, unsigned left=0);
00290
00291
00292 void set_column(unsigned i, T const * v);
00293
00294
00295 void set_column(unsigned i, T value );
00296
00297
00298 void set_column(unsigned j, vnl_vector<T> const& v);
00299
00300
00301 void set_columns(unsigned starting_column, vnl_matrix<T> const& M);
00302
00303
00304 void set_row(unsigned i, T const * v);
00305
00306
00307 void set_row(unsigned i, T value );
00308
00309
00310 void set_row(unsigned i, vnl_vector<T> const&);
00311
00312
00313
00314 vnl_matrix<T> extract(unsigned r, unsigned c,
00315 unsigned top=0, unsigned left=0) const;
00316
00317
00318
00319
00320
00321
00322 void extract ( vnl_matrix<T>& sub_matrix,
00323 unsigned top=0, unsigned left=0) const;
00324
00325
00326
00327 vnl_vector<T> get_row(unsigned r) const;
00328
00329
00330 vnl_vector<T> get_column(unsigned c) const;
00331
00332
00333 vnl_matrix<T> get_n_rows(unsigned rowstart, unsigned n) const;
00334
00335
00336 vnl_matrix<T> get_n_columns(unsigned colstart, unsigned n) const;
00337
00338
00339
00340
00341
00342 void set_identity();
00343
00344
00345 void inplace_transpose();
00346
00347
00348 void flipud();
00349
00350 void fliplr();
00351
00352
00353
00354 void normalize_rows();
00355
00356
00357
00358 void normalize_columns();
00359
00360
00361 void scale_row(unsigned row, T value);
00362
00363
00364 void scale_column(unsigned col, T value);
00365
00366
00367 void swap(vnl_matrix<T> & that);
00368
00369
00370 typedef typename vnl_c_vector<T>::abs_t abs_t;
00371
00372
00373 abs_t array_one_norm() const { return vnl_c_vector<T>::one_norm(begin(), size()); }
00374
00375
00376 abs_t array_two_norm() const { return vnl_c_vector<T>::two_norm(begin(), size()); }
00377
00378
00379 abs_t array_inf_norm() const { return vnl_c_vector<T>::inf_norm(begin(), size()); }
00380
00381
00382 abs_t absolute_value_sum() const { return array_one_norm(); }
00383
00384
00385 abs_t absolute_value_max() const { return array_inf_norm(); }
00386
00387
00388 abs_t operator_one_norm() const;
00389
00390
00391 abs_t operator_inf_norm() const;
00392
00393
00394 abs_t frobenius_norm() const { return vnl_c_vector<T>::two_norm(begin(), size()); }
00395
00396
00397 abs_t fro_norm() const { return frobenius_norm(); }
00398
00399
00400 abs_t rms() const { return vnl_c_vector<T>::rms_norm(begin(), size()); }
00401
00402
00403 T min_value() const { return vnl_c_vector<T>::min_value(begin(), size()); }
00404
00405
00406 T max_value() const { return vnl_c_vector<T>::max_value(begin(), size()); }
00407
00408
00409 unsigned arg_min() const { return vnl_c_vector<T>::arg_min(begin(), size()); }
00410
00411
00412 unsigned arg_max() const { return vnl_c_vector<T>::arg_max(begin(), size()); }
00413
00414
00415 T mean() const { return vnl_c_vector<T>::mean(begin(), size()); }
00416
00417
00418
00419
00420 bool empty() const { return !data || !num_rows || !num_cols; }
00421
00422
00423 bool is_identity() const;
00424
00425
00426 bool is_identity(double tol) const;
00427
00428
00429 bool is_zero() const;
00430
00431
00432 bool is_zero(double tol) const;
00433
00434
00435 bool is_finite() const;
00436
00437
00438 bool has_nans() const;
00439
00440
00441
00442 void assert_size(unsigned r, unsigned c) const
00443 {
00444 #ifndef NDEBUG
00445 assert_size_internal(r, c);
00446 #endif
00447 }
00448
00449
00450 void assert_finite() const
00451 {
00452 #ifndef NDEBUG
00453 assert_finite_internal();
00454 #endif
00455 }
00456
00457
00458
00459
00460 static vnl_matrix<T> read(vcl_istream& s);
00461
00462
00463 bool read_ascii(vcl_istream& s);
00464
00465
00466
00467
00468
00469 T const* data_block() const { return data[0]; }
00470
00471
00472
00473 T * data_block() { return data[0]; }
00474
00475
00476
00477 T const* const* data_array() const { return data; }
00478
00479
00480
00481 T * * data_array() { return data; }
00482
00483 typedef T element_type;
00484
00485
00486 typedef T *iterator;
00487
00488 iterator begin() { return data?data[0]:0; }
00489
00490 iterator end() { return data?data[0]+num_rows*num_cols:0; }
00491
00492
00493 typedef T const *const_iterator;
00494
00495 const_iterator begin() const { return data?data[0]:0; }
00496
00497 const_iterator end() const { return data?data[0]+num_rows*num_cols:0; }
00498
00499
00500
00501
00502
00503 vnl_matrix<T> const& as_ref() const { return *this; }
00504
00505
00506 vnl_matrix<T>& as_ref() { return *this; }
00507
00508
00509
00510
00511 bool operator_eq(vnl_matrix<T> const & rhs) const;
00512
00513
00514 bool operator==(vnl_matrix<T> const &that) const { return this->operator_eq(that); }
00515
00516
00517 bool operator!=(vnl_matrix<T> const &that) const { return !this->operator_eq(that); }
00518
00519
00520 void print(vcl_ostream& os) const;
00521
00522
00523 void clear();
00524
00525
00526
00527 bool set_size(unsigned r, unsigned c);
00528
00529
00530
00531 protected:
00532 unsigned num_rows;
00533 unsigned num_cols;
00534 T** data;
00535
00536 #if VCL_HAS_SLICED_DESTRUCTOR_BUG
00537
00538
00539 char vnl_matrix_own_data;
00540 #endif
00541
00542 void assert_size_internal(unsigned r, unsigned c) const;
00543 void assert_finite_internal() const;
00544
00545
00546 void destroy();
00547
00548 #if VCL_NEED_FRIEND_FOR_TEMPLATE_OVERLOAD
00549 # define v vnl_vector<T>
00550 # define m vnl_matrix<T>
00551 friend m operator+ VCL_NULL_TMPL_ARGS (T const&, m const&);
00552 friend m operator- VCL_NULL_TMPL_ARGS (T const&, m const&);
00553 friend m operator* VCL_NULL_TMPL_ARGS (T const&, m const&);
00554 friend m element_product VCL_NULL_TMPL_ARGS (m const&, m const&);
00555 friend m element_quotient VCL_NULL_TMPL_ARGS (m const&, m const&);
00556 friend T dot_product VCL_NULL_TMPL_ARGS (m const&, m const&);
00557 friend T inner_product VCL_NULL_TMPL_ARGS (m const&, m const&);
00558 friend T cos_angle VCL_NULL_TMPL_ARGS (m const&, m const&);
00559 friend vcl_ostream& operator<< VCL_NULL_TMPL_ARGS (vcl_ostream&, m const&);
00560 friend vcl_istream& operator>> VCL_NULL_TMPL_ARGS (vcl_istream&, m&);
00561 # undef v
00562 # undef m
00563 #endif
00564
00565
00566 static void inline_function_tickler();
00567 };
00568
00569
00570
00571
00572
00573
00574
00575
00576 template<class T>
00577 inline T vnl_matrix<T>::get(unsigned row, unsigned column) const
00578 {
00579 #ifdef ERROR_CHECKING
00580 if (row >= this->num_rows)
00581 vnl_error_matrix_row_index("get", row);
00582 if (column >= this->num_cols)
00583 vnl_error_matrix_col_index("get", column);
00584 #endif
00585 return this->data[row][column];
00586 }
00587
00588
00589
00590
00591 template<class T>
00592 inline void vnl_matrix<T>::put(unsigned row, unsigned column, T const& value)
00593 {
00594 #ifdef ERROR_CHECKING
00595 if (row >= this->num_rows)
00596 vnl_error_matrix_row_index("put", row);
00597 if (column >= this->num_cols)
00598 vnl_error_matrix_col_index("put", column);
00599 #endif
00600 this->data[row][column] = value;
00601 }
00602
00603
00604
00605
00606
00607
00608 template<class T>
00609 inline vnl_matrix<T> operator*(T const& value, vnl_matrix<T> const& m)
00610 {
00611 return vnl_matrix<T>(m, value, vnl_tag_mul());
00612 }
00613
00614
00615
00616 template<class T>
00617 inline vnl_matrix<T> operator+(T const& value, vnl_matrix<T> const& m)
00618 {
00619 return vnl_matrix<T>(m, value, vnl_tag_add());
00620 }
00621
00622
00623
00624 template<class T>
00625 inline void swap(vnl_matrix<T> &A, vnl_matrix<T> &B) { A.swap(B); }
00626
00627
00628 #endif // vnl_matrix_h_