ConjugateGradient.c 43.4 KB
Newer Older
iker_martin's avatar
iker_martin committed
1
2
3
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
4
5
6
//#include <mkl_blas.h>
//#include <mkl_spblas.h>
#include "mymkl.h"
iker_martin's avatar
iker_martin committed
7
8
9
10
#include "ScalarVectors.h"
#include "SparseMatrices.h"
#include <mpi.h>
#include <string.h>
11
#include "../malleability/MAM.h"
12
#include <sys/prctl.h>
iker_martin's avatar
iker_martin committed
13

iker_martin's avatar
iker_martin committed
14
#include <unistd.h>
15

iker_martin's avatar
iker_martin committed
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
//#define ONLY_SYM 0
#define ROOT 0
//#define DEBUG 0
#define MAX_PROCS_SET 16

typedef struct {
  double umbral, tol;
  int iter, maxiter, n;

  double beta, rho, alpha;
  double *res, *z, *d, *vec;
  SparseMatrix subm;
  double *d_full;

  int *dist_elem, *displs_elem;
  int *dist_rows, *displs_rows;
  int *vlen;
} Compute_data;

struct Dist_data {
  int ini;
  int fin;

  int tamBl; // Numero de filas
  int tot_r; // Total de filas en la matriz

  int myId;
  int numP;

  int numP_parents;

  MPI_Comm comm, comm_children, comm_parents;
  MPI_Datatype scalars, arrays;
};

51
52
typedef struct {
  SparseMatrix other_subm;
53
54
  int *recv_vlen;
  int *array_vptr, *array_vlen, *array_vpos, n, initiated;
55
56
57
58
59
60
61
  double start_time, *array_vval;
  MPI_Comm comm;
  MPI_Request reqs[2];
} user_redist_t;

static const user_redist_t empty_user_data = {
  .array_vptr = NULL,
62
63
  .array_vlen = NULL,
  .recv_vlen  = NULL,
64
65
  .array_vpos = NULL,
  .array_vval = NULL,
66
  .n = 0,
67
68
69
70
  .initiated = 0,
  .comm = MPI_COMM_NULL
};

71
void dump(Compute_data *computeData, struct Dist_data *dist_data); //FIXME Delte me
72

iker_martin's avatar
iker_martin committed
73
74
75
76
77
78
void init_app(Compute_data *computeData, struct Dist_data *dist_data, char* argv[]);
void get_mat_dist(Compute_data *computeData, struct Dist_data dist_data, SparseMatrix mat);
void get_rows_dist(Compute_data *computeData, int numP, int n);
void mat_alloc(Compute_data *computeData, SparseMatrix mat, struct Dist_data dist_data);
void computeSolution(Compute_data computeData, double **subsol, SparseMatrix mat, int myId, double **full_vec);
void pre_compute(Compute_data *computeData, struct Dist_data dist_data, double *subsol, double *full_vec);
79
80
int compute(Compute_data *computeData, struct Dist_data *dist_data, user_redist_t *user_data);
void free_computeData(Compute_data *computeData);
iker_martin's avatar
iker_martin committed
81
82
83

//===================================MALLEABILITY FUNCTIONS====================================================

84
85
void originals_set_data(struct Dist_data *dist_data, Compute_data *computeData, int num_target);
void user_func(void *args);
86
void print_counts2(struct Dist_data data_dist, int *xcounts, int *xdispls, int size, int include_zero, const char* name);
87
void targets_distribution(mam_user_reconf_t user_reconf, user_redist_t *user_data);
88
void targets_distribution_synch(mam_user_reconf_t user_reconf, user_redist_t *user_data);
89
90
void targets_update(struct Dist_data *dist_data, Compute_data *computeData, user_redist_t *user_data);
void print_global_results(double start_time);
iker_martin's avatar
iker_martin committed
91
92
//----------------------------------------------------------------------------------------------------
void get_dist(int total_r, int id, int numP, struct Dist_data *dist_data);
93
void prepare_redist_counts(int *counts, int *displs, int numP_other, int offset, struct Dist_data dist_data, int *vptr);
94
void prepare_redist_counts_vlen(int *counts, int *displs, int numP_other, int offset, struct Dist_data dist_data);
95
96
void set_counts(int id, int numP, struct Dist_data data_dist, int offset, int *sendcounts);
void getIds_intercomm(struct Dist_data dist_data, int numP_other, int *idS);
iker_martin's avatar
iker_martin committed
97
98
99
//----------------------------------------------------------------------------------------------------

int main (int argc, char *argv[]) {
100
	int init_numP;
101
	int req;
iker_martin's avatar
iker_martin committed
102
        Compute_data computeData;
103
        user_redist_t user_data;
iker_martin's avatar
iker_martin committed
104

105
        computeData.z = NULL; computeData.d_full = NULL, computeData.d = NULL;
iker_martin's avatar
iker_martin committed
106
107
108
        computeData.vec = NULL; computeData.res = NULL;
        computeData.dist_elem = NULL; computeData.displs_elem = NULL;
        computeData.dist_rows = NULL; computeData.displs_rows = NULL;
109
110
	computeData.subm.vptr = NULL;
	computeData.vlen = NULL;
iker_martin's avatar
iker_martin committed
111

112
        int num_targets = 1;
iker_martin's avatar
iker_martin committed
113
        struct Dist_data dist_data;
114
115
        if (argc >= 3) {
          num_targets = atoi(argv[2]);
116
	}
iker_martin's avatar
iker_martin committed
117

118
        MPI_Init_thread(&argc, &argv, MPI_THREAD_MULTIPLE, &req);
119
120
        MPI_Comm_size(MPI_COMM_WORLD, &dist_data.numP);
        MPI_Comm_rank(MPI_COMM_WORLD, &dist_data.myId);
iker_martin's avatar
iker_martin committed
121
	dist_data.comm = MPI_COMM_WORLD;
122
123
        user_data = empty_user_data;
	user_data.comm = dist_data.comm;
iker_martin's avatar
iker_martin committed
124

125
        prctl(PR_SET_PTRACER, PR_SET_PTRACER_ANY, 0, 0, 0);
126
        int new_group = MAM_Init(ROOT, &dist_data.comm, argv[0], user_func, (void *) &user_data);
iker_martin's avatar
iker_martin committed
127
128
129

	if( !new_group ) { //First set of processes
	  init_app(&computeData, &dist_data, argv);
130
131
          originals_set_data(&dist_data, &computeData, num_targets);

132
	  init_numP = dist_data.numP;
133
	  user_data.n = computeData.n;
134
	  user_data.array_vptr = computeData.subm.vptr;
135
	  user_data.array_vlen = computeData.vlen;
136
137
	  user_data.array_vpos = computeData.subm.vpos;
  	  user_data.array_vval = computeData.subm.vval;
138
	  MPI_Barrier(MPI_COMM_WORLD);
139
          user_data.start_time = MPI_Wtime();
iker_martin's avatar
iker_martin committed
140
        } else {
141
          targets_update(&dist_data, &computeData, &user_data);
iker_martin's avatar
iker_martin committed
142
143
	}

144
	compute(&computeData, &dist_data, &user_data);
145

146
147
148
149
150
        MPI_Barrier(dist_data.comm);
        if(dist_data.myId == ROOT) {
  	  print_global_results(user_data.start_time);
	  printf ("End(%d) --> (%d,%20.10e)\n", computeData.n, computeData.iter, computeData.tol);
	}
iker_martin's avatar
iker_martin committed
151
152

	// End of CG
153
154
        MAM_Finalize();
        free_computeData(&computeData);
155
	if(init_numP > num_targets && dist_data.myId == 0) {
iker_martin's avatar
iker_martin committed
156
157
          MPI_Abort(MPI_COMM_WORLD, -100);
	}
158
	if(dist_data.comm != MPI_COMM_WORLD && dist_data.comm != MPI_COMM_NULL) MPI_Comm_free(&dist_data.comm);
iker_martin's avatar
iker_martin committed
159
160
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
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
        MPI_Finalize();
}

