Main Page   Namespace List   Class Hierarchy   Compound List   File List   Namespace Members   Compound Members   File Members  

gsl::matrix_int Class Reference

#include <matrix_int.h>

List of all members.

Public Methods

 matrix_int ()
 matrix_int (size_t new_rows, size_t new_cols, bool clear=true)
void copy (const matrix_int &other)
template<class oclass> void copy (const oclass &other)
 matrix_int (const matrix_int &other)
template<class oclass>  matrix_int (const oclass &other)
 ~matrix_int ()
size_t get_rows () const
size_t get_cols () const
size_t size1 () const
size_t size2 () const
void dimensions (size_t *num_rows, size_t *num_cols) const
int get_element (size_t row, size_t col) const
const int& operator() (size_t row, size_t col) const
int& operator() (size_t row, size_t col)
void set_element (size_t row, size_t col, const int &v)
void set_elements (const int &new_value)
void set_all (const int &new_value)
void set_zero ()
void set_dimensions (size_t new_rows, size_t new_cols)
void load (const char *filename)
void save (const char *filename) const
int fwrite (FILE *stream) const
int fread (FILE *stream)
void load_binary (const char *filename)
void save_binary (const char *filename) const
bool operator== (const matrix_int &other) const
bool operator!= (const matrix_int &other) const
matrix_int& operator= (const matrix_int &other)
template<class omatrix> matrix_int& operator= (const omatrix &other)
 converts from any other matrix type. More...

matrix_int operator+ (const matrix_int &other) const
matrix_int operator+ (const int &f) const
matrix_int& operator+= (const int &f)
matrix_int& operator+= (const matrix_int &other)
matrix_int operator- (const matrix_int &other) const
matrix_int operator- (const int &f) const
matrix_int& operator-= (const int &f)
matrix_int& operator-= (const matrix_int &other)
matrix_int operator * (const matrix_int &other) const
matrix_int operator * (const int &f) const
matrix_int& operator *= (const int &f)
matrix_int& operator *= (const matrix_int &other)
matrix_int operator/ (const int &) const
matrix_int& operator/= (const int &)
matrix_int transpose () const
matrix_int LU_decomp (gsl::permutation *perm=NULL, int *psign=NULL) const
matrix_int LU_invert () const
int sum () const
double LU_lndet () const
vector_int_view row (size_t rowindex)
const vector_int_view row (size_t rowindex) const
vector_int_view column (size_t colindex)
const vector_int_view column (size_t colindex) const
matrix_int get_row (size_t rowindex) const
matrix_int get_col (size_t colindex) const
matrix_int row_sum () const
matrix_int column_sum () const
double trace () const
int cholesky_decomp (matrix_int &a) const
matrix_int covariance () const
bool is_square () const
void set_diagonal (int f)
void identity (size_t k)
double norm (double n) const
gsl_matrix_int* gslobj ()
const gsl_matrix_int* gslobj () const

Friends

ostream& operator<< (ostream &os, const matrix_int &m)
matrix_int operator+ (const int &f, const matrix_int &other)
matrix_int operator- (const int &f, const matrix_int &other)
matrix_int operator * (const int &f, const matrix_int &other)


Constructor & Destructor Documentation

gsl::matrix_int::matrix_int ( )
 

Definition at line 77 of file matrix_int.cc.

00077                       :m(NULL)
00078 {
00079 }

gsl::matrix_int::matrix_int ( size_t new_rows,
size_t new_cols,
bool clear = true )
 

Definition at line 81 of file matrix_int.cc.

00082 {
00083    if(clear){m = gsl_matrix_int_calloc( rows, cols );}
00084    else     {m = gsl_matrix_int_alloc ( rows, cols );}
00085 }

gsl::matrix_int::matrix_int ( const matrix_int & other ) [inline]
 

Definition at line 74 of file matrix_int.h.

00074 :m(NULL) {copy(other);}

