ConjugateGradient.c 31 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
11
12
#include "ScalarVectors.h"
#include "SparseMatrices.h"
#include <mpi.h>
#include <string.h>
#include "../malleability/malleabilityManager.h"

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

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

50
51
void dumb(Compute_data *computeData, struct Dist_data *dist_data); //FIXME Delte me

iker_martin's avatar
iker_martin committed
52
53
54
55
56
57
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);
58
int compute(Compute_data *computeData, struct Dist_data *dist_data, int sm);
59
void free_computeData(Compute_data *computeData, int terminate);
iker_martin's avatar
iker_martin committed
60
61
62
63
64

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

int n_check = 30;

65
int dist_old(struct Dist_data *dist_data, Compute_data *computeData, int num_children, int sm, int ss, int rm, int rs, int send_sync);
iker_martin's avatar
iker_martin committed
66
void dist_new(struct Dist_data *dist_data, Compute_data *computeData);
67
void update_dist_data(struct Dist_data *dist_data);
68
void print_global_results();
iker_martin's avatar
iker_martin committed
69
70
71
72
73
74
75
76
//----------------------------------------------------------------------------------------------------
void get_dist(int total_r, int id, int numP, struct Dist_data *dist_data);
void set_counts(int id, int numP, struct Dist_data data_dist, int *sendcounts);
void getIds_intercomm(struct Dist_data dist_data, int numP_other, int **idS);
//----------------------------------------------------------------------------------------------------

int main (int argc, char *argv[]) {
	int terminate;
77
	int req, num_nodes, num_cpus = 20;
78
	int sm, ss, rm, rs, send_sync;
iker_martin's avatar
iker_martin committed
79
80
81
        char *nodelist = NULL;
        Compute_data computeData;

82
        computeData.z = NULL; computeData.d_full = NULL, computeData.d = NULL;
iker_martin's avatar
iker_martin committed
83
84
85
        computeData.vec = NULL; computeData.res = NULL;
        computeData.dist_elem = NULL; computeData.displs_elem = NULL;
        computeData.dist_rows = NULL; computeData.displs_rows = NULL;
86
87
	computeData.subm.vptr = NULL;
	computeData.vlen = NULL;
iker_martin's avatar
iker_martin committed
88

89
	send_sync=1;
90
91
	sm = 1; ss = 1; rm = 0; rs = 1;

iker_martin's avatar
iker_martin committed
92
93
        int numP, myId, num_children = 0;
        struct Dist_data dist_data;
94
        if (argc >= 10) {
iker_martin's avatar
iker_martin committed
95
          num_children = atoi(argv[2]);
96
97
98
99
	  sm = atoi(argv[3]);
	  ss = atoi(argv[4]);
	  rm = atoi(argv[5]);
	  rs = atoi(argv[6]);
100
101
102
	  send_sync = atoi(argv[7]);
          nodelist = argv[8];
          num_nodes = atoi(argv[9]);
iker_martin's avatar
iker_martin committed
103
          num_cpus = num_nodes * num_cpus;
104
	}
iker_martin's avatar
iker_martin committed
105

106
        MPI_Init_thread(&argc, &argv, MPI_THREAD_MULTIPLE, &req);
iker_martin's avatar
iker_martin committed
107
108
109
110
111
        MPI_Comm_size(MPI_COMM_WORLD, &numP);
        MPI_Comm_rank(MPI_COMM_WORLD, &myId);
	dist_data.comm = MPI_COMM_WORLD;

        int new_group = init_malleability(myId, numP, ROOT, dist_data.comm, argv[0], nodelist, num_cpus, num_nodes);
112
        update_dist_data(&dist_data);
iker_martin's avatar
iker_martin committed
113
114
115

	if( !new_group ) { //First set of processes
	  init_app(&computeData, &dist_data, argv);
116
          dist_old(&dist_data, &computeData, num_children, sm, ss, rm, rs, send_sync);
117
118
	  MPI_Barrier(MPI_COMM_WORLD);
          set_global_time(MPI_Wtime());
iker_martin's avatar
iker_martin committed
119
120
121
122
        } else {
          dist_new(&dist_data, &computeData);
	}

123
//        if(computeData.iter==0)
124
	terminate = compute(&computeData, &dist_data, sm);
125

126
	if(terminate) {
iker_martin's avatar
iker_martin committed
127
          update_dist_data(&dist_data);
128
129
130
131
132
          MPI_Barrier(dist_data.comm);
          if(dist_data.myId == ROOT) {
  	    print_global_results();
	    printf ("End(%d) --> (%d,%20.10e)\n", computeData.n, computeData.iter, computeData.tol);
	  }
iker_martin's avatar
iker_martin committed
133
134
135
136
        }

	// End of CG
        free_malleability();
137
        free_computeData(&computeData, 1);
iker_martin's avatar
iker_martin committed
138
139
140
	if(sm && numP > num_children && dist_data.myId == 0) {
          MPI_Abort(MPI_COMM_WORLD, -100);
	}
iker_martin's avatar
iker_martin committed
141
142
143
144
145
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
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
        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);

