00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 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 00030 #include <iostream> 00031 #include <cstdlib> 00032 #include <string> 00033 #include <fstream> 00034 #include <sstream> 00035 00036 #include <gsl/gsl_matrix.h> 00037 #include <gsl/gsl_ieee_utils.h> 00038 00039 #include <o2scl/err_hnd.h> 00040 #include <o2scl/ovector_tlate.h> 00041 #include <o2scl/array.h> 00042 00043 #ifndef DOXYGENP 00044 namespace o2scl { 00045 #endif 00046 00047 /** 00048 \brief A matrix view of double-precision numbers 00049 00050 \future This class isn't sufficiently general for some 00051 applications, such as sub-matrices of higher-dimensional 00052 structures. It might be nice to create a more general class with 00053 a "stride" and a "tda". 00054 00055 \future The \c xmatrix class demonstrates how operator[] could 00056 return an ovector_array object and thus provide better 00057 bounds-checking. This would demand including a new parameter 00058 in omatrix_view_tlate which contains the vector type. 00059 00060 */ 00061 template<class data_t, class mparent_t, class block_t> 00062 class omatrix_view_tlate : 00063 public mparent_t { 00064 00065 public: 00066 00067 /// \name Copy constructors 00068 //@{ 00069 /// Shallow copy constructor - create a new view of the same matrix 00070 omatrix_view_tlate(const omatrix_view_tlate &v) { 00071 mparent_t::block=0; 00072 mparent_t::data=v.data; 00073 mparent_t::size1=v.size1; 00074 mparent_t::size2=v.size2; 00075 mparent_t::tda=v.tda(); 00076 mparent_t::owner=0; 00077 } 00078 00079 /// Shallow copy constructor - create a new view of the same matrix 00080 omatrix_view_tlate& operator=(const omatrix_view_tlate &v) { 00081 mparent_t::block=0; 00082 mparent_t::data=v.data; 00083 mparent_t::size1=v.size1; 00084 mparent_t::size2=v.size2; 00085 mparent_t::tda=v.tda; 00086 mparent_t::owner=0; 00087 00088 return *this; 00089 } 00090 //@} 00091 00092 ~omatrix_view_tlate() {}; 00093 00094 /// \name Get and set methods 00095 //@{ 00096 /** 00097 \brief Array-like indexing 00098 */ 00099 data_t *operator[](size_t i) { 00100 #if O2SCL_NO_RANGE_CHECK 00101 #else 00102 if (i>=mparent_t::size1) { 00103 set_err((((std::string)"Array index ")+itos(i)+" out of bounds" 00104 +" in omatrix_view_tlate::operator[]. Size: "+ 00105 itos(mparent_t::size1)+ 00106 " (index should be less than size).").c_str(),gsl_index); 00107 return (mparent_t::data); 00108 } 00109 #endif 00110 return mparent_t::data+i*mparent_t::tda; 00111 } 00112 00113 /** 00114 \brief Array-like indexing 00115 */ 00116 const data_t *operator[](size_t i) const { 00117 #if O2SCL_NO_RANGE_CHECK 00118 #else 00119 if (i>=mparent_t::size1) { 00120 set_err((((std::string)"Array index ")+itos(i)+" out of bounds" 00121 +" in omatrix_view_tlate::operator[] const. Size: "+ 00122 itos(mparent_t::size1)+ 00123 " (index should be less than size).").c_str(),gsl_index); 00124 return (mparent_t::data); 00125 } 00126 #endif 00127 return mparent_t::data+i*mparent_t::tda; 00128 } 00129 00130 /** 00131 \brief Array-like indexing 00132 */ 00133 data_t &operator()(size_t i, size_t j) { 00134 #if O2SCL_NO_RANGE_CHECK 00135 #else 00136 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00137 set_err((((std::string)"Array indices ")+itos(i)+ 00138 ", "+itos(j)+" out of bounds" 00139 +" in omatrix_cx_view_tlate::operator(). Sizes: "+ 00140 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00141 " (indices should be less than sizes).").c_str(), 00142 gsl_index); 00143 return *(mparent_t::data); 00144 } 00145 #endif 00146 return *(mparent_t::data+i*mparent_t::tda+j); 00147 } 00148 00149 /** 00150 \brief Array-like indexing 00151 */ 00152 const data_t &operator()(size_t i, size_t j) const { 00153 #if O2SCL_NO_RANGE_CHECK 00154 #else 00155 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00156 set_err((((std::string)"Array indices ")+itos(i)+ 00157 ", "+itos(j)+" out of bounds" 00158 +" in omatrix_cx_view_tlate::operator() const. Sizes: "+ 00159 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00160 " (indices should be less than sizes).").c_str(), 00161 gsl_index); 00162 return *(mparent_t::data); 00163 } 00164 #endif 00165 return *(mparent_t::data+i*mparent_t::tda+j); 00166 } 00167 00168 /** \brief Get (with optional range-checking) */ 00169 data_t get(size_t i, size_t j) const { 00170 #if O2SCL_NO_RANGE_CHECK 00171 #else 00172 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00173 set_err((((std::string)"Array indices ")+itos(i)+ 00174 ", "+itos(j)+" out of bounds" 00175 +" in omatrix_cx_view_tlate::get() const. Sizes: "+ 00176 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00177 " (indices should be less than sizes).").c_str(), 00178 gsl_index); 00179 return *(mparent_t::data); 00180 } 00181 #endif 00182 return *(mparent_t::data+i*mparent_t::tda+j); 00183 } 00184 00185 /** \brief Get pointer (with optional range-checking) */ 00186 data_t *get_ptr(size_t i, size_t j) { 00187 #if O2SCL_NO_RANGE_CHECK 00188 #else 00189 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00190 set_err((((std::string)"Array indices ")+itos(i)+ 00191 ", "+itos(j)+" out of bounds" 00192 +" in omatrix_cx_view_tlate::get_ptr(). Sizes: "+ 00193 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00194 " (indices should be less than sizes).").c_str(), 00195 gsl_index); 00196 return (mparent_t::data); 00197 } 00198 #endif 00199 return mparent_t::data+i*mparent_t::tda+j; 00200 } 00201 00202 /** \brief Get pointer (with optional range-checking) */ 00203 const data_t *get_const_ptr(size_t i, size_t j) const { 00204 #if O2SCL_NO_RANGE_CHECK 00205 #else 00206 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00207 set_err((((std::string)"Array indices ")+itos(i)+ 00208 ", "+itos(j)+" out of bounds" 00209 +" in omatrix_cx_view_tlate::get_const_ptr(). Sizes: "+ 00210 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00211 " (indices should be less than sizes).").c_str(), 00212 gsl_index); 00213 return (mparent_t::data); 00214 } 00215 #endif 00216 return (const data_t *)(mparent_t::data+i*mparent_t::tda+j); 00217 } 00218 00219 /** \brief Set (with optional range-checking) */ 00220 int set(size_t i, size_t j, data_t val) { 00221 #if O2SCL_NO_RANGE_CHECK 00222 #else 00223 if (i>=mparent_t::size1 || j>=mparent_t::size2) { 00224 set_err_ret((((std::string)"Array indices ")+itos(i)+ 00225 ", "+itos(j)+" out of bounds" 00226 +" in omatrix_cx_view_tlate::set(). Sizes: "+ 00227 itos(mparent_t::size1)+","+itos(mparent_t::size2)+ 00228 " (indices should be less than sizes).").c_str(), 00229 gsl_index); 00230 } 00231 #endif 00232 *(mparent_t::data+i*mparent_t::tda+j)=val; 00233 return 0; 00234 } 00235 00236 /** \brief Set all of the value to be the value \c val */ 00237 int set_all(double val) { 00238 for(size_t i=0;i<mparent_t::size1;i++) { 00239 for(size_t j=0;j<mparent_t::size2;j++) { 00240 *(mparent_t::data+i*mparent_t::tda+j)=val; 00241 } 00242 } 00243 return 0; 00244 } 00245 00246 /** 00247 \brief Method to return number of rows 00248 00249 If no memory has been allocated, this will quietly 00250 return zero. 00251 */ 00252 size_t rows() const { 00253 return mparent_t::size1; 00254 } 00255 00256 /** 00257 \brief Method to return number of columns 00258 00259 If no memory has been allocated, this will quietly 00260 return zero. 00261 */ 00262 size_t cols() const { 00263 return mparent_t::size2; 00264 } 00265 00266 /** 00267 \brief Method to return matrix tda 00268 00269 If no memory has been allocated, this will quietly 00270 return zero. 00271 */ 00272 size_t tda() const { 00273 return mparent_t::tda; 00274 } 00275 //@} 00276 00277 /// \name Other methods 00278 //@{ 00279 /// Return true if this object owns the data it refers to 00280 bool is_owner() const { 00281 if (mparent_t::owner==1) return true; 00282 return false; 00283 } 00284 00285 /** \brief Return a gsl matrix */ 00286 mparent_t *get_gsl_matrix() { return this; }; 00287 00288 /** \brief Return a \c const gsl matrix */ 00289 const mparent_t *get_gsl_matrix_const() const { return this; }; 00290 //@} 00291 00292 /// \name Arithmetic 00293 //@{ 00294 /** \brief operator+= */ 00295 omatrix_view_tlate<data_t,mparent_t,block_t> &operator+= 00296 (const omatrix_view_tlate<data_t,mparent_t,block_t> &x) { 00297 size_t lsize=x.size1; 00298 if (lsize>mparent_t::size1) lsize=mparent_t::size1; 00299 size_t lsize2=x.size2; 00300 if (lsize2>mparent_t::size2) lsize2=mparent_t::size2; 00301 for(size_t i=0;i<lsize;i++) { 00302 for(size_t j=0;j<lsize2;j++) { 00303 (*this)[i][j]+=x[i][j]; 00304 } 00305 } 00306 00307 return *this; 00308 } 00309 00310 /** \brief operator-= */ 00311 omatrix_view_tlate<data_t,mparent_t,block_t> &operator-= 00312 (const omatrix_view_tlate<data_t,mparent_t,block_t> &x) { 00313 size_t lsize=x.size1; 00314 if (lsize>mparent_t::size1) lsize=mparent_t::size1; 00315 size_t lsize2=x.size2; 00316 if (lsize2>mparent_t::size2) lsize2=mparent_t::size2; 00317 for(size_t i=0;i<lsize;i++) { 00318 for(size_t j=0;j<lsize2;j++) { 00319 (*this)[i][j]+=x[i][j]; 00320 } 00321 } 00322 00323 return *this; 00324 } 00325 00326 /** \brief operator+= */ 00327 omatrix_view_tlate<data_t,mparent_t,block_t> 00328 &operator+=(const data_t &y) { 00329 for(size_t i=0;i<mparent_t::size1;i++) { 00330 for(size_t j=0;j<mparent_t::size2;j++) { 00331 (*this)[i][j]+=y; 00332 } 00333 } 00334 00335 return *this; 00336 } 00337 00338 /** \brief operator-= */ 00339 omatrix_view_tlate<data_t,mparent_t,block_t> 00340 &operator-=(const data_t &y) { 00341 for(size_t i=0;i<mparent_t::size1;i++) { 00342 for(size_t j=0;j<mparent_t::size2;j++) { 00343 (*this)[i][j]-=y; 00344 } 00345 } 00346 00347 return *this; 00348 } 00349 00350 /** \brief operator*= */ 00351 omatrix_view_tlate<data_t,mparent_t,block_t> 00352 &operator*=(const data_t &y) { 00353 for(size_t i=0;i<mparent_t::size1;i++) { 00354 for(size_t j=0;j<mparent_t::size2;j++) { 00355 (*this)[i][j]*=y; 00356 } 00357 } 00358 00359 return *this; 00360 } 00361 //@} 00362 00363 #ifndef DOXYGEN_INTERNAL 00364 00365 protected: 00366 00367 /** \brief Empty constructor provided for use by 00368 omatrix_tlate(const omatrix_tlate &v) 00369 */ 00370 omatrix_view_tlate() {}; 00371 00372 #endif 00373 00374 }; 00375 00376 /** 00377 \brief A matrix of double-precision numbers 00378 00379 */ 00380 template<class data_t, class mparent_t, class vparent_t, class block_t> 00381 class omatrix_tlate : 00382 public omatrix_view_tlate<data_t,mparent_t,block_t> { 00383 public: 00384 00385 /// \name Standard constructor 00386 //@{ 00387 /** \brief Create an omatrix of size \c n with owner as \c true. 00388 */ 00389 omatrix_tlate(size_t r=0, size_t c=0) { 00390 00391 mparent_t::data=0; 00392 mparent_t::size1=0; 00393 mparent_t::size2=0; 00394 mparent_t::tda=0; 00395 00396 // This must be set to 1 even if n=0 so that future 00397 // calls to operator= work properly 00398 mparent_t::owner=1; 00399 00400 if (r>0 && c>0) { 00401 mparent_t::block=(block_t *)malloc(sizeof(block_t)); 00402 if (mparent_t::block) { 00403 mparent_t::block->data=(data_t *)malloc(r*c*sizeof(data_t)); 00404 if (mparent_t::block->data) { 00405 mparent_t::block->size=r*c; 00406 mparent_t::data=mparent_t::block->data; 00407 mparent_t::size1=r; 00408 mparent_t::size2=c; 00409 mparent_t::tda=c; 00410 } else { 00411 std::free(mparent_t::block); 00412 set_err("No memory for data in omatrix_tlate constructor", 00413 gsl_enomem); 00414 } 00415 } else { 00416 set_err("No memory for block in omatrix_tlate contructor", 00417 gsl_enomem); 00418 } 00419 } 00420 } 00421 //@} 00422 00423 /// \name Copy constructors 00424 //@{ 00425 /// Deep copy constructor, allocate new space and make a copy 00426 omatrix_tlate(const omatrix_tlate &v) : 00427 omatrix_view_tlate<data_t,mparent_t,block_t>() { 00428 size_t n=v.size1; 00429 size_t n2=v.size2; 00430 mparent_t::data=0; 00431 if (n>0 && n2>0) { 00432 mparent_t::block=(block_t *)malloc(sizeof(block_t)); 00433 if (mparent_t::block) { 00434 mparent_t::block->data=(data_t *)malloc(n*n2*sizeof(data_t)); 00435 if (mparent_t::block->data) { 00436 mparent_t::block->size=n*n2; 00437 mparent_t::data=mparent_t::block->data; 00438 mparent_t::size1=n; 00439 mparent_t::size2=n2; 00440 mparent_t::tda=n2; 00441 mparent_t::owner=1; 00442 for(size_t i=0;i<n;i++) { 00443 for(size_t j=0;j<n2;j++) { 00444 *(mparent_t::data+i*mparent_t::tda+j)=v[i][j]; 00445 } 00446 } 00447 } else { 00448 std::free(mparent_t::block); 00449 set_err("No memory for data in omatrix_tlate constructor", 00450 gsl_enomem); 00451 } 00452 } else { 00453 set_err("No memory for block in omatrix_tlate contructor", 00454 gsl_enomem); 00455 } 00456 } else { 00457 mparent_t::size1=0; 00458 mparent_t::size2=0; 00459 mparent_t::tda=0; 00460 } 00461 } 00462 00463 /// Deep copy constructor, allocate new space and make a copy 00464 omatrix_tlate 00465 (const omatrix_view_tlate<data_t,mparent_t,block_t> &v) : 00466 omatrix_view_tlate<data_t,mparent_t,block_t>() { 00467 size_t r=v.rows(); 00468 size_t c=v.cols(); 00469 mparent_t::data=0; 00470 if (r>0 && c>0) { 00471 mparent_t::block=(block_t *)malloc(sizeof(block_t)); 00472 if (mparent_t::block) { 00473 mparent_t::block->data=(data_t *)malloc(r*c*sizeof(data_t)); 00474 if (mparent_t::block->data) { 00475 mparent_t::block->size=r*c; 00476 mparent_t::data=mparent_t::block->data; 00477 mparent_t::size1=v.size1; 00478 mparent_t::size2=v.size2; 00479 mparent_t::tda=v.tda; 00480 mparent_t::owner=1; 00481 for(size_t i=0;i<r;i++) { 00482 for(size_t j=0;i<c;j++) { 00483 *(mparent_t::data+i*mparent_t::tda+j)=v[i][j]; 00484 } 00485 } 00486 } else { 00487 std::free(mparent_t::block); 00488 set_err("No memory for data in omatrix_tlate constructor", 00489 gsl_enomem); 00490 } 00491 } else { 00492 set_err("No memory for block in omatrix_tlate contructor", 00493 gsl_enomem); 00494 } 00495 } else { 00496 mparent_t::size1=0; 00497 mparent_t::size2=0; 00498 mparent_t::tda=0; 00499 } 00500 } 00501 00502 /** \brief Deep copy constructor, if owner is true, allocate 00503 space and make a new copy, otherwise, just copy into the 00504 view 00505 */ 00506 omatrix_tlate& operator=(const omatrix_tlate &v) { 00507 size_t size=v.size1; 00508 size_t size2=v.size2; 00509 mparent_t::data=0; 00510 if (mparent_t::owner) { 00511 allocate(mparent_t::size,mparent_t::size2); 00512 } else { 00513 if (mparent_t::size1!=size || mparent_t::size2!=size2) { 00514 set_err("Sizes don't match in omatrix_tlate::operator=()", 00515 gsl_ebadlen); 00516 return *this; 00517 } 00518 } 00519 for(size_t i=0;i<size;i++) { 00520 for(size_t j=0;j<size2;j++) { 00521 *(mparent_t::data+i*mparent_t::tda+j)=v[i][j]; 00522 } 00523 } 00524 return *this; 00525 } 00526 00527 /** \brief Deep copy constructor, if owner is true, allocate 00528 space and make a new copy, otherwise, just copy into the 00529 view 00530 */ 00531 omatrix_tlate& operator= 00532 (const omatrix_view_tlate<data_t,mparent_t,block_t> &v) { 00533 size_t size=v.rows(); 00534 size_t size2=v.cols(); 00535 mparent_t::data=0; 00536 if (mparent_t::owner) { 00537 allocate(size,size2); 00538 } else { 00539 if (mparent_t::size!=size || mparent_t::size2!=size2) { 00540 set_err("Sizes don't match in omatrix_tlate::operator=()", 00541 gsl_ebadlen); 00542 return *this; 00543 } 00544 } 00545 for(size_t i=0;i<size;i++) { 00546 for(size_t j=0;j<size2;j++) { 00547 *(mparent_t::data+i*mparent_t::tda+j)=v[i][j]; 00548 } 00549 } 00550 return *this; 00551 } 00552 00553 /** \brief Deep copy from an array of ovectors 00554 */ 00555 omatrix_tlate(size_t n, 00556 ovector_view_tlate<data_t,vparent_t,block_t> ova[]) { 00557 if (n>0) { 00558 size_t n2=ova[0]; 00559 if (n2>0) { 00560 allocate(n,n2); 00561 for(size_t i=0;i<n;i++) { 00562 for(size_t j=0;j<n2;j++) { 00563 (*this)[i][j]=ova[i][j]; 00564 } 00565 } 00566 } 00567 } 00568 } 00569 00570 /** \brief Deep copy from an array of uvectors 00571 */ 00572 omatrix_tlate(size_t n, uvector_view_tlate<data_t> uva[]) { 00573 if (n>0) { 00574 size_t n2=uva[0]; 00575 if (n2>0) { 00576 allocate(n,n2); 00577 for(size_t i=0;i<n;i++) { 00578 for(size_t j=0;j<n2;j++) { 00579 (*this)[i][j]=uva[i][j]; 00580 } 00581 } 00582 } 00583 } 00584 } 00585 00586 /** \brief Deep copy from a C-style 2-d array 00587 */ 00588 omatrix_tlate(size_t n, size_t n2, data_t **csa) { 00589 if (n>0 && n2>0) { 00590 allocate(n,n2); 00591 for(size_t i=0;i<n;i++) { 00592 for(size_t j=0;j<n2;j++) { 00593 (*this)[i][j]=csa[i][j]; 00594 } 00595 } 00596 } 00597 } 00598 00599 //@} 00600 00601 ~omatrix_tlate() { 00602 if (mparent_t::size1>0) { 00603 if (mparent_t::owner==1) { 00604 if (mparent_t::block->size>0) { 00605 std::free(mparent_t::block->data); 00606 } 00607 std::free(mparent_t::block); 00608 mparent_t::size1=0; 00609 mparent_t::size2=0; 00610 } 00611 } 00612 } 00613 00614 /// \name Memory allocation 00615 //@{ 00616 /** 00617 \brief Allocate memory after freeing any memory 00618 presently in use 00619 */ 00620 int allocate(size_t nrows, size_t ncols) { 00621 if (mparent_t::size1>0 || mparent_t::size2>0) free(); 00622 00623 if (nrows>0 && ncols>0) { 00624 mparent_t::block=(block_t *)malloc(sizeof(block_t)); 00625 if (mparent_t::block) { 00626 mparent_t::block->data=(data_t *) 00627 malloc(nrows*ncols*sizeof(data_t)); 00628 if (mparent_t::block->data) { 00629 mparent_t::block->size=nrows*ncols; 00630 mparent_t::data=mparent_t::block->data; 00631 mparent_t::size1=nrows; 00632 mparent_t::size2=ncols; 00633 mparent_t::tda=ncols; 00634 mparent_t::owner=1; 00635 } else { 00636 std::free(mparent_t::block); 00637 set_err_ret("No memory for data in omatrix_tlate::allocate()", 00638 gsl_enomem); 00639 } 00640 } else { 00641 set_err_ret("No memory for block in omatrix_tlate::allocate()", 00642 gsl_enomem); 00643 } 00644 } else { 00645 set_err_ret("Zero size in omatrix::allocate()",gsl_einval); 00646 } 00647 return 0; 00648 } 00649 00650 /** 00651 \brief Free the memory 00652 00653 This function will safely do nothing if used without first 00654 allocating memory or if called multiple times in succession. 00655 */ 00656 int free() { 00657 if (mparent_t::size1>0) { 00658 if (mparent_t::owner==1) { 00659 if (mparent_t::block->size>0) { 00660 std::free(mparent_t::block->data); 00661 } 00662 std::free(mparent_t::block); 00663 } 00664 mparent_t::size1=0; 00665 mparent_t::size2=0; 00666 mparent_t::tda=0; 00667 } 00668 return 0; 00669 } 00670 //@} 00671 00672 /// \name Other methods 00673 //@{ 00674 /// \brief Compute the transpose (even if matrix is not square) 00675 omatrix_tlate<data_t,mparent_t,vparent_t,block_t> transpose() { 00676 omatrix_tlate<data_t,mparent_t,vparent_t,block_t> 00677 result(mparent_t::size2,mparent_t::size1); 00678 for(size_t i=0;i<mparent_t::size1;i++) { 00679 for(size_t j=0;j<mparent_t::size2;j++) { 00680 result[j][i]=(*this)[i][j]; 00681 } 00682 } 00683 } 00684 //@} 00685 00686 }; 00687 00688 #ifdef O2SCL_NEVER_DEFINED 00689 } 00690 { 00691 #endif 00692 00693 /** 00694 \brief Create a matrix from an array 00695 */ 00696 template<class data_t, class mparent_t, class block_t> 00697 class omatrix_array_tlate : 00698 public omatrix_view_tlate<data_t,mparent_t,block_t> { 00699 00700 public: 00701 00702 /** \brief Create a vector from \c dat with size \c n */ 00703 omatrix_array_tlate(size_t tot, data_t *dat, 00704 size_t start, size_t ltda, 00705 size_t sz1, size_t sz2) { 00706 if (sz1==0 || sz2==0 || start+sz1*ltda+sz2>=tot) { 00707 mparent_t::block=0; 00708 mparent_t::data=0; 00709 mparent_t::size1=0; 00710 mparent_t::size2=0; 00711 mparent_t::tda=0; 00712 mparent_t::owner=0; 00713 set_err("Failed to create omatrix from array, insufficient size.", 00714 gsl_einval); 00715 } 00716 mparent_t::block=0; 00717 mparent_t::data=dat+start; 00718 mparent_t::size1=sz1; 00719 mparent_t::size2=sz2; 00720 mparent_t::tda=ltda; 00721 mparent_t::owner=0; 00722 } 00723 }; 00724 00725 /** \brief Create a vector from a row of a matrix 00726 */ 00727 template<class data_t, class mparent_t, class vparent_t, class block_t> 00728 class omatrix_row_tlate : 00729 public ovector_view_tlate<data_t,vparent_t,block_t> { 00730 public: 00731 /** \brief Create a vector from row \c i of matrix \c m */ 00732 omatrix_row_tlate(omatrix_view_tlate<data_t,mparent_t,block_t> &m, 00733 size_t i) { 00734 if (i<m.size1) { 00735 vparent_t::size=m.size2; 00736 vparent_t::stride=1; 00737 vparent_t::data=m.data+m.tda()*i; 00738 vparent_t::owner=0; 00739 vparent_t::block=0; 00740 } else { 00741 vparent_t::size=0; 00742 vparent_t::stride=0; 00743 vparent_t::data=0; 00744 vparent_t::owner=0; 00745 vparent_t::block=0; 00746 } 00747 } 00748 }; 00749 00750 /** \brief Create a const vector from a row of a matrix 00751 */ 00752 template<class data_t, class mparent_t, class vparent_t, class block_t> 00753 class omatrix_const_row_tlate : 00754 public ovector_view_tlate<data_t,vparent_t,block_t> { 00755 public: 00756 /** \brief Create a vector from row \c i of matrix \c m */ 00757 omatrix_const_row_tlate 00758 (const omatrix_view_tlate<data_t,mparent_t,block_t> &m, 00759 size_t i) { 00760 if (i<m.size1) { 00761 vparent_t::size=m.size2; 00762 vparent_t::stride=1; 00763 vparent_t::data=m.data+m.tda*i; 00764 vparent_t::owner=0; 00765 vparent_t::block=0; 00766 } else { 00767 vparent_t::size=0; 00768 vparent_t::stride=0; 00769 vparent_t::data=0; 00770 vparent_t::owner=0; 00771 vparent_t::block=0; 00772 } 00773 } 00774 }; 00775 00776 /** \brief Create a vector from a column of a matrix 00777 */ 00778 template<class data_t, class mparent_t, class vparent_t, class block_t> 00779 class omatrix_col_tlate : 00780 public ovector_view_tlate<data_t,vparent_t,block_t> { 00781 public: 00782 /** \brief Create a vector from col \c i of matrix \c m */ 00783 omatrix_col_tlate(omatrix_view_tlate<data_t,mparent_t,block_t> &m, 00784 size_t i) { 00785 if (i<m.size2) { 00786 vparent_t::size=m.size1; 00787 vparent_t::stride=m.tda(); 00788 vparent_t::data=m.data+i; 00789 vparent_t::owner=0; 00790 vparent_t::block=0; 00791 } else { 00792 vparent_t::size=0; 00793 vparent_t::stride=0; 00794 vparent_t::data=0; 00795 vparent_t::owner=0; 00796 vparent_t::block=0; 00797 } 00798 } 00799 }; 00800 00801 /** \brief Create a const vector from a column of a matrix 00802 */ 00803 template<class data_t, class mparent_t, class vparent_t, class block_t> 00804 class omatrix_const_col_tlate : 00805 public ovector_view_tlate<data_t,vparent_t,block_t> { 00806 public: 00807 /** \brief Create a vector from col \c i of matrix \c m */ 00808 omatrix_const_col_tlate(omatrix_view_tlate<data_t,mparent_t,block_t> &m, 00809 size_t i) { 00810 if (i<m.size2) { 00811 vparent_t::size=m.size1; 00812 vparent_t::stride=m.tda(); 00813 vparent_t::data=m.data+i; 00814 vparent_t::owner=0; 00815 vparent_t::block=0; 00816 } else { 00817 vparent_t::size=0; 00818 vparent_t::stride=0; 00819 vparent_t::data=0; 00820 vparent_t::owner=0; 00821 vparent_t::block=0; 00822 } 00823 } 00824 }; 00825 00826 /** \brief Create a vector from the main diagonal 00827 */ 00828 template<class data_t, class mparent_t, class vparent_t, class block_t> 00829 class omatrix_diag_tlate : 00830 public ovector_view_tlate<data_t,vparent_t,block_t> { 00831 public: 00832 /** \brief Create a vector of the diagonal of matrix \c m */ 00833 omatrix_diag_tlate(omatrix_view_tlate<data_t,mparent_t,block_t> &m) { 00834 00835 if (m.size2<m.size1) vparent_t::size=m.size2; 00836 else vparent_t::size=m.size1; 00837 vparent_t::stride=m.tda()+1; 00838 vparent_t::data=m.data; 00839 vparent_t::owner=0; 00840 vparent_t::block=0; 00841 } 00842 }; 00843 00844 00845 /// omatrix typedef 00846 typedef omatrix_tlate<double,gsl_matrix,gsl_vector,gsl_block> omatrix; 00847 /// omatrix_view typedef 00848 typedef omatrix_view_tlate<double,gsl_matrix,gsl_block> omatrix_view; 00849 /// omatrix_row typedef 00850 typedef omatrix_row_tlate<double,gsl_matrix,gsl_vector,gsl_block> 00851 omatrix_row; 00852 /// omatrix_col typedef 00853 typedef omatrix_col_tlate<double,gsl_matrix,gsl_vector,gsl_block> 00854 omatrix_col; 00855 /// omatrix_const_row typedef 00856 typedef omatrix_const_row_tlate<double,gsl_matrix,gsl_vector,gsl_block> 00857 omatrix_const_row; 00858 /// omatrix_const_col typedef 00859 typedef omatrix_const_col_tlate<double,gsl_matrix,gsl_vector,gsl_block> 00860 omatrix_const_col; 00861 /// omatrix_diag typedef 00862 typedef omatrix_diag_tlate<double,gsl_matrix,gsl_vector,gsl_block> 00863 omatrix_diag; 00864 /// omatrix_array typedef 00865 typedef omatrix_array_tlate<double,gsl_matrix,gsl_block> 00866 omatrix_array; 00867 00868 /// omatrix_int typedef 00869 typedef omatrix_tlate<int,gsl_matrix_int,gsl_vector_int,gsl_block_int> 00870 omatrix_int; 00871 /// omatrix_int_view typedef 00872 typedef omatrix_view_tlate<int,gsl_matrix_int,gsl_block_int> 00873 omatrix_int_view; 00874 /// omatrix_int_row typedef 00875 typedef omatrix_row_tlate<int,gsl_matrix_int,gsl_vector_int, 00876 gsl_block_int> omatrix_int_row; 00877 /// omatrix_int_col typedef 00878 typedef omatrix_col_tlate<int,gsl_matrix_int,gsl_vector_int, 00879 gsl_block_int> omatrix_int_col; 00880 /// omatrix_int_const_row typedef 00881 typedef omatrix_const_row_tlate<int,gsl_matrix_int,gsl_vector_int, 00882 gsl_block_int> omatrix_int_const_row; 00883 /// omatrix_int_const_col typedef 00884 typedef omatrix_const_col_tlate<int,gsl_matrix_int,gsl_vector_int, 00885 gsl_block_int> omatrix_int_const_col; 00886 /// omatrix_int_diag typedef 00887 typedef omatrix_diag_tlate<int,gsl_matrix_int,gsl_vector_int, 00888 gsl_block_int> omatrix_int_diag; 00889 /// omatrix_int_array typedef 00890 typedef omatrix_array_tlate<int,gsl_matrix_int,gsl_block_int> 00891 omatrix_int_array; 00892 00893 /** 00894 \brief A operator for output of omatrix objects 00895 00896 This outputs all of the matrix elements. Each row is output with 00897 an endline character at the end of each row. Positive values are 00898 preceeded by an extra space. A 2x2 example: 00899 \verbatim 00900 -3.751935e-05 -6.785864e-04 00901 -6.785864e-04 1.631984e-02 00902 \endverbatim 00903 00904 The function \c gsl_ieee_double_to_rep() is used to determine 00905 the sign of a number, so that "-0.0" as distinct from "+0.0" 00906 is handled correctly. 00907 00908 \comment 00909 I had originally thought about removing this function, since 00910 it's superceded by matrix_out, but it's kind of convenient 00911 \endcomment 00912 00913 \future This assumes that scientific mode is on and showpos 00914 is off. It'd be nice to fix this. 00915 00916 */ 00917 template<class data_t, class parent_t, class block_t> 00918 std::ostream &operator<< 00919 (std::ostream &os, 00920 const omatrix_view_tlate<data_t,parent_t,block_t> &v) { 00921 size_t i; 00922 gsl_ieee_double_rep r; 00923 for(i=0;i<v.rows()-1;i++) { 00924 for(size_t j=0;j<v.cols();j++) { 00925 gsl_ieee_double_to_rep(&(v[i][j]), &r); 00926 if (r.sign==1) os << v[i][j] << ' '; 00927 else os << ' ' << v[i][j] << ' '; 00928 } 00929 os << '\n'; 00930 } 00931 i=v.rows()-1; 00932 if (i>0) { 00933 for(size_t j=0;j<v.cols();j++) { 00934 gsl_ieee_double_to_rep(&(v[i][j]), &r); 00935 if (r.sign==1) os << v[i][j] << ' '; 00936 else os << ' ' << v[i][j] << ' '; 00937 } 00938 } 00939 return os; 00940 } 00941 00942 /** \brief A simple class to provide an \c allocate() function 00943 for \ref omatrix 00944 00945 */ 00946 class omatrix_alloc { 00947 public: 00948 /// Allocate \c v for \c i elements 00949 void allocate(omatrix &o, size_t i, size_t j) { o.allocate(i,j); } 00950 /// Free memory 00951 void free(omatrix &o, size_t i) { o.free(); } 00952 }; 00953 00954 /** \brief A matrix where the memory allocation is performed in 00955 the constructor 00956 */ 00957 #ifdef DOXYGENP 00958 template<size_t N, size_t M> class ofmatrix : 00959 public omatrix_tlate<data_t,mparent_t,vparent_t,block_t> 00960 #else 00961 template<size_t N, size_t M> class ofmatrix : 00962 public omatrix_tlate<double,gsl_matrix,gsl_vector,gsl_block> 00963 #endif 00964 { 00965 public: 00966 ofmatrix() : omatrix_tlate<double,gsl_matrix,gsl_vector, 00967 gsl_block>(N,M) { 00968 } 00969 }; 00970 00971 class xmatrix : public omatrix { 00972 public: 00973 xmatrix(size_t r=0, size_t c=0) : omatrix(r,c) { 00974 } 00975 00976 ovector_array operator[](size_t i) { 00977 #if O2SCL_NO_RANGE_CHECK 00978 #else 00979 if (i>=gsl_matrix::size1) { 00980 set_err((((std::string)"Array index ")+itos(i)+" out of bounds" 00981 +" in omatrix_view_tlate::operator[]. Size: "+ 00982 itos(gsl_matrix::size1)+ 00983 " (index should be less than size).").c_str(),gsl_index); 00984 } 00985 #endif 00986 return ovector_array(gsl_matrix::size2,gsl_matrix::data+ 00987 i*gsl_matrix::tda); 00988 } 00989 00990 const ovector_const_array operator[](size_t i) const { 00991 #if O2SCL_NO_RANGE_CHECK 00992 #else 00993 if (i>=gsl_matrix::size1) { 00994 set_err((((std::string)"Array index ")+itos(i)+" out of bounds" 00995 +" in omatrix_view_tlate::operator[]. Size: "+ 00996 itos(gsl_matrix::size1)+ 00997 " (index should be less than size).").c_str(),gsl_index); 00998 } 00999 #endif 01000 return ovector_const_array(gsl_matrix::size2,gsl_matrix::data+ 01001 i*gsl_matrix::tda); 01002 } 01003 }; 01004 01005 #ifndef DOXYGENP 01006 } 01007 #endif 01008 01009 #endif 01010
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