BiCGStab.c 32.2 KB
Newer Older
1
2
3
4
5
#include <stdio.h>
#include <stdlib.h>
#include <unistd.h>
#include <string.h>
#include <math.h>
6
7
//#include <mkl_blas.h>
#include "mymkl.h"
8
9
#include <mpi.h>
#include <hb_io.h>
10
11
//#include <vector>
#include <sys/prctl.h>
12
13
14
15
16
17
18
19

#include "reloj.h"
#include "ScalarVectors.h"
#include "SparseProduct.h"
#include "ToolsMPI.h"
#include "matrix.h"
#include "common.h"

20
21
22
#include "../malleability/MAM.h"
#include "ToolsMAM.h"

23
24
25
26
27
28
29
30
31
// ================================================================================

#define DIRECT_ERROR 1
#define PRECOND 1
// #define SPMV_OPTIMIZED 1
#ifdef SPMV_OPTIMIZED
  #define COLL_P2P_SPMV 0
#endif

32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
typedef struct {
  double tol, tol0;
  int iter, n;

  double rho;
  double *x, *b;
  double *s, *q, *r, *p, *r0, *y, *p_hat, *q_hat;
  double *aux;
  SparseMatrix matL;

#if PRECOND
  double *diags;
#endif
#if DIRECT_ERROR
  double *res_err, *x_exact;
  double direct_err;
#endif
  double t1;

  int *sizes, *dspls;
  int my_size, my_dspl;
  int *vlen;

  int myId, numP;
  MPI_Comm comm;
} Compute_data;


void originals_set_data(Compute_data *computeData, int num_target);
void targets_update(Compute_data *computeData, user_redist_t *user_data);
void user_func(void *args);
void dump(Compute_data *computeData);

