Matrix.C 100 KB
Newer Older
Francois Gygi committed
1 2
////////////////////////////////////////////////////////////////////////////////
//
Francois Gygi committed
3 4 5 6
// Copyright (c) 2008 The Regents of the University of California
//
// This file is part of Qbox
//
Francois Gygi committed
7 8
// Qbox is distributed under the terms of the GNU General Public License
// as published by the Free Software Foundation, either version 2 of
Francois Gygi committed
9 10 11 12 13 14
// the License, or (at your option) any later version.
// See the file COPYING in the root directory of this distribution
// or <http://www.gnu.org/licenses/>.
//
////////////////////////////////////////////////////////////////////////////////
//
Francois Gygi committed
15 16 17 18 19
// Matrix.C
//
////////////////////////////////////////////////////////////////////////////////

#include <cassert>
20
#include <vector>
Francois Gygi committed
21
#include <complex>
22
#include <limits>
Francois Gygi committed
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
#include <iostream>
using namespace std;

#ifdef USE_MPI
#include <mpi.h>
#endif

#include "Context.h"
#ifdef SCALAPACK
#include "blacs.h"
#endif

#include "Matrix.h"

#ifdef ADD_
#define numroc     numroc_
#define pdtran     pdtran_
#define pztranc    pztranc_
#define pdsymm     pdsymm_
#define pzsymm     pzsymm_
#define pzhemm     pzhemm_
#define pdgemm     pdgemm_
#define pzgemm     pzgemm_
#define pdsyrk     pdsyrk_
#define pzherk     pzherk_
#define pdsyr      pdsyr_
#define pdger      pdger_
50 51 52
#define pzgerc     pzgerc_
#define pzgeru     pzgeru_
#define pigemr2d   pigemr2d_
Francois Gygi committed
53 54 55 56 57 58 59 60 61 62 63 64
#define pdgemr2d   pdgemr2d_
#define pzgemr2d   pzgemr2d_
#define pdtrmm     pdtrmm_
#define pdtrsm     pdtrsm_
#define pztrsm     pztrsm_
#define pdtrtrs    pdtrtrs_
#define pdpotrf    pdpotrf_
#define pzpotrf    pzpotrf_
#define pdpotri    pdpotri_
#define pdpocon    pdpocon_
#define pdsygst    pdsygst_
#define pdsyev     pdsyev_
65 66
#define pdsyevd    pdsyevd_
#define pdsyevx    pdsyevx_
Francois Gygi committed
67
#define pzheev     pzheev_
68
#define pzheevd    pzheevd_
Francois Gygi committed
69 70 71 72
#define pdtrtri    pdtrtri_
#define pdlatra    pdlatra_
#define pdlacp2    pdlacp2_
#define pdlacp3    pdlacp3_
73
#define pdgetrf    pdgetrf_
74
#define pzgetrf    pzgetrf_
75
#define pdgetri    pdgetri_
76
#define pzgetri    pzgetri_
77 78
#define pdlapiv    pdlapiv_
#define pzlapiv    pzlapiv_
79 80
#define pdlapv2    pdlapv2_
#define pzlapv2    pzlapv2_
Francois Gygi committed
81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97

#define dscal      dscal_
#define zscal      zscal_
#define zdscal     zdscal_
#define dcopy      dcopy_
#define ddot       ddot_
#define zdotu      zdotu_
#define zdotc      zdotc_
#define daxpy      daxpy_
#define zaxpy      zaxpy_
#define dsymm      dsymm_
#define zsymm      zsymm_
#define zhemm      zhemm_
#define dgemm      dgemm_
#define zgemm      zgemm_
#define dsyr       dsyr_
#define dger       dger_
98 99
#define zgerc      zgerc_
#define zgeru      zgeru_
Francois Gygi committed
100 101 102 103 104 105 106 107 108 109 110 111 112 113 114
#define dsyrk      dsyrk_
#define zherk      zherk_
#define dtrmm      dtrmm_
#define dtrsm      dtrsm_
#define dtrtri     dtrtri_
#define ztrsm      ztrsm_
#define dtrtrs     dtrtrs_
#define dpotrf     dpotrf_
#define zpotrf     zpotrf_
#define dpotri     dpotri_
#define dpocon     dpocon_
#define dsygst     dsygst_
#define dsyev      dsyev_
#define zheev      zheev_
#define idamax     idamax_
115
#define dgetrf     dgetrf_
116
#define zgetrf     zgetrf_
117
#define dgetri     dgetri_
118
#define zgetri     zgetri_
Francois Gygi committed
119 120 121 122 123 124 125 126 127 128 129
#endif