template<class oclass>
gsl::matrix_int::matrix_int ( const oclass & other ) [inline]
 

Definition at line 77 of file matrix_int.h.

00077 :m(NULL) {copy(other);}

gsl::matrix_int::~matrix_int ( )
 

Definition at line 88 of file matrix_int.cc.

00089 {
00090    if ( m ) {gsl_matrix_int_free( m );}
00091 }


Member Function Documentation

matrix_int gsl::matrix_int::LU_decomp ( gsl::permutation * perm = NULL,
int * psign = NULL ) const
 

Definition at line 354 of file matrix_int.cc.

00355 {
00356         bool retPerm=perm!=NULL;
00357         if(!perm){perm = new permutation();}
00358         int sign;
00359         perm->resize( get_rows() );
00360         matrix result=*this;// this does conversion  if necessary
00361         gsl_linalg_LU_decomp(  result.m, perm->gsldata, &sign );
00362 
00363         if(!retPerm){delete perm;}
00364         if(psign){*psign=sign;}
00365 
00366    return result;// this does conversion  if necessary
00367 }

matrix_int gsl::matrix_int::LU_invert ( ) const
 

Definition at line 369 of file matrix_int.cc.

00370 {
00371    int i, j;
00372    permutation p;
00373    matrix a=*this;
00374    a=a.LU_decomp(&p);
00375    matrix inverse(size1(),size2());
00376    gsl_linalg_LU_invert( a.m, p.gsldata, inverse.m );
00377    return inverse;
00378 }

double gsl::matrix_int::LU_lndet ( ) const
 

returns logarithm of the determinant of the matrix.

Definition at line 397 of file matrix_int.cc.

00398 {
00399         matrix a=*this;
00400         a=a.LU_decomp();
00401 
00402         // calculate log determinant from LU decomposed matrix "a"
00403         return gsl_linalg_LU_lndet( a.m );
00404 }

int gsl::matrix_int::cholesky_decomp ( matrix_int & a ) const
 

calculates cholesky decomposition of the matrix, returning success if matrix is positive definite.

Definition at line 530 of file matrix_int.cc.

00531 {
00532         int i, j;
00533         int error;
00534         matrix result=*this;
00535         // do decomposition with call to g_s_l
00536         error = gsl_linalg_cholesky_decomp( result.m );
00537         a=result;
00538         return error;
00539 }

const vector_int_view gsl::matrix_int::column ( size_t colindex ) const
 

Definition at line 430 of file matrix_int.cc.

00431 {
00432         gsl_vector_int_view view=gsl_matrix_int_column(m, colindex);
00433         return vector_int_view::create_vector_view(view);
00434 }

vector_int_view gsl::matrix_int::column ( size_t colindex )
 

returns a vector_int_view of a single column of the matrix.

Definition at line 422 of file matrix_int.cc.

00423 {
00424         gsl_vector_int_view view=gsl_matrix_int_column(m, colindex);
00425         return vector_int_view::create_vector_view(view);
00426 }

matrix_int gsl::matrix_int::column_sum ( ) const
 

calculates sum of columns returned as a row matrix.

Definition at line 495 of file matrix_int.cc.

00496 {
00497         int     i;
00498         matrix_int sum( 1, get_cols() );
00499         
00500         sum.set_zero( );
00501         for ( i = 0; i < get_cols(); i++ ) {
00502                 sum.set_element( 0, i, get_col( i ).sum() );
00503         }
00504                 
00505         return( sum );
00506 }

template<class oclass>
void gsl::matrix_int::copy ( const oclass & other ) [inline]
 

Definition at line 63 of file matrix_int.h.

00064         {
00065                 set_dimensions( other.get_rows(), other.get_cols() );
00066                 for ( size_t i = 0; i < get_rows(); i++ ) 
00067                 {
00068                         for ( size_t j = 0; j < get_cols(); j++ ) 
00069                         {
00070                                 gsl_matrix_int_set( m, i, j, (int)other(i,j));
00071                         }
00072                 }
00073         }