void BiCGStab_init (Compute_data *computeData) {
    int size = computeData->matL.dim2, sizeR = computeData->matL.dim1; 
67
68
    int IONE = 1; 
    double DONE = 1.0, DMONE = -1.0, DZERO = 0.0;
69
70
    int n, n_dist, myId, nProcs;
    double t2;
71
#if PRECOND
72
73
74
    int i;
    int *posd = NULL;
    computeData->diags = NULL;
75
76
#endif

77
78
79
80
81
82
83
84
85
86
87
88
    computeData->s = NULL; computeData->q = NULL; computeData->r = NULL; computeData->p = NULL;
    computeData->r0 = NULL; computeData->y = NULL; computeData->p_hat = NULL; computeData->q_hat = NULL;
    computeData->aux = NULL;
    myId = computeData->myId;
    nProcs = computeData->numP;
    n = size; n_dist = sizeR;
    CreateDoubles (&computeData->s, n_dist);
    CreateDoubles (&computeData->q, n_dist);
    CreateDoubles (&computeData->r, n_dist);
    CreateDoubles (&computeData->r0, n_dist);
    CreateDoubles (&computeData->p, n_dist);
    CreateDoubles (&computeData->y, n_dist);
89
90
#if DIRECT_ERROR
    // init exact solution
91
92
93
94
    computeData->res_err = NULL; computeData->x_exact = NULL;
    CreateDoubles (&computeData->x_exact, n_dist);
    CreateDoubles (&computeData->res_err, n_dist);
    InitDoubles (computeData->x_exact, n_dist, DONE, DZERO);
95
96
97
98
#endif // DIRECT_ERROR 

#if PRECOND
    CreateInts (&posd, n_dist);
99
100
101
102
    CreateDoubles (&computeData->p_hat, n_dist);
    CreateDoubles (&computeData->q_hat, n_dist);
    CreateDoubles (&computeData->diags, n_dist);
    GetDiagonalSparseMatrix2 (computeData->matL, computeData->dspls[myId], computeData->diags, posd);
103
104
#pragma omp parallel for
    for (i=0; i<n_dist; i++) 
105
        computeData->diags[i] = DONE / computeData->diags[i];
iker_martin's avatar
iker_martin committed
106
107
    InitDoubles (computeData->p_hat, n_dist, DZERO, DZERO);
    InitDoubles (computeData->q_hat, n_dist, DZERO, DZERO);
108
#endif
109
    CreateDoubles (&computeData->aux, n); 
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128

#ifdef SPMV_OPTIMIZED
    int *permP = NULL, *ipermP = NULL;
    int *vdspP = NULL, *vdimP = NULL, *vdspR = NULL, *vdimR = NULL;
    double *vecP = NULL;
    MPI_Datatype *vectDatatypeP = NULL, *vectDatatypeR = NULL;

    CreateInts (&ipermP, size);
    CreateInts (&vdimP, nProcs); CreateInts (&vdspP, nProcs + 1);
    CreateInts (&vdimR, nProcs); CreateInts (&vdspR, nProcs + 1);
    vectDatatypeP = (MPI_Datatype *) malloc (nProcs * sizeof (MPI_Datatype));
    vectDatatypeR = (MPI_Datatype *) malloc (nProcs * sizeof (MPI_Datatype));
    createAlltoallwStruct (COLL_P2P_SPMV, MPI_COMM_WORLD, mat, sizes, dspls, vdimP, 
                vdspP, &aux, &permP, ipermP, vdimR, vdspR, vectDatatypeP, vectDatatypeR);

  // Code required before the loop  
    PermuteInts (mat.vpos, ipermP, mat.vptr[mat.dim1]);
#endif

129
    computeData->iter = 0;
130
131
132
133
134
135
#ifdef SPMV_OPTIMIZED
    joinDistributeVectorSPMV (COLL_P2P_SPMV, MPI_COMM_WORLD, x, vecP, vdimP, vdspP, 
                                vdimR, vdspR, vectDatatypeP, vectDatatypeR);
    InitDoubles (s, sizeR, DZERO, DZERO);
    ProdSparseMatrixVectorByRows (mat, 0, vecP, s);                  // s = A * x
#else
136
137
138
    MPI_Allgatherv (computeData->x, sizeR, MPI_DOUBLE, computeData->aux, computeData->sizes, computeData->dspls, MPI_DOUBLE, computeData->comm);
    InitDoubles (computeData->s, sizeR, DZERO, DZERO);
    ProdSparseMatrixVectorByRows (computeData->matL, 0, computeData->aux, computeData->s); // s = A * x
139
#endif
140
141
    rcopy (&n_dist, computeData->b, &IONE, computeData->r, &IONE);                                // r = b
    raxpy (&n_dist, &DMONE, computeData->s, &IONE, computeData->r, &IONE);           // r -= s
142

143
144
    rcopy (&n_dist, computeData->r, &IONE, computeData->p, &IONE);                                // p = r
    rcopy (&n_dist, computeData->r, &IONE, computeData->r0, &IONE);                               // r0 = r
145
    // compute tolerance and <r0,r0>
146
147
    computeData->rho = rdot (&n_dist, computeData->r, &IONE, computeData->r, &IONE);
    MPI_Allreduce (MPI_IN_PLACE, &computeData->rho, 1, MPI_DOUBLE, MPI_SUM, computeData->comm);
148

149
150
    computeData->tol0 = sqrt (computeData->rho);
    computeData->tol = computeData->tol0;
151
152
153

#if DIRECT_ERROR
    // compute direct error
154
155
    rcopy (&n_dist, computeData->x_exact, &IONE, computeData->res_err, &IONE);                    // res_err = x_exact
    raxpy (&n_dist, &DMONE, computeData->x, &IONE, computeData->res_err, &IONE);                  // res_err -= x
156
157

    // compute inf norm
158
159
    computeData->direct_err = norm_inf(n_dist, computeData->res_err);
    MPI_Allreduce(MPI_IN_PLACE, &computeData->direct_err, 1, MPI_DOUBLE, MPI_MAX, computeData->comm);
160
161

    //    // compute euclidean norm
162
    //    direct_err = rdot (&n_dist, res_err, &IONE, res_err, &IONE);            // direct_err = res_err' * res_err
163
164
165
166
    //    MPI_Allreduce(MPI_IN_PLACE, &direct_err, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
    //    direct_err = sqrt(direct_err);
#endif // DIRECT_ERROR

167
168
169
170
171
#if PRECOND
    RemoveInts (&posd);
#endif

    MPI_Barrier(computeData->comm);
172
    if (myId == 0) 
173
174
        reloj (&computeData->t1, &t2);
}
175