extern "C"
{
  int numroc(int*, int*, int*, int*, int*);
#ifdef SCALAPACK
  // PBLAS
  void pdsymm(const char*, const char*, const int*, const int*, const double*,
       const double*, const int*, const int*, const int*,
       const double*, const int*, const int*, const int*,
       const double*, double*, const int*, const int*, const int*);
130
  void pzsymm(const char*, const char*, const int*, const int*,
Francois Gygi committed
131 132 133
       const complex<double>*,
       const complex<double>*, const int*, const int*, const int*,
       const complex<double>*, const int*, const int*, const int*,
134
       const complex<double>*, complex<double>*, const int*, const int*,
Francois Gygi committed
135
       const int*);
136
  void pzhemm(const char*, const char*, const int*, const int*,
Francois Gygi committed
137 138 139
       const complex<double>*,
       const complex<double>*, const int*, const int*, const int*,
       const complex<double>*, const int*, const int*, const int*,
140
       const complex<double>*, complex<double>*, const int*, const int*,
Francois Gygi committed
141
       const int*);
142
  void pdgemm(const char*, const char*, const int*,
Francois Gygi committed
143 144 145 146
       const int*, const int*, const double*,
       const double*, const int*, const int*, const int*,
       const double*, const int*, const int*, const int*,
       const double*, double*, const int*, const int*, const int*);
147
  void pzgemm(const char*, const char*, const int*,
Francois Gygi committed
148 149 150
       const int*, const int*, const complex<double>*,
       const complex<double>*, const int*, const int*, const int*,
       const complex<double>*, const int*, const int*, const int*,
151
       const complex<double>*, complex<double>*, const int*, const int*,
Francois Gygi committed
152 153 154 155 156
       const int*);
  void pdger(const int*, const int*, const double*,
       const double*, const int*, const int*, const int*, const int*,
       const double*, const int*, const int*, const int*, const int*,
       double*, const int*, const int*, const int*);
157 158 159 160 161 162 163 164
  void pzgerc(const int*, const int*, const complex<double>*,
       const complex<double>*, const int*, const int*, const int*, const int*,
       const complex<double>*, const int*, const int*, const int*, const int*,
       complex<double>*, const int*, const int*, const int*);
  void pzgeru(const int*, const int*, const complex<double>*,
       const complex<double>*, const int*, const int*, const int*, const int*,
       const complex<double>*, const int*, const int*, const int*, const int*,
       complex<double>*, const int*, const int*, const int*);
165
  void pdsyr(const char*, const int*,
Francois Gygi committed
166 167
       const double*, const double*, const int*, const int*, const int*,
       const int*, double*, const int*, const int*, const int*);
168
  void pdsyrk(const char*, const char*, const int*, const int*,
Francois Gygi committed
169 170
       const double*, const double*, const int*, const int*, const int*,
       const double*, double*, const int*, const int*, const int*);
171
  void pzherk(const char*, const char*, const int*, const int*,
172
       const double*, const complex<double>*, const int*,
Francois Gygi committed
173
       const int*, const int*,
174
       const double*, complex<double>*, const int*,
Francois Gygi committed
175 176 177 178 179 180
       const int*, const int*);
  void pdtran(const int*,const  int*, const double*,
       const double*, const int*, const int*, const int*,
       double*, const double*, const int*, const int*, const int*);
  void pztranc(const int*, const int*, const complex<double>*,
       const complex<double>*, const int*, const int*, const int*,
181
       complex<double>*, const complex<double>*, const int*, const int*,
Francois Gygi committed
182
       const int*);
183 184
  void pdtrmm(const char*, const char*, const char*, const char*,
       const int*, const int*, const double*,
Francois Gygi committed
185 186
       const double*, const int*, const int*, const int*,
       double*, const int*, const int*, const int*);
187 188
  void pdtrsm(const char*, const char*, const char*, const char*,
       const int*, const int*, const double*,
Francois Gygi committed
189 190
       const double*, const int*, const int*, const int*,
       double*, const int*, const int*, const int*);
191 192
  void pztrsm(const char*, const char*, const char*, const char*,
       const int*, const int*, const complex<double>*,
Francois Gygi committed
193 194 195 196
       const complex<double>*, const int*, const int*, const int*,
       complex<double>*, const int*, const int*, const int*);
  double pdlatra(const int*,const double*,const int*,const int*,const int*);
  // SCALAPACK
197
  void pdtrtrs(const char*, const char*, const char*, const int*, const int*,
Francois Gygi committed
198 199
               const double*, const int*, const int*, const int*,
               double*, const int*, const int*, const int*, int*);
200 201 202
  void pigemr2d(const int*,const int*,
                const int*,const int*,const int*, const int*,
                int*,const int*,const int*,const int*,const int*);
Francois Gygi committed
203 204 205 206 207 208 209 210
  void pdgemr2d(const int*,const int*,
                const double*,const int*,const int*, const int*,
                double*,const int*,const int*,const int*,const int*);
  void pzgemr2d(const int*,const int*,
                const complex<double>*,const int*,const int*, const int*,
                complex<double>*,const int*,const int*,const int*,const int*);
  void pdpotrf(const char*, const int*, double*, const int*,
               const int*, const int*, const int*);
211
  void pzpotrf(const char*, const int*, complex<double>*, const int*,
Francois Gygi committed
212
               const int*, const int*, const int*);
213
  void pdpotri(const char*, const int*, double*, const int*,
Francois Gygi committed
214
               const int*, const int*, const int*);
215
  void pdpocon(const char*, const int*, const double*,
Francois Gygi committed
216 217
               const int*, const int*, const int*, const double*, double*,
               double*, const int*, int*, const int*, int*);
218
  void pdsygst(const int*, const char*, const int*, double*,
Francois Gygi committed
219 220
               const int*, const int*, const int*, const double*, const int*,
               const int*, const int*, double*, int*);
221
  void pdsyev(const char*, const char*, const int*,
Francois Gygi committed
222 223
              double*, const int*, const int*, const int*, double*, double*,
              const int*, const int*, const int*, double*, const int*, int*);
224
  void pdsyevd(const char*, const char*, const int*,
225 226 227 228
              double*, const int*, const int*, const int*, double*, double*,
              const int*, const int*, const int*, double*, const int*, int*,
              int*, int*);
  void pdsyevx(const char* jobz, const char* range, const char* uplo,
229
               const int* n, double* a, const int* ia, const int* ja,
230 231 232 233 234
               const int* desca, double* vl, double* vu,
               const int* il, const int* iu, double* abstol,
               int* nfound, int* nz, double* w,
               const double* orfac, double* z, const int* iz, const int* jz,
               const int* descz, double* work, const int* lwork,
235
               int* iwork, int* liwork, int* ifail,
236
               int* icluster, double* gap, int* info);
237 238
  void pzheev(const char* jobz, const char* uplo, const int* n,
              complex<double>* a, const int* ia, const int* ja,
Francois Gygi committed
239
              const int* desca, double* w, complex<double> *z,
240
              const int* iz, const int* jz, const int* descz,
241 242 243 244 245 246 247 248 249
              complex<double>* work, int* lwork,
              double* rwork, int* lrwork, int* info);
  void pzheevd(const char* jobz, const char* uplo, const int* n,
               complex<double>* a, const int* ia, const int* ja,
               const int* desca, double* w, complex<double> *z,
               const int* iz, const int* jz, const int* descz,
               complex<double>* work, int* lwork,
               double* rwork, const int* lrwork,
               int* iwork, int* liwork, int* info);
250
  void pdtrtri(const char*, const char*, const int*, double*,
Francois Gygi committed
251
               const int*, const int*, const int*, int*);
252
  void pdgetrf(const int* m, const int* n, double* val,
253
               int* ia, const int* ja, const int* desca, int* ipiv, int* info);
254 255
  void pzgetrf(const int* m, const int* n, complex<double>* val,
               int* ia, const int* ja, const int* desca, int* ipiv, int* info);
256 257
  void pdgetri(const int* n, double* val,
               const int* ia, const int* ja, int* desca, int* ipiv,
258
               double* work, int* lwork, int* iwork, int* liwork, int* info);
259 260 261
  void pzgetri(const int* n, complex<double>* val, const int* ia,
               const int* ja, int* desca, int* ipiv, complex<double>* work,
               int* lwork, int* iwork, int* liwork, int* info);
262

263 264 265 266 267 268 269 270
  void pdlapiv(const char* direc, const char* rowcol, const char* pivroc,
               const int* m, const int* n, double *a, const int* ia,
               const int* ja, const int* desca, int* ipiv, const int* ip,
               const int* jp, const int* descp, int* iwork);
  void pzlapiv(const char* direc, const char* rowcol, const char* pivroc,
               const int* m, const int* n, complex<double> *a, const int* ia,
               const int* ja, const int* desca, int* ipiv, const int* ip,
               const int* jp, const int* descp, int* iwork);
271 272 273 274 275 276 277 278
  void pdlapv2(const char* direc, const char *rowcol,
               const int* m, const int *n, double *val,
               const int *ia, const int *ja, const int* desca,
               int *ipiv, const int *ip, const int *jp, const int *descp);
  void pzlapv2(const char* direc, const char *rowcol,
               const int* m, const int *n, complex<double> *val,
               const int *ia, const int *ja, const int* desca,
               int *ipiv, const int *ip, const int *jp, const int *descp);
279

Francois Gygi committed
280 281 282 283 284
#endif
  // BLAS1
  void dscal(const int*, const double*, double*, const int*);
  void zscal(const int*, const complex<double>*, complex<double>*, const int*);
  void zdscal(const int*, const double*, complex<double>*, const int*);
285
  void daxpy(const int *, const double *, const double *, const int *,
Francois Gygi committed
286
             double *, const int *);
287
  void zaxpy(const int *, const complex<double> *, const complex<double> *,
Francois Gygi committed
288 289
             const int *, complex<double> *, const int *);
  void dcopy(const int *, const double*, const int *, double*, const int*);
290
  double ddot(const int *, const double *, const int *,
Francois Gygi committed
291
              const double *, const int *);
292
  complex<double> zdotc(const int *, const complex<double>*, const int *,
Francois Gygi committed
293
                        const complex<double>*, const int *);
294
  complex<double> zdotu(const int *, const complex<double>*, const int *,
Francois Gygi committed
295 296 297 298
                        const complex<double>*, const int *);
  int idamax(const int *, const double*, const int*);
  // BLAS3
  void dsymm(const char*, const char*, const int*, const int *,
299
             const double*, const double*, const int*,
Francois Gygi committed
300 301 302
             const double*, const int*,
             const double*, double*, const int*);
  void zsymm(const char*, const char*, const int*, const int *,
303
             const complex<double>*, const complex<double>*, const int*,
Francois Gygi committed
304 305 306
             const complex<double>*, const int*,
             const complex<double>*, complex<double>*, const int*);
  void zhemm(const char*, const char*, const int*, const int *,
307
             const complex<double>*, const complex<double>*, const int*,
Francois Gygi committed
308 309 310
             const complex<double>*, const int*,
             const complex<double>*, complex<double>*, const int*);
  void dgemm(const char*, const char*, const int*, const int *, const int*,
311
             const double*, const double*, const int*,
Francois Gygi committed
312 313 314
             const double*, const int*,
             const double*, double*, const int*);
  void zgemm(const char*, const char*, const int*, const int *, const int*,
315
             const complex<double>*, const complex<double>*, const int*,
Francois Gygi committed
316 317
             const complex<double>*, const int*,
             const complex<double>*, complex<double>*, const int*);
Francois Gygi committed
318
  void zgerc(const int*, const int *, const complex<double>*,
319 320 321
             const complex<double>*, const int*,
             const complex<double>*, const int*,
             const complex<double>*, const int*);
Francois Gygi committed
322
  void zgeru(const int*, const int *, const complex<double>*,
323 324 325
             const complex<double>*, const int*,
             const complex<double>*, const int*,
             const complex<double>*, const int*);
326
  void dger(const int *, const int*, const double *,
Francois Gygi committed
327 328
            const double *, const int *, const double *, const int *,
            double*, const int*);
329
  void dsyr(const char*, const int *, const double *,
Francois Gygi committed
330 331
            const double *, const int *, double *, const int *);
  void dsyrk(const char*, const char*, const int *, const int *,
332
             const double *, const double *, const int *,
Francois Gygi committed
333 334
             const double *, double *, const int *);
  void zherk(const char* uplo, const char* trans, const int* n, const int* k,
335
             const double* alpha, const complex<double>* a,
336
             const int*  lda,
337
             const double* beta, complex<double>* c, const int* ldc);
338 339
  void dtrmm(const char*, const char*, const char*, const char*,
             const int*, const int *, const double*, const double*,
Francois Gygi committed
340
             const int*, double*, const int*);
341 342
  void dtrsm(const char*, const char*, const char*, const char*,
             const int*, const int *, const double*, const double*,
Francois Gygi committed
343
             const int*, double*, const int*);
344 345
  void ztrsm(const char*, const char*, const char*, const char*,
             const int*, const int *, const complex<double>*,
Francois Gygi committed
346 347
             const complex<double>*, const int*, complex<double>*, const int*);
  // LAPACK
348 349
  void dtrtrs(const char*, const char*, const char*,
              const int*, const int*, const double*, const int*,
Francois Gygi committed
350 351 352 353
              double*, const int*, int*);
  void dpotrf(const char*, const int*, double*, const int*, int*);
  void zpotrf(const char*, const int*, complex<double>*, const int*, int*);
  void dpotri(const char*, const int*, double*, const int*, int*);
354
  void dpocon(const char*, const int *, const double *, const int *,
Francois Gygi committed
355
              const double *, double *, double *, const int *, int *);
356
  void dsygst(const int*, const char*, const int*,
Francois Gygi committed
357
              double*, const int*, const double*, const int*, int*);
358
  void dsyev(const char* jobz, const char* uplo, const int* n, double* a,
Francois Gygi committed
359
             const int *lda, double *w, double*work,
360
             int *lwork, int *info);
361
  void zheev(const char* jobz, const char* uplo, const int *n,
Francois Gygi committed
362
             complex<double>* a, const int *lda, double* w,
363
             complex<double>* work, int *lwork, double* rwork, int *info);
Francois Gygi committed
364
  void dtrtri(const char*, const char*, const int*, double*, const int*, int* );
365
  void dgetrf(const int* m, const int* n, double* a, const int* lda,
366
              int* ipiv, int*info);
367 368
  void zgetrf(const int* m, const int* n, complex<double>* a, const int* lda,
              int* ipiv, int*info);
369
  void dgetri(const int* m, double* val, const int* lda, int* ipiv,
370
              double* work, int* lwork, int* info);
371 372
  void zgetri(const int* m, complex<double>* val, const int* lda, int* ipiv,
              complex<double>* work, int* lwork, int* info);
Francois Gygi committed
373 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
}


