1 /* Copyright 2007,2008,2012,2014 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       : library_dgraph_stat.c                   **/
35 /**                                                        **/
36 /**   AUTHOR     : Francois PELLEGRINI                     **/
37 /**                                                        **/
38 /**   FUNCTION   : This module is the API for the source   **/
39 /**                graph handling routines of the          **/
40 /**                libSCOTCH library.                      **/
41 /**                                                        **/
42 /**   DATES      : # Version 5.0  : from : 23 jun 2007     **/
43 /**                                 to     03 apr 2008     **/
44 /**                # Version 6.0  : from : 29 nov 2012     **/
45 /**                                 to     29 oct 2014     **/
46 /**                                                        **/
47 /************************************************************/
48 
49 /*
50 **  The defines and includes.
51 */
52 
53 #define LIBRARY
54 
55 #include "module.h"
56 #include "common.h"
57 #include "dgraph.h"
58 #include "library_dgraph_stat.h"
59 #include "ptscotch.h"
60 
61 /*
62 **  The static variables.
63 */
64 
65 static int                  dgraphstatblentab[2] = { 7, 3 };
66 static MPI_Datatype         dgraphstattypetab[2] = { GNUM_MPI, MPI_DOUBLE };
67 
68 /************************************/
69 /*                                  */
70 /* These routines are the C API for */
71 /* the graph handling routines.     */
72 /*                                  */
73 /************************************/
74 
75 /* This routine is the reduction-loc operator which
76 ** returns in inout[2] the rank of the process which
77 ** holds the best partition.
78 ** It returns:
79 ** - void  : in all cases.
80 */
81 
82 static
83 void
dgraphStatReduceAll(const DgraphStatData * const in,DgraphStatData * const inout,const int * const len,const MPI_Datatype * const typedat)84 dgraphStatReduceAll (
85 const DgraphStatData * const  in,                 /* First operand                              */
86 DgraphStatData * const        inout,              /* Second and output operand                  */
87 const int * const             len,                /* Number of instances; should be 1, not used */
88 const MPI_Datatype * const    typedat)            /* MPI datatype; not used                     */
89 {
90   if (inout->velomin > in->velomin)
91     inout->velomin = in->velomin;
92   if (inout->velomax < in->velomax)
93     inout->velomax = in->velomax;
94   if (inout->degrmin > in->degrmin)
95     inout->degrmin = in->degrmin;
96   if (inout->degrmax < in->degrmax)
97     inout->degrmax = in->degrmax;
98   if (inout->edlomin > in->edlomin)
99     inout->edlomin = in->edlomin;
100   if (inout->edlomax < in->edlomax)
101     inout->edlomax = in->edlomax;
102   inout->edlosum += in->edlosum;
103   inout->velodlt += in->velodlt;
104   inout->degrdlt += in->degrdlt;
105   inout->edlodlt += in->edlodlt;
106 }
107 
108 /*+ This routine computes statistics
109 *** on the given distributed graph.
110 *** It returns:
111 *** - VOID  : in all cases.
112 +*/
113 
114 int
SCOTCH_dgraphStat(const SCOTCH_Dgraph * const grafptr,SCOTCH_Num * const velominptr,SCOTCH_Num * const velomaxptr,SCOTCH_Num * const velosumptr,double * veloavgptr,double * velodltptr,SCOTCH_Num * const degrminptr,SCOTCH_Num * const degrmaxptr,double * degravgptr,double * degrdltptr,SCOTCH_Num * const edlominptr,SCOTCH_Num * const edlomaxptr,SCOTCH_Num * const edlosumptr,double * edloavgptr,double * edlodltptr)115 SCOTCH_dgraphStat (
116 const SCOTCH_Dgraph * const grafptr,
117 SCOTCH_Num * const          velominptr,
118 SCOTCH_Num * const          velomaxptr,
119 SCOTCH_Num * const          velosumptr,
120 double *                    veloavgptr,
121 double *                    velodltptr,
122 SCOTCH_Num * const          degrminptr,
123 SCOTCH_Num * const          degrmaxptr,
124 double *                    degravgptr,
125 double *                    degrdltptr,
126 SCOTCH_Num * const          edlominptr,
127 SCOTCH_Num * const          edlomaxptr,
128 SCOTCH_Num * const          edlosumptr,
129 double *                    edloavgptr,
130 double *                    edlodltptr)
131 {
132   const Dgraph *      srcgrafptr;
133   DgraphStatData      srcgstadat;
134   DgraphStatData      srclstadat;
135   MPI_Datatype        srctypedat;
136   MPI_Aint            srcdisptab[2];
137   MPI_Op              srcoperdat;
138   Gnum                vertlocnum;
139   double              veloglbavg;
140   double              velolocdlt;
141   Gnum                degrlocmin;
142   Gnum                degrlocmax;
143   double              degrglbavg;
144   double              degrlocdlt;
145   Gnum                edloglbsum;
146   double              edloglbavg;
147   double              edlolocdlt;
148   int                 o;
149 
150   srcgrafptr = (Dgraph *) grafptr;
151 
152   velolocdlt = 0.0L;
153   if (srcgrafptr->vertglbnbr > 0) {
154     if (srcgrafptr->veloloctax != NULL) {         /* If graph has vertex loads */
155       const Gnum * restrict veloloctax;
156       Gnum                  velolocmin;
157       Gnum                  velolocmax;
158 
159       veloloctax = srcgrafptr->veloloctax;
160       velolocmin = GNUMMAX;
161       velolocmax = 0;
162       veloglbavg = (double) srcgrafptr->veloglbsum / (double) srcgrafptr->vertglbnbr;
163 
164       for (vertlocnum = srcgrafptr->baseval; vertlocnum < srcgrafptr->vertlocnnd; vertlocnum ++) {
165         Gnum                velolocval;
166 
167         velolocval = veloloctax[vertlocnum];
168         if (velolocval < velolocmin)
169           velolocmin = velolocval;
170         if (velolocval > velolocmax)
171           velolocmax = velolocval;
172         velolocdlt += fabs ((double) velolocval - veloglbavg);
173       }
174 
175       srclstadat.velomin = velolocmin;
176       srclstadat.velomax = velolocmax;
177     }
178     else {
179       srclstadat.velomin =
180       srclstadat.velomax = 1;
181       veloglbavg         = 1.0L;
182     }
183   }
184   else {
185     srclstadat.velomin =
186     srclstadat.velomax = 0;
187     veloglbavg         = 0.0L;
188   }
189   srclstadat.velodlt = velolocdlt;
190 
191   degrlocmax = 0;
192   degrlocdlt = 0.0L;
193   if (srcgrafptr->vertglbnbr > 0) {
194     degrlocmin = GNUMMAX;
195     degrglbavg = (double) srcgrafptr->edgeglbnbr / (double) srcgrafptr->vertglbnbr;
196     for (vertlocnum = srcgrafptr->baseval; vertlocnum < srcgrafptr->vertlocnnd; vertlocnum ++) {
197       Gnum                degrlocval;
198 
199       degrlocval = srcgrafptr->vendloctax[vertlocnum] - srcgrafptr->vertloctax[vertlocnum]; /* Get vertex degree */
200       if (degrlocval < degrlocmin)
201         degrlocmin = degrlocval;
202       if (degrlocval > degrlocmax)
203         degrlocmax = degrlocval;
204       degrlocdlt += fabs ((double) degrlocval - degrglbavg);
205     }
206   }
207   else {
208     degrlocmin = 0;
209     degrglbavg = 0.0L;
210   }
211   srclstadat.degrmin = degrlocmin;
212   srclstadat.degrmax = degrlocmax;
213   srclstadat.degrdlt = degrlocdlt;
214 
215   edlolocdlt = 0.0L;
216   if (srcgrafptr->edgeglbnbr > 0) {
217     if (srcgrafptr->edloloctax != NULL) {         /* If graph has edge loads */
218       Gnum                edlolocmin;
219       Gnum                edlolocmax;
220       Gnum                edlolocsum;
221 
222       edlolocmin = GNUMMAX;
223       edlolocmax = 0;
224       edlolocsum = 0;
225 
226       for (vertlocnum = srcgrafptr->baseval; vertlocnum < srcgrafptr->vertlocnnd; vertlocnum ++) {
227         Gnum                edgelocnum;
228 
229         for (edgelocnum = srcgrafptr->vertloctax[vertlocnum];
230              edgelocnum < srcgrafptr->vendloctax[vertlocnum]; edgelocnum ++) {
231           Gnum                edlolocval;
232 
233           edlolocval  = srcgrafptr->edloloctax[edgelocnum];
234           edlolocsum += edlolocval;
235           if (edlolocval < edlolocmin)            /* Account for edge load */
236             edlolocmin = edlolocval;
237           if (edlolocval > edlolocmax)
238             edlolocmax = edlolocval;
239         }
240       }
241 
242       if (MPI_Allreduce (&edlolocsum, &edloglbsum, 1, GNUM_MPI, MPI_SUM, srcgrafptr->proccomm) != MPI_SUCCESS) {
243         errorPrint ("SCOTCH_dgraphStat: communication error (1)");
244         return     (1);
245       }
246       edloglbavg = (double) edloglbsum / (double) (2 * srcgrafptr->edgeglbnbr);
247 
248       for (vertlocnum = srcgrafptr->baseval; vertlocnum < srcgrafptr->vertlocnnd; vertlocnum ++) {
249         Gnum                edgelocnum;
250 
251         for (edgelocnum = srcgrafptr->vertloctax[vertlocnum];
252              edgelocnum < srcgrafptr->vendloctax[vertlocnum]; edgelocnum ++)
253           edlolocdlt += fabs ((double) srcgrafptr->edloloctax[edgelocnum] - edloglbavg);
254       }
255     }
256     else {
257       srclstadat.edlomin =
258       srclstadat.edlomax = 1;
259       edloglbsum         = srcgrafptr->edgeglbnbr / 2;
260       edloglbavg         = 1.0L;
261     }
262   }
263   else {
264     srclstadat.edlomin =
265     srclstadat.edlomax = 0;
266     edloglbsum         = 0;
267     edloglbavg         = 0.0L;
268   }
269   srclstadat.edlodlt = edlolocdlt;
270 
271 #if ((defined COMMON_MPI_VERSION) && (COMMON_MPI_VERSION <= 100))
272   MPI_Address (&srclstadat.velomin, &srcdisptab[0]);
273   MPI_Address (&srclstadat.velodlt, &srcdisptab[1]);
274 #else /* ((defined COMMON_MPI_VERSION) && (COMMON_MPI_VERSION <= 100)) */
275   MPI_Get_address (&srclstadat.velomin, &srcdisptab[0]);
276   MPI_Get_address (&srclstadat.velodlt, &srcdisptab[1]);
277 #endif /* ((defined COMMON_MPI_VERSION) && (COMMON_MPI_VERSION <= 100)) */
278   srcdisptab[1] -= srcdisptab[0];
279   srcdisptab[0] -= srcdisptab[0];
280 
281   o = 1;                                          /* Assume something will go wrong */
282 #if ((defined COMMON_MPI_VERSION) && (COMMON_MPI_VERSION <= 100))
283   if ((MPI_Type_struct (2, dgraphstatblentab, srcdisptab, dgraphstattypetab, &srctypedat) == MPI_SUCCESS) &&
284 #else /* ((defined COMMON_MPI_VERSION) && (COMMON_MPI_VERSION <= 100)) */
285   if ((MPI_Type_create_struct (2, dgraphstatblentab, srcdisptab, dgraphstattypetab, &srctypedat) == MPI_SUCCESS) &&
286 #endif /* ((defined COMMON_MPI_VERSION) && (COMMON_MPI_VERSION <= 100)) */
287       (MPI_Type_commit (&srctypedat) == MPI_SUCCESS)) {
288     if (MPI_Op_create ((MPI_User_function *) dgraphStatReduceAll, 0, &srcoperdat) == MPI_SUCCESS) {
289       if (MPI_Allreduce (&srclstadat, &srcgstadat, 1, srctypedat, srcoperdat, srcgrafptr->proccomm) == MPI_SUCCESS)
290         o = 0;
291 
292       MPI_Op_free (&srcoperdat);
293     }
294     MPI_Type_free (&srctypedat);
295   }
296   if (o != 0) {
297     errorPrint ("SCOTCH_dgraphStat: communication error (2)");
298     return     (1);
299   }
300 
301   if (velominptr != NULL)
302     *velominptr = (SCOTCH_Num) srcgstadat.velomin;
303   if (velomaxptr != NULL)
304     *velomaxptr = (SCOTCH_Num) srcgstadat.velomax;
305   if (velosumptr != NULL)
306     *velosumptr = (SCOTCH_Num) srcgrafptr->veloglbsum;
307   if (veloavgptr != NULL)
308     *veloavgptr = (double) veloglbavg;
309   if (velodltptr != NULL)
310     *velodltptr = srcgstadat.velodlt / (double) srcgrafptr->vertglbnbr;
311 
312   if (degrminptr != NULL)
313     *degrminptr = (SCOTCH_Num) srcgstadat.degrmin;
314   if (degrmaxptr != NULL)
315     *degrmaxptr = (SCOTCH_Num) srcgstadat.degrmax;
316   if (degravgptr != NULL)
317     *degravgptr = (double) degrglbavg;
318   if (degrdltptr != NULL)
319     *degrdltptr = srcgstadat.degrdlt / (double) srcgrafptr->vertglbnbr;
320 
321   if (edlominptr != NULL)
322     *edlominptr = (SCOTCH_Num) srcgstadat.edlomin;
323   if (edlomaxptr != NULL)
324     *edlomaxptr = (SCOTCH_Num) srcgstadat.edlomax;
325   if (edlosumptr != NULL)
326     *edlosumptr = (SCOTCH_Num) edloglbsum;
327   if (edloavgptr != NULL)
328     *edloavgptr = (double) edloglbavg;
329   if (edlodltptr != NULL)
330     *edlodltptr = srcgstadat.edlodlt / (double) srcgrafptr->edgeglbnbr;
331 
332   return (0);
333 }
334