176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
void BiCGStab_compute (Compute_data *computeData, user_redist_t *user_data) {
    int size = computeData->matL.dim2, sizeR = computeData->matL.dim1; 
    int IONE = 1; 
    double DONE = 1.0, DMONE = -1.0, DZERO = 0.0;
    int n, n_dist;
    int maxiter, myId, reconfigure, rec_iter, state;
    double beta, alpha, umbral, omega, tmp;
    double t3, t4;
    double reduce[2];

    n = size; n_dist = sizeR; maxiter = 16 * size; rec_iter = maxiter / 2; umbral = 1.0e-8;
    myId = computeData->myId;
    state = -1;
    reconfigure = 0; rec_iter = -100;

    while ((computeData->iter < maxiter) && (computeData->tol > umbral)) {
192
193

#if PRECOND
194
        VvecDoubles (DONE, computeData->diags, computeData->p, DZERO, computeData->p_hat, n_dist);              // p_hat = D^-1 * p
195
#else
196
        computeData->p_hat = computeData->p;
197
#endif
iker_martin's avatar
iker_martin committed
198
     //dump(computeData);
199
200
201
202
203
204
#ifdef SPMV_OPTIMIZED
        joinDistributeVectorSPMV (COLL_P2P_SPMV, MPI_COMM_WORLD, p_hat, vecP, vdimP, 
                                    vdspP, vdimR, vdspR, vectDatatypeP, vectDatatypeR);
        InitDoubles (s, sizeR, DZERO, DZERO);
        ProdSparseMatrixVectorByRows (mat, 0, vecP, s);                   // s = A * p
#else
205
206
207
        MPI_Allgatherv (computeData->p_hat, sizeR, MPI_DOUBLE, computeData->aux, computeData->sizes, computeData->dspls, MPI_DOUBLE, computeData->comm);
        InitDoubles (computeData->s, sizeR, DZERO, DZERO);
        ProdSparseMatrixVectorByRows (computeData->matL, 0, computeData->aux, computeData->s);                   // s = A * p
208
209
210
211
#endif

        if (myId == 0) 
#if DIRECT_ERROR
212
            printf ("%d \t %g \t %g \t %g \n", computeData->iter, computeData->tol, umbral, computeData->direct_err);
213
#else        
214
        printf ("%d \t %g \n", computeData->iter, computeData->tol);
215
#endif // DIRECT_ERROR
216
217
        alpha = rdot (&n_dist, computeData->r0, &IONE, computeData->s, &IONE);
        MPI_Allreduce (MPI_IN_PLACE, &alpha, 1, MPI_DOUBLE, MPI_SUM, computeData->comm);
218

219
        alpha = computeData->rho / alpha;
220

221
        rcopy (&n_dist, computeData->r, &IONE, computeData->q, &IONE);                            // q = r
222
        tmp = -alpha;
223
        raxpy (&n_dist, &tmp, computeData->s, &IONE, computeData->q, &IONE);                      // q = r - alpha * s;
224
225
226

        // second spmv
#if PRECOND
227
        VvecDoubles (DONE, computeData->diags, computeData->q, DZERO, computeData->q_hat, n_dist);             // q_hat = D^-1 * q
228
#else
229
        computeData->q_hat = computeData->q;
230
231
232
233
234
235
236
#endif
#ifdef SPMV_OPTIMIZED
        joinDistributeVectorSPMV (COLL_P2P_SPMV, MPI_COMM_WORLD, q_hat, vecP, vdimP, 
                                  vdspP, vdimR, vdspR, vectDatatypeP, vectDatatypeR);
        InitDoubles (y, sizeR, DZERO, DZERO);
        ProdSparseMatrixVectorByRows (mat, 0, vecP, y);                // y = A * q
#else
237
238
239
        MPI_Allgatherv (computeData->q_hat, sizeR, MPI_DOUBLE, computeData->aux, computeData->sizes, computeData->dspls, MPI_DOUBLE, computeData->comm);
        InitDoubles (computeData->y, sizeR, DZERO, DZERO);
        ProdSparseMatrixVectorByRows (computeData->matL, 0, computeData->aux, computeData->y);                // y = A * q
240
241
#endif
        // omega = <q, y> / <y, y>
242
243
244
        reduce[0] = rdot (&n_dist, computeData->q, &IONE, computeData->y, &IONE);
        reduce[1] = rdot (&n_dist, computeData->y, &IONE, computeData->y, &IONE);
        MPI_Allreduce (MPI_IN_PLACE, reduce, 2, MPI_DOUBLE, MPI_SUM, computeData->comm);
245
246
247
248

        omega = reduce[0] / reduce[1];

        // x+1 = x + alpha * p + omega * q
249
250
        raxpy (&n_dist, &alpha, computeData->p_hat, &IONE, computeData->x, &IONE); 
        raxpy (&n_dist, &omega, computeData->q_hat, &IONE, computeData->x, &IONE); 
251
252

        // r+1 = q - omega * y
253
        rcopy (&n_dist, computeData->q, &IONE, computeData->r, &IONE);                            // r = q
254
        tmp = -omega;
255
        raxpy (&n_dist, &tmp, computeData->y, &IONE, computeData->r, &IONE);                      // r = q - omega * y;
256
257
        
        // rho = <r0, r+1> and tolerance
258
259
260
        reduce[0] = rdot (&n_dist, computeData->r0, &IONE, computeData->r, &IONE);
        reduce[1] = rdot (&n_dist, computeData->r, &IONE, computeData->r, &IONE);
        MPI_Allreduce (MPI_IN_PLACE, reduce, 2, MPI_DOUBLE, MPI_SUM, computeData->comm);
261
262

        tmp = reduce[0];
263
        computeData->tol = sqrt (reduce[1]) / computeData->tol0;
264
265

        // beta = (alpha / omega) * <r0, r+1> / <r0, r>
266
267
        beta = (alpha / omega) * (tmp / computeData->rho);
        computeData->rho = tmp;
268
269
270
       
        // p+1 = r+1 + beta * (p - omega * s)
        tmp = -omega; 
271
272
273
        raxpy (&n_dist, &tmp, computeData->s, &IONE, computeData->p, &IONE);                     // p -= omega * s
        rscal (&n_dist, &beta, computeData->p, &IONE);                                           // p = beta * p
        raxpy (&n_dist, &DONE, computeData->r, &IONE, computeData->p, &IONE);                    // p += r
274
275
276

#if DIRECT_ERROR
        // compute direct error
277
278
279
        rcopy (&n_dist, computeData->x_exact, &IONE, computeData->res_err, &IONE);               // res_err = x_exact
        raxpy (&n_dist, &DMONE, computeData->x, &IONE, computeData->res_err, &IONE);             // res_err -= x
  
280
        // compute inf norm
281
282
        computeData->direct_err = norm_inf(n_dist, computeData->res_err);
        MPI_Allreduce(MPI_IN_PLACE, &computeData->direct_err, 1, MPI_DOUBLE, MPI_MAX, computeData->comm);
283
284

        //        // compute euclidean norm
285
        //        direct_err = rdot (&n_dist, res_err, &IONE, res_err, &IONE);
286
287
288
289
        //        MPI_Allreduce(MPI_IN_PLACE, &direct_err, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
        //        direct_err = sqrt(direct_err);
#endif // DIRECT_ERROR

290
291
292
293
294
295
296
297
298
299
        computeData->iter++;
        if (computeData->iter == rec_iter) { reconfigure = 1;}
	if (reconfigure) {
	  MAM_Checkpoint(&state, MAM_CHECK_COMPLETION, user_func, (void *) user_data);
	  if(state == MAM_COMPLETED) {
	    reconfigure = 0; 
            //free_computeData(computeData);
            targets_update(computeData, user_data);
	    }
	}
300
301
    }

302
    MPI_Barrier(computeData->comm);
303
304
305
    if (myId == 0) 
        reloj (&t3, &t4);

306
307
308
309
310
311
    if(state == MAM_PENDING) {
      MAM_Checkpoint(&state, MAM_WAIT_COMPLETION, user_func, (void *) user_data);
      //free_computeData(computeData);
      //targets_update(dist_data, computeData, user_data);
    }

312
313
314
315
316
317
318
319
320
321
322
323
324
#ifdef SPMV_OPTIMIZED
    // Code required after the loop 
    PermuteInts (mat.vpos, permP, mat.vptr[mat.dim1]);

    // Freeing memory for Permutation
    free (vectDatatypeR); vectDatatypeR = NULL; free (vectDatatypeP); vectDatatypeP = NULL;
    RemoveDoubles (&vecP); RemoveInts (&permP);
    RemoveInts (&vdspR); RemoveInts (&vdimR); RemoveInts (&vdspP); RemoveInts (&vdimP);
    RemoveInts (&ipermP);
#endif

    if (myId == 0) {
        printf ("Size: %d \n", n);
325
326
327
328
        printf ("Iter: %d \n", computeData->iter);
        printf ("Tol: %g \n", computeData->tol);
        printf ("Time_loop: %20.10e\n", (t3-computeData->t1));
        printf ("Time_iter: %20.10e\n", (t3-computeData->t1)/computeData->iter);
329
    }
330
331
332
}