#ifndef SCALAPACK
int numroc(int* a, int* b, int* c, int* d, int* e)
{
  return *a;
}
#endif

////////////////////////////////////////////////////////////////////////////////
// reference constructor create a proxy for a ComplexMatrix rhs
DoubleMatrix::DoubleMatrix(ComplexMatrix& rhs) : ctxt_(rhs.context()),
  reference_(true)
{
  int new_m = 2 * rhs.m();
  int new_mb = 2 * rhs.mb();
  init_size(new_m,rhs.n(),new_mb,rhs.nb());
  val = (double*) rhs.valptr();
}

////////////////////////////////////////////////////////////////////////////////
// reference constructor create a proxy for a const ComplexMatrix rhs
DoubleMatrix::DoubleMatrix(const ComplexMatrix& rhs) : ctxt_(rhs.context()),
  reference_(true)
{
  int new_m = 2 * rhs.m();
  int new_mb = 2 * rhs.mb();
  init_size(new_m,rhs.n(),new_mb,rhs.nb());
  val = (double*) rhs.cvalptr();
}

////////////////////////////////////////////////////////////////////////////////
// reference constructor create a proxy for a DoubleMatrix rhs
ComplexMatrix::ComplexMatrix(DoubleMatrix& rhs) : ctxt_(rhs.context()),
  reference_(true)
{
  assert(rhs.m()%2 == 0);
  int new_m = rhs.m() / 2;
  assert(rhs.mb()%2 == 0);
  int new_mb = rhs.mb() / 2;
  init_size(new_m,rhs.n(),new_mb,rhs.nb());
  val = (complex<double>*) rhs.valptr();
}
417

Francois Gygi committed
418 419 420 421 422 423 424 425 426 427 428 429
////////////////////////////////////////////////////////////////////////////////
// reference constructor create a proxy for a const DoubleMatrix rhs
ComplexMatrix::ComplexMatrix(const DoubleMatrix& rhs) : ctxt_(rhs.context()),
  reference_(true)
{
  assert(rhs.m()%2 == 0);
  int new_m = rhs.m() / 2;
  assert(rhs.mb()%2 == 0);
  int new_mb = rhs.mb() / 2;
  init_size(new_m,rhs.n(),new_mb,rhs.nb());
  val = (complex<double>*) rhs.cvalptr();
}
430

Francois Gygi committed
431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446
////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::init_size(int m, int n, int mb, int nb)
{
  assert(m>=0);
  assert(n>=0);
  assert(mb>=0);
  assert(nb>=0);
  m_ = m;
  n_ = n;
#ifdef SCALAPACK
  mb_ = mb;
  nb_ = nb;
#else
  mb_ = m;
  nb_ = n;
#endif
447 448
  if ( mb_ == 0 ) mb_ = 1;
  if ( nb_ == 0 ) nb_ = 1;
Francois Gygi committed
449 450 451 452 453 454 455
  ictxt_ = ctxt_.ictxt();
  nprow_ = ctxt_.nprow();
  npcol_ = ctxt_.npcol();
  myrow_ = ctxt_.myrow();
  mycol_ = ctxt_.mycol();
  active_ = myrow_ >= 0;
  int isrcproc=0;
456 457
  mloc_ = 0;
  nloc_ = 0;
Francois Gygi committed
458 459 460 461 462
  if ( m_ != 0 )
    mloc_ = numroc(&m_,&mb_,&myrow_,&isrcproc,&nprow_);
  if ( n_ != 0 )
    nloc_ = numroc(&n_,&nb_,&mycol_,&isrcproc,&npcol_);
  size_ = mloc_ * nloc_;
463

Francois Gygi committed
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 495 496 497
  // set leading dimension of val array to mloc_;
  lld_ = mloc_;
  if ( lld_ == 0 ) lld_ = 1;

  // total and local number of blocks
  mblocks_ = 0;
  nblocks_ = 0;
  m_incomplete_ = false;
  n_incomplete_ = false;
  if ( active_ && mb_ > 0 && nb_ > 0 )
  {
    mblocks_ = ( mloc_ + mb_ - 1 ) / mb_;
    nblocks_ = ( nloc_ + nb_ - 1 ) / nb_;
    m_incomplete_ = mloc_ % mb_ != 0;
    n_incomplete_ = nloc_ % nb_ != 0;
  }

  if ( active_ )
  {
    desc_[0] = 1;
  }
  else
  {
    desc_[0] = -1;
  }
  desc_[1] = ictxt_;
  desc_[2] = m_;
  desc_[3] = n_;
  desc_[4] = mb_;
  desc_[5] = nb_;
  desc_[6] = 0;
  desc_[7] = 0;
  desc_[8] = lld_;
}
498

Francois Gygi committed
499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514
////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::init_size(int m, int n, int mb, int nb)
{
  assert(m>=0);
  assert(n>=0);
  assert(mb>=0);
  assert(nb>=0);
  m_ = m;
  n_ = n;
#ifdef SCALAPACK
  mb_ = mb;
  nb_ = nb;
#else
  mb_ = m;
  nb_ = n;
#endif
515 516
  if ( mb_ == 0 ) mb_ = 1;
  if ( nb_ == 0 ) nb_ = 1;
Francois Gygi committed
517 518 519 520 521 522 523 524 525
  ictxt_ = ctxt_.ictxt();
  nprow_ = ctxt_.nprow();
  npcol_ = ctxt_.npcol();
  myrow_ = ctxt_.myrow();
  mycol_ = ctxt_.mycol();
  active_ = myrow_ >= 0;
  int isrcproc=0;
  mloc_ = 0;
  nloc_ = 0;
526

Francois Gygi committed
527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565
  if ( m_ != 0 )
    mloc_ = numroc(&m_,&mb_,&myrow_,&isrcproc,&nprow_);
  if ( n_ != 0 )
    nloc_ = numroc(&n_,&nb_,&mycol_,&isrcproc,&npcol_);
  size_ = mloc_ * nloc_;

  // set leading dimension of val array to mloc_;
  lld_ = mloc_;
  if ( lld_ == 0 ) lld_ = 1;

  // total and local number of blocks
  mblocks_ = 0;
  nblocks_ = 0;
  m_incomplete_ = false;
  n_incomplete_ = false;
  if ( active_ && mb_ > 0 && nb_ > 0 )
  {
    mblocks_ = ( mloc_ + mb_ - 1 ) / mb_;
    nblocks_ = ( nloc_ + nb_ - 1 ) / nb_;
    m_incomplete_ = mloc_ % mb_ != 0;
    n_incomplete_ = nloc_ % nb_ != 0;
  }

  if ( active_ )
  {
    desc_[0] = 1;
  }
  else
  {
    desc_[0] = -1;
  }
  desc_[1] = ictxt_;
  desc_[2] = m_;
  desc_[3] = n_;
  desc_[4] = mb_;
  desc_[5] = nb_;
  desc_[6] = 0;
  desc_[7] = 0;
  desc_[8] = lld_;
566 567
}