286
287
288
        CreateInts(&(computeData->vlen), dist_data.tamBl);
        for(i=0; i<dist_data.tamBl; i++) {
          computeData->vlen[i] = computeData->subm.vptr[i+1];
289
        }
iker_martin's avatar
iker_martin committed
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
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
335
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
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
        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
/*
#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
*/
}

/*
 * 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
390
391
392
	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
393
        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
394
395
	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
396
397
398
399
400
401
402
        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
 */
403
int compute(Compute_data *computeData, struct Dist_data *dist_data, int sm) {
iker_martin's avatar
iker_martin committed
404
405
	int IZERO = 0, IONE = 1; 
	double DONE = 1.0, DMONE = -1.0, DZERO = 0.0;
406
	int state = MALL_NOT_STARTED;
iker_martin's avatar
iker_martin committed
407
        int ended_loop = 1;
408
        int cnt = 0;
409
	int reconfigure = 0, rec_iter = 10;
iker_martin's avatar
iker_martin committed
410

411
        computeData->maxiter = 110;
iker_martin's avatar
iker_martin committed
412
413
414

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

415
416
	        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
417
418
419
420
421
422
//      	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
423
        	computeData->rho = rdot (&(dist_data->tamBl), computeData->d, &IONE, computeData->z, &IONE);		// rho = (d * z)
424
	        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
425
		computeData->rho = computeData->beta / computeData->rho;                 		                // rho = beta / aux
iker_martin's avatar
iker_martin committed
426
		raxpy (&(dist_data->tamBl), &computeData->rho, computeData->d, &IONE, computeData->vec, &IONE);		// x += rho * d
iker_martin's avatar
iker_martin committed
427
		computeData->rho = -computeData->rho;
iker_martin's avatar
iker_martin committed
428
		raxpy (&(dist_data->tamBl), &computeData->rho, computeData->z, &IONE, computeData->res, &IONE);         // res -= rho * z
iker_martin's avatar
iker_martin committed
429
		computeData->alpha = computeData->beta;                                               		        // alpha = beta
iker_martin's avatar
iker_martin committed
430
		computeData->beta = rdot (&(dist_data->tamBl), computeData->res, &IONE, computeData->res, &IONE);       // beta = res' * res
431
	        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
432
		computeData->alpha = computeData->beta / computeData->alpha;                                       	// alpha = beta / alpha
iker_martin's avatar
iker_martin committed
433
434
		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
435

436
437
	        //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
438
439
440
		computeData->tol = sqrt (computeData->beta);                                          			// tol = sqrt(beta) = norm (res)
		computeData->iter++;

441
                if (computeData->iter == rec_iter) reconfigure = 1;
442
		if (reconfigure) {
443
		  state = malleability_checkpoint();
444
		  if ((state == MALL_COMPLETED && sm == 0) || state == MALL_ZOMBIE) { ended_loop = 0; break; }
445
		  else if(state == MALL_COMPLETED) {
446
		    reconfigure = 0;
447
448
449
450
451
452
453
                    free_computeData(computeData, 0);
                    update_dist_data(dist_data);
		    dist_new(dist_data, computeData);
		  }

		}

iker_martin's avatar
iker_martin committed
454
455
456
457
458
459
460
	}
#ifdef DEBUG
	if(dist_data->myId == ROOT) printf ("Ended loop\n");
#endif
	return ended_loop;
}

461
462
void dumb(Compute_data *computeData, struct Dist_data *dist_data) {
  int i;
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
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539

  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]);
      fflush(stdout);
    }
    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);
      fflush(stdout);
    }
    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]);
      fflush(stdout);
    }
    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]);
      fflush(stdout);
    }
    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]);
      fflush(stdout);
    }
    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]);
      fflush(stdout);
    }
    MPI_Barrier(dist_data->comm);
540
  }
541
  if(dist_data->myId == 0) printf("\n"); 
542
543
  fflush(stdout); MPI_Barrier(dist_data->comm);
}
iker_martin's avatar
iker_martin committed
544

