Matrix.h 18 KB
Newer Older
Francois Gygi committed
1 2 3 4 5
////////////////////////////////////////////////////////////////////////////////
//
// Matrix.h
//
////////////////////////////////////////////////////////////////////////////////
6
// $Id: Matrix.h,v 1.11 2005-02-04 21:59:30 fgygi Exp $
Francois Gygi committed
7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112

#ifndef MATRIX_H
#define MATRIX_H

#define MATRIX_DEF_BLOCK_SIZE 64

class Context;

#include <valarray>
#include <complex>
using namespace std;

class ComplexMatrix;

class DoubleMatrix
{
  private:
  
    Context ctxt_;
    int ictxt_;
    int lld_;  // leading dimension of local matrix
    int m_, n_;   // size of global matrix
    int mb_, nb_; // size of blocks
    int mloc_, nloc_;   // size of local array
    int size_; // size of local array;
    int nprow_, npcol_; // number of process rows and cols in context
    int myrow_, mycol_; // position of my process in the process grid
    int mblocks_, nblocks_; // number of local blocks
    int desc_[9];
    bool active_;
    bool m_incomplete_, n_incomplete_; // this process has an incomplete block
    bool reference_; // object was created using the copy constructor
    double* val;

  public:
  
    double* valptr(int i=0) { return &val[i]; }
    const double* cvalptr(int i=0) const { return &val[i]; }
    double& operator[] (int i) { return val[i]; }
    const double& operator[] (int i) const { return val[i]; }
    const Context& context(void) const { return ctxt_; }
    int ictxt(void) const { return ictxt_; }
    int lld(void) const { return lld_; }
    int m(void) const { return m_; } // size of global matrix
    int n(void) const { return n_; } // size of global matrix
    int mb(void) const { return mb_; } // size of blocks
    int nb(void) const { return nb_; } // size of blocks
    int mloc(void) const { return mloc_; } // size of local array
    int nloc(void) const { return nloc_; } // size of local array
    int size(void) const { return size_; } // size of local array
    int localsize(void) const { return mloc_*nloc_; } // local size of val
    double memsize(void) const { return (double)m_ * (double)n_ * sizeof(double); }
    double localmemsize(void) const
    { return (double) mloc_ * (double) nloc_ * sizeof(double); }
    const int* desc(void) const { return &desc_[0]; }
    
    // local block size of block (l,m)
    int mbs(int l) const
    {
      // if l is the last block and this process holds an incomplete
      // block, then the size is m_%mb_. Otherwise, it is mb_.
      return ( (l==mblocks_-1) && ( m_incomplete_ ) ) ? m_%mb_ : mb_;
    }
    int nbs(int m) const
    {
      // if m is the last block and this process holds an incomplete
      // block, then the size is n_%nb_. Otherwise, it is nb_.
      return ( (m==nblocks_-1) && ( n_incomplete_ ) ) ? n_%nb_ : nb_;
    }
    
    // number of local blocks
    int mblocks(void) const { return mblocks_; }     
    int nblocks(void) const { return nblocks_; }
     
    // functions to compute local indices from global indices

    // index of blocks: element (i,j) is at position (x,y)
    // in local block (l,m) of process (pr,pc)    
    int l(int i) const { return i/(nprow_*mb_); }
    int x(int i) const { return i % mb_; }
    int pr(int i) const { return (i/mb_) % nprow_; }
    
    int m(int j) const { return j/(npcol_*nb_); }
    int y(int j) const { return j % nb_; }
    int pc(int j) const { return (j/nb_) % npcol_; }
    
    // global indices:
    // (i,j) is the global index of element (x,y) of block (l,m)
    int i(int l, int x) const { return (l * nprow_ + myrow_) * mb_ + x; }     
    int j(int m, int y) const { return (m * npcol_ + mycol_) * nb_ + y; }
    
    // store element a(ii,jj) (where ii,jj are global indices)
    // in array val:
    //        int iii = l(ii) * mb_ + x(ii);
    //        int jjj = m(jj) * nb_ + y(jj);
    //        val[iii+mloc_*jjj] = a_ij;
    
    // active flag: the matrix has elements on this process
    bool active(void) const { return active_; }
    