/*
 * Init application data before
 * starting iterative computation
 */
void init_app(Compute_data *computeData, struct Dist_data *dist_data, char* argv[]) {
	SparseMatrix mat, sym;
        double *full_vec = NULL;
	double *subsol = NULL;

        if(dist_data->myId == ROOT) {
#ifdef ONLY_SYM
  	  printf ("Working with symmetric format\n");
	  CreateSparseMatrixHB (argv[1], &mat, 1);
#else
	  printf ("Working with general format\n");
	  CreateSparseMatrixHB (argv[1], &sym, 1);
	  DesymmetrizeSparseMatrices (sym, &mat);
	  RemoveSparseMatrix (&sym);
#endif
          computeData->n = mat.dim1;
        }

        // Communicate number of rows to distribute and number of elements in the matrix
        MPI_Bcast(&computeData->n, 1, MPI_INT, ROOT, MPI_COMM_WORLD);

	// Each process calcules their own distribution
        get_dist(computeData->n, dist_data->myId, dist_data->numP, dist_data);
        
        if(dist_data->myId == ROOT) { // ROOT gets rows and vpos/vval distribution
	  get_mat_dist(computeData, *dist_data, mat);
          TransformHeadertoLength(mat.vptr, computeData->n); // From vptr to vlen
        } else { // Non ROOT proceses gets row distribution
          get_rows_dist(computeData, dist_data->numP, computeData->n);
          CreateInts (&computeData->dist_elem, dist_data->numP*2);
	  InitInts (computeData->dist_elem, dist_data->numP * 2, 0.0, 0); 
          computeData->displs_elem = computeData->dist_elem + dist_data->numP;
        }
        // Allocate for each process their submatrix and get their distribution from ROOT
	mat_alloc(computeData, mat, *dist_data);

	computeSolution(*computeData, &subsol, mat, dist_data->myId, &full_vec);
	pre_compute(computeData, *dist_data, subsol, full_vec);

        //Free Initial data
	RemoveDoubles(&subsol);
	RemoveDoubles(&full_vec);
        if(dist_data->myId == ROOT) {
	  RemoveSparseMatrix(&mat);
        }
}

/*
 * MPI Dist
 * Broadcast the vptr array and each process gets the data that corresponds to itself.
 *
 * mat.vptr must be in vlen format to work correctly
 */
void get_mat_dist(Compute_data *computeData, struct Dist_data dist_data, SparseMatrix mat) {
	int i, j;
        struct Dist_data dist_data_aux;

#ifdef DEBUG
        if(dist_data.myId == ROOT) printf("Distribuyendo vptr\n");
#endif
        CreateInts (&computeData->dist_rows, dist_data.numP);
        CreateInts (&computeData->displs_rows, dist_data.numP);
        CreateInts (&computeData->dist_elem, dist_data.numP*2);
        computeData->displs_elem = computeData->dist_elem + dist_data.numP;

        InitInts (computeData->dist_rows, dist_data.numP, 0, 0);
        InitInts (computeData->displs_rows, dist_data.numP, 0, 0);
        InitInts (computeData->dist_elem, dist_data.numP*2, 0, 0);

	// Fill dist_rows and dist_elem so each process can make ScatterV or GatherV calls
        for(i=0; i<dist_data.numP; i++) {
          get_dist(computeData->n, i, dist_data.numP, &dist_data_aux);

          computeData->dist_rows[i] = dist_data_aux.tamBl;
          computeData->dist_elem[i] = mat.vptr[dist_data_aux.fin] - mat.vptr[dist_data_aux.ini];

          // Fill displacements
          if(i!=0) { 
            computeData->displs_elem[i] = computeData->displs_elem[i-1] + computeData->dist_elem[i-1];
            computeData->displs_rows[i] = computeData->displs_rows[i-1] + computeData->dist_rows[i-1];
          }
        }

#ifdef DEBUG
        printf("Proc %d almacena %d filas con %d elementos\n", dist_data.myId, computeData->dist_rows[dist_data.myId], computeData->dist_elem[dist_data.myId]);
        fflush(stdout);
#endif
}

/*
 * MPI Dist
 * Get the rows distribution of n rows in a given number of processes
 */
void get_rows_dist(Compute_data *computeData, int numP, int n) {
	int i, j;
        struct Dist_data dist_data;

        CreateInts (&(computeData->dist_rows), numP);
        CreateInts (&(computeData->displs_rows), numP);

        InitInts (computeData->dist_rows, numP, 0, 0);
        InitInts (computeData->displs_rows, numP, 0, 0);

	// Fill dist_rows and dist_elem so each process can make ScatterV or GatherV calls
        for(i=0; i<numP; i++) {
          get_dist(n, i, numP, &dist_data);
          computeData->dist_rows[i] = dist_data.tamBl;

          // Fill displacements
          if(i!=0) { 
            computeData->displs_rows[i] = computeData->displs_rows[i-1] + computeData->dist_rows[i-1];
          }
        }
}

/*
 * Matrix allocation
 *
 * The matrix that each process will use is allocated and
 * their vptr array initialised.
 *
 * MPI Dist
 * Distribute vpos and vvalues data among processes
 * Both arrays have the same distribution
 */
