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