Francois Gygi committed
568 569 570
////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::clear(void)
{
571
  assert(val!=0||size_==0);
Francois Gygi committed
572 573 574 575 576 577
  memset(val,0,size_*sizeof(double));
}

////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::clear(void)
{
578
  assert(val!=0||size_==0);
Francois Gygi committed
579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664
  memset(val,0,size_*sizeof(complex<double>));
}

////////////////////////////////////////////////////////////////////////////////
// real identity: initialize matrix to identity
////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::identity(void)
{
  clear();
  set('d',1.0);
}

////////////////////////////////////////////////////////////////////////////////
// complex identity: initialize matrix to identity
////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::identity(void)
{
  clear();
  set('d',complex<double>(1.0,0.0));
}

////////////////////////////////////////////////////////////////////////////////
// set value of diagonal or off-diagonal elements to a constant
// uplo=='u': set strictly upper part to x
// uplo=='l': set strictly lower part to x
// uplo=='d': set diagonal to x
////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::set(char uplo, double xx)
{
  if ( active_ )
  {
    if ( uplo=='l' || uplo=='L' )
    {
      // initialize strictly lower part
      for (int li=0; li < mblocks_;li++)
      {
        for (int lj=0; lj < nblocks_;lj++)
        {
          for (int ii=0; ii < mbs(li); ii++)
          {
            for (int jj=0; jj < nbs(lj);jj++)
            {
              if ( i(li,ii) > j(lj,jj) )
                val[ (ii+li*mb_)+(jj+lj*nb_)*mloc_ ] = xx;
            }
          }
        }
      }
    }
    else if ( uplo=='u' || uplo=='U' )
    {
      // initialize strictly upper part
      for ( int li=0; li < mblocks_; li++ )
      {
        for ( int lj=0; lj < nblocks_; lj++ )
        {
          for ( int ii=0; ii < mbs(li); ii++ )
          {
            for ( int jj=0; jj < nbs(lj); jj++ )
            {
              if ( i(li,ii) < j(lj,jj) )
                val[ (ii+li*mb_)+(jj+lj*nb_)*mloc_ ] = xx;
            }
          }
        }
      }
    }
    else if ( uplo=='d' || uplo=='D' )
    {
      // initialize diagonal elements
      if ( active() )
      {
        // loop through all local blocks (ll,mm)
        for ( int ll = 0; ll < mblocks(); ll++)
        {
          for ( int mm = 0; mm < nblocks(); mm++)
          {
            // check if block (ll,mm) has diagonal elements
            int imin = i(ll,0);
            int imax = imin + mbs(ll)-1;
            int jmin = j(mm,0);
            int jmax = jmin + nbs(mm)-1;
            // cout << " process (" << myrow_ << "," << mycol_ << ")"
            // << " block (" << ll << "," << mm << ")"
            // << " imin/imax=" << imin << "/" << imax
            // << " jmin/jmax=" << jmin << "/" << jmax << endl;
665

Francois Gygi committed
666 667 668 669 670
            if ((imin <= jmax) && (imax >= jmin))
            {
              // block (ll,mm) holds diagonal elements
              int idiagmin = max(imin,jmin);
              int idiagmax = min(imax,jmax);
671

Francois Gygi committed
672 673 674
              // cout << " process (" << myrow_ << "," << mycol_ << ")"
              // << " holds diagonal elements " << idiagmin << " to " <<
              // idiagmax << " in block (" << ll << "," << mm << ")" << endl;
675

Francois Gygi committed
676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760
              for ( int ii = idiagmin; ii <= idiagmax; ii++ )
              {
                // access element (ii,ii)
                int jj = ii;
                int iii = ll * mb_ + x(ii);
                int jjj = mm * nb_ + y(jj);
                val[iii+mloc_*jjj] = xx;
              }
            }
          }
        }
      }
    }
    else
    {
      cout << " DoubleMatrix::set: invalid argument" << endl;
#ifdef USE_MPI
      MPI_Abort(MPI_COMM_WORLD,2);
#else
      exit(2);
#endif
    }
  }
}

////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::set(char uplo, complex<double> xx)
{
  if ( active_ )
  {
    if ( uplo=='l' || uplo=='L' )
    {
      // initialize strictly lower part
      for (int li=0; li < mblocks_;li++)
      {
        for (int lj=0; lj < nblocks_;lj++)
        {
          for (int ii=0; ii < mbs(li); ii++)
          {
            for (int jj=0; jj < nbs(lj);jj++)
            {
              if ( i(li,ii) > j(lj,jj) )
                val[ (ii+li*mb_)+(jj+lj*nb_)*mloc_ ] = xx;
            }
          }
        }
      }
    }
    else if ( uplo=='u' || uplo=='U' )
    {
      // initialize strictly upper part
      for ( int li=0; li < mblocks_; li++ )
      {
        for ( int lj=0; lj < nblocks_; lj++ )
        {
          for ( int ii=0; ii < mbs(li); ii++ )
          {
            for ( int jj=0; jj < nbs(lj); jj++ )
            {
              if ( i(li,ii) < j(lj,jj) )
                val[ (ii+li*mb_)+(jj+lj*nb_)*mloc_ ] = xx;
            }
          }
        }
      }
    }
    else if ( uplo=='d' || uplo=='D' )
    {
      // initialize diagonal elements
      if ( active() )
      {
        // loop through all local blocks (ll,mm)
        for ( int ll = 0; ll < mblocks(); ll++)
        {
          for ( int mm = 0; mm < nblocks(); mm++)
          {
            // check if block (ll,mm) has diagonal elements
            int imin = i(ll,0);
            int imax = imin + mbs(ll)-1;
            int jmin = j(mm,0);
            int jmax = jmin + nbs(mm)-1;
            // cout << " process (" << myrow_ << "," << mycol_ << ")"
            // << " block (" << ll << "," << mm << ")"
            // << " imin/imax=" << imin << "/" << imax
            // << " jmin/jmax=" << jmin << "/" << jmax << endl;
761

Francois Gygi committed
762 763 764 765 766
            if ((imin <= jmax) && (imax >= jmin))
            {
              // block (ll,mm) holds diagonal elements
              int idiagmin = max(imin,jmin);
              int idiagmax = min(imax,jmax);
767

Francois Gygi committed
768 769 770
              // cout << " process (" << myrow_ << "," << mycol_ << ")"
              // << " holds diagonal elements " << idiagmin << " to " <<
              // idiagmax << " in block (" << ll << "," << mm << ")" << endl;
771

Francois Gygi committed
772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868
              for ( int ii = idiagmin; ii <= idiagmax; ii++ )
              {
                // access element (ii,ii)
                int jj = ii;
                int iii = ll * mb_ + x(ii);
                int jjj = mm * nb_ + y(jj);
                val[iii+mloc_*jjj] = xx;
              }
            }
          }
        }
      }
    }
    else
    {
      cout << " DoubleMatrix::set: invalid argument" << endl;
#ifdef USE_MPI
      MPI_Abort(MPI_COMM_WORLD,2);
#else
      exit(2);
#endif
    }
  }
}

////////////////////////////////////////////////////////////////////////////////
// initialize *this using a replicated matrix a
void DoubleMatrix::init(const double* const a, int lda)
{
  if ( active_ )
  {
    for ( int li=0; li < mblocks_; li++ )
    {
      for ( int lj=0; lj < nblocks_; lj++ )
      {
        for ( int ii=0; ii < mbs(li); ii++ )
        {
          for ( int jj=0; jj < nbs(lj); jj++ )
          {
            val[ (ii+li*mb_)+(jj+lj*nb_)*mloc_ ]
                = a[ i(li,ii) + j(lj,jj)*lda ];
          }
        }
      }
    }
  }
}

////////////////////////////////////////////////////////////////////////////////
double DoubleMatrix::dot(const DoubleMatrix &x) const
{
  assert( ictxt_ == x.ictxt() );
  double  sum=0.;
  double  tsum=0.;
  if ( active_ )
  {
    assert( m_ == x.m() );
    assert( n_ == x.n() );
    assert( mb_ == x.mb() );
    assert( nb_ == x.nb() );
    assert( mloc_ == x.mloc() );
    assert( nloc_ == x.nloc() );
    assert(size_==x.size());
    int ione=1;
    tsum=ddot(&size_, val, &ione, x.val, &ione);
  }
#ifdef SCALAPACK
  if ( active_ )
    MPI_Allreduce(&tsum, &sum, 1, MPI_DOUBLE, MPI_SUM, ctxt_.comm() );
#else
  sum=tsum;
#endif
  return sum;
}