void mat_alloc(Compute_data *computeData, SparseMatrix mat, struct Dist_data dist_data) {
	int i;
	int elems; // Number of elements this process has
#ifdef DEBUG
        if(dist_data.myId == ROOT) printf("Distribuyendo vpos y vvalue\n");
#endif

	// dist_rows[myId] is the number of rows, n the number of columns, and dist_elem[myId] is the number of elements this process will have in the matrix
        CreateSparseMatrixVptr(&(computeData->subm), dist_data.tamBl, computeData->n, 0);
        computeData->subm.vptr[0] = 0;

        MPI_Scatterv((mat.vptr)+1, computeData->dist_rows, computeData->displs_rows, MPI_INT, (computeData->subm.vptr)+1, dist_data.tamBl, MPI_INT, ROOT, MPI_COMM_WORLD);

304
305
306
        CreateInts(&(computeData->vlen), dist_data.tamBl);
        for(i=0; i<dist_data.tamBl; i++) {
          computeData->vlen[i] = computeData->subm.vptr[i+1];
307
        }
iker_martin's avatar
iker_martin committed
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
        TransformLengthtoHeader(computeData->subm.vptr, computeData->subm.dim1); // The array is converted from vlen to vptr
        elems = computeData->subm.vptr[dist_data.tamBl];
        CreateSparseMatrixValues(&(computeData->subm), dist_data.tamBl, computeData->n, elems, 0);

        MPI_Scatterv(mat.vpos, computeData->dist_elem, computeData->displs_elem, MPI_INT,    computeData->subm.vpos, elems, MPI_INT,    ROOT, MPI_COMM_WORLD);
        MPI_Scatterv(mat.vval, computeData->dist_elem, computeData->displs_elem, MPI_DOUBLE, computeData->subm.vval, elems, MPI_DOUBLE, ROOT, MPI_COMM_WORLD);

	// Free elem arrays, as they are not going to be used again
        RemoveInts (&computeData->dist_elem);
}

/*
 * Compute solution
 */
void computeSolution(Compute_data computeData, double **subsol, SparseMatrix mat, int myId, double **full_vec) {
        
	CreateDoubles (subsol, computeData.dist_rows[myId]);
	InitDoubles (*subsol, computeData.dist_rows[myId], 0.0, 0.0);
	CreateDoubles(full_vec, computeData.n);
	InitDoubles (*full_vec, computeData.n, 1.0, 0.0);

//Compute SOLUTION
#ifdef ONLY_SYM
	ProdSymSparseMatrixVector (computeData.subm, *full_vec, *subsol);                  // sol += A * x
#else
	ProdSparseMatrixVector (computeData.subm, *full_vec, *subsol);                    	// sol += A * x
#endif
335

iker_martin's avatar
iker_martin committed
336
337
338
339
340
341
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
370
371
#ifdef DEBUG
	int aux, i;
	double *solD = NULL, *sol = NULL;
	if(myId == ROOT) {
          printf("Computing solution\n");
	  CreateDoubles (&sol, computeData.n);
	  CreateDoubles (&solD, computeData.n);
	  InitDoubles (sol, computeData.n, 0.0, 0.0);
	  InitDoubles (solD, computeData.n, 0.0, 0.0);

          TransformLengthtoHeader(mat.vptr, mat.dim1); // vlen to vptr (At mat_alloc was needed as vlen)
        }

	MPI_Gatherv(*subsol, computeData.dist_rows[myId], MPI_DOUBLE, sol, computeData.dist_rows, computeData.displs_rows, MPI_DOUBLE, ROOT, MPI_COMM_WORLD);

        if(myId == ROOT) {

#ifdef ONLY_SYM
	  ProdSymSparseMatrixVector (mat, *full_vec, solD);                   // solD += A * x
#else
	  ProdSparseMatrixVector (mat, *full_vec, solD);                      // solD += A * x
#endif // ONLY_SIM
          aux = 1;
          printf("Checking sol array is ok\n");
          for(i=0; i<mat.dim1; i++) {
            if(sol[i] != solD[i]) {
              printf("[%d]Expected %lf - Result %lf\n", i, solD[i],sol[i]);
              aux = 0;
            }
          }
          if(aux) printf("sol array is correct\n");
          
        }
	RemoveDoubles (&sol);
	RemoveDoubles (&solD);
#endif // DEBUG
372

iker_martin's avatar
iker_martin 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
}

/*
 * Realiza los preparativos para pasar al bucle de computo principal
 * inicializando los datos y realizando una primera iteración
 */