    void init_size(int m, int n, int mb = MATRIX_DEF_BLOCK_SIZE, 
      int nb = MATRIX_DEF_BLOCK_SIZE);
          
    void resize(int m, int n, int mb = MATRIX_DEF_BLOCK_SIZE, 
      int nb = MATRIX_DEF_BLOCK_SIZE)
    {
113 114 115
      const int old_size = size_;
      init_size(m,n,mb,nb);
      if ( size_ == old_size ) return;
Francois Gygi committed
116 117 118
      init_size(m,n,mb,nb);
      delete[] val;
      val = new double[size_];
119
      clear();
Francois Gygi committed
120 121 122 123 124
    }    
    
    void print(ostream& os) const;
    
    explicit DoubleMatrix(const Context& ctxt) : ctxt_(ctxt),
125
        m_(0), n_(0), mb_(0), nb_(0), size_(0), reference_(false), val(0) {}
Francois Gygi committed
126 127 128 129
        
    // Construct a DoubleMatrix of dimensions m,n
    explicit DoubleMatrix(const Context& ctxt, int m, int n,
        int mb=MATRIX_DEF_BLOCK_SIZE, 
130
        int nb=MATRIX_DEF_BLOCK_SIZE) : ctxt_(ctxt),
131
        m_(0), n_(0), mb_(0), nb_(0), size_(0), reference_(false), val(0)
Francois Gygi committed
132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215
    {
      resize(m,n,mb,nb);
    }
    
    // copy constructor: create a separate copy of rhs
    explicit DoubleMatrix(const DoubleMatrix& rhs) : ctxt_(rhs.context()),
      reference_(false)
    {
      init_size(rhs.m(),rhs.n(),rhs.mb(),rhs.nb());
      val = new double[size_];
      memcpy(val, rhs.val, size_*sizeof(double));
    }
    
    // reference constructor create a proxy for a ComplexMatrix rhs
    explicit DoubleMatrix(ComplexMatrix& rhs);
    explicit DoubleMatrix(const ComplexMatrix& rhs);
        
    ~DoubleMatrix(void) 
    { 
      if ( !reference_ ) delete[] val;
    }
                  
    DoubleMatrix& operator=(const DoubleMatrix& a);
    
    DoubleMatrix& operator+=(const DoubleMatrix& a);
    DoubleMatrix& operator-=(const DoubleMatrix& a);
    
    DoubleMatrix& operator*=(double a);
    
    void matgather(double *a, int lda) const;

    void initdiag(const double* const a); // initialize diagonal from a[]
    void init(const double* const a, int lda);
    void set(char uplo, double x);
    void identity(void);
    void clear(void);
    
    double dot(const DoubleMatrix &a) const;
    void axpy(double alpha, const DoubleMatrix &a);
    void scal(double alpha);
    
    double nrm2(void) const;
    double asum(void) const;
    double amax(void) const;
    double trace(void) const;
    void symmetrize(char uplo);
    
    // rank-1 update: *this += alpha * x(kx) * (y(ky))^T
    // where x(kx) is row kx of x, and y(ky) is row ky of y
    void ger(double alpha, const DoubleMatrix& x, int kx,
                           const DoubleMatrix& y, int ky);

    // symmetric rank-1 update
    void syr(char uplo, double alpha, const DoubleMatrix& x,
             int ix, char rowcol);
    
    // get submatrix A(ia:ia+m,ja:ja+n) of A
    void getsub(const DoubleMatrix& a,int m,int n,int ia,int ja);
        
    // matrix * matrix
    // this = alpha*op(A)*op(B)+beta*this
    void gemm(char transa, char transb,
         double alpha, const DoubleMatrix& a,
         const DoubleMatrix& b, double beta);
    
    // symmetric_matrix * matrix
    // *this = alpha * A * B + beta * this
    void symm(char side, char uplo,
         double alpha, const DoubleMatrix& a,
         const DoubleMatrix& b, double beta);
    
    // symmetric rank k update
    // this = beta * this + alpha * A * A^T  (trans=='n')
    // this = beta * this + alpha * A^T * A  (trans=='t')
    void syrk(char uplo, char trans,
         double alpha, const DoubleMatrix& a, double beta);

