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