void BiCGStab_free (Compute_data *computeData) {
333

334
335
336
    RemoveDoubles (&computeData->aux); RemoveDoubles (&computeData->s); 
    RemoveDoubles (&computeData->q); RemoveDoubles (&computeData->r); 
    RemoveDoubles (&computeData->p); RemoveDoubles (&computeData->r0); RemoveDoubles (&computeData->y);
337
#if PRECOND
338
339
    RemoveDoubles (&computeData->diags);
    RemoveDoubles(&computeData->p_hat); RemoveDoubles (&computeData->q_hat); 
340
#endif
iker_martin's avatar
iker_martin committed
341
342
343
#if DIRECT_ERROR
    RemoveDoubles (&computeData->res_err); RemoveDoubles (&computeData->x_exact); 
#endif
344
345
346
347

    RemoveDoubles (&computeData->x); 
    RemoveDoubles (&computeData->b);
    RemoveInts (&computeData->sizes); RemoveInts (&computeData->dspls); 
iker_martin's avatar
iker_martin committed
348
    RemoveInts (&computeData->vlen);
349
350
351
352
353
    RemoveSparseMatrix (&computeData->matL);
}

void originals_free() {

354
355
356
357
358
359
360
361
362
363
}

/*********************************************************************************/

int main (int argc, char **argv) {
    int dim; 
    double *sol1 = NULL, *sol2 = NULL;
    int index = 0, indexL = 0;
    SparseMatrix mat  = {0, 0, NULL, NULL, NULL}, sym = {0, 0, NULL, NULL, NULL};

364
    int root = 0, myId, nProcs, isTarget, numTarget, req;
365
366
367
    int dimL, dspL, *vdimL = NULL, *vdspL = NULL;
    SparseMatrix matL = {0, 0, NULL, NULL, NULL};
    double *sol1L = NULL, *sol2L = NULL;
368
    double beta;
369

370
371
    int IONE = 1;
    double DMONE = -1.0;
372

373
374
375
    int mat_from_file, nodes, size_param, stencil_points;
    Compute_data computeData;
    user_redist_t user_data;
376
377
378

    /***************************************/

379
    MPI_Init_thread (&argc, &argv, MPI_THREAD_MULTIPLE, &req);
380
381
382
383
384
385

    // Definition of the variables nProcs and myId
    MPI_Comm_size(MPI_COMM_WORLD, &nProcs);
    MPI_Comm_rank(MPI_COMM_WORLD, &myId);
    root = nProcs-1;
    root = 0;
386
387
388
389
390
391
392
393
    computeData.myId = myId;
    computeData.numP = nProcs;
    computeData.comm = MPI_COMM_WORLD;
    user_data = empty_user_data;
    user_data.comm = computeData.comm;

    prctl(PR_SET_PTRACER, PR_SET_PTRACER_ANY, 0, 0, 0);
    isTarget = MAM_Init(root, &computeData.comm, argv[0], user_func, (void *) &user_data);
394

395
396
397
398
399
400
    if(isTarget) {
      targets_update(&computeData, &user_data);
    } else {
    /***************************************/
      if (argc == 4) {
          mat_from_file = atoi(argv[2]);
iker_martin's avatar
iker_martin committed
401
	  numTarget = atoi(argv[3]);
402
403
404
405
406
407
      } else {
          mat_from_file = atoi(argv[2]);
          nodes = atoi(argv[3]);
          size_param = atoi(argv[4]);
          stencil_points = atoi(argv[5]);
      }
408
409
    /***************************************/

410
411
412
413
414
415
416
417
418
      printf ("A\n");
      CreateInts (&vdimL, nProcs); CreateInts (&vdspL, nProcs); 
      if(mat_from_file) {
          if (myId == root) {
              // Creating the matrix
              ReadMatrixHB (argv[1], &sym);
              TransposeSparseMatrices (sym, 0, &mat, 0);
              dim = mat.dim1;
          }
419
420

        // Distributing the matrix
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
          dim = DistributeMatrix (mat, index, &matL, indexL, vdimL, vdspL, root, MPI_COMM_WORLD);
          dimL = vdimL[myId]; dspL = vdspL[myId];
          if (myId == root) {
            RemoveSparseMatrix (&mat);
            RemoveSparseMatrix (&sym);
          } 
          printf ("B\n");
      }
      else {
          dim = size_param * size_param * size_param;
          int divL, rstL, i;
          divL = (dim / nProcs); rstL = (dim % nProcs);
          for (i=0; i<nProcs; i++) vdimL[i] = divL + (i < rstL);
          vdspL[0] = 0; for (i=1; i<nProcs; i++) vdspL[i] = vdspL[i-1] + vdimL[i-1];
          dimL = vdimL[myId]; dspL = vdspL[myId];
          int band_width = size_param * (size_param + 1) + 1;
          band_width = 100 * nodes;
          long nnz_here = ((long) (stencil_points + 2 * band_width)) * dimL;
          printf ("dimL: %d, nodes: %d, size_param: %d, band_width: %d, stencil_points: %d, nnz_here: %ld\n",
                  dimL, nodes, size_param, band_width, stencil_points, nnz_here);
          allocate_matrix(dimL, dim, nnz_here, &matL);
          generate_Poisson3D_filled(&matL, size_param, stencil_points, band_width, dspL, dimL, dim);

          // To generate ill-conditioned matrices
  //        double factor = 1.0e6;
  //        ScaleFirstRowCol(matL, dspL, dimL, myId, root, factor);
      }
      MPI_Barrier(MPI_COMM_WORLD);

      // Creating the vectors
      CreateDoubles (&sol1, dim);
//      CreateDoubles (&sol2, dim);
      CreateDoubles (&sol1L, dimL);
      CreateDoubles (&sol2L, dimL);

//      InitDoubles (sol2, dim, 0.0, 0.0);
      InitDoubles (sol1L, dimL, 0.0, 0.0);
      InitDoubles (sol2L, dimL, 0.0, 0.0);
459
460
461

    /***************************************/

462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
      printf ("C\n");

      beta = 1.0 / sqrt(dim);
      if(mat_from_file) {
          // compute b = A * x_c, x_c = 1/sqrt(nbrows)
          InitDoubles (sol1, dim, 1.0, 0.0);
          ProdSparseMatrixVectorByRows (matL, 0, sol1, sol1L);                  // s = A * x
          rscal (&dimL, &beta, sol1L, &IONE);                                         // s = beta * s
      } else {
          InitDoubles (sol1, dim, 0.0, 0.0);

          int k=0;
          int *vptrM = matL.vptr;
          for (int i=0; i < matL.dim1; i++) {
              for(int j=vptrM[i]; j<vptrM[i+1]; j++) {
                  sol1L[k] += matL.vval[j];
              }
          }
      }

      printf ("D\n");

//      MPI_Scatterv (sol2, vdimL, vdspL, MPI_DOUBLE, sol2L, dimL, MPI_DOUBLE, root, MPI_COMM_WORLD); //FIXME It does not seem to do anything

      printf ("E\n");
      computeData.sizes = vdimL;
      computeData.my_size = dimL;
      computeData.dspls = vdspL;
      computeData.my_dspl = dspL;
      computeData.b = sol1L;
      computeData.x = sol2L;
      computeData.matL = matL;
      computeData.n = computeData.matL.dim2;
      RemoveDoubles (&sol1); 
      BiCGStab_init (&computeData);
      originals_set_data(&computeData, numTarget);
498
499
500
    }


501
    BiCGStab_compute (&computeData, &user_data);
502
503
504
505
506

    printf ("F\n");

    // Error computation ||b-Ax||
//    if(mat_from_file) {
507
508
509
510
511
512
513
514
515
        dim = matL.dim2;
        CreateDoubles (&sol2, dim);
        InitDoubles (sol2, dim, 0.0, 0.0);
        MPI_Allgatherv (computeData.x, computeData.my_size, MPI_DOUBLE, sol2, computeData.sizes, computeData.dspls, MPI_DOUBLE, computeData.comm);
        InitDoubles (computeData.x, computeData.my_size, 0, 0);
        ProdSparseMatrixVectorByRows (computeData.matL, 0, sol2, computeData.x);
        raxpy (&dimL, &DMONE, computeData.x, &IONE, computeData.b, &IONE);          
        beta = rdot (&computeData.my_size, computeData.b, &IONE, computeData.b, &IONE);
        MPI_Allreduce (MPI_IN_PLACE, &beta, 1, MPI_DOUBLE, MPI_SUM, computeData.comm);
516
517
518
519
520
        
//    } else {
//        // case with x_exact = {1.0}
//        for (int i=0; i<dimL; i++)
//            sol2L[i] -= 1.0;
521
//        beta = rdot (&dimL, sol2L, &IONE, sol2L, &IONE);            
522
523
524
525
526
527
528
529
//    } 

    beta = sqrt(beta);
    if (myId == 0) 
        printf ("Error: %20.10e\n", beta);

    /***************************************/
    // Freeing memory
530
    BiCGStab_free (&computeData);
531
532
    RemoveDoubles (&sol2); 

533
    MAM_Finalize ();
534
535
536
537
538
    MPI_Finalize ();

    return 0;
}

