testjacobi.C 5.26 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 16 17 18 19
////////////////////////////////////////////////////////////////////////////////
//
// testjacobi.C
//
////////////////////////////////////////////////////////////////////////////////
//
// Test the Jacobi implementation of the Matrix class
Francois Gygi committed
20
//
21
// use: testjacobi nprow npcol n nb
Francois Gygi committed
22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42
//

#include <cassert>
#include <cstdlib>
#include <cmath>
#include <iostream>
#include <iomanip>
#include <string>
#include <fstream>
#include <valarray>
#include <algorithm>

#include "Timer.h"

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

#include "Context.h"
#include "Matrix.h"
#include "jacobi.h"
Francois Gygi committed
43
using namespace std;
Francois Gygi committed
44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67

double aa(int i, int j) { return 1.0/(i+1)+2.0*i/(j+1); }
double frank(int n, int i, int j) { return n - max(i,j); }
//double frank(int n, int i, int j) { return n - fabs(i-j); }
double bb(int i, int j) { return i-j-3; }

int main(int argc, char **argv)
{
  // use: testjacobi nprow npcol n nb
  int mype;
  int npes;
#ifdef USE_MPI
  MPI_Init(&argc,&argv);
  MPI_Comm_size(MPI_COMM_WORLD, &npes);
  MPI_Comm_rank(MPI_COMM_WORLD, &mype);
#else
  npes=1;
  mype=0;
#endif

  //char* infilename = argv[1];
  //ifstream infile(infilename);

  Timer tm;
Francois Gygi committed
68 69 70 71 72
  if ( argc != 5 )
  {
    cout << "use: testjacobi nprow npcol n nb" << endl;
    return 1;
  }
Francois Gygi committed
73 74 75 76 77 78 79 80 81 82 83 84 85 86
  int nprow=atoi(argv[1]);
  int npcol=atoi(argv[2]);
  int m_a=atoi(argv[3]);
  int mb_a=atoi(argv[4]);
  int n_a = m_a;
  int nb_a = mb_a;
  if(mype == 0)
  {
    //infile >> nprow >> npcol;
    cout<<"nprow="<<nprow<<", npcol="<<npcol<<endl;
    //infile >> m_a >> n_a >> mb_a >> nb_a >> ta;
    cout<<"m_a="<<m_a<<", n_a="<<n_a<<endl;
  }
#ifdef USE_MPI
87 88 89 90 91 92
  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);
Francois Gygi committed
93
#endif
94
  {
95
    Context ctxt(MPI_COMM_WORLD,nprow,npcol);
Francois Gygi committed
96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115

    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 u(ctxt,m_a,n_a,mb_a,nb_a);
    vector<double> e(n_a);

    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() << endl;
    }

    for ( int m = 0; m < a.nblocks(); m++ )
      for ( int l = 0; l < a.mblocks(); l++ )
116
        for ( int y = 0; y < a.nbs(m); y++ )
Francois Gygi committed
117 118 119 120 121 122 123 124 125 126 127 128
          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 = frank(a.n(),i,j);
            int iii = x + l*a.mb();
            int jjj = y + m*a.nb();
            int ival = iii + jjj * a.mloc();
            a[ival] = aij;
          }
    //cout << a;
129

Francois Gygi committed
130 131
    const int maxsweep = 30;
    tm.start();
Francois Gygi committed
132 133
    double tol = 1.e-10;
    int nsweep = jacobi(maxsweep,tol,a,u,e);
Francois Gygi committed
134 135 136
    tm.stop();
    if ( ctxt.onpe0() ) cout << " nsweep: " << nsweep << endl;
    if (ctxt.mype() == 0) cout << "Jacobi time: " << tm.real() << endl;
137

Francois Gygi committed
138
    sort(e.begin(),e.end());
139

Francois Gygi committed
140 141 142 143 144 145
    vector<double> e_exact(e.size());
    for ( int i = 0; i < e_exact.size(); i++ )
    {
      e_exact[i] = 1.0 / ( 2.0 * ( 1.0 - cos( ((2*i+1)*M_PI)/(2*a.n()+1)) ) );
    }
    sort(e_exact.begin(),e_exact.end());
146

Francois Gygi committed
147 148 149 150 151 152 153 154 155 156 157 158 159
    double asum = 0.0;
    for ( int i = 0; i < e.size(); i++ )
    {
      //cout << ctxt.mype() << ": eig[" << i << "]= " << e[i]
      //     << "  " << e_exact[i] << endl;
      asum += fabs(e[i]-e_exact[i]);
    }
    cout << " abs eigenvalue error: " << asum << endl;
    //cout << " u=" << endl;
    //cout << u;
    //cout << " au=" << endl;
    //cout << a;

160
#if 1
Francois Gygi committed
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
    // Eigenvalues and eigenvectors of c if c is square
    if ( a.m() == a.n() && a.mb() == a.nb() )
    {
      for ( int m = 0; m < a.nblocks(); m++ )
        for ( int l = 0; l < a.mblocks(); l++ )
          for ( int y = 0; y < a.nbs(m); y++ )
            for ( int x = 0; x < a.mbs(l); x++ )
            {
              int i = a.i(l,x);
              int j = a.j(m,y);
              int iii = x + l*a.mb();
              int jjj = y + m*a.nb();
              int ival = iii + jjj * a.mloc();
              if ( i == j )
                a[ival] = i + 1.e-5*drand48();
              else
                a[ival] = 1.e-5*drand48();
            }
      tm.reset();
      tm.start();
      if (mype == 0) cout << "Eigenproblem... ";
      DoubleMatrix z(a.context(),a.n(),a.n(),a.nb(),a.nb());
      valarray<double> w(a.m());
      a.syev('l',w,z);
      //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;
    }
#endif
  }

#ifdef USE_MPI
  MPI_Finalize();
#endif
}