    // matrix transpose
    // this = alpha * transpose(A) + beta * this
    void transpose(double alpha, const DoubleMatrix& a, double beta);
    void transpose(const DoubleMatrix& a);

    void trmm(char side, char uplo, char trans, char diag, 
              double alpha, const DoubleMatrix& a);
216 217
    
    // solve triangular system
Francois Gygi committed
218 219 220 221 222 223 224 225 226
    void trsm(char side, char uplo, char trans, char diag, 
              double alpha, const DoubleMatrix& a);
    void trtrs(char uplo, char trans, char diag, DoubleMatrix& b) const;
    
    // Cholesky decomposition of a symmetric matrix
    void potrf(char uplo);
    // Inverse of a symmetric matrix from Cholesky factor
    void potri(char uplo);
    
227 228 229 230 231 232
    // LU decomposition
    void lu(valarray<int>& ipiv);
    
    // compute inverse of a square matrix
    void inverse(void);
    
Francois Gygi committed
233 234 235 236 237 238 239 240 241 242 243
    // Inverse of triangular matrix
    void trtri(char uplo,char diag);
    
    double pocon(char) const;
    void sygst(int, char, const DoubleMatrix&);
    
    // compute eigenvalues and eigenvectors of symmetric matrix *this
    void syev(char uplo, valarray<double>& w, DoubleMatrix& z);
    
    // compute eigenvalues (only) of symmetric matrix *this
    void syev(char uplo, valarray<double>& w);
244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259

    // compute eigenvalues and eigenvectors of symmetric matrix *this
    // using the divide and conquer method of Tisseur and Dongarra
    void syevd(char uplo, valarray<double>& w, DoubleMatrix& z);
    
    // compute eigenvalues (only) of symmetric matrix *this
    // using the divide and conquer method of Tisseur and Dongarra
    void syevd(char uplo, valarray<double>& w);

    // compute eigenvalues and eigenvectors of symmetric matrix *this
    // using the expert driver
    void syevx(char uplo, valarray<double>& w, DoubleMatrix& z, double abstol);
    
    // compute eigenvalues (only) of symmetric matrix *this
    // using the divide and conquer method of Tisseur and Dongarra
    //void syevx(char uplo, valarray<double>& w);
Francois Gygi committed
260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355
};
ostream& operator << ( ostream& os, const DoubleMatrix& a );

class ComplexMatrix
{
  private:
  
    Context ctxt_;
    int ictxt_;
    int lld_;  // leading dimension of local matrix
    int m_, n_;   // size of global matrix
    int mb_, nb_; // size of blocks
    int mloc_, nloc_;   // size of local array
    int size_; // size of local array;
    int nprow_, npcol_; // number of process rows and cols in  context
    int myrow_, mycol_; // position of my process in the process grid
    int mblocks_, nblocks_; // number of local blocks
    int desc_[9];
    bool active_;
    bool m_incomplete_, n_incomplete_; // this process has an incomplete block
    bool reference_; // object was created using the copy constructor
    complex<double>* val;

  public:
  
    complex<double>* valptr(int i=0) { return &val[i]; }
    const complex<double>* cvalptr(int i=0) const { return &val[i]; }
    complex<double>& operator[] (int i) { return val[i]; }
    const complex<double>& operator[] (int i) const { return val[i]; }
    const Context& context(void) const { return ctxt_; }
    int ictxt(void) const { return ictxt_; }
    int lld(void) const { return lld_; }
    int m(void) const { return m_; } // size of global matrix
    int n(void) const { return n_; } // size of global matrix
    int mb(void) const { return mb_; } // size of blocks
    int nb(void) const { return nb_; } // size of blocks
    int mloc(void) const { return mloc_; } // size of local array
    int nloc(void) const { return nloc_; } // size of local array
    int size(void) const { return size_; } // size of local array
    int localsize(void) const { return mloc_*nloc_; } // local size of val
    double memsize(void) const 
    { return (double) m_ * (double) n_ * sizeof(complex<double>); }
    double localmemsize(void) const
    { return (double) mloc_ * (double) nloc_ * sizeof(complex<double>); }
    const int* desc(void) const { return &desc_[0]; }
    