void pre_compute(Compute_data *computeData, struct Dist_data dist_data, double *subsol, double *full_vec) {

	int IZERO = 0, IONE = 1; 
	double DONE = 1.0, DMONE = -1.0, DZERO = 0.0;

        if(dist_data.myId == ROOT) {
	  printf("Start CG\n");
        }

        computeData->res = NULL; computeData->z = NULL; computeData->d = NULL;
	computeData->umbral = 1.0e-8;

	CreateDoubles(&computeData->res, dist_data.tamBl); 
	CreateDoubles(&computeData->z, dist_data.tamBl); 
	CreateDoubles(&computeData->d, dist_data.tamBl);
	CreateDoubles (&computeData->vec, dist_data.tamBl);
	CreateDoubles (&computeData->d_full, computeData->n);

	InitDoubles (computeData->vec, dist_data.tamBl, DZERO, DZERO); // x = 0
	InitDoubles (full_vec, computeData->n, DZERO, DZERO); // full_x = 0
	
	computeData->iter = 0;

#ifdef ONLY_SYM
	ProdSymSparseMatrixVector (computeData->subm, full_vec, computeData->z);                     				// z += A * full_x
//	mkl_dcsrsymv ("U", &n, mat.vval, mat.vptr, mat.vpos, vec, z); 			   // z = A * full_x
#else
	ProdSparseMatrixVector (computeData->subm, full_vec, computeData->z);                       				// z += A * full_x
#endif
iker_martin's avatar
iker_martin committed
408
409
410
	rcopy (&(dist_data.tamBl), subsol, &IONE, computeData->res, &IONE);             					// res = b
	raxpy (&(dist_data.tamBl), &DMONE, computeData->z, &IONE, computeData->res, &IONE);           				// res -= z
	//rcopy (&(computeData.subm.dim1), computeData.res, &IONE, &(computeData.d+computeData.displs_rows[myId]), &IONE);      // d_full = res
iker_martin's avatar
iker_martin committed
411
        MPI_Allgatherv(computeData->res, dist_data.tamBl, MPI_DOUBLE, computeData->d_full, computeData->dist_rows, computeData->displs_rows, MPI_DOUBLE, MPI_COMM_WORLD);
iker_martin's avatar
iker_martin committed
412
413
	rcopy (&(dist_data.tamBl), &(computeData->d_full[dist_data.ini]), &IONE, computeData->d, &IONE);             		// d = d_full[ini] to d_full[ini+tamBl]
	computeData->beta = rdot (&(dist_data.tamBl), computeData->res, &IONE, computeData->res, &IONE);      			// beta = res' * res
iker_martin's avatar
iker_martin committed
414
415
416
417
418
419
420
        MPI_Allreduce(MPI_IN_PLACE, &computeData->beta, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
	computeData->tol = sqrt (computeData->beta);                                          			   		// tol = sqrt(beta) = norm (res)
}

/*
 * Bucle de computo principal
 */
421
int compute(Compute_data *computeData, struct Dist_data *dist_data, user_redist_t *user_data) {
iker_martin's avatar
iker_martin committed
422
423
	int IZERO = 0, IONE = 1; 
	double DONE = 1.0, DMONE = -1.0, DZERO = 0.0;
424
	int state = MAM_UNRESERVED;
iker_martin's avatar
iker_martin committed
425
        int ended_loop = 1;
426
	int reconfigure = 0, rec_iter = 500;
iker_martin's avatar
iker_martin committed
427

428
        computeData->maxiter = 1000;
iker_martin's avatar
iker_martin committed
429
430
431

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

432
433
	        MPI_Allgatherv(computeData->d, dist_data->tamBl, MPI_DOUBLE, computeData->d_full, 
				computeData->dist_rows, computeData->displs_rows, MPI_DOUBLE, dist_data->comm);		// d_full = Gather(d)
iker_martin's avatar
iker_martin committed
434
435
436
437
438
439
//      	COMPUTATION
#ifdef ONLY_SYM
		ProdSymSparseMatrixVector (computeData->subm, computeData->d_full, computeData->z);                     // z += A * d_full
#else
		ProdSparseMatrixVector (computeData->subm, computeData->d_full, computeData->z);                    	// z += A * d_full
#endif
iker_martin's avatar
iker_martin committed
440
        	computeData->rho = rdot (&(dist_data->tamBl), computeData->d, &IONE, computeData->z, &IONE);		// rho = (d * z)
441
	        MPI_Allreduce(MPI_IN_PLACE, &computeData->rho, 1, MPI_DOUBLE, MPI_SUM, dist_data->comm);		// Reduce(rho, SUM)
iker_martin's avatar
iker_martin committed
442
		computeData->rho = computeData->beta / computeData->rho;                 		                // rho = beta / aux
iker_martin's avatar
iker_martin committed
443
		raxpy (&(dist_data->tamBl), &computeData->rho, computeData->d, &IONE, computeData->vec, &IONE);		// x += rho * d
iker_martin's avatar
iker_martin committed
444
		computeData->rho = -computeData->rho;
iker_martin's avatar
iker_martin committed
445
		raxpy (&(dist_data->tamBl), &computeData->rho, computeData->z, &IONE, computeData->res, &IONE);         // res -= rho * z
iker_martin's avatar
iker_martin committed
446
		computeData->alpha = computeData->beta;                                               		        // alpha = beta
iker_martin's avatar
iker_martin committed
447
		computeData->beta = rdot (&(dist_data->tamBl), computeData->res, &IONE, computeData->res, &IONE);       // beta = res' * res
448
	        MPI_Allreduce(MPI_IN_PLACE, &computeData->beta, 1, MPI_DOUBLE, MPI_SUM, dist_data->comm);		// Reduce(beta, SUM)
iker_martin's avatar
iker_martin committed
449
		computeData->alpha = computeData->beta / computeData->alpha;                                       	// alpha = beta / alpha
iker_martin's avatar
iker_martin committed
450
451
		rscal (&(dist_data->tamBl), &computeData->alpha, computeData->d, &IONE);                   		// d = alpha * d
		raxpy (&(dist_data->tamBl), &DONE, computeData->res, &IONE, computeData->d, &IONE);        		// d += res
452

iker_martin's avatar
iker_martin committed
453
454
455
		computeData->tol = sqrt (computeData->beta);                                          			// tol = sqrt(beta) = norm (res)
		computeData->iter++;

456
                if (computeData->iter == rec_iter) { reconfigure = 1;}
457
		if (reconfigure) {
458
459
460
461
462
		  MAM_Checkpoint(&state, MAM_CHECK_COMPLETION, user_func, (void *) user_data);
		  if(state == MAM_COMPLETED) {
		    reconfigure = 0; 
                    free_computeData(computeData);
                    targets_update(dist_data, computeData, user_data);
463
464
465
		  }

		}
466
	
iker_martin's avatar
iker_martin committed
467
	}
468
469
470
471
472
473

	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);
	}
iker_martin's avatar
iker_martin committed
474
475
476
477
478
479
#ifdef DEBUG
	if(dist_data->myId == ROOT) printf ("Ended loop\n");
#endif
	return ended_loop;
}

480
void dump(Compute_data *computeData, struct Dist_data *dist_data) {
481
  int i;
482

483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
  if(dist_data->myId == 0) printf("TamBL="); 
  fflush(stdout); MPI_Barrier(dist_data->comm);
  for(i=0; i<dist_data->numP; i++) {
    if(dist_data->myId == i) {
      printf("%d, ", dist_data->tamBl);
    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(dist_data->comm);
  }
  if(dist_data->myId == 0) printf("\n"); 
  fflush(stdout); MPI_Barrier(dist_data->comm);

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

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

    }
    fflush(stdout);
    sleep(1);
    MPI_Barrier(dist_data->comm);
  }
  if(dist_data->myId == 0) printf("\n"); 
  fflush(stdout); MPI_Barrier(dist_data->comm);

513
514
515
516
517
518
519
  if(dist_data->myId == 0) printf("Vptr="); 
  fflush(stdout); MPI_Barrier(dist_data->comm);
  for(i=0; i<dist_data->numP; i++) {
    if(dist_data->myId == i) {

      printf("%d, ", computeData->subm.vptr[dist_data->tamBl]);
    }
520
521
    fflush(stdout);
    sleep(1);
522
523
524
525
526
527
528
529
530
531
532
533
534
    MPI_Barrier(dist_data->comm);
  }
  if(dist_data->myId == 0) printf("\n"); 
  fflush(stdout); MPI_Barrier(dist_data->comm);


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

      printf("%lf, ", computeData->tol);
    }
535
536
    fflush(stdout);
    sleep(1);