void gsl::matrix_int::copy ( const matrix_int & other ) [inline]
 

Definition at line 57 of file matrix_int.h.

Referenced by matrix_int(), and operator=().

00058         {
00059                 set_dimensions(other.size1(),other.size2());
00060                 gsl_matrix_int_memcpy( m, other.m );
00061         }

matrix_int gsl::matrix_int::covariance ( ) const
 

calculates covariance of the matrix columns.

void gsl::matrix_int::dimensions ( size_t * num_rows,
size_t * num_cols ) const
 

Definition at line 98 of file matrix_int.cc.

00099 {
00100    *num_rows = m->size1;
00101    *num_cols = m->size2;
00102 }

int gsl::matrix_int::fread ( FILE * stream ) [inline]
 

Definition at line 117 of file matrix_int.h.

00117 {return gsl_matrix_int_fread (stream, m);}

int gsl::matrix_int::fwrite ( FILE * stream ) const [inline]
 

Definition at line 114 of file matrix_int.h.

00114 {return gsl_matrix_int_fwrite (stream, m);}

matrix_int gsl::matrix_int::get_col ( size_t colindex ) const
 

returns a column matrix containing a single column of the matrix.

Definition at line 458 of file matrix_int.cc.

Referenced by column_sum().

00459 {
00460         matrix_int columnmatrix( get_rows(), 1 );
00461         gsl_vector_int *tempvector = gsl_vector_int_calloc( get_rows() );
00462         
00463         if ( colindex < 0 || colindex >= get_cols() )
00464         {
00465                 cerr << "column index must be in range 0 to " << get_cols() - 1 << endl;
00466                 exit( 1 );
00467         }
00468         
00469         gsl_matrix_int_get_col( tempvector, m, colindex );
00470         gsl_matrix_int_set_col( columnmatrix.m, 0, tempvector );
00471         for ( int i = 0; i < get_rows(); i++ )
00472                 cout << gsl_vector_int_get( tempvector, i ) << endl;
00473 
00474         // tidy up
00475         gsl_vector_int_free( tempvector );
00476         
00477         return( columnmatrix );
00478 }

size_t gsl::matrix_int::get_cols ( ) const [inline]
 

Definition at line 86 of file matrix_int.h.

Referenced by column_sum(), copy(), get_col(), get_row(), is_square(), norm(), set_diagonal(), set_dimensions(), sum(), trace(), and transpose().

00086 {return m->size2;}

int gsl::matrix_int::get_element ( size_t row,
size_t col ) const [inline]
 

Definition at line 96 of file matrix_int.h.

Referenced by norm(), operator==(), and trace().

00096 {return  gsl_matrix_int_get( m, row, col ) ;}

matrix_int gsl::matrix_int::get_row ( size_t rowindex ) const
 

returns a column matrix containing a single row of the matrix.

Definition at line 437 of file matrix_int.cc.

Referenced by row_sum().

00438 {
00439         matrix_int rowmatrix( 1, get_cols() );
00440         gsl_vector_int *tempvector = gsl_vector_int_calloc( get_cols() );
00441         
00442         if ( rowindex < 0 || rowindex >= get_rows() )
00443         {
00444                 cerr << "row index must be in range 0 to " << get_rows() - 1 << endl;
00445                 exit( 1 );
00446         }
00447 
00448         gsl_matrix_int_get_row( tempvector, m, rowindex );
00449         gsl_matrix_int_set_row( rowmatrix.m, 0, tempvector );
00450 
00451         // tidy up
00452         gsl_vector_int_free( tempvector );
00453         
00454         return( rowmatrix );
00455 }

size_t gsl::matrix_int::get_rows ( ) const [inline]
 

Definition at line 84 of file matrix_int.h.

Referenced by LU_decomp(), copy(), get_col(), get_row(), is_square(), norm(), operator *(), row_sum(), set_diagonal(), set_dimensions(), sum(), trace(), and transpose().