////////////////////////////////////////////////////////////////////////////////
complex<double> ComplexMatrix::dot(const ComplexMatrix &x) const
{
  assert( ictxt_ == x.ictxt() );
  complex<double>  sum=0.0;
  complex<double>  tsum=0.0;
  if ( active_ )
  {
    assert( m_ == x.m() );
    assert( n_ == x.n() );
    assert( mb_ == x.mb() );
    assert( nb_ == x.nb() );
    assert( mloc_ == x.mloc() );
    assert( nloc_ == x.nloc() );
    assert(size_==x.size());
    //int ione=1;
    //tsum=zdotc(&size_, val, &ione, x.val, &ione);
    for ( int i = 0; i < size_; i++ )
      tsum += conj(val[i]) * x.val[i];
  }
#ifdef SCALAPACK
  if ( active_ )
869
    MPI_Allreduce((double*)&tsum, (double*)&sum, 2,
Francois Gygi committed
870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898
                  MPI_DOUBLE, MPI_SUM, ctxt_.comm() );
#else
  sum=tsum;
#endif
  return sum;
}

////////////////////////////////////////////////////////////////////////////////
complex<double> ComplexMatrix::dotu(const ComplexMatrix &x) const
{
  assert( ictxt_ == x.ictxt() );
  complex<double>  sum=0.0;
  complex<double>  tsum=0.0;
  if ( active_ )
  {
    assert( m_ == x.m() );
    assert( n_ == x.n() );
    assert( mb_ == x.mb() );
    assert( nb_ == x.nb() );
    assert( mloc_ == x.mloc() );
    assert( nloc_ == x.nloc() );
    assert(size_==x.size());
    //int ione=1;
    //tsum=zdotu(&size_, val, &ione, x.val, &ione);
    for ( int i = 0; i < size_; i++ )
      tsum += val[i] * x.val[i];
  }
#ifdef SCALAPACK
  if ( active_ )
899
    MPI_Allreduce((double*)&tsum, (double*)&sum, 2,
Francois Gygi committed
900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990
                  MPI_DOUBLE, MPI_SUM, ctxt_.comm() );
#else
  sum=tsum;
#endif
  return sum;
}

////////////////////////////////////////////////////////////////////////////////
double DoubleMatrix::amax(void) const
{
  double am = 0.0, tam = 0.0;
  if ( active_ )
  {
    int ione=1;
    tam = val[idamax(&size_,val,&ione) - 1];
  }
#ifdef SCALAPACK
  if ( active_ )
    MPI_Allreduce(&tam, &am, 1, MPI_DOUBLE, MPI_MAX, ctxt_.comm() );
#else
  am=tam;
#endif
  return am;
}

////////////////////////////////////////////////////////////////////////////////
// axpy: *this = *this + alpha * x
void DoubleMatrix::axpy(double alpha, const DoubleMatrix &x)
{
  assert( ictxt_ == x.ictxt() );
  int ione=1;
  assert(m_==x.m());
  assert(n_==x.n());
  assert(mloc_==x.mloc());
  assert(nloc_==x.nloc());
  if( active_ )
    daxpy(&size_, &alpha, x.val, &ione, val, &ione);
}

////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::axpy(complex<double> alpha, const ComplexMatrix &x)
{
  assert( ictxt_ == x.ictxt() );
  int ione=1;
  assert(m_==x.m());
  assert(n_==x.n());
  assert(mloc_==x.mloc());
  assert(nloc_==x.nloc());
  if( active_ )
    zaxpy(&size_, &alpha, x.val, &ione, val, &ione);
}

////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::axpy(double alpha, const ComplexMatrix &x)
{
  assert( ictxt_ == x.ictxt() );
  int ione=1;
  assert(m_==x.m());
  assert(n_==x.n());
  assert(mloc_==x.mloc());
  assert(nloc_==x.nloc());
  int len = 2 * size_;
  if( active_ )
    daxpy(&len, &alpha, (double*) x.val, &ione, (double*) val, &ione);
}

////////////////////////////////////////////////////////////////////////////////
// real getsub: *this = sub(A)
// copy submatrix A(ia:ia+m, ja:ja+n) into *this;
// *this and A may live in different contexts
void DoubleMatrix::getsub(const DoubleMatrix &a,
  int m, int n, int ia, int ja)
{
#if SCALAPACK
  int iap=ia+1;
  int jap=ja+1;
  assert(n<=n_);
  assert(n<=a.n());
  assert(m<=m_);
  assert(m<=a.m());
  int ione = 1;
  int gictxt;
  Cblacs_get( 0, 0, &gictxt );
  pdgemr2d(&m,&n,a.val,&iap,&jap,a.desc_,val,&ione,&ione,desc_,&gictxt);
#else
  for ( int j = 0; j < n; j++ )
    for ( int i = 0; i < m; i++ )
      val[i+j*m_] = a.val[(i+ia) + (j+ja)*a.m()];
#endif
}

991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016
////////////////////////////////////////////////////////////////////////////////
// real getsub: *this = sub(A)
// copy submatrix A(ia:ia+m, ja:ja+n) into *this(idest:idest+m,jdest:jdest+n)
// *this and A may live in different contexts
void DoubleMatrix::getsub(const DoubleMatrix &a,
  int m, int n, int isrc, int jsrc, int idest, int jdest)
{
#if SCALAPACK
  int iap=isrc+1;
  int jap=jsrc+1;
  int idp=idest+1;
  int jdp=jdest+1;
  assert(n<=n_);
  assert(n<=a.n());
  assert(m<=m_);
  assert(m<=a.m());
  int gictxt;
  Cblacs_get( 0, 0, &gictxt );
  pdgemr2d(&m,&n,a.val,&iap,&jap,a.desc_,val,&idp,&jdp,desc_,&gictxt);
#else
  for ( int j = 0; j < n; j++ )
    for ( int i = 0; i < m; i++ )
      val[(idest+i)+(jdest+j)*m_] = a.val[(i+isrc) + (j+jsrc)*a.m()];
#endif
}

Francois Gygi committed
1017 1018
////////////////////////////////////////////////////////////////////////////////
// complex getsub: *this = sub(A)
1019
// copy submatrix A(ia:ia+m, ja:ja+n) into *this
Francois Gygi committed
1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041
// *this and A may live in different contexts
void ComplexMatrix::getsub(const ComplexMatrix &a,
  int m, int n, int ia, int ja)
{
#if SCALAPACK
  int iap=ia+1;
  int jap=ja+1;
  assert(n<=n_);
  assert(n<=a.n());
  assert(m<=m_);
  assert(m<=a.m());
  int ione = 1;
  int gictxt;
  Cblacs_get( 0, 0, &gictxt );
  pzgemr2d(&m,&n,a.val,&iap,&jap,a.desc_,val,&ione,&ione,desc_,&gictxt);
#else
  for ( int j = 0; j < n; j++ )
    for ( int i = 0; i < m; i++ )
      val[i+j*m_] = a.val[(i+ia) + (j+ja)*a.m()];
#endif
}

1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067
////////////////////////////////////////////////////////////////////////////////
// complex getsub: *this = sub(A)
// copy submatrix A(ia:ia+m, ja:ja+n) into *this(idest:idest+m,jdest:jdest+n)
// *this and A may live in different contexts
void ComplexMatrix::getsub(const ComplexMatrix &a,
  int m, int n, int isrc, int jsrc, int idest, int jdest)
{
#if SCALAPACK
  int iap=isrc+1;
  int jap=jsrc+1;
  int idp=idest+1;
  int jdp=jdest+1;
  assert(n<=n_);
  assert(n<=a.n());
  assert(m<=m_);
  assert(m<=a.m());
  int gictxt;
  Cblacs_get( 0, 0, &gictxt );
  pzgemr2d(&m,&n,a.val,&iap,&jap,a.desc_,val,&idp,&jdp,desc_,&gictxt);
#else
  for ( int j = 0; j < n; j++ )
    for ( int i = 0; i < m; i++ )
      val[(idest+i)+(jdest+j)*m_] = a.val[(i+isrc) + (j+jsrc)*a.m()];
#endif
}