537
538
539
540
541
542
543
544
545
546
547
548
549
    MPI_Barrier(dist_data->comm);
  }
  if(dist_data->myId == 0) printf("\n"); 
  fflush(stdout); MPI_Barrier(dist_data->comm);


  if(dist_data->myId == 0) printf("Z[last]="); 
  fflush(stdout); MPI_Barrier(dist_data->comm);
  for(i=0; i<dist_data->numP; i++) {
    if(dist_data->myId == i) {

      printf("%lf, ", computeData->z[dist_data->tamBl-1]);
    }
550
551
    fflush(stdout);
    sleep(1);
552
553
554
555
556
557
558
559
560
561
562
563
    MPI_Barrier(dist_data->comm);
  }
  if(dist_data->myId == 0) printf("\n"); 
  fflush(stdout); MPI_Barrier(dist_data->comm);

  if(dist_data->myId == 0) printf("D[last]="); 
  fflush(stdout); MPI_Barrier(dist_data->comm);
  for(i=0; i<dist_data->numP; i++) {
    if(dist_data->myId == i) {

      printf("%lf, ", computeData->d[dist_data->tamBl-1]);
    }
564
565
    fflush(stdout);
    sleep(1);
566
567
568
569
570
571
572
573
574
575
576
577
    MPI_Barrier(dist_data->comm);
  }
  if(dist_data->myId == 0) printf("\n"); 
  fflush(stdout); MPI_Barrier(dist_data->comm);

  if(dist_data->myId == 0) printf("res[last]="); 
  fflush(stdout); MPI_Barrier(dist_data->comm);
  for(i=0; i<dist_data->numP; i++) {
    if(dist_data->myId == i) {

      printf("%lf, ", computeData->res[dist_data->tamBl-1]);
    }
578
579
    fflush(stdout);
    sleep(1);
580
581
582
583
584
585
586
587
588
589
590
    MPI_Barrier(dist_data->comm);
  }
  if(dist_data->myId == 0) printf("\n"); 
  fflush(stdout); MPI_Barrier(dist_data->comm);
  if(dist_data->myId == 0) printf("Vec[last]="); 
  fflush(stdout); MPI_Barrier(dist_data->comm);
  for(i=0; i<dist_data->numP; i++) {
    if(dist_data->myId == i) {

      printf("%lf, ", computeData->vec[dist_data->tamBl-1]);
    }
591
592
    fflush(stdout);
    sleep(1);
593
    MPI_Barrier(dist_data->comm);
594
  }
595
  if(dist_data->myId == 0) printf("\n"); 
596
597
  fflush(stdout); MPI_Barrier(dist_data->comm);
}
iker_martin's avatar
iker_martin committed
598

599
void free_computeData(Compute_data *computeData) {
600
	if(computeData->res != NULL) {
iker_martin's avatar
iker_martin committed
601
	RemoveDoubles (&computeData->res); 
602
603
	}
	if(computeData->z != NULL) {
iker_martin's avatar
iker_martin committed
604
        RemoveDoubles (&computeData->z); 
605
606
	}
	if(computeData->d != NULL) {
iker_martin's avatar
iker_martin committed
607
        RemoveDoubles (&computeData->d);
608
609
	}
	if(computeData->vec != NULL) {
iker_martin's avatar
iker_martin committed
610
	RemoveDoubles (&computeData->vec);
611
	}
iker_martin's avatar
iker_martin committed
612

613
	if(computeData->d_full != NULL) {
614
615
        RemoveDoubles (&computeData->d_full);
	}
616
	
617
	if(computeData->subm.vptr != NULL) {
iker_martin's avatar
iker_martin committed
618
	RemoveSparseMatrix2 (&computeData->subm);
619
620
621
	}

	if(computeData->dist_rows != NULL) {
iker_martin's avatar
iker_martin committed
622
        RemoveInts (&computeData->dist_rows);
623
624
	}
	if(computeData->displs_rows != NULL) {
iker_martin's avatar
iker_martin committed
625
        RemoveInts (&computeData->displs_rows);
626
	}
627
	
628
	if(computeData->vlen != NULL) {
iker_martin's avatar
iker_martin committed
629
        RemoveInts (&computeData->vlen);
630
	}
iker_martin's avatar
iker_martin committed
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
}

/*
 *  _____________________________________________________________________________________
 * ||                                                                                   ||
 * ||                                                                                   ||
 * ||                            DISTRIBUTION FUNCTIONS                                 ||
 * ||                                                                                   ||
 * ||                                                                                   ||
 * \_____________________________________________________________________________________/
*/

/*
 * Las siguientes funciones están todas relacionadas con la distribución de los datos
 * o procesos.
 */

/*
 * ========================================================================================
 * ========================================================================================
 * ========================PARENTS COMMUNICATION FUNCTIONS=================================
 * ========================================================================================
 * ========================================================================================
*/

/*
657
 * Función para declarar los datos a comunicar por parte de MAM
iker_martin's avatar
iker_martin committed
658
 */
659
660
void originals_set_data(struct Dist_data *dist_data, Compute_data *computeData, int num_target) {
    size_t index;
661

662
663
664
665
    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->umbral), NULL, 1, MPI_DOUBLE, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);
666
667
    //MAM_Data_add(&(computeData->n), NULL, 1, MPI_INT, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
    //MAM_Data_add(&(computeData->umbral), NULL, 1, MPI_DOUBLE, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
668
669
670
671
672
673
674
675
676
677
678
    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->beta), NULL, 1, MPI_DOUBLE, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);

    MAM_Data_add(computeData->d, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);

    MAM_Data_add(computeData->vec, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    MAM_Data_add(computeData->res, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    MAM_Data_add(computeData->z, NULL, computeData->n, MPI_DOUBLE, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);

    MAM_Data_add(computeData->vlen, NULL, computeData->n, MPI_INT, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
iker_martin's avatar
iker_martin committed
679
680
681
682
683
684
685
686
687
688
689
}

/*
 * ========================================================================================
 * ========================================================================================
 * ========================CHILDREN COMMUNICATION FUNCTIONS================================
 * ========================================================================================
 * ========================================================================================
*/

/*
690
 * Función llamada por MAM como callback.
iker_martin's avatar
iker_martin committed
691
 *
692
693
694
695
696
697
698
699
700
701
702
703
704
 * La misma realiza la redistribucion de datos por parte del usuario.
 * Como se usan comunicaciones no bloqueantes, primero se inicia
 * la comunicación y en las siguientes llamadas se comprueba si
 * la misma ha terminado.
 */
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) {
      MPI_Bcast(&user_data->start_time, 1, MPI_DOUBLE, 0, user_reconf.comm);
705
706
707
708

      //targets_distribution_synch(user_reconf, user_data);
      //flag = 1;

709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
      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);
}

/*
 * Funcion encargada de realizar la redistribucion de datos
727
 * asíncrona del usuario.
iker_martin's avatar
iker_martin committed
728
 *
729
730
 * Calcula el total de elementos a enviar/recibir por cada proceso
 * y tras ello llama a la funcion Ialltoallv dos veces.
iker_martin's avatar
iker_martin committed
731
 *
732
733
 * Además inicializa la memoria para aquellos procesos que vayan
 * a recibir datos.
iker_martin's avatar
iker_martin committed
734
 */