00084 {return m->size1;}

const gsl_matrix_int * gsl::matrix_int::gslobj ( ) const [inline]
 

Definition at line 217 of file matrix_int.h.

00217 {if (!m){cout << "matrix_int::gslobj ERROR, data not initialized!! " << endl; exit(-1);}return m;}

gsl_matrix_int * gsl::matrix_int::gslobj ( ) [inline]
 

for interfacing with gsl c

Definition at line 216 of file matrix_int.h.

00216 {if (!m){cout << "matrix_int::gslobj ERROR, data not initialized!! " << endl; exit(-1);}return m;}

void gsl::matrix_int::identity ( size_t k )
 

sets matrix to a k dimensional unit matrix.

Definition at line 542 of file matrix_int.cc.

00543 {
00544         set_dimensions( k, k );
00545         set_zero();
00546         set_diagonal( 1 );
00547 }

bool gsl::matrix_int::is_square ( ) const
 

returns 1 if matrix is square, 0 otherwise.

Definition at line 559 of file matrix_int.cc.

00560 {
00561         if ( get_rows() == get_cols() )
00562                 return true;
00563         return false;
00564 }

void gsl::matrix_int::load ( const char * filename )
 

Definition at line 122 of file matrix_int.cc.

00123 {
00124    FILE * f = fopen( filename, "r" ) ;
00125    matrix temp;
00126    gsl_matrix_fread ( f, temp.m );
00127    fclose (f);
00128    *this = temp;
00129 }

void gsl::matrix_int::load_binary ( const char * filename )
 

Definition at line 165 of file matrix_int.cc.

00166 {
00167    ;
00168 }

double gsl::matrix_int::norm ( double n ) const
 

returns sum of nth power of all elements.

Definition at line 567 of file matrix_int.cc.

00568 {
00569         int i, j;
00570         double sum = 0.0;
00571         
00572         for ( i = 0; i < get_rows(); i++ ) {
00573                 for ( j = 0; j < get_cols(); j++ ) {
00574                         sum += pow( get_element( i, j ), n );
00575                 }
00576         }
00577         
00578         return sum;
00579 }

matrix_int gsl::matrix_int::operator * ( const int & f ) const
 

Definition at line 285 of file matrix_int.cc.

00286 {
00287    matrix_int result( *this );
00288    gsl_matrix_int_scale( result.m, f );
00289    return( result );
00290 }

matrix_int gsl::matrix_int::operator * ( const matrix_int & other ) const
 

Definition at line 272 of file matrix_int.cc.

00273 {
00274         matrix result( get_rows(), other.get_cols() );
00275 #ifdef type_is_double
00276         gsl_linalg_matmult(m,other.m,result.m);
00277 #else //type_is_double
00278         matrix a=*this;
00279         matrix b=other;
00280         gsl_linalg_matmult(a.m,b.m,result.m);
00281 #endif //type_is_double
00282         return result ;
00283 }

matrix_int & gsl::matrix_int::operator *= ( const matrix_int & other )
 

Definition at line 307 of file matrix_int.cc.

00308 {
00309    *this=(*this)*other;
00310    return( *this );
00311 }

matrix_int & gsl::matrix_int::operator *= ( const int & f )
 

Definition at line 300 of file matrix_int.cc.

00301 {
00302    gsl_matrix_int_scale( m, f );
00303 
00304    return( *this );
00305 }

bool gsl::matrix_int::operator!= ( const matrix_int & other ) const [inline]
 

Definition at line 125 of file matrix_int.h.

00125 {return !((*this)==other);}

int & gsl::matrix_int::operator() ( size_t row,
size_t col ) [inline]
 

Definition at line 98 of file matrix_int.h.

00098 {return *gsl_matrix_int_ptr( m, row, col ) ;}

const int & gsl::matrix_int::operator() ( size_t row,
size_t col ) const [inline]
 

Definition at line 97 of file matrix_int.h.