Francois Gygi committed
1068 1069 1070 1071 1072 1073 1074 1075
////////////////////////////////////////////////////////////////////////////////
// real matrix transpose
// this = alpha * transpose(A) + beta * this
////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::transpose(double alpha, const DoubleMatrix& a, double beta)
{
  assert(this != &a);
  assert( ictxt_ == a.ictxt() );
1076

Francois Gygi committed
1077 1078 1079 1080
  if ( active() )
  {
    assert(a.m() == n_);
    assert(a.n() == m_);
1081

Francois Gygi committed
1082 1083 1084 1085 1086 1087 1088 1089 1090
#ifdef SCALAPACK
    int ione = 1;
    pdtran(&m_, &n_, &alpha,
         a.val, &ione, &ione, a.desc_,
         &beta, val, &ione, &ione, desc_);
#else
    scal(beta);
    for ( int i=0; i<m_; i++ )
      for ( int j=0; j<i; j++ )
1091
      {
Francois Gygi committed
1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119
        val[i*m_+j] += alpha * a.val[j*m_+i];
        val[j*m_+i] += alpha * a.val[i*m_+j];
      }
    for ( int i=0; i<m_; i++ )
      val[i*m_+i] += alpha * a.val[i*m_+i];
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// real matrix transpose
// *this = transpose(a)
////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::transpose(const DoubleMatrix& a)
{
  assert(this != &a);
  transpose(1.0,a,0.0);
}

////////////////////////////////////////////////////////////////////////////////
// complex hermitian transpose
// this = alpha * A^H + beta * this
////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::transpose(complex<double> alpha, const ComplexMatrix& a,
  complex<double> beta)
{
  assert(this != &a);
  assert( ictxt_ == a.ictxt() );
1120

Francois Gygi committed
1121 1122 1123 1124
  if ( active() )
  {
    assert(a.m() == n_);
    assert(a.n() == m_);
1125

Francois Gygi committed
1126 1127 1128 1129 1130 1131 1132 1133 1134
#ifdef SCALAPACK
    int ione = 1;
    pztranc(&m_, &n_, &alpha,
         a.val, &ione, &ione, a.desc_,
         &beta, val, &ione, &ione, desc_);
#else
    scal(beta);
    for ( int i=0; i<m_; i++ )
      for ( int j=0; j<i; j++ )
1135
      {
Francois Gygi committed
1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161
        val[i*m_+j] += alpha * conj(a.val[j*m_+i]);
        val[j*m_+i] += alpha * conj(a.val[i*m_+j]);
      }
    for ( int i=0; i<m_; i++ )
      val[i*m_+i] += alpha * a.val[i*m_+i];
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// complex matrix transpose
// *this = transpose(a)
////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::transpose(const ComplexMatrix& a)
{
  assert(this != &a);
  transpose(complex<double>(1.0,0.0),a,complex<double>(0.0,0.0));
}

////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::symmetrize(char uplo)
{
  // symmetrize
  // if uplo == 'l' : copy strictly lower triangle to strictly upper triangle
  // if uplo == 'u' : copy strictly upper triangle to strictly lower triangle
  // if uplo == 'n' : A = 0.5 * ( A^T + A )
1162

Francois Gygi committed
1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196
  if ( uplo == 'n' )
  {
    DoubleMatrix tmp(*this);
    transpose(0.5,tmp,0.5);
  }
  else if ( uplo == 'l' )
  {
    set('u',0.0);
    DoubleMatrix tmp(*this);
    tmp.set('d',0.0);
    transpose(1.0,tmp,1.0);
  }
  else if ( uplo == 'u' )
  {
    set('l',0.0);
    DoubleMatrix tmp(*this);
    tmp.set('d',0.0);
    transpose(1.0,tmp,1.0);
  }
  else
  {
    cout << " DoubleMatrix::symmetrize: invalid argument" << endl;
#ifdef USE_MPI
      MPI_Abort(MPI_COMM_WORLD, 2);
#else
      exit(2);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::symmetrize(char uplo)
{
  // symmetrize
1197 1198
  // uplo == 'l' : copy conjugate of strictly lower triangle to strictly upper
  // uplo == 'u' : copy conjugate of strictly upper triangle to strictly lower
Francois Gygi committed
1199
  // uplo == 'n' : A = 0.5 * ( A^H + A )
1200

Francois Gygi committed
1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252
  if ( uplo == 'n' )
  {
    ComplexMatrix tmp(*this);
    transpose(complex<double>(0.5,0.0),tmp,complex<double>(0.5,0.0));
  }
  else if ( uplo == 'l' )
  {
    set('u',complex<double>(0.0,0.0));
    ComplexMatrix tmp(*this);
    tmp.set('d',complex<double>(0.0,0.0));
    transpose(complex<double>(1.0,0.0),tmp,complex<double>(1.0,0.0));
  }
  else if ( uplo == 'u' )
  {
    set('l',complex<double>(0.0,0.0));
    ComplexMatrix tmp(*this);
    tmp.set('d',complex<double>(0.0,0.0));
    transpose(complex<double>(1.0,0.0),tmp,complex<double>(1.0,0.0));
  }
  else
  {
    cout << " ComplexMatrix::symmetrize: invalid argument" << endl;
#ifdef USE_MPI
      MPI_Abort(MPI_COMM_WORLD, 2);
#else
      exit(2);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
double DoubleMatrix::nrm2(void) const
{
  return sqrt(dot(*this));
}

////////////////////////////////////////////////////////////////////////////////
double ComplexMatrix::nrm2(void) const
{
  return abs(dot(*this));
}

////////////////////////////////////////////////////////////////////////////////
// rank-1 update using row kx of x and (row ky of y)^T
// *this = *this + alpha * x(kx) * y(ky)^T
void DoubleMatrix::ger(double alpha,
  const DoubleMatrix& x, int kx, const DoubleMatrix& y, int ky)
{
  assert(x.n()==m_);
  assert(y.n()==n_);
#if SCALAPACK
  int ione=1;
1253

Francois Gygi committed
1254 1255 1256
  int ix = kx+1;
  int jx = 1;
  int incx = x.m();
1257

Francois Gygi committed
1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271
  int iy = ky+1;
  int jy = 1;
  int incy = y.m();
  pdger(&m_,&n_,&alpha,x.val,&ix,&jx,x.desc_,&incx,
                       y.val,&iy,&jy,y.desc_,&incy,
                       val,&ione,&ione,desc_);
#else
  int incx = x.m();
  int incy = y.m();
  dger(&m_,&n_,&alpha,&x.val[kx*x.m()],&incx,
                      &y.val[ky*y.m()],&incy,val,&m_);
#endif
}

1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329
////////////////////////////////////////////////////////////////////////////////
// rank-1 update using row kx of x and conj(row ky of y)^T
// *this = *this + alpha * x(kx) * conj(y(ky))^T
void ComplexMatrix::gerc(complex<double> alpha,
  const ComplexMatrix& x, int kx, const ComplexMatrix& y, int ky)
{
  assert(x.n()==m_);
  assert(y.n()==n_);
#if SCALAPACK
  int ione=1;

  int ix = kx+1;
  int jx = 1;
  int incx = x.m();

  int iy = ky+1;
  int jy = 1;
  int incy = y.m();
  pzgerc(&m_,&n_,&alpha,x.val,&ix,&jx,x.desc_,&incx,
                       y.val,&iy,&jy,y.desc_,&incy,
                       val,&ione,&ione,desc_);
#else
  int incx = x.m();
  int incy = y.m();
  zgerc(&m_,&n_,&alpha,&x.val[kx*x.m()],&incx,
                      &y.val[ky*y.m()],&incy,val,&m_);
#endif
}

////////////////////////////////////////////////////////////////////////////////
// rank-1 update using row kx of x and conj(row ky of y)^T
// *this = *this + alpha * x(kx) * y(ky)^T
void ComplexMatrix::geru(complex<double> alpha,
  const ComplexMatrix& x, int kx, const ComplexMatrix& y, int ky)
{
  assert(x.n()==m_);
  assert(y.n()==n_);
#if SCALAPACK
  int ione=1;

  int ix = kx+1;
  int jx = 1;
  int incx = x.m();

  int iy = ky+1;
  int jy = 1;
  int incy = y.m();
  pzgeru(&m_,&n_,&alpha,x.val,&ix,&jx,x.desc_,&incx,
                       y.val,&iy,&jy,y.desc_,&incy,
                       val,&ione,&ione,desc_);
#else
  int incx = x.m();
  int incy = y.m();
  zgeru(&m_,&n_,&alpha,&x.val[kx*x.m()],&incx,
                      &y.val[ky*y.m()],&incy,val,&m_);
#endif
}

Francois Gygi committed
1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389
////////////////////////////////////////////////////////////////////////////////
// symmetric rank-1 update using a row or a column of a Matrix x
void DoubleMatrix::syr(char uplo, double alpha,
  const DoubleMatrix& x, int k, char rowcol)
{
  assert(n_==m_);
#if SCALAPACK
  int ix,jx,incx,ione=1;
  if ( rowcol == 'c' )
  {
    // use column k of matrix x
    assert(x.m()==n_);
    ix = 1;
    jx = k+1;
    incx = 1;
  }
  else if ( rowcol == 'r' )
  {
    // use row k of matrix x
    assert(x.n()==n_);
    ix = k+1;
    jx = 1;
    incx = x.m();
  }
  else
  {
    cout << " DoubleMatrix::syr: invalid argument rowcol" << endl;
    MPI_Abort(MPI_COMM_WORLD,2);
  }
  pdsyr(&uplo,&n_,&alpha,x.val,&ix,&jx,x.desc_,&incx,
        val,&ione,&ione,desc_);
#else
  if ( rowcol == 'c' )
  {
    // use column k of matrix x
    assert(x.m()==n_);
    int incx = 1;
    dsyr(&uplo,&n_,&alpha,&x.val[k*x.m()],&incx,val,&m_);
  }
  else if ( rowcol == 'r' )
  {
    // use row k of matrix x
    assert(x.n()==n_);
    int incx = x.m();
    dsyr(&uplo,&n_,&alpha,&x.val[k],&incx,val,&m_);
  }
  else
  {
    cout << " DoubleMatrix::syr: invalid argument rowcol" << endl;
    exit(2);
  }
#endif
}

////////////////////////////////////////////////////////////////////////////////
DoubleMatrix& DoubleMatrix::operator=(const DoubleMatrix& a)
{
  if ( this == &a ) return *this;

  // operator= works only for matrices having same distribution on same context
1390
  assert( a.ictxt() == ictxt_ && a.m() == m_ && a.mb() == mb_ &&
Francois Gygi committed
1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407
          a.n() == n_ && a.nb() == nb_ );
  if ( active() )
  {
    for ( int i = 0; i < 9; i++ )
    {
      assert( desc_[i] == a.desc_[i] );
    }
    memcpy(val, a.val, mloc_*nloc_*sizeof(double));
  }
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
ComplexMatrix& ComplexMatrix::operator=(const ComplexMatrix& a)
{
  if ( this == &a ) return *this;

1408
  assert( a.ictxt() == ictxt_ && a.m() == m_ && a.mb() == mb_ &&
Francois Gygi committed
1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546
          a.n() == n_ && a.nb() == nb_ );
  if ( active() )
  {
    for ( int i = 0; i < 9; i++ )
    {
      assert( desc_[i] == a.desc_[i] );
    }
    memcpy(val, a.val, mloc_*nloc_*sizeof(complex<double>));
  }
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
// operator+=
DoubleMatrix& DoubleMatrix::operator+=(const DoubleMatrix &x)
{
  assert( ictxt_ == x.ictxt() );
  int ione=1;
  assert(m_==x.m());
  assert(n_==x.n());
  assert(mloc_==x.mloc());
  assert(nloc_==x.nloc());
  double alpha = 1.0;
  if( active_ )
    daxpy(&size_, &alpha, x.val, &ione, val, &ione);
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
// operator-=
DoubleMatrix& DoubleMatrix::operator-=(const DoubleMatrix &x)
{
  assert( ictxt_ == x.ictxt() );
  int ione=1;
  assert(m_==x.m());
  assert(n_==x.n());
  assert(mloc_==x.mloc());
  assert(nloc_==x.nloc());
  double alpha = -1.0;
  if( active_ )
    daxpy(&size_, &alpha, x.val, &ione, val, &ione);
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
// operator+=
ComplexMatrix& ComplexMatrix::operator+=(const ComplexMatrix& x)
{
  assert( ictxt_ == x.ictxt() );
  int ione=1;
  assert(m_==x.m());
  assert(n_==x.n());
  assert(mloc_==x.mloc());
  assert(nloc_==x.nloc());
  double alpha = 1.0;
  int two_size = 2 * size_;
  if( active_ )
    daxpy(&two_size, &alpha, (double*) x.val, &ione, (double*) val, &ione);
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
// operator-=
ComplexMatrix& ComplexMatrix::operator-=(const ComplexMatrix& x)
{
  assert( ictxt_ == x.ictxt() );
  int ione=1;
  assert(m_==x.m());
  assert(n_==x.n());
  assert(mloc_==x.mloc());
  assert(nloc_==x.nloc());
  double alpha = -1.0;
  int two_size = 2 * size_;
  if( active_ )
    daxpy(&two_size, &alpha, (double*) x.val, &ione, (double*) val, &ione);
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
// operator*=
DoubleMatrix& DoubleMatrix::operator*=(double alpha)
{
  int ione=1;
  if( active_ )
    dscal(&size_, &alpha, val, &ione);
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
ComplexMatrix& ComplexMatrix::operator*=(double alpha)
{
  int ione=1;
  if( active_ )
    zdscal(&size_, &alpha, val, &ione);
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
// operator*=
ComplexMatrix& ComplexMatrix::operator*=(complex<double> alpha)
{
  int ione=1;
  if( active_ )
    zscal(&size_, &alpha, val, &ione);
  return *this;
}

////////////////////////////////////////////////////////////////////////////////
// scal
void DoubleMatrix::scal(double alpha)
{
  *this *= alpha;
}

////////////////////////////////////////////////////////////////////////////////
// scal
void ComplexMatrix::scal(double alpha)
{
  *this *= alpha;
}

////////////////////////////////////////////////////////////////////////////////
// scal
void ComplexMatrix::scal(complex<double> alpha)
{
  *this *= alpha;
}

////////////////////////////////////////////////////////////////////////////////
// matrix multiplication
// this = alpha*op(A)*op(B)+beta*this
////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::gemm(char transa, char transb,
                double alpha, const DoubleMatrix& a,
                const DoubleMatrix& b, double beta)
{
  assert( ictxt_ == a.ictxt() );
  assert( ictxt_ == b.ictxt() );
1547

Francois Gygi committed
1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572
  if ( active() )
  {
    int m, n, k;
    if ( transa == 'N' || transa == 'n' )
    {
      m = a.m();
      k = a.n();
      assert(a.m()==m_);
    }
    else
    {
      m = a.n();
      k = a.m();
      assert(a.n()==m_);
    }
    if ( transb == 'N' || transb == 'n' )
    {
      n = b.n();
      assert(k==b.m());
    }
    else
    {
      n = b.m();
      assert(k==b.n());
    }
1573

Francois Gygi committed
1574 1575 1576 1577 1578 1579 1580
#ifdef SCALAPACK
    int ione=1;
    pdgemm(&transa, &transb, &m, &n, &k, &alpha,
         a.val, &ione, &ione, a.desc_,
         b.val, &ione, &ione, b.desc_,
         &beta, val, &ione, &ione, desc_);
#else
1581
    dgemm(&transa, &transb, &m, &n, &k, &alpha, a.val, &a.lld_,
Francois Gygi committed
1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596
          b.val, &b.lld_, &beta, val, &lld_);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// complex matrix multiplication
// this = alpha*op(A)*op(B)+beta*this
////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::gemm(char transa, char transb,
                complex<double> alpha, const ComplexMatrix& a,
                const ComplexMatrix& b, complex<double> beta)
{
  assert( ictxt_ == a.ictxt() );
  assert( ictxt_ == b.ictxt() );
1597

Francois Gygi committed
1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622
  if ( active() )
  {
    int m, n, k;
    if ( transa == 'N' || transa == 'n' )
    {
      m = a.m();
      k = a.n();
      assert(a.m()==m_);
    }
    else
    {
      m = a.n();
      k = a.m();
      assert(a.n()==m_);
    }
    if ( transb == 'N' || transb == 'n' )
    {
      n = b.n();
      assert(k==b.m());
    }
    else
    {
      n = b.m();
      assert(k==b.n());
    }
1623

Francois Gygi committed
1624 1625 1626 1627 1628 1629 1630
#ifdef SCALAPACK
    int ione=1;
    pzgemm(&transa, &transb, &m, &n, &k, &alpha,
         a.val, &ione, &ione, a.desc_,
         b.val, &ione, &ione, b.desc_,
         &beta, val, &ione, &ione, desc_);
#else
1631
    zgemm(&transa, &transb, &m, &n, &k, &alpha, a.val, &a.lld_,
Francois Gygi committed
1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647
          b.val, &b.lld_, &beta, val, &lld_);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// symmetric_matrix * matrix multiplication
// this = beta * this + alpha * a * b
// this = beta * this + alpha * b * a
////////////////////////////////////////////////////////////////////////////////
void DoubleMatrix::symm(char side, char uplo,
                double alpha, const DoubleMatrix& a,
                const DoubleMatrix& b, double beta)
{
  assert( ictxt_ == a.ictxt() );
  assert( ictxt_ == b.ictxt() );
1648

Francois Gygi committed
1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659
  if ( active() )
  {
    assert(a.n()==a.m());
    if ( side == 'L' || side == 'l' )
    {
      assert(a.n()==b.m());
    }
    else
    {
      assert(a.m()==b.n());
    }
1660

Francois Gygi committed
1661 1662 1663 1664 1665 1666 1667
#ifdef SCALAPACK
    int ione=1;
    pdsymm(&side, &uplo, &m_, &n_, &alpha,
         a.val, &ione, &ione, a.desc_,
         b.val, &ione, &ione, b.desc_,
         &beta, val, &ione, &ione, desc_);
#else
1668
    dsymm(&side, &uplo, &m_, &n_, &alpha, a.val, &a.lld_,
Francois Gygi committed
1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684
          b.val, &b.lld_, &beta, val, &lld_);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// hermitian_matrix * matrix multiplication
// this = beta * this + alpha * a * b
// this = beta * this + alpha * b * a
////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::hemm(char side, char uplo,
                complex<double> alpha, const ComplexMatrix& a,
                const ComplexMatrix& b, complex<double> beta)
{
  assert( ictxt_ == a.ictxt() );
  assert( ictxt_ == b.ictxt() );
1685

Francois Gygi committed
1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696
  if ( active() )
  {
    assert(a.n()==a.m());
    if ( side == 'L' || side == 'l' )
    {
      assert(a.n()==b.m());
    }
    else
    {
      assert(a.m()==b.n());
    }
1697

Francois Gygi committed
1698 1699 1700 1701 1702 1703 1704
#ifdef SCALAPACK
    int ione=1;
    pzhemm(&side, &uplo, &m_, &n_, &alpha,
         a.val, &ione, &ione, a.desc_,
         b.val, &ione, &ione, b.desc_,
         &beta, val, &ione, &ione, desc_);
#else
1705
    zhemm(&side, &uplo, &m_, &n_, &alpha, a.val, &a.lld_,
Francois Gygi committed
1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721
          b.val, &b.lld_, &beta, val, &lld_);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// complex_symmetric_matrix * matrix multiplication
// this = beta * this + alpha * a * b
// this = beta * this + alpha * b * a
////////////////////////////////////////////////////////////////////////////////
void ComplexMatrix::symm(char side, char uplo,
                complex<double> alpha, const ComplexMatrix& a,
                const ComplexMatrix& b, complex<double> beta)
{
  assert( ictxt_ == a.ictxt() );
  assert( ictxt_ == b.ictxt() );
1722

Francois Gygi committed
1723 1724 1725 1726 1727 1728 1729 1730 1731 1732 1733
  if ( active() )
  {
    assert(a.n()==a.m());
    if ( side == 'L' || side == 'l' )
    {
      assert(a.n()==b.m());
    }
    else
    {
      assert(a.m()==b.n());
    }
1734

Francois Gygi committed
1735 1736 1737 1738 1739 1740 1741
#ifdef SCALAPACK
    int ione=1;
    pzsymm(&side, &uplo, &m_, &n_, &alpha,
         a.val, &ione, &ione, a.desc_,
         b.val, &ione, &ione, b.desc_,
         &beta, val, &ione, &ione, desc_);
#else
1742
    zsymm(&side, &uplo, &m_, &n_, &alpha, a.val, &a.lld_,
Francois Gygi committed
1743 1744 1745 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756
          b.val, &b.lld_, &beta, val, &lld_);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// Compute a matrix-matrix product for a real triangular
// matrix or its transpose.
// *this = alpha op(A) * *this    if side=='l'
// *this = alpha * *this * op(A)  if side=='r'
// where op(A) = A or trans(A)
// alpha is a scalar, *this is an m by n matrix, and A is a unit or non-unit,
// upper- or lower-triangular matrix.
////////////////////////////////////////////////////////////////////////////////
1757
void DoubleMatrix::trmm(char side, char uplo, char trans, char diag,
Francois Gygi committed
1758 1759 1760 1761 1762 1763 1764 1765 1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 1776
                        double alpha, const DoubleMatrix& a)
{
  if ( active() )
  {
    assert(a.m_==a.n_);
    if ( side=='L' || side=='l' )
    {
      assert(a.n_==m_);
    }
    else
    {
      assert(a.n_==n_);
    }
#ifdef SCALAPACK
    int ione=1;
    pdtrmm(&side, &uplo, &trans, &diag, &m_, &n_,
           &alpha, a.val, &ione, &ione, a.desc_,
           val, &ione, &ione, desc_);
#else
1777
    dtrmm(&side, &uplo, &trans, &diag,
Francois Gygi committed
1778 1779 1780 1781 1782 1783 1784 1785 1786 1787 1788 1789
          &m_, &n_, &alpha, a.val, &a.m_, val, &m_);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// Solve op(A) * X = alpha * *this  (if side=='l')
// or    X * op(A) = alpha * *this  (if side=='r')
// where op(A) = A or trans(A)
// alpha is a scalar, *this is an m by n matrix, and A is a unit or non-unit,
// upper- or lower-triangular matrix.
////////////////////////////////////////////////////////////////////////////////
1790
void DoubleMatrix::trsm(char side, char uplo, char trans, char diag,
Francois Gygi committed
1791 1792 1793 1794 1795 1796 1797 1798 1799 1800 1801 1802 1803 1804 1805 1806 1807 1808 1809
                        double alpha, const DoubleMatrix& a)
{
  if ( active() )
  {
    assert(a.m_==a.n_);
    if ( side=='L' || side=='l' )
    {
      assert(a.n_==m_);
    }
    else
    {
      assert(a.n_==n_);
    }
#ifdef SCALAPACK
    int ione=1;
    pdtrsm(&side, &uplo, &trans, &diag, &m_, &n_,
           &alpha, a.val, &ione, &ione, a.desc_,
           val, &ione, &ione, desc_);
#else
1810
    dtrsm(&side, &uplo, &trans, &diag,
Francois Gygi committed
1811 1812 1813 1814 1815 1816 1817 1818 1819 1820 1821 1822
          &m_, &n_, &alpha, a.val, &a.m_, val, &m_);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// Solve op(A) * X = alpha * *this  (if side=='l')
// or    X * op(A) = alpha * *this  (if side=='r')
// where op(A) = A or trans(A)
// alpha is a scalar, *this is an m by n matrix, and A is a unit or non-unit,
// upper- or lower-triangular matrix.
////////////////////////////////////////////////////////////////////////////////
1823
void ComplexMatrix::trsm(char side, char uplo, char trans,
Francois Gygi committed
1824 1825 1826 1827 1828 1829 1830 1831 1832 1833 1834 1835 1836 1837 1838 1839 1840 1841 1842
  char diag, complex<double> alpha, const ComplexMatrix& a)
{
  if ( active() )
  {
    assert(a.m_==a.n_);
    if ( side=='L' || side=='l' )
    {
      assert(a.n_==m_);
    }
    else
    {
      assert(a.n_==n_);
    }
#ifdef SCALAPACK
    int ione=1;
    pztrsm(&side, &uplo, &trans, &diag, &m_, &n_,
           &alpha, a.val, &ione, &ione, a.desc_,
           val, &ione, &ione, desc_);
#else
1843
    ztrsm(&side, &uplo, &trans, &diag,
Francois Gygi committed
1844 1845 1846 1847 1848 1849 1850 1851 1852 1853 1854
          &m_, &n_, &alpha, a.val, &a.m_, val, &m_);
#endif
  }
}

////////////////////////////////////////////////////////////////////////////////
// Solves a triangular system of the form A * X = B or
// A**T * X = B, where A is a triangular matrix of  order  N,
// and  B  is an N-by-NRHS matrix.
// Output in B.
////////////////////////////////////////////////////////////////////////////////
1855
void DoubleMatrix::trtrs(char uplo, char trans, char diag,
Francois Gygi committed
1856 1857 1858 1859 1860 1861 1862 1863 1864
                         DoubleMatrix& b) const
{
  int info;
  if ( active() )
  {
    assert(m_==n_);

#ifdef SCALAPACK
    int ione=1;
1865
    pdtrtrs(&uplo, &trans, &diag, &m_, &b.n_,
Francois Gygi committed
1866 1867 1868
    val, &ione, &ione, desc_,
    b.val, &ione, &ione, b.desc_, &info);
#else
1869
    dtrtrs(&uplo, &trans, &diag, &m_, &b.n_, val, &m_,
Francois Gygi committed
1870 1871 1872 1873 1874 1875 1876 1877 1878 1879 1880 1881 1882 1883
           b.val, &b.m_, &info);
#endif
    if(info!=0)
    {
      cout <<" Matrix::trtrs, info=" << info << endl;
#ifdef USE_MPI
      MPI_Abort(MPI_COMM_WORLD, 2);
#else
      exit(2);
#endif
    }
  }
}

1884
////////////////////////////////////////////////////////////////////////////////
1885
// LU decomposition of a double matrix
1886 1887 1888 1889 1890 1891 1892 1893 1894 1895 1896 1897 1898 1899 1900 1901 1902 1903 1904 1905 1906 1907 1908 1909 1910 1911 1912 1913