735
736
737
738
void targets_distribution(mam_user_reconf_t user_reconf, user_redist_t *user_data) {
    int i, n, offset, elems, numP, *vlen, *rank_states;
    int *scounts, *rcounts, *sdispls, *rdispls;
    size_t total_qty;
iker_martin's avatar
iker_martin committed
739
    void *value = NULL;
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
    struct Dist_data dist_data;
    MPI_Datatype type;

    int aux_int;
    int *recv_vpos = &aux_int;
    double aux_double;
    double *recv_vval = &aux_double;

    MPI_Comm_size(user_reconf.comm, &numP);
    scounts = calloc(numP, sizeof(int));
    sdispls = calloc(numP, sizeof(int));
    rcounts = calloc(numP, sizeof(int));
    rdispls = calloc(numP, sizeof(int));
    offset = 0;
	  
    rank_states = (int *) malloc(numP * sizeof(int));
    MPI_Allgather(&user_reconf.rank_state, 1, MPI_INT, rank_states, 1, MPI_INT, user_reconf.comm);

    MAM_Data_get_pointer(&value, 0, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
    vlen = ((int *)value);
    n = (int) total_qty;

    if(user_reconf.rank_state != MAM_PROC_ZOMBIE) {
      MPI_Comm_rank(user_data->comm, &dist_data.myId);
      dist_data.numP = user_reconf.numT;
      if(user_reconf.rank_state == MAM_PROC_NEW_RANK) {
	user_data->array_vpos = &aux_int;
	user_data->array_vval = &aux_double;
	for(i=0; i<user_reconf.numS; i++) {
          if(rank_states[i] == MAM_PROC_CONTINUE) {
            dist_data.myId += user_reconf.numS;
	    break;
	  }
	}
      }
      get_dist(n, dist_data.myId, dist_data.numP, &dist_data);
    
      CreateSparseMatrixVptr(&user_data->other_subm, dist_data.tamBl, n, 0);
      user_data->other_subm.vptr[0] = 0;
      //memcpy(user_data->other_subm.vptr+1, vlen, dist_data.tamBl * sizeof(int));
      for(i=0; i<dist_data.tamBl; i++) {
        user_data->other_subm.vptr[i+1] = vlen[i];
      }
      TransformLengthtoHeader(user_data->other_subm.vptr, user_data->other_subm.dim1); // The array is converted from vlen to vptr
      elems = user_data->other_subm.vptr[dist_data.tamBl];
      CreateSparseMatrixValues(&user_data->other_subm, dist_data.tamBl, n, elems, 0);
      recv_vpos = user_data->other_subm.vpos;
      recv_vval = user_data->other_subm.vval;

      prepare_redist_counts(rcounts, rdispls, user_reconf.numS, offset, dist_data, user_data->other_subm.vptr);
790
    //  print_counts2(dist_data, rcounts, rdispls, numP, 0, "TARGETS");
791
792
793
794
795
796
797
798
799
    } 

    if(user_reconf.rank_state != MAM_PROC_NEW_RANK) {
      MPI_Comm_rank(user_data->comm, &dist_data.myId);
      dist_data.numP = user_reconf.numS;
      get_dist(n, dist_data.myId, dist_data.numP, &dist_data);
      offset = (user_reconf.numS + user_reconf.numT) == numP ? 
	      user_reconf.numS : 0; 
      prepare_redist_counts(scounts, sdispls, user_reconf.numT, offset, dist_data, user_data->array_vptr);
800
    //  print_counts2(dist_data, scounts, sdispls, numP, 0, "SOURCES");
801
802
803
804
805
806
807
808
809
810
    }

    // COMUNICACION DE DATOS //
    MPI_Ialltoallv(user_data->array_vpos, scounts, sdispls, MPI_INT,    recv_vpos, rcounts, rdispls, MPI_INT,    user_reconf.comm, &user_data->reqs[0]);
    MPI_Ialltoallv(user_data->array_vval, scounts, sdispls, MPI_DOUBLE, recv_vval, rcounts, rdispls, MPI_DOUBLE, user_reconf.comm, &user_data->reqs[1]);

    free(rank_states);
    free(scounts); free(sdispls); free(rcounts); free(rdispls);
}

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
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
899
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
/*
 * Funcion encargada de realizar la redistribucion de datos
 * síncrona del usuario.
 *
 * Calcula el total de elementos a enviar/recibir por cada proceso
 * y tras ello llama a la funcion Ialltoallv dos veces.
 *
 * Además inicializa la memoria para aquellos procesos que vayan
 * a recibir datos.
 */
void targets_distribution_synch(mam_user_reconf_t user_reconf, user_redist_t *user_data) {
    int i, n, offset, elems, rank, numP, *vlen, *rank_states;
    int *scounts, *rcounts, *sdispls, *rdispls;
    size_t total_qty;
    void *value = NULL;
    struct Dist_data dist_data;
    MPI_Datatype type;

    int aux_int;
    int *recv_vpos = &aux_int;
    double aux_double;
    double *recv_vval = &aux_double;
    user_data->recv_vlen = &aux_int;


    MPI_Comm_rank(user_reconf.comm, &rank);
    MPI_Comm_size(user_reconf.comm, &numP);
    scounts = calloc(numP, sizeof(int));
    sdispls = calloc(numP, sizeof(int));
    rcounts = calloc(numP, sizeof(int));
    rdispls = calloc(numP, sizeof(int));
    offset = 0;
	  
    rank_states = (int *) malloc(numP * sizeof(int));
    MPI_Allgather(&user_reconf.rank_state, 1, MPI_INT, rank_states, 1, MPI_INT, user_reconf.comm);

    if(rank == 0) n = user_data->n;
    MPI_Bcast(&n, 1, MPI_INT, 0, user_reconf.comm);

    if(user_reconf.rank_state != MAM_PROC_ZOMBIE) {
      MPI_Comm_rank(user_data->comm, &dist_data.myId);
      dist_data.numP = user_reconf.numT;
      if(user_reconf.rank_state == MAM_PROC_NEW_RANK) {
	user_data->array_vlen = &aux_int;
	for(i=0; i<user_reconf.numS; i++) {
          if(rank_states[i] == MAM_PROC_CONTINUE) {
            dist_data.myId += user_reconf.numS;
	    break;
	  }
	}
      }
      get_dist(n, dist_data.myId, dist_data.numP, &dist_data);
      CreateInts(&user_data->recv_vlen, dist_data.tamBl);

      prepare_redist_counts_vlen(rcounts, rdispls, user_reconf.numS, offset, dist_data);
    //  print_counts2(dist_data, rcounts, rdispls, numP, 0, "TARGETS");
    } 

    if(user_reconf.rank_state != MAM_PROC_NEW_RANK) {
      MPI_Comm_rank(user_data->comm, &dist_data.myId);
      dist_data.numP = user_reconf.numS;
      get_dist(n, dist_data.myId, dist_data.numP, &dist_data);
      offset = (user_reconf.numS + user_reconf.numT) == numP ? 
	      user_reconf.numS : 0; 
      prepare_redist_counts_vlen(scounts, sdispls, user_reconf.numT, offset, dist_data);
    //  print_counts2(dist_data, scounts, sdispls, numP, 0, "SOURCES");
    }

    // COMUNICACION DE DATOS //
    MPI_Alltoallv(user_data->array_vlen, scounts, sdispls, MPI_INT, user_data->recv_vlen, rcounts, rdispls, MPI_INT, user_reconf.comm);
    free(scounts); free(sdispls); free(rcounts); free(rdispls);
    scounts = calloc(numP, sizeof(int));
    sdispls = calloc(numP, sizeof(int));
    rcounts = calloc(numP, sizeof(int));
    rdispls = calloc(numP, sizeof(int));
    offset = 0;

    if(user_reconf.rank_state != MAM_PROC_ZOMBIE) {
      MPI_Comm_rank(user_data->comm, &dist_data.myId);
      dist_data.numP = user_reconf.numT;
      if(user_reconf.rank_state == MAM_PROC_NEW_RANK) {
	user_data->array_vlen = &aux_int;
	for(i=0; i<user_reconf.numS; i++) {
          if(rank_states[i] == MAM_PROC_CONTINUE) {
            dist_data.myId += user_reconf.numS;
	    break;
	  }
	}
      }
      get_dist(n, dist_data.myId, dist_data.numP, &dist_data);
      CreateSparseMatrixVptr(&user_data->other_subm, dist_data.tamBl, n, 0);
      user_data->other_subm.vptr[0] = 0;
      //memcpy(user_data->other_subm.vptr+1, vlen, dist_data.tamBl * sizeof(int));
      for(i=0; i<dist_data.tamBl; i++) {
        user_data->other_subm.vptr[i+1] = user_data->recv_vlen[i];
      }
      TransformLengthtoHeader(user_data->other_subm.vptr, user_data->other_subm.dim1); // The array is converted from vlen to vptr
      elems = user_data->other_subm.vptr[dist_data.tamBl];
      CreateSparseMatrixValues(&user_data->other_subm, dist_data.tamBl, n, elems, 0);
      recv_vpos = user_data->other_subm.vpos;
      recv_vval = user_data->other_subm.vval;

      prepare_redist_counts(rcounts, rdispls, user_reconf.numS, offset, dist_data, user_data->other_subm.vptr);
    }
    if(user_reconf.rank_state != MAM_PROC_NEW_RANK) {
      MPI_Comm_rank(user_data->comm, &dist_data.myId);
      dist_data.numP = user_reconf.numS;
      get_dist(n, dist_data.myId, dist_data.numP, &dist_data);
      offset = (user_reconf.numS + user_reconf.numT) == numP ? 
	      user_reconf.numS : 0; 
      prepare_redist_counts(scounts, sdispls, user_reconf.numT, offset, dist_data, user_data->array_vptr);
    }
    // COMUNICACION DE DATOS //
    MPI_Alltoallv(user_data->array_vpos, scounts, sdispls, MPI_INT,    recv_vpos, rcounts, rdispls, MPI_INT,    user_reconf.comm);
    MPI_Alltoallv(user_data->array_vval, scounts, sdispls, MPI_DOUBLE, recv_vval, rcounts, rdispls, MPI_DOUBLE, user_reconf.comm);

    free(rank_states);
    free(scounts); free(sdispls); free(rcounts); free(rdispls);
}

931
932
933
934
935
936
void targets_update(struct Dist_data *dist_data, Compute_data *computeData, user_redist_t *user_data) {
    int IONE = 1, i;
    size_t entry, total_qty;
    double start_time;
    void *value = NULL;
    MPI_Datatype type;
937

938
939
940
941
942
    MPI_Comm_size(dist_data->comm, &dist_data->numP);
    MPI_Comm_rank(dist_data->comm, &dist_data->myId);

    entry = 0;
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);
943
    //MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
944
    computeData->n = *((int *)value);
945
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_CONSTANT);
946
    //MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
