1 /* Copyright 2007-2010,2012 IPB, Universite de Bordeaux, 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       : parmetis_dgraph_order.c                 **/
35 /**                                                        **/
36 /**   AUTHOR     : Francois PELLEGRINI                     **/
37 /**                                                        **/
38 /**   FUNCTION   : This module is the compatibility        **/
39 /**                library for the ParMeTiS ordering       **/
40 /**                routines.                               **/
41 /**                                                        **/
42 /**   DATES      : # Version 5.0  : from : 17 oct 2007     **/
43 /**                                 to     07 dec 2007     **/
44 /**                # Version 5.1  : from : 18 mar 2009     **/
45 /**                                 to     30 jun 2010     **/
46 /**                # Version 6.0  : from : 13 sep 2012     **/
47 /**                                 to     13 sep 2012     **/
48 /**                                                        **/
49 /************************************************************/
50 
51 /*
52 **  The defines and includes.
53 */
54 
55 #define LIBRARY
56 
57 #include "module.h"
58 #include "common.h"
59 #include "ptscotch.h"
60 #include "parmetis.h"                             /* Our "parmetis.h" file */
61 
62 /************************************/
63 /*                                  */
64 /* These routines are the C API for */
65 /* the ParMeTiS graph ordering      */
66 /* routine.                         */
67 /*                                  */
68 /************************************/
69 
70 static
71 void
_SCOTCH_ParMETIS_V3_NodeNDTree(SCOTCH_Num * const sizeglbtnd,SCOTCH_Num * const sizeglbtab,SCOTCH_Num * const sepaglbtab,const SCOTCH_Num levlmax,const SCOTCH_Num levlnum,const SCOTCH_Num cblknum,SCOTCH_Num cblkidx)72 _SCOTCH_ParMETIS_V3_NodeNDTree (
73 SCOTCH_Num * const          sizeglbtnd,
74 SCOTCH_Num * const          sizeglbtab,
75 SCOTCH_Num * const          sepaglbtab,
76 const SCOTCH_Num            levlmax,
77 const SCOTCH_Num            levlnum,
78 const SCOTCH_Num            cblknum,
79 SCOTCH_Num                  cblkidx)
80 {
81   SCOTCH_Num          sizeval;
82 
83   sizeval = sizeglbtab[cblknum];                  /* Assume node is terminal or has no nested dissection sons */
84   if (levlnum < levlmax) {
85     if ((sepaglbtab[3 * cblknum]     >= 0) &&     /* If node has at least two sons, assume is is a nested dissection node */
86         (sepaglbtab[3 * cblknum + 1] >= 0)) {
87       _SCOTCH_ParMETIS_V3_NodeNDTree (sizeglbtnd, sizeglbtab, sepaglbtab, levlmax, levlnum + 1, sepaglbtab[3 * cblknum],     (cblkidx << 1) + 1);
88       _SCOTCH_ParMETIS_V3_NodeNDTree (sizeglbtnd, sizeglbtab, sepaglbtab, levlmax, levlnum + 1, sepaglbtab[3 * cblknum + 1], (cblkidx << 1));
89       sizeval = (sepaglbtab[3 * cblknum + 2] < 0) ? 0 : sizeglbtab[sepaglbtab[3 * cblknum + 2]]; /* Get size of separator, if any */
90     }
91   }
92 
93   sizeglbtnd[- cblkidx] = sizeval;                /* Set size of current column block */
94 }
95 
96 /*
97 **
98 */
99 
100 void
METISNAMEU(ParMETIS_V3_NodeND)101 METISNAMEU(ParMETIS_V3_NodeND) (
102 const SCOTCH_Num * const    vtxdist,
103 SCOTCH_Num * const          xadj,
104 SCOTCH_Num * const          adjncy,
105 const SCOTCH_Num * const    numflag,
106 const SCOTCH_Num * const    options,              /* Not used */
107 SCOTCH_Num * const          order,
108 SCOTCH_Num * const          sizes,                /* Of size twice the number of processors ; not used */
109 MPI_Comm *                  comm)
110 {
111   MPI_Comm            proccomm;
112   int                 procglbnbr;
113   int                 proclocnum;
114   SCOTCH_Num          baseval;
115   SCOTCH_Dgraph       grafdat;                    /* Scotch distributed graph object to interface with libScotch    */
116   SCOTCH_Dordering    ordedat;                    /* Scotch distributed ordering object to interface with libScotch */
117   SCOTCH_Strat        stradat;
118   SCOTCH_Num          vertlocnbr;
119   SCOTCH_Num          edgelocnbr;
120 
121   proccomm = *comm;
122   if (SCOTCH_dgraphInit (&grafdat, proccomm) != 0)
123     return;
124 
125   MPI_Comm_size (proccomm, &procglbnbr);
126   MPI_Comm_rank (proccomm, &proclocnum);
127   baseval    = *numflag;
128   vertlocnbr = vtxdist[proclocnum + 1] - vtxdist[proclocnum];
129   edgelocnbr = xadj[vertlocnbr] - baseval;
130 
131   if (sizes != NULL)
132     memSet (sizes, ~0, (2 * procglbnbr - 1) * sizeof (SCOTCH_Num)); /* Array not used if procglbnbr is not a power of 2 or if error */
133 
134   if (SCOTCH_dgraphBuild (&grafdat, baseval,
135                           vertlocnbr, vertlocnbr, xadj, xadj + 1, NULL, NULL,
136                           edgelocnbr, edgelocnbr, adjncy, NULL, NULL) == 0) {
137     SCOTCH_stratInit (&stradat);
138 #ifdef SCOTCH_DEBUG_ALL
139     if (SCOTCH_dgraphCheck (&grafdat) == 0)       /* TRICK: next instruction called only if graph is consistent */
140 #endif /* SCOTCH_DEBUG_ALL */
141     {
142       if (SCOTCH_dgraphOrderInit (&grafdat, &ordedat) == 0) {
143         SCOTCH_Num          levlmax;
144         SCOTCH_Num          bitsnbr;
145         SCOTCH_Num          proctmp;
146 
147         SCOTCH_dgraphOrderCompute (&grafdat, &ordedat, &stradat);
148         SCOTCH_dgraphOrderPerm    (&grafdat, &ordedat, order);
149 
150         for (levlmax = -1, bitsnbr = 0, proctmp = procglbnbr; /* Count number of bits set to 1 in procglbnbr */
151              proctmp != 0; levlmax ++, proctmp >>= 1)
152           bitsnbr += proctmp & 1;
153 
154         if (bitsnbr == 1) {
155           SCOTCH_Num          cblkglbnbr;
156 
157           if ((cblkglbnbr = SCOTCH_dgraphOrderCblkDist (&grafdat, &ordedat)) >= 0) {
158             SCOTCH_Num *        treeglbtab;
159             SCOTCH_Num *        sizeglbtab;
160             SCOTCH_Num *        sepaglbtab;
161 
162             if (memAllocGroup ((void **) (void *)
163                                &treeglbtab, (size_t) (cblkglbnbr * sizeof (SCOTCH_Num)),
164                                &sizeglbtab, (size_t) (cblkglbnbr * sizeof (SCOTCH_Num)),
165                                &sepaglbtab, (size_t) (cblkglbnbr * sizeof (SCOTCH_Num) * 3), NULL) != NULL) {
166               if (SCOTCH_dgraphOrderTreeDist (&grafdat, &ordedat, treeglbtab, sizeglbtab) == 0) {
167                 SCOTCH_Num          rootnum;
168                 SCOTCH_Num          cblknum;
169 
170                 memSet (sepaglbtab, ~0, cblkglbnbr * sizeof (SCOTCH_Num) * 3);
171 
172                 for (rootnum = -1, cblknum = 0; cblknum < cblkglbnbr; cblknum ++) {
173                   SCOTCH_Num          fathnum;
174 
175                   fathnum = treeglbtab[cblknum] - baseval; /* Use un-based indices  */
176                   if (fathnum < 0) {              /* If father index indicates root */
177                     if (rootnum != -1) {          /* If another root already found  */
178                       rootnum = -1;               /* Indicate an error              */
179                       break;
180                     }
181                     rootnum = cblknum;            /* Record index of root node */
182                   }
183                   else {
184                     SCOTCH_Num          i;
185 
186                     for (i = 0; i < 3; i ++) {
187                       SCOTCH_Num          j;
188 
189                       j = 3 * fathnum + i;        /* Slot number of prospective son  */
190                       if (sepaglbtab[j] < 0) {    /* If potentially empty slot found */
191                         if (sepaglbtab[j] == -1)  /* If we don't have too many sons  */
192                           sepaglbtab[j] = cblknum; /* Add link to son in slot        */
193                         break;
194                       }
195                     }
196                     if (i == 3) {                 /* If no empty slot found             */
197                       sepaglbtab[3 * fathnum] = -2; /* Indicate there are too many sons */
198                       break;
199                     }
200                   }
201                 }
202 
203                 if ((rootnum >= 0) && (sizes != NULL)) { /* If no error above, go on processing separator tree         */
204                   memSet (sizes, 0, (2 * procglbnbr - 1) * sizeof (SCOTCH_Num)); /* Set array of sizes to 0 by default */
205                   _SCOTCH_ParMETIS_V3_NodeNDTree (sizes + (2 * procglbnbr - 1), sizeglbtab, sepaglbtab, levlmax, 0, rootnum, 1);
206                 }
207               }
208 
209               memFree (treeglbtab);               /* Free group leader */
210             }
211           }
212         }
213 
214         SCOTCH_dgraphOrderExit (&grafdat, &ordedat);
215       }
216     }
217     SCOTCH_stratExit (&stradat);
218   }
219   SCOTCH_dgraphExit (&grafdat);
220 }
221