539
540
541
542
543
544
545
546
547
548

/* MAM New functions */

/*
 * Función para declarar los datos a comunicar por parte de MAM
 */
void originals_set_data(Compute_data *computeData, int num_target) {

    TransformHeadertoLength (computeData->matL.vptr, computeData->n);
    CreateInts (&computeData->vlen, computeData->n); 
iker_martin's avatar
iker_martin committed
549
    CopyInts (computeData->matL.vptr+1, computeData->vlen, computeData->n); 
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
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
665
666
667
668
669
670
671
672
673
674
675
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
    TransformLengthtoHeader (computeData->matL.vptr, computeData->n);

    MAM_Set_target_number(num_target);

    MAM_Data_add(&(computeData->n), NULL, 1, MPI_INT, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);
    MAM_Data_add(&(computeData->tol0), NULL, 1, MPI_DOUBLE, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);
    MAM_Data_add(&(computeData->t1), NULL, 1, MPI_DOUBLE, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);

    MAM_Data_add(&(computeData->iter), NULL, 1, MPI_INT, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
    MAM_Data_add(&(computeData->tol), NULL, 1, MPI_DOUBLE, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
    MAM_Data_add(&(computeData->rho), NULL, 1, MPI_DOUBLE, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
#if DIRECT_ERROR
    MAM_Data_add(&(computeData->direct_err), NULL, 1, MPI_DOUBLE, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
#endif

    MAM_Data_add(computeData->vlen, NULL, computeData->n, MPI_INT, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT); //TODO Calcular vlen
    MAM_Data_add(computeData->r0, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
    MAM_Data_add(computeData->b, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
#if PRECOND
    MAM_Data_add(computeData->diags, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
#endif
#if DIRECT_ERROR
    MAM_Data_add(computeData->x_exact, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
#endif

    MAM_Data_add(computeData->p, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    MAM_Data_add(computeData->r, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    MAM_Data_add(computeData->x, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);

}


void targets_update(Compute_data *computeData, user_redist_t *user_data) {
    size_t entry, total_qty;
    void *value = NULL;
    MPI_Datatype type;

    MPI_Comm_size(computeData->comm, &computeData->numP);
    MPI_Comm_rank(computeData->comm, &computeData->myId);

    entry = 0;
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);
    computeData->n = *((int *)value);
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);
    computeData->tol0 = *((double *)value);
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);
    computeData->t1 = *((double *)value);

    entry = 0;
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
    computeData->iter = *((int *)value); 
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
    computeData->tol = *((double *)value); 
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
    computeData->rho = *((double *)value); 
#if DIRECT_ERROR
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
    computeData->direct_err = *((double *)value); 
#endif

    entry = 0;
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
    computeData->vlen = ((int *)value);
    //computeData->vlen = user_data->recv_vlen;
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
    computeData->r0 = ((double *)value);
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
    computeData->b = ((double *)value);
#if PRECOND
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
    computeData->diags = ((double *)value);
#endif
#if DIRECT_ERROR
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
    computeData->x_exact = ((double *)value);
#endif

    entry = 0;
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    computeData->p = ((double *)value);
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    computeData->r = ((double *)value);
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    computeData->x = ((double *)value);
  
    int n_dist = computeData->matL.dim1;
    int n = computeData->n;
    CreateInts (&computeData->sizes, computeData->numP); 
    CreateInts (&computeData->dspls, computeData->numP); 
    CreateDoubles (&computeData->s, n_dist);
    CreateDoubles (&computeData->q, n_dist);
    CreateDoubles (&computeData->y, n_dist);
    CreateDoubles (&computeData->aux, n); 
#if PRECOND
    CreateDoubles (&computeData->p_hat, n_dist);
    CreateDoubles (&computeData->q_hat, n_dist);
#endif
#if DIRECT_ERROR
    CreateDoubles (&computeData->res_err, n_dist);
#endif
    ComputeMatrixSizes (n, computeData->sizes, computeData->dspls, computeData->comm);
    computeData->my_size = computeData->sizes[computeData->myId];
    computeData->my_dspl = computeData->dspls[computeData->myId];

    computeData->matL = user_data->other_subm;
    *user_data = empty_user_data;
    user_data->array_vptr = computeData->matL.vptr;
    user_data->array_vlen = computeData->vlen;
    user_data->array_vpos = computeData->matL.vpos;
    user_data->array_vval = computeData->matL.vval;
    user_data->comm = computeData->comm;
}


void user_func(void *args) {
    int local_flag, flag = 0;
    mam_user_reconf_t user_reconf;

    MAM_Get_Reconf_Info(&user_reconf);
    user_redist_t *user_data = (user_redist_t *) args;
    if(!user_data->initiated) {
      //targets_distribution_synch(user_reconf, user_data);
      //flag = 1;

      targets_distribution(user_reconf, user_data);
      user_data->initiated = 1;

      if(user_reconf.rank_state == MAM_PROC_NEW_RANK) {
        MPI_Waitall(2, user_data->reqs, MPI_STATUSES_IGNORE);

	flag = 1;
      }
    } else {
      MPI_Testall(2, user_data->reqs, &local_flag, MPI_STATUSES_IGNORE);
      MPI_Allreduce(&local_flag, &flag, 1, MPI_INT, MPI_MIN, user_data->comm);
    }

    if(flag) MAM_Resume_redistribution(NULL);
}


void dump(Compute_data *computeData) {
  int i;

  if(computeData->myId == 0) printf("TamBL="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {
      printf("%d, ", computeData->my_size);
    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }
  if(computeData->myId == 0) printf("\n"); 
  fflush(stdout); MPI_Barrier(computeData->comm);

  if(computeData->myId == 0) printf("Vlen="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%d, ", computeData->vlen[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }
iker_martin's avatar
iker_martin committed
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
761
762
763
764
765
766
767
768
769
770
771
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
869
870
871
872

  if(computeData->myId == 0) printf("Vptr="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size+1; j++) {
        printf("%d, ", computeData->matL.vptr[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("Vpos="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->matL.vptr[computeData->my_size]; j++) {
        printf("%d, ", computeData->matL.vpos[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("Vval="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->matL.vptr[computeData->my_size]; j++) {
        printf("%lf, ", computeData->matL.vval[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\nX="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->x[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\nB="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->b[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\nr="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->r[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\nr0="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->r0[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\np="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->p[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\ndiags="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->diags[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\np_hat="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->p_hat[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }
/*
  if(computeData->myId == 0) printf("\nq_hat="); 
873
  fflush(stdout); MPI_Barrier(computeData->comm);
iker_martin's avatar
iker_martin committed
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
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->q_hat[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\ns="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->s[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }

  if(computeData->myId == 0) printf("\ny="); 
  fflush(stdout); MPI_Barrier(computeData->comm);
  for(i=0; i<computeData->numP; i++) {
    if(computeData->myId == i) {

      for(int j=0; j<computeData->my_size; j++) {
        printf("%lf, ", computeData->y[j]);
      }

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(computeData->comm);
  }
  */
917
}
iker_martin's avatar
iker_martin committed
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

/*
typedef struct {
  double tol, tol0;
  int iter, n;

  double rho;
  double *x, *b;
  double *s, *q, *r, *p, *r0, *y, *p_hat, *q_hat;
  double *aux;
  SparseMatrix matL;

#if PRECOND
  double *diags;
#endif
#if DIRECT_ERROR
  double *res_err, *x_exact;
  double direct_err;
#endif
  double t1;

  int *sizes, *dspls;
  int my_size, my_dspl;
  int *vlen;

  int myId, numP;
  MPI_Comm comm;
} Compute_data;
*/