    // local block size of block (l,m)
    int mbs(int l) const
    {
      // if l is the last block and this process holds an incomplete
      // block, then the size is m_%mb_. Otherwise, it is mb_.
      return ( (l==mblocks_-1) && ( m_incomplete_ ) ) ? m_%mb_ : mb_;
    }
    int nbs(int m) const
    {
      // if m is the last block and this process holds an incomplete
      // block, then the size is n_%nb_. Otherwise, it is nb_.
      return ( (m==nblocks_-1) && ( n_incomplete_ ) ) ? n_%nb_ : nb_;
    }
    
    // number of local blocks
    int mblocks(void) const { return mblocks_; }     
    int nblocks(void) const { return nblocks_; }
     
    // functions to compute local indices from global indices

    // index of blocks: element (i,j) is at position (x,y)
    // in local block (l,m) of process (pr,pc)    
    int l(int i) const { return i/(nprow_*mb_); }
    int x(int i) const { return i % mb_; }
    int pr(int i) const { return (i/mb_) % nprow_; }
    
    int m(int j) const { return j/(npcol_*nb_); }
    int y(int j) const { return j % nb_; }
    int pc(int j) const { return (j/nb_) % npcol_; }
    
    // global indices:
    // (i,j) is the global index of element (x,y) of block (l,m)
    int i(int l, int x) const { return (l * nprow_ + myrow_) * mb_ + x; }     
    int j(int m, int y) const { return (m * npcol_ + mycol_) * nb_ + y; }
    
    // store element a(ii,jj) (where ii,jj are global indices)
    // in array val:
    //        int iii = l(ii) * mb_ + x(ii);
    //        int jjj = m(jj) * nb_ + y(jj);
    //        val[iii+mloc_*jjj] = a_ij;
    
    // active flag: the matrix has elements on this process
    bool active(void) const { return active_; }
    
    void init_size(int m, int n, int mb = MATRIX_DEF_BLOCK_SIZE, 
      int nb = MATRIX_DEF_BLOCK_SIZE);
          
    void resize(int m, int n, int mb = MATRIX_DEF_BLOCK_SIZE, 
      int nb = MATRIX_DEF_BLOCK_SIZE)
    {
356
      const int old_size = size_;
Francois Gygi committed
357
      init_size(m,n,mb,nb);
358
      if ( size_ == old_size ) return;
Francois Gygi committed
359 360
      delete[] val;
      val = new complex<double>[size_];
361
      clear();
Francois Gygi committed
362 363 364 365 366
    }    
    
    void print(ostream& os) const;
    
    explicit ComplexMatrix(const Context& ctxt) : ctxt_(ctxt),
367
        m_(0), n_(0), mb_(0), nb_(0), size_(0), reference_(false), val(0) {}
Francois Gygi committed
368 369 370 371
        
    // Construct a ComplexMatrix of dimensions m,n
    explicit ComplexMatrix(const Context& ctxt, int m, int n,
        int mb=MATRIX_DEF_BLOCK_SIZE, 
372
        int nb=MATRIX_DEF_BLOCK_SIZE) : ctxt_(ctxt), 
373
        m_(0), n_(0), mb_(0), nb_(0), size_(0), reference_(false), val(0)
Francois Gygi committed
374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494
    {
      resize(m,n,mb,nb);
    }
    
    // copy constructor: create a separate copy of rhs
    explicit ComplexMatrix(const ComplexMatrix& rhs) : ctxt_(rhs.context()),
      reference_(false)
    {
      init_size(rhs.m(),rhs.n(),rhs.mb(),rhs.nb());
      val = new complex<double>[size_];
      memcpy(val, rhs.val, size_*sizeof(complex<double>));
    }
    
    // reference constructor: create a proxy for a DoubleMatrix rhs
    explicit ComplexMatrix(DoubleMatrix& rhs);
    explicit ComplexMatrix(const DoubleMatrix& rhs);
    
    ~ComplexMatrix(void)
    { 
      if ( !reference_ ) delete[] val;
    }
                  