947
948
949
950
951
952
953
    computeData->umbral = *((double *)value);

    get_dist(computeData->n, dist_data->myId, dist_data->numP, dist_data);
    get_rows_dist(computeData, dist_data->numP, computeData->n);

    entry = 0;
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
iker_martin's avatar
iker_martin committed
954
    computeData->iter = *((int *)value);
955
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
956
    computeData->tol = *((double *)value);
957
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_REPLICATED, MAM_DATA_VARIABLE);
958
959
    computeData->beta = *((double *)value);

960
961
    entry = 0;
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
962
    computeData->d = ((double *)value);
963

964
965
    CreateDoubles(&computeData->d_full, computeData->n);
    rcopy (&(dist_data->tamBl), computeData->d, &IONE, &(computeData->d_full[dist_data->ini]), &IONE); // d_full[ini] to d_full[ini+tamBl] = d 
iker_martin's avatar
iker_martin committed
966

967
968
969
970
971
972
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    computeData->vec = ((double *)value);
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    computeData->res = ((double *)value);
    MAM_Data_get_pointer(&value, entry++, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_VARIABLE);
    computeData->z = ((double *)value);
973

974
    MAM_Data_get_pointer(&value, 0, &total_qty, &type, MAM_DATA_DISTRIBUTED, MAM_DATA_CONSTANT);
975
    computeData->vlen = ((int *)value);
976
    //computeData->vlen = user_data->recv_vlen;
977

978
979
980
981
982
    start_time = user_data->start_time;
    computeData->subm = user_data->other_subm;
    *user_data = empty_user_data;
    user_data->start_time = start_time;
    user_data->array_vptr = computeData->subm.vptr;
983
    user_data->array_vlen = computeData->vlen;
984
985
986
    user_data->array_vpos = computeData->subm.vpos;
    user_data->array_vval = computeData->subm.vval;
    user_data->comm = dist_data->comm;
iker_martin's avatar
iker_martin committed
987
988
989
990
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
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
}

/*
 * ========================================================================================
 * ========================================================================================
 * ================================DISTRIBUTION FUNCTIONS==================================
 * ========================================================================================
 * ========================================================================================
*/

/*
 * Obtiene para el Id que se pasa junto a su
 * numero de procesos total, con cuantas filas (tamBl),
 * elementos por fila, y total de filas (fin - ini)
 * con las que va a trabajar el proceso
 */
