00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 2009, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef O2SCL_OMATRIX_TLATE_H 00024 #define O2SCL_OMATRIX_TLATE_H 00025 00026 /** \file omatrix_tlate.h 00027 \brief File for definitions of matrices 00028 00029 \future The \c xmatrix class demonstrates how operator[] could 00030 return an ovector_array object and thus provide more 00031 bounds-checking. This would demand including a new parameter in 00032 omatrix_view_tlate which contains the vector type. 00033 */ 00034 00035 #include <iostream> 00036 #include <cstdlib> 00037 #include <string> 00038 #include <fstream> 00039 #include <sstream> 00040 00041 #include <gsl/gsl_matrix.h> 00042 #include <gsl/gsl_ieee_utils.h> 00043 00044 #include <o2scl/err_hnd.h> 00045 #include <o2scl/ovector_tlate.h> 00046 #include <o2scl/array.h> 00047 00048 #ifndef DOXYGENP 00049 namespace o2scl { 00050 #endif 00051 00052 /** 00053 \brief A const matrix view of omatrix objects 00054 */ 00055 template<class data_t, class mparent_t, class block_t> 00056 class omatrix_const_view_tlate : 00057 public mparent_t { 00058 00059 public: 00060 00061 /// \name Copy constructors 00062 //@{ 00063 /// Shallow copy constructor - create a new view of the same matrix 00064 omatrix_const_view_tlate(const omatrix_const_view_tlate &v) { 00065 mparent_t::block=0; 00066 mparent_t::data=v.data; 00067 mparent_t::size1=v.size1; 00068 mparent_t::size2=v.size2; 00069 mparent_t::tda=v.tda(); 00070 mparent_t::owner=0; 00071 } 00072 00073 /// Shallow copy constructor - create a new view of the same matrix 00074 omatrix_const_view_tlate& operator=(const omatrix_const_view_tlate &v) { 00075 mparent_t::block=0; 00076 mparent_t::data=v.data; 00077 mparent_t::size1=v.size1; 00078 mparent_t::size2=v.size2; 00079 mparent_t::tda=v.tda; 00080 mparent_t::owner=0; 00081 00082 return *this; 00083 } 00084 //@} 00085 00086 ~omatrix_const_view_tlate() {}; 00087 00088 /// \name Get and set methods 00089 //@{ 00090 /** 00091 \brief Array-like indexing 00092 */ 00093 const data_t *operator[](size_t i) const { 00094 #if O2SCL_NO_RANGE_CHECK 00095 #else 00096 if (i>=mparent_t::size1) { 00097 O2SCL_ERR((((std::string)"Array index ")+itos(i)+" out of bounds" 00098 +" in omatrix_const_view_tlate::operator[] const. Size: "+ 00099 itos(mparent_t::size1)+ 00100 " (index should be less than size).").c_str(),gsl_eindex); 00101 return (mparent_t::data); 00102 } 00103 #endif 00104 return mparent_t::data+i*mparent_t::tda; 00105 } 00106 00107 /** 00108 \brief Array-like indexing 00109 */ 00110 const data_t &operator()(size_t i, size_t j) const { 00111 #if O2SCL_NO_RANGE_CHECK 00112 #else 00113 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00114 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00115 ", "+itos(j)+" out of bounds" 00116 +" in omatrix_cx_view_tlate::operator() const. Sizes: "+ 00117 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00118 " (indices should be less than sizes).").c_str(), 00119 gsl_eindex); 00120 return *(mparent_t::data); 00121 } 00122 #endif 00123 return *(mparent_t::data+i*mparent_t::tda+j); 00124 } 00125 00126 /** \brief Get (with optional range-checking) */ 00127 data_t get(size_t i, size_t j) const { 00128 #if O2SCL_NO_RANGE_CHECK 00129 #else 00130 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00131 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00132 ", "+itos(j)+" out of bounds" 00133 +" in omatrix_cx_view_tlate::get() const. Sizes: "+ 00134 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00135 " (indices should be less than sizes).").c_str(), 00136 gsl_eindex); 00137 return *(mparent_t::data); 00138 } 00139 #endif 00140 return *(mparent_t::data+i*mparent_t::tda+j); 00141 } 00142 00143 /** \brief Get pointer (with optional range-checking) */ 00144 data_t *get_ptr(size_t i, size_t j) { 00145 #if O2SCL_NO_RANGE_CHECK 00146 #else 00147 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00148 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00149 ", "+itos(j)+" out of bounds" 00150 +" in omatrix_cx_view_tlate::get_ptr(). Sizes: "+ 00151 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00152 " (indices should be less than sizes).").c_str(), 00153 gsl_eindex); 00154 return (mparent_t::data); 00155 } 00156 #endif 00157 return mparent_t::data+i*mparent_t::tda+j; 00158 } 00159 00160 /** \brief Get pointer (with optional range-checking) */ 00161 const data_t *get_const_ptr(size_t i, size_t j) const { 00162 #if O2SCL_NO_RANGE_CHECK 00163 #else 00164 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00165 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00166 ", "+itos(j)+" out of bounds" 00167 +" in omatrix_cx_view_tlate::get_const_ptr(). Sizes: "+ 00168 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00169 " (indices should be less than sizes).").c_str(), 00170 gsl_eindex); 00171 return (mparent_t::data); 00172 } 00173 #endif 00174 return (const data_t *)(mparent_t::data+i*mparent_t::tda+j); 00175 } 00176 00177 /** 00178 \brief Method to return number of rows 00179 00180 If no memory has been allocated, this will quietly 00181 return zero. 00182 */ 00183 size_t rows() const { 00184 return mparent_t::size1; 00185 } 00186 00187 /** 00188 \brief Method to return number of columns 00189 00190 If no memory has been allocated, this will quietly 00191 return zero. 00192 */ 00193 size_t cols() const { 00194 return mparent_t::size2; 00195 } 00196 00197 /** 00198 \brief Method to return matrix tda 00199 00200 If no memory has been allocated, this will quietly 00201 return zero. 00202 */ 00203 size_t tda() const { 00204 return mparent_t::tda; 00205 } 00206 //@} 00207 00208 /// \name Other methods 00209 //@{ 00210 /** \brief Return true if this object owns the data it refers to 00211 00212 This can be used to determine if an object is a "matrix_view", 00213 or a "matrix". If is_owner() is true, then it is an \ref 00214 omatrix_tlate object. 00215 */ 00216 bool is_owner() const { 00217 if (mparent_t::owner==1) return true; 00218 return false; 00219 } 00220 00221 /// Desc 00222 data_t max() const { 00223 data_t mx=*mparent_t::data; 00224 for(size_t i=0;i<mparent_t::size1;i++) { 00225 for(size_t j=0;j<mparent_t::size2;j++) { 00226 if (*(mparent_t::data+i*mparent_t::tda+j)>mx) { 00227 mx=*(mparent_t::data+i*mparent_t::tda+j); 00228 } 00229 } 00230 } 00231 return mx; 00232 } 00233 00234 /// Desc 00235 data_t min() const { 00236 data_t mx=*mparent_t::data; 00237 for(size_t i=0;i<mparent_t::size1;i++) { 00238 for(size_t j=0;j<mparent_t::size2;j++) { 00239 if (*(mparent_t::data+i*mparent_t::tda+j)<mx) { 00240 mx=*(mparent_t::data+i*mparent_t::tda+j); 00241 } 00242 } 00243 } 00244 return mx; 00245 } 00246 00247 /** \brief Return a \c const gsl matrix */ 00248 const mparent_t *get_gsl_matrix_const() const { return this; }; 00249 //@} 00250 00251 #ifndef DOXYGEN_INTERNAL 00252 00253 protected: 00254 00255 /** \brief Desc 00256 */ 00257 omatrix_const_view_tlate() {}; 00258 00259 #endif 00260 00261 }; 00262 00263 /** 00264 \brief A base class for omatrix and omatrix_view 00265 */ 00266 template<class data_t, class mparent_t, class block_t> 00267 class omatrix_base_tlate : 00268 public omatrix_const_view_tlate<data_t,mparent_t,block_t> { 00269 00270 public: 00271 00272 /// \name Copy constructors 00273 //@{ 00274 /// Shallow copy constructor - create a new view of the same matrix 00275 omatrix_base_tlate(const omatrix_base_tlate &v) { 00276 mparent_t::block=0; 00277 mparent_t::data=v.data; 00278 mparent_t::size1=v.size1; 00279 mparent_t::size2=v.size2; 00280 mparent_t::tda=v.tda(); 00281 mparent_t::owner=0; 00282 } 00283 00284 /// Shallow copy constructor - create a new view of the same matrix 00285 omatrix_base_tlate& operator=(const omatrix_base_tlate &v) { 00286 mparent_t::block=0; 00287 mparent_t::data=v.data; 00288 mparent_t::size1=v.size1; 00289 mparent_t::size2=v.size2; 00290 mparent_t::tda=v.tda; 00291 mparent_t::owner=0; 00292 00293 return *this; 00294 } 00295 //@} 00296 00297 ~omatrix_base_tlate() {}; 00298 00299 /// \name Get and set methods 00300 //@{ 00301 /** 00302 \brief Array-like indexing 00303 */ 00304 data_t *operator[](size_t i) { 00305 #if O2SCL_NO_RANGE_CHECK 00306 #else 00307 if (i>=mparent_t::size1) { 00308 O2SCL_ERR((((std::string)"Array index ")+itos(i)+" out of bounds" 00309 +" in omatrix_base_tlate::operator[]. Size: "+ 00310 itos(mparent_t::size1)+ 00311 " (index should be less than size).").c_str(),gsl_eindex); 00312 return (mparent_t::data); 00313 } 00314 #endif 00315 return mparent_t::data+i*mparent_t::tda; 00316 } 00317 00318 /** 00319 \brief Array-like indexing 00320 */ 00321 const data_t *operator[](size_t i) const { 00322 #if O2SCL_NO_RANGE_CHECK 00323 #else 00324 if (i>=mparent_t::size1) { 00325 O2SCL_ERR((((std::string)"Array index ")+itos(i)+" out of bounds" 00326 +" in omatrix_base_tlate::operator[] const. Size: "+ 00327 itos(mparent_t::size1)+ 00328 " (index should be less than size).").c_str(),gsl_eindex); 00329 return (mparent_t::data); 00330 } 00331 #endif 00332 return mparent_t::data+i*mparent_t::tda; 00333 } 00334 00335 /** 00336 \brief Array-like indexing 00337 */ 00338 data_t &operator()(size_t i, size_t j) { 00339 #if O2SCL_NO_RANGE_CHECK 00340 #else 00341 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00342 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00343 ", "+itos(j)+" out of bounds" 00344 +" in omatrix_cx_view_tlate::operator(). Sizes: "+ 00345 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00346 " (indices should be less than sizes).").c_str(), 00347 gsl_eindex); 00348 return *(mparent_t::data); 00349 } 00350 #endif 00351 return *(mparent_t::data+i*mparent_t::tda+j); 00352 } 00353 00354 /** 00355 \brief Array-like indexing 00356 */ 00357 const data_t &operator()(size_t i, size_t j) const { 00358 #if O2SCL_NO_RANGE_CHECK 00359 #else 00360 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00361 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00362 ", "+itos(j)+" out of bounds" 00363 +" in omatrix_cx_view_tlate::operator() const. Sizes: "+ 00364 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00365 " (indices should be less than sizes).").c_str(), 00366 gsl_eindex); 00367 return *(mparent_t::data); 00368 } 00369 #endif 00370 return *(mparent_t::data+i*mparent_t::tda+j); 00371 } 00372 00373 /** \brief Get pointer (with optional range-checking) */ 00374 data_t *get_ptr(size_t i, size_t j) { 00375 #if O2SCL_NO_RANGE_CHECK 00376 #else 00377 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00378 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00379 ", "+itos(j)+" out of bounds" 00380 +" in omatrix_cx_view_tlate::get_ptr(). Sizes: "+ 00381 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00382 " (indices should be less than sizes).").c_str(), 00383 gsl_eindex); 00384 return (mparent_t::data); 00385 } 00386 #endif 00387 return mparent_t::data+i*mparent_t::tda+j; 00388 } 00389 00390 /** \brief Set (with optional range-checking) */ 00391 int set(size_t i, size_t j, data_t val) { 00392 #if O2SCL_NO_RANGE_CHECK 00393 #else 00394 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00395 O2SCL_ERR_RET((((std::string)"Array indices ")+itos(i)+ 00396 ", "+itos(j)+" out of bounds" 00397 +" in omatrix_cx_view_tlate::set(). Sizes: "+ 00398 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00399 " (indices should be less than sizes).").c_str(), 00400 gsl_eindex); 00401 } 00402 #endif 00403 *(mparent_t::data+i*mparent_t::tda+j)=val; 00404 return 0; 00405 } 00406 00407 /** \brief Set all of the value to be the value \c val */ 00408 int set_all(double val) { 00409 for(size_t i=0;i<mparent_t::size1;i++) { 00410 for(size_t j=0;j<mparent_t::size2;j++) { 00411 *(mparent_t::data+i*mparent_t::tda+j)=val; 00412 } 00413 } 00414 return 0; 00415 } 00416 //@} 00417 00418 /// \name Other methods 00419 //@{ 00420 /** \brief Return a gsl matrix */ 00421 mparent_t *get_gsl_matrix() { return this; }; 00422 //@} 00423 00424 /// \name Arithmetic 00425 //@{ 00426 /** \brief operator+= */ 00427 omatrix_base_tlate<data_t,mparent_t,block_t> &operator+= 00428 (const omatrix_base_tlate<data_t,mparent_t,block_t> &x) { 00429 size_t lsize=x.size1; 00430 if (lsize>mparent_t::size1) lsize=mparent_t::size1; 00431 size_t lsize2=x.size2; 00432 if (lsize2>mparent_t::size2) lsize2=mparent_t::size2; 00433 for(size_t i=0;i<lsize;i++) { 00434 for(size_t j=0;j<lsize2;j++) { 00435 (*this)[i][j]+=x[i][j]; 00436 } 00437 } 00438 00439 return *this; 00440 } 00441 00442 /** \brief operator-= */ 00443 omatrix_base_tlate<data_t,mparent_t,block_t> &operator-= 00444 (const omatrix_base_tlate<data_t,mparent_t,block_t> &x) { 00445 size_t lsize=x.size1; 00446 if (lsize>mparent_t::size1) lsize=mparent_t::size1; 00447 size_t lsize2=x.size2; 00448 if (lsize2>mparent_t::size2) lsize2=mparent_t::size2; 00449 for(size_t i=0;i<lsize;i++) { 00450 for(size_t j=0;j<lsize2;j++) { 00451 (*this)[i][j]+=x[i][j]; 00452 } 00453 } 00454 00455 return *this; 00456 } 00457 00458 /** \brief operator+= */ 00459 omatrix_base_tlate<data_t,mparent_t,block_t> 00460 &operator+=(const data_t &y) { 00461 for(size_t i=0;i<mparent_t::size1;i++) { 00462 for(size_t j=0;j<mparent_t::size2;j++) { 00463 (*this)[i][j]+=y; 00464 } 00465 } 00466 00467 return *this; 00468 } 00469 00470 /** \brief operator-= */ 00471 omatrix_base_tlate<data_t,mparent_t,block_t> 00472 &operator-=(const data_t &y) { 00473 for(size_t i=0;i<mparent_t::size1;i++) { 00474 for(size_t j=0;j<mparent_t::size2;j++) { 00475 (*this)[i][j]-=y; 00476 } 00477 } 00478 00479 return *this; 00480 } 00481 00482 /** \brief operator*= */ 00483 omatrix_base_tlate<data_t,mparent_t,block_t> 00484 &operator*=(const data_t &y) { 00485 for(size_t i=0;i<mparent_t::size1;i++) { 00486 for(size_t j=0;j<mparent_t::size2;j++) { 00487 (*this)[i][j]*=y; 00488 } 00489 } 00490 00491 return *this; 00492 } 00493 //@} 00494 00495 #ifndef DOXYGEN_INTERNAL 00496 00497 protected: 00498 00499 /** \brief Desc 00500 */ 00501 omatrix_base_tlate() {}; 00502 00503 #endif 00504 00505 }; 00506 00507 /** 00508 \brief A matrix view of double-precision numbers 00509 00510 This is a matrix view, which views a matrix stored somewhere 00511 else. For a full matrix template class with its own memory, see 00512 \ref omatrix_tlate . The basic matrix classes are built upon 00513 these templates. A matrix of double-precision numbers is an 00514 object of type \ref omatrix . See \ref vecmat_section in the 00515 User's Guide for more information. 00516 */ 00517 template<class data_t, class mparent_t, class block_t> 00518 class omatrix_view_tlate : 00519 public omatrix_base_tlate<data_t,mparent_t,block_t> { 00520 00521 public: 00522 00523 /// \name Copy constructors 00524 //@{ 00525 /// Shallow copy constructor - create a new view of the same matrix 00526 omatrix_view_tlate(const omatrix_view_tlate &v) : 00527 omatrix_base_tlate<data_t,mparent_t,block_t>() { 00528 mparent_t::block=0; 00529 mparent_t::data=v.data; 00530 mparent_t::size1=v.size1; 00531 mparent_t::size2=v.size2; 00532 mparent_t::tda=v.tda(); 00533 mparent_t::owner=0; 00534 } 00535 00536 /// Shallow copy constructor - create a new view of the same matrix 00537 omatrix_view_tlate& operator=(const omatrix_view_tlate &v) { 00538 mparent_t::block=0; 00539 mparent_t::data=v.data; 00540 mparent_t::size1=v.size1; 00541 mparent_t::size2=v.size2; 00542 mparent_t::tda=v.tda; 00543 mparent_t::owner=0; 00544 00545 return *this; 00546 } 00547 00548 /// Shallow copy constructor - create a new view of the same matrix 00549 omatrix_view_tlate(omatrix_base_tlate<data_t,mparent_t,block_t> &v) : 00550 omatrix_base_tlate<data_t,mparent_t,block_t>() { 00551 mparent_t::block=0; 00552 mparent_t::data=v.data; 00553 mparent_t::size1=v.size1; 00554 mparent_t::size2=v.size2; 00555 mparent_t::tda=v.tda(); 00556 mparent_t::owner=0; 00557 } 00558 00559 /// Shallow copy constructor - create a new view of the same matrix 00560 omatrix_view_tlate& operator= 00561 (omatrix_base_tlate<data_t,mparent_t,block_t> &v) { 00562 mparent_t::block=0; 00563 mparent_t::data=v.data; 00564 mparent_t::size1=v.size1; 00565 mparent_t::size2=v.size2; 00566 mparent_t::tda=v.tda; 00567 mparent_t::owner=0; 00568 00569 return *this; 00570 } 00571 //@} 00572 00573 ~omatrix_view_tlate() {}; 00574 00575 /// \name Get and set methods 00576 //@{ 00577 /** 00578 \brief Array-like indexing 00579 */ 00580 data_t *operator[](size_t i) const { 00581 #if O2SCL_NO_RANGE_CHECK 00582 #else 00583 if (i>=mparent_t::size1) { 00584 O2SCL_ERR((((std::string)"Array index ")+itos(i)+" out of bounds" 00585 +" in omatrix_view_tlate::operator[]. Size: "+ 00586 itos(mparent_t::size1)+ 00587 " (index should be less than size).").c_str(),gsl_eindex); 00588 return (mparent_t::data); 00589 } 00590 #endif 00591 return mparent_t::data+i*mparent_t::tda; 00592 } 00593 00594 /** 00595 \brief Array-like indexing 00596 */ 00597 data_t &operator()(size_t i, size_t j) const { 00598 #if O2SCL_NO_RANGE_CHECK 00599 #else 00600 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00601 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00602 ", "+itos(j)+" out of bounds" 00603 +" in omatrix_cx_view_tlate::operator(). Sizes: "+ 00604 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00605 " (indices should be less than sizes).").c_str(), 00606 gsl_eindex); 00607 return *(mparent_t::data); 00608 } 00609 #endif 00610 return *(mparent_t::data+i*mparent_t::tda+j); 00611 } 00612 00613 /** \brief Get pointer (with optional range-checking) */ 00614 data_t *get_ptr(size_t i, size_t j) const { 00615 #if O2SCL_NO_RANGE_CHECK 00616 #else 00617 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00618 O2SCL_ERR((((std::string)"Array indices ")+itos(i)+ 00619 ", "+itos(j)+" out of bounds" 00620 +" in omatrix_cx_view_tlate::get_ptr(). Sizes: "+ 00621 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00622 " (indices should be less than sizes).").c_str(), 00623 gsl_eindex); 00624 return (mparent_t::data); 00625 } 00626 #endif 00627 return mparent_t::data+i*mparent_t::tda+j; 00628 } 00629 00630 /** \brief Set (with optional range-checking) */ 00631 int set(size_t i, size_t j, data_t val) const { 00632 #if O2SCL_NO_RANGE_CHECK 00633 #else 00634 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00635 O2SCL_ERR_RET((((std::string)"Array indices ")+itos(i)+ 00636 ", "+itos(j)+" out of bounds" 00637 +" in omatrix_cx_view_tlate::set(). Sizes: "+ 00638 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00639 " (indices should be less than sizes).").c_str(), 00640 gsl_eindex); 00641 } 00642 #endif 00643 *(mparent_t::data+i*mparent_t::tda+j)=val; 00644 return 0; 00645 } 00646 00647 /** \brief Set all of the value to be the value \c val */ 00648 int set_all(double val) const { 00649 for(size_t i=0;i<mparent_t::size1;i++) { 00650 for(size_t j=0;j<mparent_t::size2;j++) { 00651 *(mparent_t::data+i*mparent_t::tda+j)=val; 00652 } 00653 } 00654 return 0; 00655 } 00656 //@} 00657 00658 #ifndef DOXYGEN_INTERNAL 00659 00660 protected: 00661 00662 /** \brief Desc 00663 */ 00664 omatrix_view_tlate() {}; 00665 00666 #endif 00667 00668 }; 00669 00670 /** 00671 \brief A matrix of double-precision numbers 00672 00673 The basic matrix classes are built upon this template. A matrix 00674 of double-precision numbers is an object of type \ref omatrix , 00675 which is just a <tt>typedef</tt> defined using this class 00676 template. See \ref vecmat_section in the User's Guide for more 00677 information. 00678 */ 00679 template<class data_t, class mparent_t, class vparent_t, class block_t> 00680 class omatrix_tlate : 00681 public omatrix_base_tlate<data_t,mparent_t,block_t> { 00682 public: 00683 00684 /// \name Standard constructor 00685 //@{ 00686 /** \brief Create an omatrix of size \c n with owner as \c true. 00687 */ 00688 omatrix_tlate(size_t r=0, size_t c=0) { 00689 00690 mparent_t::block=0; 00691 mparent_t::data=0; 00692 mparent_t::size1=0; 00693 mparent_t::size2=0; 00694 mparent_t::tda=0; 00695 00696 // This must be set to 1 even if n=0 so that future 00697 // calls to operator= work properly 00698 mparent_t::owner=1; 00699 00700 if (r>0 && c>0) { 00701 mparent_t::block=(block_t *)malloc(sizeof(block_t)); 00702 if (mparent_t::block) { 00703 mparent_t::block->data=(data_t *)malloc(r*c*sizeof(data_t)); 00704 if (mparent_t::block->data) { 00705 mparent_t::block->size=r*c; 00706 mparent_t::data=mparent_t::block->data; 00707 mparent_t::size1=r; 00708 mparent_t::size2=c; 00709 mparent_t::tda=c; 00710 } else { 00711 std::free(mparent_t::block); 00712 O2SCL_ERR("No memory for data in omatrix_tlate constructor", 00713 gsl_enomem); 00714 } 00715 } else { 00716 O2SCL_ERR("No memory for block in omatrix_tlate contructor", 00717 gsl_enomem); 00718 } 00719 } 00720 } 00721 //@} 00722 00723 /// \name Copy constructors 00724 //@{ 00725 /// Deep copy constructor, allocate new space and make a copy 00726 omatrix_tlate(const omatrix_tlate &v) : 00727 omatrix_base_tlate<data_t,mparent_t,block_t>() { 00728 00729 mparent_t::block=0; 00730 mparent_t::data=0; 00731 mparent_t::size1=0; 00732 mparent_t::size2=0; 00733 mparent_t::tda=0; 00734 mparent_t::owner=1; 00735 00736 size_t n=v.size1; 00737 size_t n2=v.size2; 00738 if (n>0 && n2>0) { 00739 mparent_t::block=(block_t *)malloc(sizeof(block_t)); 00740 if (mparent_t::block) { 00741 mparent_t::block->data=(data_t *)malloc(n*n2*sizeof(data_t)); 00742 if (mparent_t::block->data) { 00743 mparent_t::block->size=n*n2; 00744 mparent_t::data=mparent_t::block->data; 00745 mparent_t::size1=n; 00746 mparent_t::size2=n2; 00747 mparent_t::tda=n2; 00748 mparent_t::owner=1; 00749 for(size_t i=0;i<n;i++) { 00750 for(size_t j=0;j<n2;j++) { 00751 *(mparent_t::data+i*mparent_t::tda+j)=v[i][j]; 00752 } 00753 } 00754 } else { 00755 std::free(mparent_t::block); 00756 O2SCL_ERR("No memory for data in omatrix_tlate constructor", 00757 gsl_enomem); 00758 } 00759 } else { 00760 O2SCL_ERR("No memory for block in omatrix_tlate contructor", 00761 gsl_enomem); 00762 } 00763 } else { 00764 mparent_t::size1=0; 00765 mparent_t::size2=0; 00766 mparent_t::tda=0; 00767 } 00768 } 00769 00770 /// Deep copy constructor, allocate new space and make a copy 00771 omatrix_tlate 00772 (const omatrix_const_view_tlate<data_t,mparent_t,block_t> &v) : 00773 omatrix_view_tlate<data_t,mparent_t,block_t>() { 00774 00775 mparent_t::block=0; 00776 mparent_t::data=0; 00777 mparent_t::size1=0; 00778 mparent_t::size2=0; 00779 mparent_t::tda=0; 00780 mparent_t::owner=1; 00781 00782 size_t r=v.rows(); 00783 size_t c=v.cols(); 00784 if (r>0 && c>0) { 00785 mparent_t::block=(block_t *)malloc(sizeof(block_t)); 00786 if (mparent_t::block) { 00787 mparent_t::block->data=(data_t *)malloc(r*c*sizeof(data_t)); 00788 if (mparent_t::block->data) { 00789 mparent_t::block->size=r*c; 00790 mparent_t::data=mparent_t::block->data; 00791 mparent_t::size1=v.size1; 00792 mparent_t::size2=v.size2; 00793 mparent_t::tda=v.tda; 00794 mparent_t::owner=1; 00795 for(size_t i=0;i<r;i++) { 00796 for(size_t j=0;i<c;j++) { 00797 *(mparent_t::data+i*mparent_t::tda+j)=v[i][j]; 00798 } 00799 } 00800 } else { 00801 std::free(mparent_t::block); 00802 O2SCL_ERR("No memory for data in omatrix_tlate constructor", 00803 gsl_enomem); 00804 } 00805 } else { 00806 O2SCL_ERR("No memory for block in omatrix_tlate contructor", 00807 gsl_enomem); 00808 } 00809 } else { 00810 mparent_t::size1=0; 00811 mparent_t::size2=0; 00812 mparent_t::tda=0; 00813 } 00814 } 00815 00816 /** \brief Deep copy constructor, if owner is true, allocate 00817 space and make a new copy, otherwise, just copy into the 00818 view 00819 */ 00820 omatrix_tlate& operator=(const omatrix_tlate &v) { 00821 00822 // Check for self-assignment 00823 if (this==&v) return *this; 00824 00825 mparent_t::block=0; 00826 mparent_t::data=0; 00827 mparent_t::size1=0; 00828 mparent_t::size2=0; 00829 mparent_t::tda=0; 00830 mparent_t::owner=1; 00831 00832 size_t sz=v.rows(); 00833 size_t sz2=v.cols(); 00834 mparent_t::data=0; 00835 if (mparent_t::owner) { 00836 allocate(mparent_t::size1,mparent_t::size2); 00837 } else { 00838 if (mparent_t::size1!=sz || mparent_t::size2!=sz2) { 00839 O2SCL_ERR("Sizes don't match in omatrix_tlate::operator=()", 00840 gsl_ebadlen); 00841 return *this; 00842 } 00843 } 00844 for(size_t i=0;i<sz;i++) { 00845 for(size_t j=0;j<sz2;j++) { 00846 *(mparent_t::data+i*mparent_t::tda+j)=v[i][j]; 00847 } 00848 } 00849 return *this; 00850 } 00851 00852 /** \brief Deep copy constructor, if owner is true, allocate 00853 space and make a new copy, otherwise, just copy into the 00854 view 00855 */ 00856 omatrix_tlate& operator= 00857 (const omatrix_const_view_tlate<data_t,mparent_t,block_t> &v) { 00858 00859 mparent_t::block=0; 00860 mparent_t::data=0; 00861 mparent_t::size1=0; 00862 mparent_t::size2=0; 00863 mparent_t::tda=0; 00864 mparent_t::owner=1; 00865 00866 // Check for self-assignment 00867 if (this==&v) return *this; 00868 00869 size_t size=v.rows(); 00870 size_t size2=v.cols(); 00871 mparent_t::data=0; 00872 if (mparent_t::owner) { 00873 allocate(size,size2); 00874 } else { 00875 if (mparent_t::size!=size || mparent_t::size2!=size2) { 00876 O2SCL_ERR("Sizes don't match in omatrix_tlate::operator=()", 00877 gsl_ebadlen); 00878 return *this; 00879 } 00880 } 00881 for(size_t i=0;i<size;i++) { 00882 for(size_t j=0;j<size2;j++) { 00883 *(mparent_t::data+i*mparent_t::tda+j)=v[i][j]; 00884 } 00885 } 00886 return *this; 00887 } 00888 00889 /** \brief Deep copy from an array of ovectors 00890 */ 00891 omatrix_tlate(size_t n, 00892 ovector_base_tlate<data_t,vparent_t,block_t> ova[]) { 00893 if (n>0) { 00894 size_t n2=ova[0]; 00895 if (n2>0) { 00896 allocate(n,n2); 00897 for(size_t i=0;i<n;i++) { 00898 for(size_t j=0;j<n2;j++) { 00899 (*this)[i][j]=ova[i][j]; 00900 } 00901 } 00902 } 00903 } 00904 } 00905 00906 /** \brief Deep copy from an array of uvectors 00907 */ 00908 omatrix_tlate(size_t n, uvector_base_tlate<data_t> uva[]) { 00909 if (n>0) { 00910 size_t n2=uva[0]; 00911 if (n2>0) { 00912 allocate(n,n2); 00913 for(size_t i=0;i<n;i++) { 00914 for(size_t j=0;j<n2;j++) { 00915 (*this)[i][j]=uva[i][j]; 00916 } 00917 } 00918 } 00919 } 00920 } 00921 00922 /** \brief Deep copy from a C-style 2-d array 00923 */ 00924 omatrix_tlate(size_t n, size_t n2, data_t **csa) { 00925 if (n>0 && n2>0) { 00926 allocate(n,n2); 00927 for(size_t i=0;i<n;i++) { 00928 for(size_t j=0;j<n2;j++) { 00929 (*this)[i][j]=csa[i][j]; 00930 } 00931 } 00932 } 00933 } 00934 00935 //@} 00936 00937 ~omatrix_tlate() { 00938 if (mparent_t::size1>0) { 00939 if (mparent_t::owner==1) { 00940 if (mparent_t::block->size>0) { 00941 std::free(mparent_t::block->data); 00942 } 00943 std::free(mparent_t::block); 00944 mparent_t::size1=0; 00945 mparent_t::size2=0; 00946 } 00947 } 00948 } 00949 00950 /// \name Memory allocation 00951 //@{ 00952 /** 00953 \brief Allocate memory after freeing any memory 00954 presently in use 00955 */ 00956 int allocate(size_t nrows, size_t ncols) { 00957 if (mparent_t::size1>0 || mparent_t::size2>0) free(); 00958 00959 if (nrows>0 && ncols>0) { 00960 mparent_t::block=(block_t *)malloc(sizeof(block_t)); 00961 if (mparent_t::block) { 00962 mparent_t::block->data=(data_t *) 00963 malloc(nrows*ncols*sizeof(data_t)); 00964 if (mparent_t::block->data) { 00965 mparent_t::block->size=nrows*ncols; 00966 mparent_t::data=mparent_t::block->data; 00967 mparent_t::size1=nrows; 00968 mparent_t::size2=ncols; 00969 mparent_t::tda=ncols; 00970 mparent_t::owner=1; 00971 } else { 00972 std::free(mparent_t::block); 00973 O2SCL_ERR_RET("No memory for data in omatrix_tlate::allocate()", 00974 gsl_enomem); 00975 } 00976 } else { 00977 O2SCL_ERR_RET("No memory for block in omatrix_tlate::allocate()", 00978 gsl_enomem); 00979 } 00980 } else { 00981 O2SCL_ERR_RET("Zero size in omatrix::allocate()",gsl_einval); 00982 } 00983 return 0; 00984 } 00985 00986 /** 00987 \brief Free the memory 00988 00989 This function will safely do nothing if used without first 00990 allocating memory or if called multiple times in succession. 00991 */ 00992 int free() { 00993 if (mparent_t::size1>0) { 00994 if (mparent_t::owner==1) { 00995 if (mparent_t::block->size>0) { 00996 std::free(mparent_t::block->data); 00997 } 00998 std::free(mparent_t::block); 00999 } 01000 mparent_t::size1=0; 01001 mparent_t::size2=0; 01002 mparent_t::tda=0; 01003 } 01004 return 0; 01005 } 01006 //@} 01007 01008 /// \name Other methods 01009 //@{ 01010 /// \brief Compute the transpose (even if matrix is not square) 01011 omatrix_tlate<data_t,mparent_t,vparent_t,block_t> transpose() { 01012 omatrix_tlate<data_t,mparent_t,vparent_t,block_t> 01013 result(mparent_t::size2,mparent_t::size1); 01014 for(size_t i=0;i<mparent_t::size1;i++) { 01015 for(size_t j=0;j<mparent_t::size2;j++) { 01016 result[j][i]=(*this)[i][j]; 01017 } 01018 } 01019 } 01020 //@} 01021 01022 }; 01023 01024 #ifdef O2SCL_NEVER_DEFINED 01025 } 01026 { 01027 #endif 01028 01029 /** 01030 \brief Create a matrix from an array 01031 */ 01032 template<class data_t, class mparent_t, class block_t> 01033 class omatrix_array_tlate : 01034 public omatrix_view_tlate<data_t,mparent_t,block_t> { 01035 01036 public: 01037 01038 /** \brief Create a vector from \c dat with size \c n */ 01039 omatrix_array_tlate(data_t *dat, size_t ltda, 01040 size_t sz1, size_t sz2) { 01041 if (sz1==0 || sz2==0) { 01042 O2SCL_ERR("Failed to create omatrix from array, insufficient size.", 01043 gsl_einval); 01044 mparent_t::block=0; 01045 mparent_t::data=0; 01046 mparent_t::size1=0; 01047 mparent_t::size2=0; 01048 mparent_t::tda=0; 01049 mparent_t::owner=0; 01050 } 01051 mparent_t::block=0; 01052 mparent_t::data=dat; 01053 mparent_t::size1=sz1; 01054 mparent_t::size2=sz2; 01055 mparent_t::tda=ltda; 01056 mparent_t::owner=0; 01057 } 01058 }; 01059 01060 /** \brief Create a vector from a row of a matrix 01061 */ 01062 template<class data_t, class mparent_t, class vparent_t, class block_t> 01063 class omatrix_row_tlate : 01064 public ovector_view_tlate<data_t,vparent_t,block_t> { 01065 public: 01066 /** \brief Create a vector from row \c i of matrix \c m */ 01067 omatrix_row_tlate 01068 (omatrix_base_tlate<data_t,mparent_t,block_t> &m, 01069 size_t i) { 01070 if (i<m.size1) { 01071 vparent_t::size=m.size2; 01072 vparent_t::stride=1; 01073 vparent_t::data=m.data+m.tda()*i; 01074 vparent_t::owner=0; 01075 vparent_t::block=0; 01076 } else { 01077 O2SCL_ERR("Invalid row in omatrix_row_tlate constructor.", 01078 gsl_einval); 01079 vparent_t::size=0; 01080 vparent_t::stride=0; 01081 vparent_t::data=0; 01082 vparent_t::owner=0; 01083 vparent_t::block=0; 01084 } 01085 } 01086 }; 01087 01088 /** \brief Create a const vector from a row of a matrix 01089 */ 01090 template<class data_t, class mparent_t, class vparent_t, class block_t> 01091 class omatrix_const_row_tlate : 01092 public ovector_const_view_tlate<data_t,vparent_t,block_t> { 01093 public: 01094 /** \brief Create a vector from row \c i of matrix \c m */ 01095 omatrix_const_row_tlate 01096 (const omatrix_base_tlate<data_t,mparent_t,block_t> &m, 01097 size_t i) { 01098 if (i<m.size1) { 01099 vparent_t::size=m.size2; 01100 vparent_t::stride=1; 01101 vparent_t::data=m.data+m.tda*i; 01102 vparent_t::owner=0; 01103 vparent_t::block=0; 01104 } else { 01105 O2SCL_ERR("Invalid row in omatrix_const_row_tlate constructor.", 01106 gsl_einval); 01107 vparent_t::size=0; 01108 vparent_t::stride=0; 01109 vparent_t::data=0; 01110 vparent_t::owner=0; 01111 vparent_t::block=0; 01112 } 01113 } 01114 }; 01115 01116 /** \brief Create a vector from a column of a matrix 01117 01118 \future Also do a umatrix constructor since we 01119 can't do that with uvectors 01120 */ 01121 template<class data_t, class mparent_t, class vparent_t, class block_t> 01122 class omatrix_col_tlate : 01123 public ovector_view_tlate<data_t,vparent_t,block_t> { 01124 public: 01125 /** \brief Create a vector from col \c i of matrix \c m */ 01126 omatrix_col_tlate(omatrix_base_tlate<data_t,mparent_t,block_t> &m, 01127 size_t i) { 01128 if (i<m.size2) { 01129 vparent_t::size=m.size1; 01130 vparent_t::stride=m.tda(); 01131 vparent_t::data=m.data+i; 01132 vparent_t::owner=0; 01133 vparent_t::block=0; 01134 } else { 01135 O2SCL_ERR("Invalid column in omatrix_col_tlate constructor.", 01136 gsl_einval); 01137 vparent_t::size=0; 01138 vparent_t::stride=0; 01139 vparent_t::data=0; 01140 vparent_t::owner=0; 01141 vparent_t::block=0; 01142 } 01143 } 01144 }; 01145 01146 /** \brief Create a const vector from a column of a matrix 01147 */ 01148 template<class data_t, class mparent_t, class vparent_t, class block_t> 01149 class omatrix_const_col_tlate : 01150 public ovector_const_view_tlate<data_t,vparent_t,block_t> { 01151 public: 01152 /** \brief Create a vector from col \c i of matrix \c m */ 01153 omatrix_const_col_tlate 01154 (const omatrix_base_tlate<data_t,mparent_t,block_t> &m, size_t i) { 01155 if (i<m.size2) { 01156 vparent_t::size=m.size1; 01157 vparent_t::stride=m.tda(); 01158 vparent_t::data=m.data+i; 01159 vparent_t::owner=0; 01160 vparent_t::block=0; 01161 } else { 01162 O2SCL_ERR("Invalid column in omatrix_const_col_tlate constructor.", 01163 gsl_einval); 01164 vparent_t::size=0; 01165 vparent_t::stride=0; 01166 vparent_t::data=0; 01167 vparent_t::owner=0; 01168 vparent_t::block=0; 01169 } 01170 } 01171 }; 01172 01173 /** \brief Create a vector from the main diagonal 01174 */ 01175 template<class data_t, class mparent_t, class vparent_t, class block_t> 01176 class omatrix_diag_tlate : 01177 public ovector_view_tlate<data_t,vparent_t,block_t> { 01178 public: 01179 /** \brief Create a vector of the diagonal of matrix \c m */ 01180 omatrix_diag_tlate(omatrix_view_tlate<data_t,mparent_t,block_t> &m) { 01181 01182 if (m.size2<m.size1) vparent_t::size=m.size2; 01183 else vparent_t::size=m.size1; 01184 vparent_t::stride=m.tda()+1; 01185 vparent_t::data=m.data; 01186 vparent_t::owner=0; 01187 vparent_t::block=0; 01188 } 01189 }; 01190 01191 /** \brief Create a vector from the main diagonal 01192 */ 01193 template<class data_t, class mparent_t, class vparent_t, class block_t> 01194 class omatrix_const_diag_tlate : 01195 public ovector_const_view_tlate<data_t,vparent_t,block_t> { 01196 public: 01197 /** \brief Create a vector of the diagonal of matrix \c m */ 01198 omatrix_const_diag_tlate 01199 (const omatrix_view_tlate<data_t,mparent_t,block_t> &m) { 01200 01201 if (m.size2<m.size1) vparent_t::size=m.size2; 01202 else vparent_t::size=m.size1; 01203 vparent_t::stride=m.tda()+1; 01204 vparent_t::data=m.data; 01205 vparent_t::owner=0; 01206 vparent_t::block=0; 01207 } 01208 }; 01209 01210 01211 /// omatrix typedef 01212 typedef omatrix_tlate<double,gsl_matrix,gsl_vector_norm,gsl_block> omatrix; 01213 /// omatrix_view typedef 01214 typedef omatrix_view_tlate<double,gsl_matrix,gsl_block> omatrix_view; 01215 /// omatrix_const_view typedef 01216 typedef omatrix_const_view_tlate<double,gsl_matrix,gsl_block> 01217 omatrix_const_view; 01218 /// omatrix_base typedef 01219 typedef omatrix_base_tlate<double,gsl_matrix,gsl_block> omatrix_base; 01220 /// omatrix_row typedef 01221 typedef omatrix_row_tlate<double,gsl_matrix,gsl_vector_norm,gsl_block> 01222 omatrix_row; 01223 /// omatrix_col typedef 01224 typedef omatrix_col_tlate<double,gsl_matrix,gsl_vector_norm,gsl_block> 01225 omatrix_col; 01226 /// omatrix_const_row typedef 01227 typedef omatrix_const_row_tlate<double,gsl_matrix,gsl_vector_norm,gsl_block> 01228 omatrix_const_row; 01229 /// omatrix_const_col typedef 01230 typedef omatrix_const_col_tlate<double,gsl_matrix,gsl_vector_norm,gsl_block> 01231 omatrix_const_col; 01232 /// omatrix_diag typedef 01233 typedef omatrix_diag_tlate<double,gsl_matrix,gsl_vector_norm,gsl_block> 01234 omatrix_diag; 01235 /// omatrix_const_diag typedef 01236 typedef omatrix_const_diag_tlate<double,gsl_matrix,gsl_vector_norm,gsl_block> 01237 omatrix_const_diag; 01238 /// omatrix_array typedef 01239 typedef omatrix_array_tlate<double,gsl_matrix,gsl_block> 01240 omatrix_array; 01241 01242 /// omatrix_int typedef 01243 typedef omatrix_tlate<int,gsl_matrix_int,gsl_vector_int,gsl_block_int> 01244 omatrix_int; 01245 /// omatrix_int_view typedef 01246 typedef omatrix_view_tlate<int,gsl_matrix_int,gsl_block_int> 01247 omatrix_int_view; 01248 /// omatrix_int_base typedef 01249 typedef omatrix_base_tlate<int,gsl_matrix_int,gsl_block_int> 01250 omatrix_int_base; 01251 /// omatrix_int_const_view typedef 01252 typedef omatrix_const_view_tlate<int,gsl_matrix_int,gsl_block_int> 01253 omatrix_int_const_view; 01254 /// omatrix_int_row typedef 01255 typedef omatrix_row_tlate<int,gsl_matrix_int,gsl_vector_int, 01256 gsl_block_int> omatrix_int_row; 01257 /// omatrix_int_col typedef 01258 typedef omatrix_col_tlate<int,gsl_matrix_int,gsl_vector_int, 01259 gsl_block_int> omatrix_int_col; 01260 /// omatrix_int_const_row typedef 01261 typedef omatrix_const_row_tlate<int,gsl_matrix_int,gsl_vector_int, 01262 gsl_block_int> omatrix_int_const_row; 01263 /// omatrix_int_const_col typedef 01264 typedef omatrix_const_col_tlate<int,gsl_matrix_int,gsl_vector_int, 01265 gsl_block_int> omatrix_int_const_col; 01266 /// omatrix_int_diag typedef 01267 typedef omatrix_diag_tlate<int,gsl_matrix_int,gsl_vector_int, 01268 gsl_block_int> omatrix_int_diag; 01269 /// omatrix_int_const_diag typedef 01270 typedef omatrix_const_diag_tlate<int,gsl_matrix_int,gsl_vector_int, 01271 gsl_block_int> omatrix_int_const_diag; 01272 /// omatrix_int_array typedef 01273 typedef omatrix_array_tlate<int,gsl_matrix_int,gsl_block_int> 01274 omatrix_int_array; 01275 01276 /** 01277 \brief A operator for output of omatrix objects 01278 01279 This outputs all of the matrix elements. Each row is output with 01280 an endline character at the end of each row. Positive values are 01281 preceeded by an extra space. A 2x2 example: 01282 \verbatim 01283 -3.751935e-05 -6.785864e-04 01284 -6.785864e-04 1.631984e-02 01285 \endverbatim 01286 01287 The function \c gsl_ieee_double_to_rep() is used to determine 01288 the sign of a number, so that "-0.0" as distinct from "+0.0" 01289 is handled correctly. 01290 01291 \comment 01292 I had originally thought about removing this function, since 01293 it's superceded by matrix_out, but it's simpler to be able 01294 to just use the << operator so I keep it for now. 01295 \endcomment 01296 01297 \future This assumes that scientific mode is on and showpos 01298 is off. It'd be nice to fix this. 01299 01300 */ 01301 template<class data_t, class parent_t, class block_t> 01302 std::ostream &operator<< 01303 (std::ostream &os, 01304 const omatrix_const_view_tlate<data_t,parent_t,block_t> &v) { 01305 size_t i; 01306 gsl_ieee_double_rep r; 01307 for(i=0;i<v.rows()-1;i++) { 01308 for(size_t j=0;j<v.cols();j++) { 01309 gsl_ieee_double_to_rep(&(v[i][j]), &r); 01310 if (r.sign==1) os << v[i][j] << ' '; 01311 else os << ' ' << v[i][j] << ' '; 01312 } 01313 os << '\n'; 01314 } 01315 i=v.rows()-1; 01316 if (i>0) { 01317 for(size_t j=0;j<v.cols();j++) { 01318 gsl_ieee_double_to_rep(&(v[i][j]), &r); 01319 if (r.sign==1) os << v[i][j] << ' '; 01320 else os << ' ' << v[i][j] << ' '; 01321 } 01322 } 01323 return os; 01324 } 01325 01326 /** \brief A simple class to provide an \c allocate() function 01327 for \ref omatrix 01328 01329 */ 01330 class omatrix_alloc { 01331 public: 01332 /// Allocate \c v for \c i elements 01333 void allocate(omatrix &o, size_t i, size_t j) { o.allocate(i,j); } 01334 /// Free memory 01335 void free(omatrix &o, size_t i) { o.free(); } 01336 }; 01337 01338 #ifdef DOXGYENP2 01339 01340 /** \brief A matrix where the memory allocation is performed in 01341 the constructor 01342 */ 01343 template<size_t N, size_t M> class ofmatrix : public omatrix_tlate { 01344 01345 public: 01346 01347 /// Create a matrix 01348 ofmatrix(); 01349 }; 01350 01351 #else 01352 01353 template<size_t N, size_t M> class ofmatrix : 01354 public omatrix_tlate<double,gsl_matrix,gsl_vector,gsl_block> 01355 { 01356 public: 01357 ofmatrix() : omatrix_tlate<double,gsl_matrix,gsl_vector, 01358 gsl_block>(N,M) { 01359 } 01360 }; 01361 01362 #endif 01363 01364 class xmatrix : public omatrix { 01365 public: 01366 xmatrix(size_t r=0, size_t c=0) : omatrix(r,c) { 01367 } 01368 01369 ovector_array operator[](size_t i) { 01370 #if O2SCL_NO_RANGE_CHECK 01371 #else 01372 if (i>=gsl_matrix::size1) { 01373 O2SCL_ERR((((std::string)"Array index ")+itos(i)+" out of bounds" 01374 +" in omatrix_view_tlate::operator[]. Size: "+ 01375 itos(gsl_matrix::size1)+ 01376 " (index should be less than size).").c_str(),gsl_eindex); 01377 } 01378 #endif 01379 return ovector_array(gsl_matrix::size2,gsl_matrix::data+ 01380 i*gsl_matrix::tda); 01381 } 01382 01383 const ovector_const_array operator[](size_t i) const { 01384 #if O2SCL_NO_RANGE_CHECK 01385 #else 01386 if (i>=gsl_matrix::size1) { 01387 O2SCL_ERR((((std::string)"Array index ")+itos(i)+" out of bounds" 01388 +" in omatrix_view_tlate::operator[]. Size: "+ 01389 itos(gsl_matrix::size1)+ 01390 " (index should be less than size).").c_str(),gsl_eindex); 01391 } 01392 #endif 01393 return ovector_const_array(gsl_matrix::size2,gsl_matrix::data+ 01394 i*gsl_matrix::tda); 01395 } 01396 }; 01397 01398 #ifndef DOXYGENP 01399 } 01400 #endif 01401 01402 #endif 01403
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page