testjade.C 6.22 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 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
// testjade.C
//
// use: testjade nprow npcol m mb
//
////////////////////////////////////////////////////////////////////////////////

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

#include "Timer.h"

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

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

const bool print = false;

enum MatrixType { FRANK, SCALED_FRANK, LAPLACIAN, SMALL };
const MatrixType mtype = SCALED_FRANK;

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; }

double frank(int n, int i, int j) { return n - max(i,j); }
double aijf(int n, int k, int i, int j)
52
{
Francois Gygi committed
53 54 55 56 57 58 59
  switch ( mtype )
  {
    case FRANK :
      if ( k == 0 )
        return frank(n,i,j);
      else
        return frank(n,n-i,n-j);
60

Francois Gygi committed
61 62 63 64 65
    case SCALED_FRANK :
      if ( k == 0 )
        return frank(n,i,j)/(n*n);
      else
        return frank(n,n-i,n-j)/(n*n);
66

Francois Gygi committed
67 68 69 70 71 72 73 74 75
    case LAPLACIAN :
      if ( i == j )
      {
        return 2.0 + k;
      }
      else if ( (i-j)*(i-j) == 1 )
      {
        return 1.0;
      }
76

Francois Gygi committed
77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100
    case SMALL :
    {
      if ( k == 0 )
      {
        if ( i == j && i == 0 )
          return 1.0;
        else if ( i == j && i > 0 )
          return 4.0;
        else
          return 0.0;
      }
      else
      {
        if ( i == j && i == 0 )
          return 1.0;
        else if ( i == j && i > 0 )
          return 1.0;
        else
          return 0.5;
      }
    }
  }
  return 0.0;
}
101

Francois Gygi committed
102 103 104 105 106 107
int main(int argc, char **argv)
{
  // use: testjade nprow npcol n nb
  const int maxsweep = 30;
  const double tol = 1.e-8;
  const int nmat = 2;
108

Francois Gygi committed
109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136
  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

  Timer tm;
  assert(argc==5);
  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 << "nmat=" << nmat << " m_a=" << m_a << ", n_a=" << n_a << endl;
    cout << "tol=" << tol << endl;
  }
#ifdef USE_MPI
137 138 139 140 141 142
  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
143
#endif
144
  {
145
    Context ctxt(MPI_COMM_WORLD,nprow,npcol);
Francois Gygi committed
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

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

    vector<DoubleMatrix*> a(nmat);
    vector<vector<double> > adiag(nmat);
    for ( int k = 0; k < nmat; k++ )
    {
      a[k] = new DoubleMatrix(ctxt,m_a,n_a,mb_a,nb_a);
      adiag[k].resize(n_a);
    }
    DoubleMatrix u(ctxt,m_a,n_a,mb_a,nb_a);

    if ( mype == 0 )
    {
      cout << " m_a x n_a / mb_a x nb_a / ta = "
           << a[0]->m() << "x" << a[0]->n() << " / "
           << a[0]->mb() << "x" << a[0]->nb() << endl;
    }

    for ( int k = 0; k < nmat; k++ )
    {
      for ( int m = 0; m < a[k]->nblocks(); m++ )
        for ( int l = 0; l < a[k]->mblocks(); l++ )
173
          for ( int y = 0; y < a[k]->nbs(m); y++ )
Francois Gygi committed
174 175 176 177 178 179 180 181 182 183 184 185 186 187
            for ( int x = 0; x < a[k]->mbs(l); x++ )
            {
              int i = a[k]->i(l,x);
              int j = a[k]->j(m,y);
              //double aij = a.i(l,x) * 10 + a.j(m,y);
              double aij = aijf(a[k]->n(),k,i,j);
              int iii = x + l*a[k]->mb();
              int jjj = y + m*a[k]->nb();
              int ival = iii + jjj * a[k]->mloc();
              (*a[k])[ival] = aij;
            }
      //cout << " a[" << k << "]=" << endl;
      //cout << (*a[k]);
    }
188

Francois Gygi committed
189 190 191
    tm.start();
    int nsweep = jade(maxsweep,tol,a,u,adiag);
    tm.stop();
192
    if ( ctxt.onpe0() )
Francois Gygi committed
193
    {
194
      cout << " m=" << m_a << " mb=" << mb_a
Francois Gygi committed
195 196 197
           << " ctxt: " << ctxt.nprow() << "x" << ctxt.npcol()
           << " nsweep=" << nsweep << " time: " << tm.real() << endl;
    }
198

Francois Gygi committed
199 200
    for ( int k = 0; k < nmat; k++ )
      sort(adiag[k].begin(),adiag[k].end());
201

Francois Gygi committed
202 203 204 205 206 207 208 209 210 211 212 213 214
    if ( nmat == 1 && (mtype == FRANK || mtype == SCALED_FRANK) )
    {
      vector<double> e_exact(adiag[0].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*n_a+1)) ) );
        if ( mtype == SCALED_FRANK )
        {
          int n = a[0]->n();
          e_exact[i] /= (n*n);
        }
      }
      sort(e_exact.begin(),e_exact.end());
215

Francois Gygi committed
216 217 218 219 220 221 222 223 224 225 226 227
      if ( mype == 0 )
      {
        for ( int k = 0; k < nmat; k++ )
        {
          double asum = 0.0;
          for ( int i = 0; i < e_exact.size(); i++ )
          {
            //cout << ctxt.mype() << ": eig[" << k << "][" << i << "]= "
            //     << adiag[k][i]
            //     << "  " << e_exact[i] << endl;
            asum += fabs(adiag[k][i]-e_exact[i]);
          }
228
          cout << "a[" << k << "] sum of abs eigenvalue errors: "
Francois Gygi committed
229 230 231 232 233
               << asum << endl;
        }
      }

    }
234

Francois Gygi committed
235 236 237 238 239
    if ( print )
    {
      // a[k] contains AU.
      // compute the product C = U^T A U
      DoubleMatrix c(*a[0]);
240 241
      for ( int k = 0; k < nmat; k++ )
      {
Francois Gygi committed
242
        c.gemm('t','n',1.0,u,(*a[k]),0.0);
243 244
        if ( mype == 0 )
        {
Francois Gygi committed
245
          cout << " a[" << k << "]=" << endl;
246 247
          cout << c;
        }
Francois Gygi committed
248
      }
249 250 251
      cout << " u=" << endl;
      cout << u;
    }
Francois Gygi committed
252 253 254 255 256 257
  }

#ifdef USE_MPI
  MPI_Finalize();
#endif
}