void get_dist(int total_r, int id, int numP, struct Dist_data *dist_data) {
  int rem;

  dist_data->tot_r = total_r;
  dist_data->tamBl = total_r / numP;
  rem = total_r % numP;

  if(id < rem) { // First subgroup
    dist_data->ini = id * dist_data->tamBl + id;
    dist_data->fin = (id+1) * dist_data->tamBl + (id+1);
  } else { // Second subgroup
    dist_data->ini = id * dist_data->tamBl + rem;
    dist_data->fin = (id+1) * dist_data->tamBl + rem;
  }
  
  if(dist_data->fin > total_r) {
    dist_data->fin = total_r;
  }
  if(dist_data->ini > dist_data->fin) {
    dist_data->ini = dist_data->fin;
  }

  dist_data->tamBl = dist_data->fin - dist_data->ini;
}

1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
void prepare_redist_counts(int *counts, int *displs, int numP_other, int offset, struct Dist_data dist_data, int *vptr) {
  int idS[2], i, idS_zero; 
  int last_index, first_index;

  getIds_intercomm(dist_data, numP_other, idS);
  idS[0] += offset;
  idS[1] += offset;
  idS_zero = 0;

  if(!idS[0]) {
    set_counts(0, numP_other, dist_data, offset, counts);
    idS_zero = 1;
  }
  for(i=idS[0] + idS_zero; i<idS[1]; i++) {
    set_counts(i, numP_other, dist_data, offset, counts);
    displs[i] = displs[i-1] + counts[i-1];
  }

  if(!idS[0]) {
    last_index = counts[0];
    first_index = 0;
    counts[0] = vptr[last_index] - vptr[first_index];
  }
  for(i=idS[0] + idS_zero; i<idS[1]; i++) {
    last_index = displs[i] + counts[i];
    first_index = displs[i];
    counts[i] = vptr[last_index] - vptr[first_index];
    displs[i] = displs[i-1] + counts[i-1];
  }
}

1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
void prepare_redist_counts_vlen(int *counts, int *displs, int numP_other, int offset, struct Dist_data dist_data) {
  int idS[2], i, idS_zero; 
  int last_index, first_index;

  getIds_intercomm(dist_data, numP_other, idS);
  idS[0] += offset;
  idS[1] += offset;
  idS_zero = 0;

  if(!idS[0]) {
    set_counts(0, numP_other, dist_data, offset, counts);
    idS_zero = 1;
  }
  for(i=idS[0] + idS_zero; i<idS[1]; i++) {
    set_counts(i, numP_other, dist_data, offset, counts);
    displs[i] = displs[i-1] + counts[i-1];
  }
}

iker_martin's avatar
iker_martin committed
1078
1079
1080
1081
/*
 * Obtiene para un Id de proceso, cuantos elementos va 
 * a enviar/recibir el proceso myId
 */
1082
void set_counts(int id, int numP, struct Dist_data data_dist, int offset, int *sendcounts) {
iker_martin's avatar
iker_martin committed
1083
1084
1085
  struct Dist_data other;
  int biggest_ini, smallest_end, tot_rows;

1086
  get_dist(data_dist.tot_r, id-offset, numP, &other);
iker_martin's avatar
iker_martin committed
1087
1088
1089
1090
1091
1092
1093

  // Si el rango de valores no coincide, se pasa al siguiente proceso
  if(data_dist.ini >= other.fin || data_dist.fin <= other.ini) {
    return;
  }

  // Obtiene el proceso con mayor ini entre los dos procesos
1094
  biggest_ini = (data_dist.ini > other.ini) ? data_dist.ini : other.ini;
iker_martin's avatar
iker_martin committed
1095
  // Obtiene el proceso con menor fin entre los dos procesos
1096
1097
  smallest_end = (data_dist.fin < other.fin) ? data_dist.fin : other.fin;

iker_martin's avatar
iker_martin committed
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
  sendcounts[id] = smallest_end - biggest_ini; // Numero de elementos a enviar/recibir del proceso Id
}

/*
 * Obtiene para un proceso de un grupo a que rango procesos de 
 * otro grupo tiene que enviar o recibir datos.
 *
 * Devuelve el primer identificador y el último (Excluido) con el que
 * comunicarse.
 */
1108
void getIds_intercomm(struct Dist_data dist_data, int numP_other, int *idS) {
iker_martin's avatar
iker_martin committed
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
    int idI, idE;
    int tamOther = dist_data.tot_r / numP_other;
    int remOther = dist_data.tot_r % numP_other;
    int middle = (tamOther + 1) * remOther;

    if(middle > dist_data.ini) { // First subgroup
      idI = dist_data.ini / (tamOther + 1);
    } else { // Second subgroup
      idI = ((dist_data.ini - middle) / tamOther) + remOther;
    }

    if(middle >= dist_data.fin) { // First subgroup
      idE = dist_data.fin / (tamOther + 1);
      idE = (dist_data.fin % (tamOther + 1) > 0 && idE+1 <= numP_other) ? idE+1 : idE;
    } else { // Second subgroup
      idE = ((dist_data.fin - middle) / tamOther) + remOther;
      idE = ((dist_data.fin - middle) % tamOther > 0 && idE+1 <= numP_other) ? idE+1 : idE;
    }

1128
1129
    idS[0] = idI;
    idS[1] = idE;
iker_martin's avatar
iker_martin committed
1130
1131
}

1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
void print_counts2(struct Dist_data data_dist, int *xcounts, int *xdispls, int size, int include_zero, const char* name) {
  int i;

  for(i=0; i < size; i++) {
    if(xcounts[i] != 0 || include_zero) {
      printf("P%d of %d | %scounts[%d]=%d disp=%d\n", data_dist.myId, data_dist.numP, name, i, xcounts[i], xdispls[i]);
    }
  }
}

1142
void print_global_results(double start_time) {
1143
1144
1145
  size_t i;
  double sp_time, sy_time, asy_time, mall_time, global_time;

1146
1147
  MAM_Retrieve_times(&sp_time, &sy_time, &asy_time, &mall_time);
  global_time = MPI_Wtime() - start_time;
1148
1149
1150
1151
1152
1153
1154
1155
  printf("T_spawn: %lf", sp_time);
  printf("\nT_SR: %lf", sy_time);
  printf("\nT_AR: %lf", asy_time);
  printf("\nT_Malleability: %lf", mall_time);
  printf("\nT_total: %lf\n", global_time);
}


iker_martin's avatar
iker_martin committed
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
/*
 
	  double starttime, endtime, total, res;
          MPI_Barrier(MPI_COMM_WORLD);
	  starttime = MPI_Wtime();
	  endtime = MPI_Wtime();
          total = endtime - starttime;
          MPI_Reduce(&total, &res, 1, MPI_DOUBLE, MPI_MAX, ROOT, MPI_COMM_WORLD);
          if(dist_data.myId == ROOT) {printf("Tiempo BCAST PADRE %f\n", total); fflush(stdout);}
 */