Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
* minimum of datatype handling for alltoall tuned
[simgrid.git] / src / smpi / smpi_coll.c
1 /* $Id$tag */
2
3 /* smpi_coll.c -- various optimized routing for collectives                   */
4
5 /* Copyright (c) 2009 Stephane Genaud.                                        */
6 /* All rights reserved.                                                       */
7
8 /* This program is free software; you can redistribute it and/or modify it
9  *  * under the terms of the license (GNU LGPL) which comes with this package. */
10
11
12 #include <stdio.h>
13 #include <stdlib.h>
14 #include <string.h>
15
16 #include "private.h"
17 #include "smpi_coll_private.h"
18
19 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_coll, smpi,
20                                 "Logging specific to SMPI (coll)");
21
22 /* proc_tree taken and translated from P2P-MPI */
23
24 struct proc_tree {
25         int PROCTREE_A;
26         int numChildren;
27         int * child;
28         int parent;
29         int me;
30         int root;
31         int isRoot;
32 };
33 typedef struct proc_tree * proc_tree_t;
34
35
36
37 /* prototypes */
38 void build_tree( int index, int extent, proc_tree_t *tree);
39 proc_tree_t alloc_tree( int arity );
40 void free_tree( proc_tree_t tree);
41 void print_tree(proc_tree_t tree);
42
43
44
45
46 /**
47  * alloc and init
48  **/
49 proc_tree_t alloc_tree( int arity ) {
50         proc_tree_t tree = malloc(1*sizeof(struct proc_tree));
51         int i;
52
53         tree->PROCTREE_A = arity;
54         tree->isRoot = 0; 
55         tree->numChildren = 0;
56         tree->child = malloc(arity*sizeof(int));
57         for (i=0; i < arity; i++) {
58                 tree->child[i] = -1;
59         }
60         tree->root = -1;
61         tree->parent = -1;
62         return( tree );
63 }
64
65 /**
66  * free
67  **/
68 void free_tree( proc_tree_t tree) {
69         free (tree->child );
70         free(tree);
71 }
72
73
74
75 /**
76  * Build the tree depending on a process rank (index) and the group size (extent)
77  * @param index the rank of the calling process
78  * @param extent the total number of processes
79  **/
80 void build_tree( int index, int extent, proc_tree_t *tree) {
81         int places = (*tree)->PROCTREE_A * index;
82         int i;
83         int ch;
84         int pr;
85
86         (*tree)->me = index;
87         (*tree)->root = 0 ;
88
89         for (i = 1; i <= (*tree)->PROCTREE_A; i++) {
90                 ++places;
91                 ch = (*tree)->PROCTREE_A * index + i + (*tree)->root;
92                 //printf("places %d\n",places);
93                 ch %= extent;
94                 if (places < extent) {
95                          //printf("ch <%d> = <%d>\n",i,ch);
96                          //printf("adding to the tree at index <%d>\n\n",i-1);
97                         (*tree)->child[i - 1] = ch;
98                         (*tree)->numChildren++;
99                 }
100                 else {
101                        //fprintf(stderr,"not adding to the tree\n");
102                 }
103         }
104         //fprintf(stderr,"procTree.numChildren <%d>\n",(*tree)->numChildren);
105
106         if (index == (*tree)->root) {
107                 (*tree)->isRoot = 1;
108         }
109         else {
110                 (*tree)->isRoot = 0;
111                 pr = (index - 1) / (*tree)->PROCTREE_A;
112                 (*tree)->parent = pr;
113         }
114 }
115
116 /**
117  * prints the tree as a graphical representation
118  **/
119 void print_tree(proc_tree_t tree) {
120       int i;
121       char *spacer;
122       if (-1 != tree->parent ) {
123             printf("[%d]\n   +---",tree->parent);
124             spacer= strdup("     ");
125       }
126       else {
127             spacer=strdup("");
128       }
129       printf("<%d>\n",tree->me);
130       for (i=0;i < tree->numChildren; i++) {
131               printf("%s   +--- %d\n", spacer,tree->child[i]);
132       }
133       free(spacer);
134 }
135       
136 /**
137  * bcast
138  **/
139 int tree_bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm, proc_tree_t tree);
140 int tree_bcast( void *buf, int count, MPI_Datatype datatype, int root,
141                 MPI_Comm comm, proc_tree_t tree) 
142 {
143         int system_tag = 999;  // used negative int but smpi_create_request() declares this illegal (to be checked)
144         int rank;
145         int retval = MPI_SUCCESS;
146         int i;
147         smpi_mpi_request_t request;
148         smpi_mpi_request_t * requests;
149
150         rank = smpi_mpi_comm_rank(comm);
151
152         /* wait for data from my parent in the tree */
153         if (!tree->isRoot) {
154                 DEBUG4("[%d] recv(%d  from %d, tag=%d)\n",rank,rank, tree->parent, system_tag+rank);
155                 retval = smpi_create_request(buf, count, datatype, 
156                                 tree->parent, rank, 
157                                 system_tag + rank, 
158                                 comm, &request);
159                 if (MPI_SUCCESS != retval) {
160                         printf("** internal error: smpi_create_request() rank=%d returned retval=%d, %s:%d\n",
161                                         rank,retval,__FILE__,__LINE__);
162                 }
163                 smpi_mpi_irecv(request);
164                 DEBUG2("[%d] waiting on irecv from %d\n",rank , tree->parent);
165                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
166                 xbt_mallocator_release(smpi_global->request_mallocator, request);
167         }
168
169         requests = xbt_malloc( tree->numChildren * sizeof(smpi_mpi_request_t));
170         DEBUG2("[%d] creates %d requests\n",rank,tree->numChildren);
171
172         /* iniates sends to ranks lower in the tree */
173         for (i=0; i < tree->numChildren; i++) {
174                 if (tree->child[i] != -1) {
175                         DEBUG4("[%d] send(%d->%d, tag=%d)\n",rank,rank, tree->child[i], system_tag+tree->child[i]);
176                         retval = smpi_create_request(buf, count, datatype, 
177                                         rank, tree->child[i], 
178                                         system_tag + tree->child[i], 
179                                         comm, &(requests[i]));
180                         DEBUG5("[%d] after create req[%d]=%p req->(src=%d,dst=%d)\n",rank,i,requests[i],requests[i]->src,requests[i]->dst );
181                         if (MPI_SUCCESS != retval) {
182                               printf("** internal error: smpi_create_request() rank=%d returned retval=%d, %s:%d\n",
183                                               rank,retval,__FILE__,__LINE__);
184                         }
185                         smpi_mpi_isend(requests[i]);
186                         /* FIXME : we should not wait immediately here. See next FIXME. */
187                         smpi_mpi_wait( requests[i], MPI_STATUS_IGNORE);
188                         xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
189                 }
190         }
191         /* FIXME : normally, we sould wait only once all isend have been issued:
192          * this is the following commented code. It deadlocks, probably because
193          * of a bug in the sender process */
194         
195         /* wait for completion of sends */
196         /* printf("[%d] wait for %d send completions\n",rank,tree->numChildren);
197           smpi_mpi_waitall( tree->numChildren, requests, MPI_STATUS_IGNORE);
198          printf("[%d] reqs completed\n)",rank);
199            */
200         
201         xbt_free(requests);
202         return(retval);
203         /* checked ok with valgrind --leak-check=full*/
204 }
205
206
207 /**
208  * anti-bcast
209  **/
210 int tree_antibcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm, proc_tree_t tree);
211 int tree_antibcast( void *buf, int count, MPI_Datatype datatype, int root,
212                 MPI_Comm comm, proc_tree_t tree) 
213 {
214         int system_tag = 999;  // used negative int but smpi_create_request() declares this illegal (to be checked)
215         int rank;
216         int retval = MPI_SUCCESS;
217         int i;
218         smpi_mpi_request_t request;
219         smpi_mpi_request_t * requests;
220
221         rank = smpi_mpi_comm_rank(comm);
222
223           //------------------anti-bcast-------------------
224         
225         // everyone sends to its parent, except root.
226         if (!tree->isRoot) {
227                 retval = smpi_create_request(buf, count, datatype,
228                                 rank,tree->parent,  
229                                 system_tag + rank, 
230                                 comm, &request);
231                 if (MPI_SUCCESS != retval) {
232                         printf("** internal error: smpi_create_request() rank=%d returned retval=%d, %s:%d\n",
233                                         rank,retval,__FILE__,__LINE__);
234                 }
235                 smpi_mpi_isend(request);
236                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
237                 xbt_mallocator_release(smpi_global->request_mallocator, request);
238         }
239
240         //every one receives as many messages as it has children
241         requests = xbt_malloc( tree->numChildren * sizeof(smpi_mpi_request_t));
242         for (i=0; i < tree->numChildren; i++) {
243                 if (tree->child[i] != -1) {
244                         retval = smpi_create_request(buf, count, datatype, 
245                                         tree->child[i], rank, 
246                                         system_tag + tree->child[i], 
247                                         comm, &(requests[i]));
248                         if (MPI_SUCCESS != retval) {
249                                 printf("** internal error: smpi_create_request() rank=%d returned retval=%d, %s:%d\n",
250                                                 rank,retval,__FILE__,__LINE__);
251                         }
252                         smpi_mpi_irecv(requests[i]);
253                         /* FIXME : we should not wait immediately here. See next FIXME. */
254                         smpi_mpi_wait( requests[i], MPI_STATUS_IGNORE);
255                         xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
256                 }
257         }
258         xbt_free(requests);
259         return(retval);
260
261         /* checked ok with valgrind --leak-check=full*/
262
263
264 /**
265  * bcast with a binary, ternary, or whatever tree .. 
266  **/
267 int nary_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root,
268                 MPI_Comm comm, int arity)
269 {
270 int rank;
271 int retval;
272
273         rank = smpi_mpi_comm_rank( comm );
274         // arity=2: a binary tree, arity=4 seem to be a good setting (see P2P-MPI))
275         proc_tree_t tree = alloc_tree( arity ); 
276         build_tree( rank, comm->size, &tree );
277
278         retval = tree_bcast( buf, count, datatype, root, comm, tree );
279
280         free_tree( tree );
281         return( retval );
282 }
283
284
285 /**
286  * Barrier
287  **/
288 int nary_tree_barrier( MPI_Comm comm , int arity)
289 {
290         int rank;
291         int retval = MPI_SUCCESS;
292         char dummy='$';
293
294         rank = smpi_mpi_comm_rank(comm);
295         // arity=2: a binary tree, arity=4 seem to be a good setting (see P2P-MPI))
296         proc_tree_t tree = alloc_tree( arity ); 
297
298         build_tree( rank, comm->size, &tree );
299
300         retval = tree_antibcast( &dummy, 1, MPI_CHAR, 0, comm, tree);
301         if (MPI_SUCCESS != retval) {
302                 printf("[%s:%d] ** Error: tree_antibcast() returned retval=%d\n",__FILE__,__LINE__,retval);
303         }
304         else {
305             retval = tree_bcast(&dummy,  1, MPI_CHAR, 0, comm, tree);
306         }
307
308         free_tree( tree );
309         return(retval);
310
311         /* checked ok with valgrind --leak-check=full*/
312 }
313
314
315 /**
316  * Alltoall pairwise
317  *
318  * this algorithm performs size steps (1<=s<=size) and
319  * at each step s, a process p sends iand receive to.from a unique distinct remote process
320  * size=5 : s=1:  4->0->1, 0->1->2, 1->2->3, ... 
321  *          s=2:  3->0->2, 4->1->3, 0->2->4, 1->3->0 , 2->4->1
322  *          .... 
323  * Openmpi calls this routine when the message size sent to each rank is greater than 3000 bytes
324  **/
325
326 int smpi_coll_tuned_alltoall_pairwise (void *sendbuf, int sendcount, MPI_Datatype datatype,
327                     void* recvbuf, int recvcount, MPI_Datatype recvdatatype, MPI_Comm comm)
328 {
329         int retval = MPI_SUCCESS;
330           int rank;
331           int size = comm->size;
332           int step;
333           int sendto, recvfrom;
334           int tag_alltoall=999;
335           void * tmpsend, *tmprecv;
336
337           rank = smpi_mpi_comm_rank(comm);
338         INFO1("[%d] algorithm alltoall_pairwise() called.\n",rank);
339
340
341           /* Perform pairwise exchange - starting from 1 so the local copy is last */
342           for (step = 1; step < size+1; step++) {
343
344                     /* who do we talk to in this step? */
345                     sendto  = (rank+step)%size;
346                     recvfrom = (rank+size-step)%size;
347
348                     /* where from are we sending and where from are we receiving actual data ? */
349                     tmpsend = (char*)sendbuf+sendto*datatype->size*sendcount;
350                     tmprecv = (char*)recvbuf+recvfrom*recvdatatype->size*recvcount;
351
352                     /* send and receive */
353                     /* in OpenMPI, they use :
354                          err = ompi_coll_tuned_sendrecv( tmpsend, scount, sdtype, sendto, MCA_COLL_BASE_TAG_ALLTOALL,
355                          tmprecv, rcount, rdtype, recvfrom, MCA_COLL_BASE_TAG_ALLTOALL,
356                          comm, MPI_STATUS_IGNORE, rank);
357                      */
358                     retval = MPI_Sendrecv(tmpsend, sendcount, datatype, sendto, tag_alltoall,
359                                                 tmprecv, recvcount, recvdatatype, recvfrom, tag_alltoall,
360                                                   comm, MPI_STATUS_IGNORE);
361           }
362           return(retval);
363 }
364
365 /**
366  * helper: copy a datatype into another (in the simple case dt1=dt2) 
367 **/
368 int copy_dt( void *sbuf, int scount, const MPI_Datatype sdtype, void *rbuf, int rcount, const MPI_Datatype rdtype);
369 int copy_dt( void *sbuf, int scount, const MPI_Datatype sdtype,
370              void *rbuf, int rcount, const MPI_Datatype rdtype)
371 {
372     /* First check if we really have something to do */
373     if (0 == rcount) {
374         return ((0 == scount) ? MPI_SUCCESS : MPI_ERR_TRUNCATE);
375     }
376    /* If same datatypes used, just copy. */
377    if (sdtype == rdtype) {
378       int count = ( scount < rcount ? scount : rcount );
379       memcpy( rbuf, sbuf, sdtype->size*count);
380       return ((scount > rcount) ? MPI_ERR_TRUNCATE : MPI_SUCCESS);
381    }
382    /* FIXME:  cases 
383     * - If receive packed. 
384     * - If send packed
385     * to be treated once we have the MPI_Pack things ...
386     **/
387      return( MPI_SUCCESS );
388 }
389
390 /**
391  *
392  **/
393 int smpi_coll_tuned_alltoall_basic_linear(void *sbuf, int scount, MPI_Datatype sdtype,
394                     void* rbuf, int rcount, MPI_Datatype rdtype, MPI_Comm comm)
395 {
396           int i;
397         int system_alltoall_tag = 888;
398           int rank;
399           int size = comm->size;
400           int err;
401           char *psnd;
402           char *prcv;
403         int nreq = 0;
404           MPI_Aint lb;
405           MPI_Aint sndinc;
406           MPI_Aint rcvinc;
407           MPI_Request *reqs;
408
409           /* Initialize. */
410           rank = smpi_mpi_comm_rank(comm);
411         INFO1("[%d] algorithm alltoall_basic_linear() called.\n",rank);
412
413         err = smpi_mpi_type_get_extent(sdtype, &lb, &sndinc);
414         err = smpi_mpi_type_get_extent(rdtype, &lb, &rcvinc);
415           /* simple optimization */
416           psnd = ((char *) sbuf) + (rank * sndinc);
417           prcv = ((char *) rbuf) + (rank * rcvinc);
418
419           err = copy_dt( psnd, scount, sdtype, prcv, rcount, rdtype );
420           if (MPI_SUCCESS != err) {
421                     return err;
422           }
423
424           /* If only one process, we're done. */
425           if (1 == size) {
426                     return MPI_SUCCESS;
427           }
428
429           /* Initiate all send/recv to/from others. */
430           reqs =  xbt_malloc(2*(size-1) * sizeof(smpi_mpi_request_t));
431
432           prcv = (char *) rbuf;
433           psnd = (char *) sbuf;
434
435           /* Post all receives first -- a simple optimization */
436           for (i = (rank + 1) % size; i != rank; i = (i + 1) % size) {
437                     err = smpi_create_request( prcv + (i * rcvinc), rcount, rdtype,
438                                         i, rank,
439                                         system_alltoall_tag,
440                                         comm, &(reqs[nreq]));
441                 if (MPI_SUCCESS != err) {
442                         DEBUG2("[%d] failed to create request for rank %d\n",rank,i);
443                         for (i=0;i< nreq;i++) 
444                                 xbt_mallocator_release(smpi_global->request_mallocator, reqs[i]);
445                         return err;
446                 }
447                 nreq++;
448           }
449           /* Now post all sends in reverse order 
450            *        - We would like to minimize the search time through message queue
451            *                 when messages actually arrive in the order in which they were posted.
452            *                      */
453           for (i = (rank + size - 1) % size; i != rank; i = (i + size - 1) % size ) {
454                     err = smpi_create_request (psnd + (i * sndinc), scount, sdtype, 
455                                          rank, i,
456                                                      system_alltoall_tag, 
457                                          comm, &(reqs[nreq]));
458                 if (MPI_SUCCESS != err) {
459                         DEBUG2("[%d] failed to create request for rank %d\n",rank,i);
460                         for (i=0;i< nreq;i++) 
461                                 xbt_mallocator_release(smpi_global->request_mallocator, reqs[i]);
462                         return err;
463                 }
464                 nreq++;
465         }
466
467         /* Start your engines.  This will never return an error. */
468         for ( i=0; i< nreq/2; i++ ) {
469             smpi_mpi_irecv( reqs[i] );
470         }
471         for ( i= nreq/2; i<nreq; i++ ) {
472             smpi_mpi_isend( reqs[i] );
473         }
474
475
476         /* Wait for them all.  If there's an error, note that we don't
477          * care what the error was -- just that there *was* an error.  The
478          * PML will finish all requests, even if one or more of them fail.
479          * i.e., by the end of this call, all the requests are free-able.
480          * So free them anyway -- even if there was an error, and return
481          * the error after we free everything. */
482
483           err = smpi_mpi_waitall(nreq, reqs, MPI_STATUS_IGNORE);
484
485           /* Free the reqs */
486         for (i=0;i< 2*(size-1);i++) {
487             xbt_mallocator_release(smpi_global->request_mallocator, reqs[i]);
488         }
489           return err;
490 }
491
492
493 /**
494  * Alltoall Bruck 
495  *
496  * Openmpi calls this routine when the message size sent to each rank < 2000 bytes and size < 12
497  **/
498
499
500 int smpi_coll_tuned_alltoall_bruck(void *sendbuf, int sendcount, MPI_Datatype sdtype,
501                     void* recvbuf, int recvcount, MPI_Datatype rdtype, MPI_Comm comm)
502 {
503 /*        int size = comm->size;
504           int i, k, line = -1;
505           int sendto, recvfrom, distance, *displs=NULL, *blen=NULL;
506           int maxpacksize, packsize, position;
507           char * tmpbuf=NULL, *packbuf=NULL;
508           ptrdiff_t lb, sext, rext;
509           int err = 0;
510           int weallocated = 0;
511           MPI_Datatype iddt;
512
513           rank = smpi_mpi_comm_rank(comm);
514 */
515           INFO0("coll:tuned:alltoall_intra_bruck ** NOT IMPLEMENTED YET**");
516 /*
517           displs = xbt_malloc(size*sizeof(int));
518           blen =   xbt_malloc(size*sizeof(int));
519
520           weallocated = 1;
521 */
522           /* Prepare for packing data */
523 /*
524           err = MPI_Pack_size( scount*size, sdtype, comm, &maxpacksize );
525           if (err != MPI_SUCCESS) {  }
526 */
527           /* pack buffer allocation */
528 /*        packbuf = (char*) malloc((unsigned) maxpacksize);
529           if (packbuf == NULL) { }
530 */
531           /* tmp buffer allocation for message data */
532 /*        tmpbuf = xbt_malloc(scount*size*sext);
533           if (tmpbuf == NULL) {  }
534 */
535
536           /* Step 1 - local rotation - shift up by rank */
537 /*        err = ompi_ddt_copy_content_same_ddt (sdtype, (int32_t) ((size-rank)*scount),
538                                 tmpbuf, ((char*)sbuf)+rank*scount*sext);
539           if (err<0) {
540                     line = __LINE__; err = -1; goto err_hndl;
541           }
542
543           if (rank != 0) {
544                     err = ompi_ddt_copy_content_same_ddt (sdtype, (int32_t) (rank*scount),
545                                           tmpbuf+(size-rank)*scount*sext, (char*)sbuf);
546                     if (err<0) {
547                                 line = __LINE__; err = -1; goto err_hndl;
548                     }
549           }
550 */
551           /* perform communication step */
552 /*        for (distance = 1; distance < size; distance<<=1) {
553 */
554                     /* send data to "sendto" */
555 /*                  sendto = (rank+distance)%size;
556                     recvfrom = (rank-distance+size)%size;
557                     packsize = 0;
558                     k = 0;
559 */
560                     /* create indexed datatype */
561 //                  for (i = 1; i < size; i++) {
562 //                              if ((i&distance) == distance) {
563 //                                        displs[k] = i*scount; blen[k] = scount;
564 //                                        k++;
565 //                              }
566 //                  }
567                     /* Set indexes and displacements */
568 //                  err = MPI_Type_indexed(k, blen, displs, sdtype, &iddt);
569 //                  if (err != MPI_SUCCESS) { line = __LINE__; goto err_hndl;  }
570 //                  /* Commit the new datatype */
571 ///                 err = MPI_Type_commit(&iddt);
572 //                  if (err != MPI_SUCCESS) { line = __LINE__; goto err_hndl;  }
573
574                     /* have the new distribution ddt, pack and exchange data */
575 //                  err = MPI_Pack(tmpbuf, 1, iddt, packbuf, maxpacksize, &packsize, comm);
576 //                  if (err != MPI_SUCCESS) { line = __LINE__; goto err_hndl;  }
577
578                     /* Sendreceive */
579 //                  err = ompi_coll_tuned_sendrecv ( packbuf, packsize, MPI_PACKED, sendto,
580 //                                        MCA_COLL_BASE_TAG_ALLTOALL,
581 //                                        rbuf, packsize, MPI_PACKED, recvfrom,
582 //                                        MCA_COLL_BASE_TAG_ALLTOALL,
583 //                                        comm, MPI_STATUS_IGNORE, rank);
584 //                  if (err != MPI_SUCCESS) { line = __LINE__; goto err_hndl; }
585
586                     /* Unpack data from rbuf to tmpbuf */
587 //                  position = 0;
588 //          err = MPI_Unpack(rbuf, packsize, &position,
589 //                                        tmpbuf, 1, iddt, comm);
590 //                  if (err != MPI_SUCCESS) { line = __LINE__; goto err_hndl; }
591
592                     /* free ddt */
593 //                  err = MPI_Type_free(&iddt);
594 //                  if (err != MPI_SUCCESS) { line = __LINE__; goto err_hndl;  }
595 //        } /* end of for (distance = 1... */
596
597           /* Step 3 - local rotation - */
598 //        for (i = 0; i < size; i++) {
599
600 //                  err = ompi_ddt_copy_content_same_ddt (rdtype, (int32_t) rcount,
601 //                                        ((char*)rbuf)+(((rank-i+size)%size)*rcount*rext),
602 //                                        tmpbuf+i*rcount*rext);
603 //
604 //        if (err<0) {
605 //                              line = __LINE__; err = -1; goto err_hndl;
606 //                  }
607 //        }
608
609           /* Step 4 - clean up */
610 /*        if (tmpbuf != NULL) free(tmpbuf);
611           if (packbuf != NULL) free(packbuf);
612           if (weallocated) {
613                     if (displs != NULL) free(displs);
614                     if (blen != NULL) free(blen);
615           }
616           return OMPI_SUCCESS;
617
618 err_hndl:
619           OPAL_OUTPUT((ompi_coll_tuned_stream,"%s:%4d\tError occurred %d, rank %2d", __FILE__,line,err,rank));
620           if (tmpbuf != NULL) free(tmpbuf);
621           if (packbuf != NULL) free(packbuf);
622           if (weallocated) {
623                     if (displs != NULL) free(displs);
624                     if (blen != NULL) free(blen);
625           }
626           return err;
627           */
628           return -1; /* FIXME: to be changed*/
629 }
630
631
632 /**
633  * -----------------------------------------------------------------------------------------------------
634  * example usage
635  * -----------------------------------------------------------------------------------------------------
636  **/
637 /*
638  * int main() {
639
640       int rank; 
641       int size=12;
642
643       proc_tree_t tree;
644       for (rank=0;rank<size;rank++) {
645             printf("--------------tree for rank %d ----------\n",rank);
646             tree = alloc_tree( 2 );
647             build_tree( rank, size, &tree );
648             print_tree( tree );
649             free_tree( tree );
650    
651       }
652       printf("-------------- bcast ----------\n");
653       for (rank=0;rank<size;rank++) {
654               bcast( rank, size );
655       }
656
657
658 }
659 */
660
661                 
662