00097 {return *gsl_matrix_int_ptr( m, row, col ) ;}

matrix_int gsl::matrix_int::operator+ ( const int & f ) const
 

Definition at line 203 of file matrix_int.cc.

00204 {
00205    matrix_int result( *this );
00206    gsl_matrix_int_add_constant( result.m, f );
00207 
00208    return( result );
00209 }

matrix_int gsl::matrix_int::operator+ ( const matrix_int & other ) const
 

Definition at line 196 of file matrix_int.cc.

00197 {
00198         matrix_int result(*this);
00199         gsl_matrix_int_add( result.m, other.m );
00200         return result;
00201 }

matrix_int & gsl::matrix_int::operator+= ( const matrix_int & other )
 

Definition at line 226 of file matrix_int.cc.

00227 {
00228    gsl_matrix_int_add( m, other.m );
00229 
00230    return( *this );
00231 }

matrix_int & gsl::matrix_int::operator+= ( const int & f )
 

Definition at line 219 of file matrix_int.cc.

00220 {
00221    gsl_matrix_int_add_constant( m, f );
00222 
00223    return( *this );
00224 }

matrix_int gsl::matrix_int::operator- ( const int & f ) const
 

Definition at line 241 of file matrix_int.cc.

00242 {
00243    matrix_int result( *this );
00244    gsl_matrix_int_add_constant( result.m, -f );
00245 
00246    return( result );
00247 }

matrix_int gsl::matrix_int::operator- ( const matrix_int & other ) const
 

Definition at line 233 of file matrix_int.cc.

00234 {
00235    matrix_int result( *this );
00236    gsl_matrix_int_sub( result.m, other.m );
00237 
00238    return( result );
00239 }

matrix_int & gsl::matrix_int::operator-= ( const matrix_int & other )
 

Definition at line 264 of file matrix_int.cc.

00265 {
00266    gsl_matrix_int_sub( m, other.m );
00267 
00268    return( *this );
00269 }

matrix_int & gsl::matrix_int::operator-= ( const int & f )
 

Definition at line 257 of file matrix_int.cc.

00258 {
00259    gsl_matrix_int_add_constant( m, -f );
00260 
00261    return( *this );
00262 }

matrix_int gsl::matrix_int::operator/ ( const int & f ) const
 

Definition at line 313 of file matrix_int.cc.

00314 {
00315    matrix_int result( *this );
00316 
00317    if ( f != 0.0 ) {
00318       gsl_matrix_int_scale( result.m, 1.0 / f );
00319    } else {
00320       cout << "e_r_r_o_r: division by zero." << endl;
00321       return( result );
00322    }
00323 
00324    return( result );
00325 }

matrix_int & gsl::matrix_int::operator/= ( const int & f )
 

Definition at line 327 of file matrix_int.cc.

00328 {
00329    if ( f != 0.0 ) {
00330       gsl_matrix_int_scale( m, 1.0 / f );
00331    } else {
00332       cout << "e_rr_or: division by zero." << endl;
00333       return( *this );
00334    }
00335 
00336    return( *this );
00337 }

template<class omatrix>
matrix_int & gsl::matrix_int::operator= ( const omatrix & other ) [inline]
 

converts from any other matrix type.

Definition at line 130 of file matrix_int.h.

00130 {copy(other);return *this;}

matrix_int & gsl::matrix_int::operator= ( const matrix_int & other ) [inline]
 

Definition at line 127 of file matrix_int.h.

00127 {copy(other);return *this;}

bool gsl::matrix_int::operator== ( const matrix_int & other ) const
 

Definition at line 175 of file matrix_int.cc.

00176 {
00177    size_t i, j;
00178    
00179    // first check for same dimensions
00180    if ( size1() != other.size1() || size2() != other.size2() )
00181    {
00182       return false;
00183    }
00184 
00185    for ( i = 0; i < size1(); i++ ) {
00186       for ( j = 0; j < size2(); j++ ) {
00187          if ( this->get_element( i, j ) != other.get_element( i, j ) ) {
00188             return false;
00189          }
00190       }
00191    }
00192 
00193    return true;
00194 }