545
void free_computeData(Compute_data *computeData, int terminate) {
546
	if(computeData->res != NULL) {
iker_martin's avatar
iker_martin committed
547
	RemoveDoubles (&computeData->res); 
548
549
	}
	if(computeData->z != NULL) {
iker_martin's avatar
iker_martin committed
550
        RemoveDoubles (&computeData->z); 
551
552
	}
	if(computeData->d != NULL) {
iker_martin's avatar
iker_martin committed
553
        RemoveDoubles (&computeData->d);
554
555
	}
	if(computeData->vec != NULL) {
iker_martin's avatar
iker_martin committed
556
	RemoveDoubles (&computeData->vec);
557
	}
iker_martin's avatar
iker_martin committed
558
559


560
	if(computeData->d_full != NULL && terminate) {
561
562
563
        RemoveDoubles (&computeData->d_full);
	}
	if(computeData->subm.vptr != NULL) {
iker_martin's avatar
iker_martin committed
564
	RemoveSparseMatrix2 (&computeData->subm);
565
566
567
	}

	if(computeData->dist_rows != NULL) {
iker_martin's avatar
iker_martin committed
568
        RemoveInts (&computeData->dist_rows);
569
570
	}
	if(computeData->displs_rows != NULL) {
iker_martin's avatar
iker_martin committed
571
        RemoveInts (&computeData->displs_rows);
572
573
	}
	if(computeData->vlen != NULL) {
iker_martin's avatar
iker_martin committed
574
        RemoveInts (&computeData->vlen);
575
	}
iker_martin's avatar
iker_martin committed
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
}

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

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

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

/*
 */
603
int dist_old(struct Dist_data *dist_data, Compute_data *computeData, int num_children, int sm, int ss, int rm, int rs, int send_sync) {
604
    int phy_dist = 2;
605

606
607
608
    set_malleability_configuration(sm, ss, phy_dist, rm, rs);
    set_children_number(num_children);
    
609
610
611
612
613
    malleability_add_data(&(computeData->n), 1, MAL_INT, MAL_DATA_ALONE, 1, 1);
    malleability_add_data(&(computeData->iter), 1, MAL_INT, MAL_DATA_ALONE, 1, 1);
    malleability_add_data(&(computeData->tol), 1, MAL_DOUBLE, MAL_DATA_ALONE, 1, 1);
    malleability_add_data(&(computeData->beta), 1, MAL_DOUBLE, MAL_DATA_ALONE, 1, 1);
    malleability_add_data(&(computeData->umbral), 1, MAL_DOUBLE, MAL_DATA_ALONE, 1, 1);
iker_martin's avatar
iker_martin committed
614

615
616
    //malleability_add_data(computeData->d_full, computeData->n, MAL_DOUBLE, MAL_DATA_ALONE, 1, 1);
    malleability_add_data(computeData->d, computeData->n, MAL_DOUBLE, MAL_DATA_ALONE, 0, 1);
iker_martin's avatar
iker_martin committed
617

618
619
620
    malleability_add_data(computeData->vec, computeData->n, MAL_DOUBLE, MAL_DATA_ALONE, 0, 1);
    malleability_add_data(computeData->res, computeData->n, MAL_DOUBLE, MAL_DATA_ALONE, 0, 1);
    malleability_add_data(computeData->z, computeData->n, MAL_DOUBLE, MAL_DATA_ALONE, 0, 1);
iker_martin's avatar
iker_martin committed
621

622
    //FIXME SIguientes valores pueden ser asincronos
623
624
625
    malleability_add_data(computeData->vlen, computeData->n, MAL_INT, 1+MAL_DATA_INDEPENDENT, 0, send_sync);
    malleability_add_data(computeData->subm.vpos, computeData->n, MAL_INT, 1+MAL_DATA_DEPENDENT, 0, send_sync);
    malleability_add_data(computeData->subm.vval, computeData->n, MAL_DOUBLE, 1+MAL_DATA_DEPENDENT, 0, send_sync);
iker_martin's avatar
iker_martin committed
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
}

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

/*
 * Función llamada por un set de procesos hijos.
 *
 * Primero los hijos obtienen de los padres una información iniciar 
 * con la que conocer el tamaño de sus vectores y matriz, como asi 
 * tambien cuantos datos van a recibir de cada padre.
 *
 * Tras esto se preparan para recibir los datos de los padres.
 *
 */
