testMatrix.C 12.4 KB
Newer Older
Francois Gygi committed
1 2 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
// 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
13 14 15
////////////////////////////////////////////////////////////////////////////////
//
// testMatrix.C
Francois Gygi committed
16
//
Francois Gygi committed
17
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
18 19 20 21
//
// multiply a matrix a(m,k) by b(k,n) to get c(m,n)
// using blocks of size (mb,nb) on a process grid (nprow,npcol)
//
22
// use: testMatrix input_file [-check]
Francois Gygi committed
23
// input_file:
24
// nprow npcol
Francois Gygi committed
25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48
// m_a n_a mb_a nb_a transa
// m_b n_b mb_b nb_b transb
// m_c n_c mb_c nb_c
//

#include <cassert>
#include <cstdlib>
#include <cmath>
#include <iostream>
#include <iomanip>
#include <string>
#include <fstream>
#include <valarray>
using namespace std;

#include "Timer.h"

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

#include "Context.h"
#include "Matrix.h"

49 50
double aa(int i, int j) { return 1.0/(i+1)+2.0*i/(j+1); }
double bb(int i, int j) { return i-j-3; }
Francois Gygi committed
51 52 53 54 55 56 57 58 59

int main(int argc, char **argv)
{
  int mype;
  int npes;
  MPI_Init(&argc,&argv);
  MPI_Comm_size(MPI_COMM_WORLD, &npes);
  MPI_Comm_rank(MPI_COMM_WORLD, &mype);

60
  char* infilename = argv[1];
61 62
  ifstream infile(infilename);

63
  assert(argc == 2 || argc == 3);
Francois Gygi committed
64
  bool tcheck = false;
65
  if ( argc == 3 )
Francois Gygi committed
66
  {
67
    if ( !strcmp(argv[2],"-check") )
Francois Gygi committed
68 69 70
      tcheck = true;
    else
    {
71
      cerr << " invalid argv[2]" << endl;
Francois Gygi committed
72 73 74 75
      MPI_Abort(MPI_COMM_WORLD,2);
    }
  }
  Timer tm;
76
  int nprow, npcol;
Francois Gygi committed
77 78 79 80 81 82
  int m_a, n_a, mb_a, nb_a;
  int m_b, n_b, mb_b, nb_b;
  int m_c, n_c, mb_c, nb_c;
  char ta, tb;
  if(mype == 0)
  {
83 84
    infile >> nprow >> npcol;
    cout<<"nprow="<<nprow<<", npcol="<<npcol<<endl;
85
    infile >> m_a >> n_a >> mb_a >> nb_a >> ta;
Francois Gygi committed
86
    cout<<"m_a="<<m_a<<", n_a="<<n_a<<endl;
87
    infile >> m_b >> n_b >> mb_b >> nb_b >> tb;
88
    cout<<"m_b="<<m_b<<", n_b="<<n_b<<endl;
89
    infile >> m_c >> n_c >> mb_c >> nb_c;
Francois Gygi committed
90 91
    cout<<"m_c="<<m_c<<", n_c="<<n_c<<endl;
  }
92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108
  MPI_Bcast(&nprow, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&npcol, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&m_a, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&n_a, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&mb_a, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&nb_a, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&m_b, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&n_b, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&mb_b, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&nb_b, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&m_c, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&n_c, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&mb_c, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&nb_c, 1, MPI_INT, 0, MPI_COMM_WORLD);
  MPI_Bcast(&ta, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
  MPI_Bcast(&tb, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
  {
Francois Gygi committed
109 110 111
    if ( ta == 'N' ) ta = 'n';
    if ( tb == 'N' ) tb = 'n';

112
    Context ctxt(MPI_COMM_WORLD,nprow,npcol);
Francois Gygi committed
113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138

    if ( mype == 0 )
    {
      cout << " Context " << ctxt.ictxt()
           << ": " << ctxt.nprow() << "x" << ctxt.npcol() << endl;
    }

    DoubleMatrix a(ctxt,m_a,n_a,mb_a,nb_a);
    DoubleMatrix b(ctxt,m_b,n_b,mb_b,nb_b);
    DoubleMatrix c(ctxt,m_c,n_c,mb_c,nb_c);

    if ( mype == 0 )
    {
      cout << " m_a x n_a / mb_a x nb_a / ta = "
           << a.m() << "x" << a.n() << " / "
           << a.mb() << "x" << a.nb() << " / " << ta << endl;
      cout << " m_b x n_b / mb_b x nb_b / tb = "
           << b.m() << "x" << b.n() << " / "
           << b.mb() << "x" << b.nb() << " / " << tb << endl;
      cout << " m_c x n_c / mb_c x nb_c      = "
           << c.m() << "x" << c.n() << " / "
           << c.mb() << "x" << c.nb() << endl;
    }

    for ( int m = 0; m < a.nblocks(); m++ )
      for ( int l = 0; l < a.mblocks(); l++ )
139
        for ( int y = 0; y < a.nbs(m); y++ )
Francois Gygi committed
140 141 142 143 144 145 146 147 148 149 150 151 152 153
          for ( int x = 0; x < a.mbs(l); x++ )
          {
            int i = a.i(l,x);
            int j = a.j(m,y);
            // double aij = a.i(l,x) * 10 + a.j(m,y);
            double aij = aa(i,j);
            int iii = x + l*a.mb();
            int jjj = y + m*a.nb();
            int ival = iii + jjj * a.mloc();
            a[ival] = aij;
          }

    for ( int m = 0; m < b.nblocks(); m++ )
      for ( int l = 0; l < b.mblocks(); l++ )
154
        for ( int y = 0; y < b.nbs(m); y++ )
Francois Gygi committed
155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176
          for ( int x = 0; x < b.mbs(l); x++ )
          {
            int i = b.i(l,x);
            int j = b.j(m,y);
            // double bij = b.i(l,x) * 10 + b.j(m,y);
            double bij = bb(i,j);
            int iii = x + l*b.mb();
            int jjj = y + m*b.nb();
            int ival = iii + jjj * b.mloc();
            b[ival] = bij;
          }

    tm.start();
    c.gemm(ta,tb,1.0,a,b,0.0);
    tm.stop();


    if ( tcheck )
    {
    cout << " checking results..." << endl;
    for ( int m = 0; m < c.nblocks(); m++ )
      for ( int l = 0; l < c.mblocks(); l++ )
177
        for ( int y = 0; y < c.nbs(m); y++ )
Francois Gygi committed
178 179 180 181 182 183
          for ( int x = 0; x < c.mbs(l); x++ )
          {
            int i = c.i(l,x);
            int j = c.j(m,y);
            double sum = 0.0;
            int kmax = ( ta == 'n' ) ? a.n() : a.m();
184

Francois Gygi committed
185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204
            if ( ( ta == 'n' ) && ( tb == 'n' ) )
            {
              for ( int k = 0; k < kmax; k++ )
                sum += aa(i,k) * bb(k,j);
            }
            else if ( ( ta != 'n' ) && ( tb == 'n' ) )
            {
              for ( int k = 0; k < kmax; k++ )
                sum += aa(k,i) * bb(k,j);
            }
            else if ( ( ta == 'n' ) && ( tb != 'n' ) )
            {
              for ( int k = 0; k < kmax; k++ )
                sum += aa(i,k) * bb(j,k);
            }
            else if ( ( ta != 'n' ) && ( tb != 'n' ) )
            {
              for ( int k = 0; k < kmax; k++ )
                sum += aa(k,i) * bb(j,k);
            }
205

Francois Gygi committed
206 207 208 209 210 211 212 213 214 215
            int iii = x + l*c.mb();
            int jjj = y + m*c.nb();
            int ival = iii + jjj * c.mloc();
            if ( fabs( c[ival] - sum ) > 1.e-8 )
            {
              cout << " error at element (" << i << "," << j << ") "
                   << c[ival] << " " << sum << endl;
              exit(1);
            }
          }
216 217


218
       cout << " results checked" << endl;
Francois Gygi committed
219 220
    }

221
    cout << " CPU/Real: " << setw(8) << tm.cpu()
Francois Gygi committed
222 223 224 225
         << " / " << setw(8) << tm.real();
    if ( tm.real() > 0.0 )
    {
      int kmax = ( ta == 'n' ) ? a.n() : a.m();
226
      cout << "  MFlops: "
Francois Gygi committed
227 228
           << (2.0e-6*m_c*n_c*kmax) / tm.real() << endl;
    }
229
#if 1
Francois Gygi committed
230 231
    double norma=a.nrm2();
    if(mype == 0)cout<<"Norm(a)="<<norma<<endl;
232

233 234
    double norm;
#if 0
Francois Gygi committed
235 236 237 238 239
    if(mype == 0)cout<<"DoubleMatrix::matgather..."<<endl;
    double*  aa=new double[a.m()*a.n()];
    a.matgather(aa, a.m());
    if(mype == 0)cout<<"DoubleMatrix::init..."<<endl;
    b.init(aa, a.m());
240
    norm=b.nrm2();
Francois Gygi committed
241 242 243
    if ( mype == 0 ) cout << "Norm(b)=" << norm << endl;
    if ( fabs(norm-norma)>0.000001 )
       cout << "DoubleMatrix: problem with matgather/init" << endl;
244
#endif
Francois Gygi committed
245 246 247 248

    if ( c.n() == b.m() && c.m() == b.n() )
    {
      if(mype == 0)cout<<"DoubleMatrix::transpose..."<<endl;
249 250
      tm.reset();
      tm.start();
Francois Gygi committed
251
      c.transpose(1.0,b,0.0);
252 253
      tm.stop();
      if ( mype == 0 ) cout << " transpose time: " << tm.real() << endl;
Francois Gygi committed
254 255 256
      norm=c.nrm2();
      if(mype == 0)cout<<"Norm(c)="<<norm<<endl;
    }
257

Francois Gygi committed
258 259
    if(mype == 0)cout<<"DoubleMatrix::scal..."<<endl;
    c.scal(0.5);
260

Francois Gygi committed
261 262 263 264 265
    if ( a.m() == b.m() && a.n() == b.n() )
    {
      if(mype == 0)cout<<"DoubleMatrix::axpy..."<<endl;
      a.axpy(-2., b);
    }
266

Francois Gygi committed
267
    if ( a.m()==c.m() && a.n()==c.n() && a.mb()==c.mb() && a.nb()==c.nb() )
Francois Gygi committed
268 269 270 271
    {
      if(mype == 0)cout<<"DoubleMatrix::operator=..."<<endl;
      c=a;
    }
272

Francois Gygi committed
273 274 275
    if(mype == 0)cout<<"DoubleMatrix::nrm2..."<<endl;
    norm=c.nrm2();
    if (mype == 0) cout<<"Norm="<<norm<<endl;
276

277
#if 1
Francois Gygi committed
278 279 280 281 282
    a.identity();
    DoubleMatrix a2(a);
    a -= a2;
    norm = a.nrm2();
    if (mype == 0) cout << "Norm(a)=" << norm << endl;
283
#endif
284

Francois Gygi committed
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
#if 1
    // Inverse of c if c is square
    if ( c.m() == c.n() && c.mb() == c.nb() )
    {
      for ( int m = 0; m < c.nblocks(); m++ )
        for ( int l = 0; l < c.mblocks(); l++ )
          for ( int y = 0; y < c.nbs(m); y++ )
            for ( int x = 0; x < c.mbs(l); x++ )
            {
              int i = c.i(l,x);
              int j = c.j(m,y);
              int iii = x + l*c.mb();
              int jjj = y + m*c.nb();
              int ival = iii + jjj * c.mloc();
              if ( i == j )
                c[ival] = i + 1.e-5*drand48();
              else
                c[ival] = 1.e-5*drand48();
            }
      tm.reset();
      tm.start();
      if (mype == 0) cout << "Inverse ... ";
      c.inverse();
      if (mype == 0) cout << " done" << endl;
      tm.stop();
      if (mype == 0) cout << "Inverse time: " << tm.real() << endl;
    }
#endif
313

314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336
    // Eigenvalues and eigenvectors of c if c is square
    if ( c.m() == c.n() && c.mb() == c.nb() )
    {
      for ( int m = 0; m < c.nblocks(); m++ )
        for ( int l = 0; l < c.mblocks(); l++ )
          for ( int y = 0; y < c.nbs(m); y++ )
            for ( int x = 0; x < c.mbs(l); x++ )
            {
              int i = c.i(l,x);
              int j = c.j(m,y);
              int iii = x + l*c.mb();
              int jjj = y + m*c.nb();
              int ival = iii + jjj * c.mloc();
              if ( i == j )
                c[ival] = i + 1.e-5*drand48();
              else
                c[ival] = 1.e-5*drand48();
            }
      tm.reset();
      tm.start();
      if (mype == 0) cout << "Eigenproblem... ";
      DoubleMatrix z(c.context(),c.n(),c.n(),c.nb(),c.nb());
      valarray<double> w(c.m());
Francois Gygi committed
337
      c.syevd('l',w,z);
338 339 340 341
      //c.syevx('l',w,z,1.e-5);
      if (mype == 0) cout << " done" << endl;
      tm.stop();
      if (mype == 0) cout << "Eigenproblem time: " << tm.real() << endl;
342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369

      // complex eigenvalue problem
      ComplexMatrix cc(c.context(),m_c,n_c,mb_c,nb_c);
      for ( int m = 0; m < cc.nblocks(); m++ )
        for ( int l = 0; l < cc.mblocks(); l++ )
          for ( int y = 0; y < cc.nbs(m); y++ )
            for ( int x = 0; x < cc.mbs(l); x++ )
            {
              int i = cc.i(l,x);
              int j = cc.j(m,y);
              int iii = x + l*cc.mb();
              int jjj = y + m*cc.nb();
              int ival = iii + jjj * cc.mloc();
              if ( i == j )
                cc[ival] = i + 1.e-5*drand48();
              else
                cc[ival] = complex<double>(1.e-5*drand48(), 1.e-5*drand48());
            }
      tm.reset();
      tm.start();
      if (mype == 0) cout << "Complex Eigenproblem... ";
      ComplexMatrix zz(cc.context(),cc.n(),cc.n(),cc.nb(),cc.nb());
      valarray<double> ww(cc.m());
      cc.heev('l',ww,zz);
      //c.syevx('l',w,z,1.e-5);
      if (mype == 0) cout << " done" << endl;
      tm.stop();
      if (mype == 0) cout << "Complex Eigenproblem time: " << tm.real() << endl;
Francois Gygi committed
370

371
    }
372

373 374 375
//  Gram-Schmidt orthogonalization of matrix a
    for ( int m = 0; m < a.nblocks(); m++ )
      for ( int l = 0; l < a.mblocks(); l++ )
376
        for ( int y = 0; y < a.nbs(m); y++ )
377 378 379 380
          for ( int x = 0; x < a.mbs(l); x++ )
          {
            int i = a.i(l,x);
            int j = a.j(m,y);
381
            //double aij = aa(i,j);
382 383 384
            int iii = x + l*a.mb();
            int jjj = y + m*a.nb();
            int ival = iii + jjj * a.mloc();
385
            if ( i == j )
386 387 388 389 390 391 392 393 394
              a[ival] = i + 1.e-6*drand48();
            else
              a[ival] = 1.e-6*drand48();
          }
    DoubleMatrix s(a.context(),a.n(),a.n(),a.nb(),a.nb());
    tm.reset();
    tm.start();
    s.syrk('l','t',2.0,a,0.0);
    if (mype == 0) cout << "Gram syrk time: " << tm.real() << endl;
395

396 397 398 399 400
    tm.reset();
    tm.start();
    s.syr('l',-1.0,a,0,'r');
    tm.stop();
    if (mype == 0) cout << "Gram syr time: " << tm.real() << endl;
401

402 403 404 405 406 407
    // Cholesky decomposition
    tm.reset();
    tm.start();
    s.potrf('l'); // Cholesky decomposition: S = L * L^T
    tm.stop();
    if (mype == 0) cout << "Gram Cholesky time: " << tm.real() << endl;
408

409 410 411 412 413 414 415 416
    // Triangular solve
    tm.reset();
    tm.start();
    // solve triangular system X * L^T = C
    a.trsm('r','l','t','n',1.0,s);
    tm.stop();
    if (mype == 0) cout << "Gram triangular solve time: " << tm.real() << endl;
#endif
Francois Gygi committed
417 418 419 420
  }

  MPI_Finalize();
}