const vector_int_view gsl::matrix_int::row ( size_t rowindex ) const
 

Definition at line 415 of file matrix_int.cc.

00416 {
00417         gsl_vector_int_view view=gsl_matrix_int_row(m, rowindex);
00418         return vector_int_view::create_vector_view(view);
00419 }

vector_int_view gsl::matrix_int::row ( size_t rowindex )
 

returns a vector_int_view of a single row of the matrix.

Definition at line 407 of file matrix_int.cc.

00408 {
00409         gsl_vector_int_view view=gsl_matrix_int_row(m, rowindex);
00410         return vector_int_view::create_vector_view(view);
00411 }

matrix_int gsl::matrix_int::row_sum ( ) const
 

calculates sum of rows returned as a column matrix.

Definition at line 481 of file matrix_int.cc.

00482 {
00483         int     i;
00484         matrix_int sum( get_rows(), 1 );
00485         
00486         sum.set_zero();
00487         for ( i = 0; i < get_rows(); i++ ) {
00488                 sum.set_element( i, 0, get_row( i ).sum() );
00489         }
00490         
00491         return( sum );
00492 }

void gsl::matrix_int::save ( const char * filename ) const
 

Definition at line 131 of file matrix_int.cc.

00132 {
00133    FILE * f = fopen( filename, "w" ) ;
00134    matrix temp = *this;
00135    gsl_matrix_fwrite ( f, temp.m );
00136    fclose ( f );
00137 }

void gsl::matrix_int::save_binary ( const char * filename ) const
 

Definition at line 170 of file matrix_int.cc.

00171 {
00172    ;
00173 }

void gsl::matrix_int::set_all ( const int & new_value ) [inline]
 

Definition at line 103 of file matrix_int.h.

00103 {gsl_matrix_int_set_all ( m, new_value );}

void gsl::matrix_int::set_diagonal ( int f )
 

set diagonal elements of a square matrix to f.

Definition at line 550 of file matrix_int.cc.

Referenced by identity().

00551 {
00552         size_t i;
00553         int mn=(get_rows()<get_cols() ? get_rows() : get_cols());
00554         for ( i = 0; i < mn; i++ )
00555                         set_element( i, i, f );
00556 }

void gsl::matrix_int::set_dimensions ( size_t new_rows,
size_t new_cols )
 

Definition at line 110 of file matrix_int.cc.

Referenced by copy(), and identity().

00111 {
00112    // if dimensions have changed re-allocate matrix
00113    if ( m && (get_rows() != new_rows || get_cols() != new_cols )) {
00114               gsl_matrix_int_free( m );
00115    }
00116 
00117    // allocate
00118    m = gsl_matrix_int_calloc( new_rows, new_cols );
00119 }

void gsl::matrix_int::set_element ( size_t row,
size_t col,
const int & v ) [inline]
 

Definition at line 100 of file matrix_int.h.

Referenced by set_diagonal().

00100 { gsl_matrix_int_set( m, row, col, v );}

void gsl::matrix_int::set_elements ( const int & new_value )
 

Definition at line 105 of file matrix_int.cc.

00106 {
00107    gsl_matrix_int_set_all( m, new_value );
00108 }

void gsl::matrix_int::set_zero ( ) [inline]
 

Definition at line 104 of file matrix_int.h.

Referenced by identity().

00104 {gsl_matrix_int_set_zero( m );}

size_t gsl::matrix_int::size1 ( ) const [inline]
 

Definition at line 88 of file matrix_int.h.

Referenced by LU_invert(), and operator==().

00088 {return m->size1;}

size_t gsl::matrix_int::size2 ( ) const [inline]
 

Definition at line 90 of file matrix_int.h.

Referenced by LU_invert(), and operator==().