void dist_new(struct Dist_data *dist_data, Compute_data *computeData) {
647
648
    int IONE = 1, i, is_synch;
    size_t entry, entries;
iker_martin's avatar
iker_martin committed
649
    void *value = NULL;
650
    is_synch = 1;
651
    entry = 0;
652

iker_martin's avatar
iker_martin committed
653
    malleability_get_data(&value, 0, 1, 1);
654
655
    computeData->n = *((int *)value);
    malleability_get_data(&value, 1, 1, 1);
iker_martin's avatar
iker_martin committed
656
657
    computeData->iter = *((int *)value);
    malleability_get_data(&value, 2, 1, 1);
658
    computeData->tol = *((double *)value);
iker_martin's avatar
iker_martin committed
659
    malleability_get_data(&value, 3, 1, 1);
660
661
    computeData->beta = *((double *)value);
    malleability_get_data(&value, 4, 1, 1);
iker_martin's avatar
iker_martin committed
662
    computeData->umbral = *((double *)value);
663

664
665
666
667
    //malleability_get_data(&value, 5, 1, 1);
    //computeData->d_full = ((double *)value);
    malleability_get_data(&value, entry++, 0, 1);
    computeData->d = ((double *)value);
668

669
    malleability_get_data(&value, entry++, 0, 1);
iker_martin's avatar
iker_martin committed
670
    computeData->vec = ((double *)value);
671
    malleability_get_data(&value, entry++, 0, 1);
iker_martin's avatar
iker_martin committed
672
    computeData->res = ((double *)value);
673
    malleability_get_data(&value, entry++, 0, 1);
iker_martin's avatar
iker_martin committed
674
675
    computeData->z = ((double *)value);

676
    get_dist(computeData->n, dist_data->myId, dist_data->numP, dist_data);
iker_martin's avatar
iker_martin committed
677
    get_rows_dist(computeData, dist_data->numP, computeData->n);
678
679
680
681
    //CreateDoubles(&computeData->d, dist_data->tamBl);
    //rcopy (&(dist_data->tamBl), &(computeData->d_full[dist_data->ini]), &IONE, computeData->d, &IONE);  // d = d_full[ini] to d_full[ini+tamBl]
    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
682

683
    malleability_get_entries(&entries, 0, 0); //Get if there is any asynch data to recover
684
    if(entries) { is_synch=0; entry=0; }
685

686
    malleability_get_data(&value, entry++, 0, is_synch);
687
    computeData->vlen = ((int *)value);
688
    
689
690
691
692
    CreateSparseMatrixVptr(&(computeData->subm), dist_data->tamBl, computeData->n, 0);
    computeData->subm.vptr[0] = 0;
    for(i=0; i<dist_data->tamBl; i++) {
      computeData->subm.vptr[i+1] = computeData->vlen[i];
iker_martin's avatar
iker_martin committed
693
    }
694
    TransformLengthtoHeader(computeData->subm.vptr, computeData->subm.dim1); // The array is converted from vlen to vptr
695

696
    malleability_get_data(&value, entry++, 0, is_synch);
697
    computeData->subm.vpos = ((int *)value);
698
    malleability_get_data(&value, entry++, 0, is_synch);
699
700
    computeData->subm.vval = ((double *)value);
}
iker_martin's avatar
iker_martin committed
701

702
703
704
705
706
707
708
void update_dist_data(struct Dist_data *dist_data) {
  int myId, numP;
  get_malleability_user_comm(&(dist_data->comm));
  MPI_Comm_size(dist_data->comm, &numP);
  MPI_Comm_rank(dist_data->comm, &myId);
  dist_data->myId = myId;
  dist_data->numP = numP;
iker_martin's avatar
iker_martin committed
709
710
711
712
713
714
715
716
717
718
719
720
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
}

/*
 * ========================================================================================
 * ========================================================================================
 * ================================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;
}

/*
 * Obtiene para un Id de proceso, cuantos elementos va 
 * a enviar/recibir el proceso myId
 */
void set_counts(int id, int numP, struct Dist_data data_dist, int *sendcounts) {
  struct Dist_data other;
  int biggest_ini, smallest_end, tot_rows;

  get_dist(data_dist.tot_r, id, numP, &other);

  // 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
  if(data_dist.ini > other.ini) { 
    biggest_ini = data_dist.ini;
  } else {
    biggest_ini = other.ini;
  }

  // Obtiene el proceso con menor fin entre los dos procesos
  if(data_dist.fin < other.fin) {
    smallest_end = data_dist.fin;
  } else {
    smallest_end = other.fin;
  }
  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.
 */
void getIds_intercomm(struct Dist_data dist_data, int numP_other, int **idS) {
    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;
    }

    //free(*idS);
    CreateInts(idS, 2);
    (*idS)[0] = idI;
    (*idS)[1] = idE;
}

814
815
816
817
818
819
820
821
822
823
824
825
826
827
void print_global_results() {
  size_t i;
  double sp_time, sy_time, asy_time, mall_time, global_time;

  retrieve_results(&sp_time, &sy_time, &asy_time, &mall_time, &global_time);
  global_time = MPI_Wtime() - global_time;
  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
828
829
830
831
832
833
834
835
836
837
/*
 
	  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);}
 */