    ComplexMatrix& operator=(const ComplexMatrix& a);
    
    ComplexMatrix& operator+=(const ComplexMatrix& a);
    ComplexMatrix& operator-=(const ComplexMatrix& a);
    
    ComplexMatrix& operator*=(double a);
    ComplexMatrix& operator*=(complex<double> a);
    
    void matgather(complex<double> *a, int lda) const;

    void initdiag(const complex<double>* const a); // initialize diagonal 
    void init(const complex<double>* const a, int lda);
    void set(char uplo, complex<double> x);
    void identity(void);
    void clear(void);
    
    complex<double> dot(const ComplexMatrix &a) const;  // tr A^H * A
    complex<double> dotu(const ComplexMatrix &a) const; // tr A^T * A
    void axpy(complex<double> alpha, const ComplexMatrix &a);
    void axpy(double alpha, const ComplexMatrix &a);
    void scal(complex<double> alpha);
    void scal(double alpha);
    
    double nrm2(void) const;
    double asum(void) const;
    double amax(void) const;
    complex<double> trace(void) const;
    void symmetrize(char uplo);
    
    // rank-1 update: *this += alpha * x(kx) * (y(ky))^T
    // where x(kx) is row kx of x, and y(ky) is row ky of y
    void ger(complex<double> alpha, const ComplexMatrix& x,int kx,
                                    const ComplexMatrix& y,int ky);
    void geru(complex<double> alpha, const ComplexMatrix& x,int kx,
                                     const ComplexMatrix& y,int ky);
    void gerc(complex<double> alpha, const ComplexMatrix& x,int kx,
                                     const ComplexMatrix& y,int ky);

    // symmetric rank-1 update
    void her(char uplo, complex<double> alpha, 
             const ComplexMatrix& x, int ix, char rowcol);
    
    // get submatrix A(ia:ia+m,ja:ja+n) of A
    void getsub(const ComplexMatrix& a, int m, int n, int ia, int ja);
        
    // matrix * matrix
    // this = alpha*op(A)*op(B)+beta*this
    void gemm(char transa, char transb,
         complex<double> alpha, const ComplexMatrix& a,
         const ComplexMatrix& b, complex<double> beta);
    
    // hermitian_matrix * matrix
    // *this = alpha * A * B + beta * this
    void hemm(char side, char uplo,
         complex<double> alpha, const ComplexMatrix& a,
         const ComplexMatrix& b, complex<double> beta);
    
    // complex_symmetric_matrix * matrix
    // *this = alpha * A * B + beta * this
    void symm(char side, char uplo,
         complex<double> alpha, const ComplexMatrix& a,
         const ComplexMatrix& b, complex<double> beta);
    
    // hermitian rank k update
    // this = beta * this + alpha * A * A^T  (trans=='n')
    // this = beta * this + alpha * A^T * A  (trans=='t')
    void herk(char uplo, char trans,
         complex<double> alpha, const ComplexMatrix& a, 
         complex<double> beta);

    // matrix transpose
    // this = alpha * transpose(A) + beta * this
    void transpose(complex<double> alpha, const ComplexMatrix& a, 
                   complex<double> beta);
    void transpose(const ComplexMatrix& a);

    void trmm(char side, char uplo, char trans, char diag, 
              complex<double> alpha, const ComplexMatrix& a);
    void trsm(char side, char uplo, char trans, char diag, 
              complex<double> alpha, const ComplexMatrix& a);
    void trtrs(char uplo, char trans, char diag, ComplexMatrix& b) const;
    
    // Cholesky decomposition of a hermitian matrix
    void potrf(char uplo);
    // Inverse of a symmetric matrix from Cholesky factor
    void potri(char uplo);
    
    // Inverse of triangular matrix
    void trtri(char uplo,char diag);
    
    double pocon(char) const;
    
    // compute eigenvalues and eigenvectors of hermitian matrix *this
    void heev(char uplo, valarray<double>& w, ComplexMatrix& z);
    // compute eigenvalues (only) of hermitian matrix *this
    void heev(char uplo, valarray<double>& w);
};
ostream& operator << ( ostream& os, const ComplexMatrix& a );
#endif