00090 {return m->size2;}

int gsl::matrix_int::sum ( ) const
 

returns sum of all the matrix elements.

Definition at line 382 of file matrix_int.cc.

Referenced by column_sum(), and row_sum().

00383 {
00384         int i, j;
00385         int sum = 0;
00386 
00387         for ( i = 0; i < get_rows(); i++ ) {
00388                 for ( j = 0; j < get_cols(); j++ ) {
00389                         sum += gsl_matrix_int_get( m, i, j );
00390                 }
00391         }
00392 
00393         return( sum );
00394 }

double gsl::matrix_int::trace ( ) const
 

returns trace (diagonal sum) of a square matrix.

Definition at line 509 of file matrix_int.cc.

00510 {
00511         int     i;
00512         double sum = 0.0;
00513         
00514         if ( get_rows() != get_cols() ) {
00515                 cerr << "e_r_r_o_r: cannot calculate trace of non-square matrix.";
00516                 cerr << endl;
00517                 exit( 1 );
00518         }
00519 
00520         // calculate sum of diagonal elements
00521         for ( i = 0; i < get_rows(); i++ ) {
00522                 sum += get_element( i, i );
00523         }
00524 
00525         return( sum );
00526 }

matrix_int gsl::matrix_int::transpose ( ) const
 

Definition at line 339 of file matrix_int.cc.

00340 {
00341    int i, j;
00342    matrix_int result( get_cols(), get_rows() );
00343 
00344    for ( i = 0; i < get_rows(); i++ ) {
00345       for ( j = 0; j < get_cols(); j++ ) {
00346          gsl_matrix_int_set( result.m, j, i, gsl_matrix_int_get( m, i, j ) );
00347       }
00348    }
00349    
00350    return( result );
00351 }


Friends And Related Function Documentation

matrix_int operator * ( const int & f,
const matrix_int & other ) [friend]
 

Definition at line 292 of file matrix_int.cc.

00293 {
00294    matrix_int result( other );
00295    gsl_matrix_int_scale( result.m, f );
00296 
00297    return( result );
00298 }

matrix_int operator+ ( const int & f,
const matrix_int & other ) [friend]
 

Definition at line 211 of file matrix_int.cc.

00212 {
00213    matrix_int result( other );
00214    gsl_matrix_int_add_constant( result.m, f );
00215 
00216    return( result );
00217 }

matrix_int operator- ( const int & f,
const matrix_int & other ) [friend]
 

Definition at line 249 of file matrix_int.cc.

00250 {
00251    matrix_int result( -1 * other );
00252    gsl_matrix_int_add_constant( result.m, f );
00253 
00254    return( result );
00255 }

ostream & operator<< ( ostream & os,
const matrix_int & m ) [friend]
 

Definition at line 139 of file matrix_int.cc.

00140 {
00141    size_t i, j;
00142 
00143    os.setf( ios::fixed );
00144 
00145 //FIXME for aCC (doesn't find correct outstream function
00146 //     for ( i = 0; i < m.get_rows(); i++ ) {
00147 //         for ( j = 0; j < m.get_cols() - 1; j++ ) {
00148 //                 os << setprecision( 6 ) << setw( 11 ) ;//<< m.get_element( i, j ) << " ";
00149 //         }
00150 //         os << setprecision( 6 ) << setw( 11 ) ;//<< m.get_element( i, j ) << endl;
00151 //     }
00152 
00153    for ( i = 0; i < m.get_rows(); i++ ) {
00154            for ( j = 0; j < m.get_cols() - 1; j++ ) {
00155                    os << m.get_element( i, j ) << " ";
00156            }
00157            os << m.get_element( i, j ) << endl;
00158    }
00159 
00160    return os;
00161 }


The documentation for this class was generated from the following files:
Generated at Sun Dec 16 23:44:45 2001 for gslwrap by doxygen1.2.8.1 written by Dimitri van Heesch, © 1997-2001