1 /* Copyright 2004,2007,2010 ENSEIRB, INRIA & CNRS
2 **
3 ** This file is part of the Scotch software package for static mapping,
4 ** graph partitioning and sparse matrix ordering.
5 **
6 ** This software is governed by the CeCILL-C license under French law
7 ** and abiding by the rules of distribution of free software. You can
8 ** use, modify and/or redistribute the software under the terms of the
9 ** CeCILL-C license as circulated by CEA, CNRS and INRIA at the following
10 ** URL: "http://www.cecill.info".
11 **
12 ** As a counterpart to the access to the source code and rights to copy,
13 ** modify and redistribute granted by the license, users are provided
14 ** only with a limited warranty and the software's author, the holder of
15 ** the economic rights, and the successive licensors have only limited
16 ** liability.
17 **
18 ** In this respect, the user's attention is drawn to the risks associated
19 ** with loading, using, modifying and/or developing or reproducing the
20 ** software by the user in light of its specific status of free software,
21 ** that may mean that it is complicated to manipulate, and that also
22 ** therefore means that it is reserved for developers and experienced
23 ** professionals having in-depth computer knowledge. Users are therefore
24 ** encouraged to load and test the software's suitability as regards
25 ** their requirements in conditions enabling the security of their
26 ** systems and/or data to be ensured and, more generally, to use and
27 ** operate it in the same conditions as regards security.
28 **
29 ** The fact that you are presently reading this means that you have had
30 ** knowledge of the CeCILL-C license and that you accept its terms.
31 */
32 /************************************************************/
33 /**                                                        **/
34 /**   NAME       : dorder_io.c                             **/
35 /**                                                        **/
36 /**   AUTHOR     : Francois PELLEGRINI                     **/
37 /**                                                        **/
38 /**   FUNCTION   : This module handles distributed         **/
39 /**                orderings.                              **/
40 /**                                                        **/
41 /**   DATES      : # Version 5.0  : from : 27 apr 2006     **/
42 /**                                 to     13 jun 2008     **/
43 /**                # Version 5.1  : from : 30 jul 2010     **/
44 /**                                 to     11 aug 2010     **/
45 /**                                                        **/
46 /************************************************************/
47 
48 /*
49 **  The defines and includes.
50 */
51 
52 #define DORDER
53 
54 #include "module.h"
55 #include "common.h"
56 #include "comm.h"
57 #include "dgraph.h"
58 #include "dorder.h"
59 #include "order.h"
60 
61 /************************************/
62 /*                                  */
63 /* These routines handle orderings. */
64 /*                                  */
65 /************************************/
66 
67 /* This routine saves a distributed ordering.
68 ** The distributed graph structure is provided
69 ** to access the distribution of vertex labels,
70 ** whenever present.
71 ** It returns:
72 ** - 0   : on success.
73 ** - !0  : on error.
74 */
75 
76 int
dorderSave(const Dorder * restrict const ordeptr,const Dgraph * restrict const grafptr,FILE * restrict const stream)77 dorderSave (
78 const Dorder * restrict const ordeptr,
79 const Dgraph * restrict const grafptr,
80 FILE * restrict const         stream)
81 {
82   Gnum * restrict       peritab;
83   Gnum * restrict       permtab;
84   Gnum * restrict       vlbltax;
85   int                   procglbnbr;
86   int                   reduloctab[3];
87   int                   reduglbtab[3];
88   int                   protnum;
89 
90   if (stream != NULL) {                           /* If file provided         */
91     reduloctab[0] = 1;                            /* This process is the root */
92     reduloctab[1] = ordeptr->proclocnum;          /* Get its rank             */
93   }
94   else {
95     reduloctab[0] =                               /* This process is not the root */
96     reduloctab[1] = 0;
97   }
98   reduloctab[2] = (grafptr->vlblloctax != NULL) ? 1 : 0; /* See if vertex labels provided */
99   if (MPI_Allreduce (reduloctab, reduglbtab, 3, MPI_INT, MPI_SUM, ordeptr->proccomm) != MPI_SUCCESS) {
100     errorPrint ("dorderSave: communication error (1)");
101     return     (1);
102   }
103   if (reduglbtab[0] != 1) {
104     errorPrint ("dorderSave: should have only one root");
105     return     (1);
106   }
107   MPI_Comm_size (ordeptr->proccomm, &procglbnbr);
108   if ((reduglbtab[2] != 0) && (reduglbtab[2] != procglbnbr)) {
109     errorPrint ("dorderSave: inconsistent parameters");
110     return     (1);
111   }
112   protnum = (int) reduglbtab[1];                  /* Get rank of root process */
113 
114   reduloctab[0] = 0;
115   permtab       = NULL;
116   if (protnum == ordeptr->proclocnum) {
117     Gnum                  vlblnbr;
118 
119     vlblnbr = (grafptr->vlblloctax != NULL) ? ordeptr->vnodglbnbr : 0;
120     if (memAllocGroup ((void **) (void *)
121                        &permtab, (size_t) (ordeptr->vnodglbnbr * sizeof (Gnum)),
122                        &peritab, (size_t) (ordeptr->vnodglbnbr * sizeof (Gnum)),
123                        &vlbltax, (size_t) (vlblnbr             * sizeof (Gnum)), NULL) == NULL) {
124       errorPrint ("dorderSave: out of memory");
125 #ifdef SCOTCH_DEBUG_DORDER1
126       reduloctab[0] = 1;
127 #else /* SCOTCH_DEBUG_DORDER1 */
128       return (1);
129 #endif /* SCOTCH_DEBUG_DORDER1 */
130     }
131   }
132 #ifdef SCOTCH_DEBUG_DORDER1                       /* This communication cannot be covered by a useful one */
133   if (MPI_Bcast (&reduloctab[0], 1, MPI_INT, protnum, ordeptr->proccomm) != MPI_SUCCESS) {
134     errorPrint ("dorderSave: communication error (2)");
135     if (permtab != NULL);
136       memFree (permtab);                          /* Free group leader */
137     return (1);
138   }
139   if (reduloctab[0] != 0)
140     return (1);
141 #endif /* SCOTCH_DEBUG_DORDER1 */
142 
143   if (grafptr->vlblloctax != NULL) {
144     if (commGatherv (grafptr->vlblloctax + grafptr->baseval, grafptr->vertlocnbr, GNUM_MPI,
145                      vlbltax, grafptr->proccnttab, grafptr->procdsptab, GNUM_MPI, protnum, grafptr->proccomm) != MPI_SUCCESS) {
146       errorPrint ("dorderSave: communication error (3)");
147       return     (1);
148     }
149   }
150 
151   if (protnum == ordeptr->proclocnum) {
152     Gnum                  vertnum;
153 
154     for (vertnum = 0; vertnum < ordeptr->vnodglbnbr; ) { /* Till all inverse permutation indices collected */
155       const DorderLink *    linkptr;
156 
157       for (linkptr = ordeptr->linkdat.nextptr; linkptr != &ordeptr->linkdat; linkptr = linkptr->nextptr) {
158         const DorderCblk *    cblkptr;
159 
160         cblkptr = (DorderCblk *) linkptr;         /* TRICK: FIRST */
161         if (((cblkptr->typeval & DORDERCBLKLEAF) != 0) &&
162             (cblkptr->data.leaf.ordelocval == vertnum) && /* If column block fragment starts at proper index       */
163             (cblkptr->data.leaf.vnodlocnbr > 0)) { /* And is not an empty local block with relevent data elsewhere */
164           memCpy (peritab + vertnum, cblkptr->data.leaf.periloctab, cblkptr->data.leaf.vnodlocnbr * sizeof (Gnum));
165           vertnum += cblkptr->data.leaf.vnodlocnbr;
166           break;
167         }
168       }
169       if (linkptr == &ordeptr->linkdat) {         /* If fragment not found locally */
170         MPI_Status            statdat;
171         int                   recvnbr;
172 
173         if (MPI_Bcast (&vertnum, 1, GNUM_MPI, protnum, ordeptr->proccomm) != MPI_SUCCESS) { /* Broadcast missing fragment */
174           errorPrint ("dorderSave: communication error (4)");
175           memFree    (permtab);                       /* Free group leader */
176           return     (1);
177         }
178         if (MPI_Recv (peritab + vertnum, ordeptr->vnodglbnbr - vertnum, GNUM_MPI,
179                       MPI_ANY_SOURCE, DORDERTAGPERI, ordeptr->proccomm, &statdat) != MPI_SUCCESS) {
180           errorPrint ("dorderSave: communication error (5)");
181           return     (1);
182         }
183         MPI_Get_count (&statdat, GNUM_MPI, &recvnbr);
184         vertnum += recvnbr;
185       }
186     }
187     vertnum = -1;                                 /* Indicate termination */
188     if (MPI_Bcast (&vertnum, 1, GNUM_MPI, protnum, ordeptr->proccomm) != MPI_SUCCESS) {
189       errorPrint ("dorderSave: communication error (6)");
190       memFree    (permtab);                       /* Free group leader */
191       return     (1);
192     }
193 
194     if (fprintf (stream, GNUMSTRING "\n", (Gnum) ordeptr->vnodglbnbr) == EOF) {
195       errorPrint ("dorderSave: bad output (1)");
196       memFree    (permtab);
197       return     (1);
198     }
199 
200     orderPeri (peritab, ordeptr->baseval, ordeptr->vnodglbnbr, permtab, ordeptr->baseval); /* Compute direct permutation */
201 
202     if (grafptr->vlblloctax != NULL) {            /* If ordering has label array */
203       vlbltax -= ordeptr->baseval;                /* Base label array            */
204 
205       for (vertnum = 0; vertnum < ordeptr->vnodglbnbr; vertnum ++) {
206         if (fprintf (stream, GNUMSTRING "\t" GNUMSTRING "\n",
207                      (Gnum) vlbltax[vertnum + ordeptr->baseval],
208                      (Gnum) vlbltax[permtab[vertnum]]) == EOF) {
209           errorPrint ("dorderSave: bad output (2)");
210           memFree    (permtab);
211           return     (1);
212         }
213       }
214     }
215     else {
216       for (vertnum = 0; vertnum < ordeptr->vnodglbnbr; vertnum ++) {
217         if (fprintf (stream, GNUMSTRING "\t" GNUMSTRING "\n",
218                      (Gnum) (vertnum + ordeptr->baseval),
219                      (Gnum) permtab[vertnum]) == EOF) {
220           errorPrint ("dorderSave: bad output (3)");
221           memFree    (permtab);
222           return     (1);
223         }
224       }
225     }
226 
227     memFree (permtab);                            /* Free group leader */
228   }
229   else {
230     while (1) {
231       const DorderLink *    linkptr;
232       Gnum                  vertnum;
233 
234       if (MPI_Bcast (&vertnum, 1, GNUM_MPI, protnum, ordeptr->proccomm) != MPI_SUCCESS) {
235         errorPrint ("dorderSave: communication error (7)");
236         return     (1);
237       }
238       if (vertnum == -1)                          /* If asked to quit */
239         break;                                    /* Finish           */
240 
241       for (linkptr = ordeptr->linkdat.nextptr; linkptr != &ordeptr->linkdat; linkptr = linkptr->nextptr) {
242         const DorderCblk *    cblkptr;
243 
244         cblkptr = (DorderCblk *) linkptr;         /* TRICK: FIRST */
245 
246         if (((cblkptr->typeval & DORDERCBLKLEAF) != 0) && /* If matching column block fragment found */
247             (cblkptr->data.leaf.ordelocval == vertnum) &&
248             (cblkptr->data.leaf.vnodlocnbr > 0)) { /* And is not an empty local block with relevent data elsewhere */
249           if (MPI_Send (cblkptr->data.leaf.periloctab, cblkptr->data.leaf.vnodlocnbr,
250                         GNUM_MPI, protnum, DORDERTAGPERI, ordeptr->proccomm) != MPI_SUCCESS) {
251             errorPrint ("dorderSave: communication error (8)");
252             return     (1);
253           }
254           break;
255         }
256       }
257     }
258   }
259 
260   return  (0);
261 }
262