コード例 #1
0
int
SCOTCH_graphPart (
SCOTCH_Graph * const        grafptr,              /*+ Graph to map     +*/
const SCOTCH_Num            partnbr,              /*+ Number of parts  +*/
SCOTCH_Strat * const        stratptr,             /*+ Mapping strategy +*/
SCOTCH_Num * const          maptab)               /*+ Mapping array    +*/
{
  SCOTCH_Arch         archdat;
  int                 o;

  SCOTCH_archInit  (&archdat);
  SCOTCH_archCmplt (&archdat, partnbr);
  o = SCOTCH_graphMap (grafptr, &archdat, stratptr, maptab);
  SCOTCH_archExit  (&archdat);

  return (o);
}
コード例 #2
0
int
SCOTCH_graphRepartFixed (
SCOTCH_Graph * const        grafptr,              /*+ Graph to map                +*/
const SCOTCH_Num            partnbr,              /*+ Number of parts             +*/
SCOTCH_Num * const          parotab,              /*+ Old partition array         +*/
const double                emraval,              /*+ Edge migration ratio        +*/
const SCOTCH_Num *          vmlotab,              /*+ Vertex migration cost array +*/
SCOTCH_Strat * const        straptr,              /*+ Mapping strategy            +*/
SCOTCH_Num * const          parttab)              /*+ Partition array             +*/
{
  SCOTCH_Arch         archdat;
  int                 o;

  SCOTCH_archInit  (&archdat);
  SCOTCH_archCmplt (&archdat, partnbr);
  o = SCOTCH_graphRemapFixed (grafptr, &archdat, parotab, emraval, vmlotab, straptr, parttab);
  SCOTCH_archExit (&archdat);

  return (o);
}
コード例 #3
0
// Call scotch with options from dictionary.
Foam::label Foam::scotchDecomp::decompose
(
    const List<int>& adjncy,
    const List<int>& xadj,
    const scalarField& cWeights,

    List<int>& finalDecomp
)
{
    // Dump graph
    if (decompositionDict_.found("scotchCoeffs"))
    {
        const dictionary& scotchCoeffs =
            decompositionDict_.subDict("scotchCoeffs");

        if (scotchCoeffs.found("writeGraph"))
        {
            Switch writeGraph(scotchCoeffs.lookup("writeGraph"));

            if (writeGraph)
            {
                OFstream str(mesh_.time().path() / mesh_.name() + ".grf");

                Info<< "Dumping Scotch graph file to " << str.name() << endl
                    << "Use this in combination with gpart." << endl;

                label version = 0;
                str << version << nl;
                // Numer of vertices
                str << xadj.size()-1 << ' ' << adjncy.size() << nl;
                // Numbering starts from 0
                label baseval = 0;
                // Has weights?
                label hasEdgeWeights = 0;
                label hasVertexWeights = 0;
                label numericflag = 10*hasEdgeWeights+hasVertexWeights;
                str << baseval << ' ' << numericflag << nl;
                for (label cellI = 0; cellI < xadj.size()-1; cellI++)
                {
                    label start = xadj[cellI];
                    label end = xadj[cellI+1];
                    str << end-start;

                    for (label i = start; i < end; i++)
                    {
                        str << ' ' << adjncy[i];
                    }
                    str << nl;
                }
            }
        }
    }


    // Strategy
    // ~~~~~~~~

    // Default.
    SCOTCH_Strat stradat;
    check(SCOTCH_stratInit(&stradat), "SCOTCH_stratInit");

    if (decompositionDict_.found("scotchCoeffs"))
    {
        const dictionary& scotchCoeffs =
            decompositionDict_.subDict("scotchCoeffs");


        string strategy;
        if (scotchCoeffs.readIfPresent("strategy", strategy))
        {
            if (debug)
            {
                Info<< "scotchDecomp : Using strategy " << strategy << endl;
            }
            SCOTCH_stratGraphMap(&stradat, strategy.c_str());
            //fprintf(stdout, "S\tStrat=");
            //SCOTCH_stratSave(&stradat, stdout);
            //fprintf(stdout, "\n");
        }
    }


    // Graph
    // ~~~~~

    List<int> velotab;


    // Check for externally provided cellweights and if so initialise weights
    scalar minWeights = gMin(cWeights);
    if (cWeights.size() > 0)
    {
        if (minWeights <= 0)
        {
            WarningIn
            (
                "scotchDecomp::decompose"
                "(const pointField&, const scalarField&)"
            )   << "Illegal minimum weight " << minWeights
                << endl;
        }

        if (cWeights.size() != xadj.size()-1)
        {
            FatalErrorIn
            (
                "scotchDecomp::decompose"
                "(const pointField&, const scalarField&)"
            )   << "Number of cell weights " << cWeights.size()
                << " does not equal number of cells " << xadj.size()-1
                << exit(FatalError);
        }

        // Convert to integers.
        velotab.setSize(cWeights.size());
        forAll(velotab, i)
        {
            velotab[i] = int(cWeights[i]/minWeights);
        }
    }



    SCOTCH_Graph grafdat;
    check(SCOTCH_graphInit(&grafdat), "SCOTCH_graphInit");
    check
    (
        SCOTCH_graphBuild
        (
            &grafdat,
            0,                      // baseval, c-style numbering
            xadj.size()-1,          // vertnbr, nCells
            xadj.begin(),           // verttab, start index per cell into adjncy
            &xadj[1],               // vendtab, end index  ,,
            velotab.begin(),        // velotab, vertex weights
            NULL,                   // vlbltab
            adjncy.size(),          // edgenbr, number of arcs
            adjncy.begin(),         // edgetab
            NULL                    // edlotab, edge weights
        ),
        "SCOTCH_graphBuild"
    );
    check(SCOTCH_graphCheck(&grafdat), "SCOTCH_graphCheck");


    // Architecture
    // ~~~~~~~~~~~~
    // (fully connected network topology since using switch)

    SCOTCH_Arch archdat;
    check(SCOTCH_archInit(&archdat), "SCOTCH_archInit");

    List<label> processorWeights;
    if (decompositionDict_.found("scotchCoeffs"))
    {
        const dictionary& scotchCoeffs =
            decompositionDict_.subDict("scotchCoeffs");

        scotchCoeffs.readIfPresent("processorWeights", processorWeights);
    }
    if (processorWeights.size())
    {
        if (debug)
        {
            Info<< "scotchDecomp : Using procesor weights " << processorWeights
                << endl;
        }
        check
        (
            SCOTCH_archCmpltw(&archdat, nProcessors_, processorWeights.begin()),
            "SCOTCH_archCmpltw"
        );
    }
    else
    {
        check
        (
            SCOTCH_archCmplt(&archdat, nProcessors_),
            "SCOTCH_archCmplt"
        );
    }


    //SCOTCH_Mapping mapdat;
    //SCOTCH_graphMapInit(&grafdat, &mapdat, &archdat, NULL);
    //SCOTCH_graphMapCompute(&grafdat, &mapdat, &stradat); /* Perform mapping */
    //SCOTCH_graphMapExit(&grafdat, &mapdat);


    // Hack:switch off fpu error trapping
#   ifdef LINUX_GNUC
    int oldExcepts = fedisableexcept
    (
        FE_DIVBYZERO
      | FE_INVALID
      | FE_OVERFLOW
    );
#   endif

    finalDecomp.setSize(xadj.size()-1);
    finalDecomp = 0;
    check
    (
        SCOTCH_graphMap
        (
            &grafdat,
            &archdat,
            &stradat,           // const SCOTCH_Strat *
            finalDecomp.begin() // parttab
        ),
        "SCOTCH_graphMap"
    );

#   ifdef LINUX_GNUC
    feenableexcept(oldExcepts);
#   endif



    //finalDecomp.setSize(xadj.size()-1);
    //check
    //(
    //    SCOTCH_graphPart
    //    (
    //        &grafdat,
    //        nProcessors_,       // partnbr
    //        &stradat,           // const SCOTCH_Strat *
    //        finalDecomp.begin() // parttab
    //    ),
    //    "SCOTCH_graphPart"
    //);

    // Release storage for graph
    SCOTCH_graphExit(&grafdat);
    // Release storage for strategy
    SCOTCH_stratExit(&stradat);
    // Release storage for network topology
    SCOTCH_archExit(&archdat);

    return 0;
}
コード例 #4
0
ファイル: scotch.c プロジェクト: firedrakeproject/petsc
static PetscErrorCode MatPartitioningApply_PTScotch_Private(MatPartitioning part, PetscBool useND, IS *partitioning)
{
  MPI_Comm                 pcomm,comm;
  MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
  PetscErrorCode           ierr;
  PetscMPIInt              rank;
  Mat                      mat  = part->adj;
  Mat_MPIAdj               *adj = (Mat_MPIAdj*)mat->data;
  PetscBool                flg,distributed;
  PetscBool                proc_weight_flg;
  PetscInt                 i,j,p,bs=1,nold;
  PetscInt                 *NDorder = NULL;
  PetscReal                *vwgttab,deltval;
  SCOTCH_Num               *locals,*velotab,*veloloctab,*edloloctab,vertlocnbr,edgelocnbr,nparts=part->n;

  PetscFunctionBegin;
  ierr = PetscObjectGetComm((PetscObject)part,&pcomm);CHKERRQ(ierr);
  /* Duplicate the communicator to be sure that PTSCOTCH attribute caching does not interfere with PETSc. */
  ierr = MPI_Comm_dup(pcomm,&comm);CHKERRQ(ierr);
  ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
  ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);CHKERRQ(ierr);
  if (!flg) {
    /* bs indicates if the converted matrix is "reduced" from the original and hence the
       resulting partition results need to be stretched to match the original matrix */
    nold = mat->rmap->n;
    ierr = MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&mat);CHKERRQ(ierr);
    if (mat->rmap->n > 0) bs = nold/mat->rmap->n;
    adj  = (Mat_MPIAdj*)mat->data;
  }

  proc_weight_flg = PETSC_TRUE;
  ierr = PetscOptionsGetBool(NULL, NULL, "-mat_partitioning_ptscotch_proc_weight", &proc_weight_flg, NULL);CHKERRQ(ierr);

  ierr = PetscMalloc1(mat->rmap->n+1,&locals);CHKERRQ(ierr);

  if (useND) {
#if defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
    PetscInt    *sizes, *seps, log2size, subd, *level, base = 0;
    PetscMPIInt size;

    ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
    log2size = PetscLog2Real(size);
    subd = PetscPowInt(2,log2size);
    if (subd != size) SETERRQ(comm,PETSC_ERR_SUP,"Only power of 2 communicator sizes");
    ierr = PetscMalloc1(mat->rmap->n,&NDorder);CHKERRQ(ierr);
    ierr = PetscMalloc3(2*size,&sizes,4*size,&seps,size,&level);CHKERRQ(ierr);
    SCOTCH_ParMETIS_V3_NodeND(mat->rmap->range,adj->i,adj->j,&base,NULL,NDorder,sizes,&comm);
    ierr = MatPartitioningSizesToSep_Private(subd,sizes,seps,level);CHKERRQ(ierr);
    for (i=0;i<mat->rmap->n;i++) {
      PetscInt loc;

      ierr = PetscFindInt(NDorder[i],2*subd,seps,&loc);CHKERRQ(ierr);
      if (loc < 0) {
        loc = -(loc+1);
        if (loc%2) { /* part of subdomain */
          locals[i] = loc/2;
        } else {
          ierr = PetscFindInt(NDorder[i],2*(subd-1),seps+2*subd,&loc);CHKERRQ(ierr);
          loc = loc < 0 ? -(loc+1)/2 : loc/2;
          locals[i] = level[loc];
        }
      } else locals[i] = loc/2;
    }
    ierr = PetscFree3(sizes,seps,level);CHKERRQ(ierr);
#else
    SETERRQ(pcomm,PETSC_ERR_SUP,"Need libptscotchparmetis.a compiled with -DSCOTCH_METIS_PREFIX");
#endif
  } else {
    velotab = NULL;
    if (proc_weight_flg) {
      ierr = PetscMalloc1(nparts,&vwgttab);CHKERRQ(ierr);
      ierr = PetscMalloc1(nparts,&velotab);CHKERRQ(ierr);
      for (j=0; j<nparts; j++) {
        if (part->part_weights) vwgttab[j] = part->part_weights[j]*nparts;
        else vwgttab[j] = 1.0;
      }
      for (i=0; i<nparts; i++) {
        deltval = PetscAbsReal(vwgttab[i]-PetscFloorReal(vwgttab[i]+0.5));
        if (deltval>0.01) {
          for (j=0; j<nparts; j++) vwgttab[j] /= deltval;
        }
      }
      for (i=0; i<nparts; i++) velotab[i] = (SCOTCH_Num)(vwgttab[i] + 0.5);
      ierr = PetscFree(vwgttab);CHKERRQ(ierr);
    }

    vertlocnbr = mat->rmap->range[rank+1] - mat->rmap->range[rank];
    edgelocnbr = adj->i[vertlocnbr];
    veloloctab = part->vertex_weights;
    edloloctab = adj->values;

    /* detect whether all vertices are located at the same process in original graph */
    for (p = 0; !mat->rmap->range[p+1] && p < nparts; ++p);
    distributed = (mat->rmap->range[p+1] == mat->rmap->N) ? PETSC_FALSE : PETSC_TRUE;

    if (distributed) {
      SCOTCH_Arch              archdat;
      SCOTCH_Dgraph            grafdat;
      SCOTCH_Dmapping          mappdat;
      SCOTCH_Strat             stradat;

      ierr = SCOTCH_dgraphInit(&grafdat,comm);CHKERRQ(ierr);
      ierr = SCOTCH_dgraphBuild(&grafdat,0,vertlocnbr,vertlocnbr,adj->i,adj->i+1,veloloctab,
                                NULL,edgelocnbr,edgelocnbr,adj->j,NULL,edloloctab);CHKERRQ(ierr);

#if defined(PETSC_USE_DEBUG)
      ierr = SCOTCH_dgraphCheck(&grafdat);CHKERRQ(ierr);
#endif

      ierr = SCOTCH_archInit(&archdat);CHKERRQ(ierr);
      ierr = SCOTCH_stratInit(&stradat);CHKERRQ(ierr);
      ierr = SCOTCH_stratDgraphMapBuild(&stradat,scotch->strategy,nparts,nparts,scotch->imbalance);CHKERRQ(ierr);

      if (velotab) {
        ierr = SCOTCH_archCmpltw(&archdat,nparts,velotab);CHKERRQ(ierr);
      } else {
        ierr = SCOTCH_archCmplt( &archdat,nparts);CHKERRQ(ierr);
      }
      ierr = SCOTCH_dgraphMapInit(&grafdat,&mappdat,&archdat,locals);CHKERRQ(ierr);
      ierr = SCOTCH_dgraphMapCompute(&grafdat,&mappdat,&stradat);CHKERRQ(ierr);

      SCOTCH_dgraphMapExit(&grafdat,&mappdat);
      SCOTCH_archExit(&archdat);
      SCOTCH_stratExit(&stradat);
      SCOTCH_dgraphExit(&grafdat);

    } else if (rank == p) {
      SCOTCH_Graph   grafdat;
      SCOTCH_Strat   stradat;

      ierr = SCOTCH_graphInit(&grafdat);CHKERRQ(ierr);
      ierr = SCOTCH_graphBuild(&grafdat,0,vertlocnbr,adj->i,adj->i+1,veloloctab,NULL,edgelocnbr,adj->j,edloloctab);CHKERRQ(ierr);

#if defined(PETSC_USE_DEBUG)
      ierr = SCOTCH_graphCheck(&grafdat);CHKERRQ(ierr);
#endif

      ierr = SCOTCH_stratInit(&stradat);CHKERRQ(ierr);
      ierr = SCOTCH_stratGraphMapBuild(&stradat,scotch->strategy,nparts,scotch->imbalance);CHKERRQ(ierr);

      ierr = SCOTCH_graphPart(&grafdat,nparts,&stradat,locals);CHKERRQ(ierr);

      SCOTCH_stratExit(&stradat);
      SCOTCH_graphExit(&grafdat);
    }

    ierr = PetscFree(velotab);CHKERRQ(ierr);
  }
  ierr = MPI_Comm_free(&comm);CHKERRQ(ierr);

  if (bs > 1) {
    PetscInt *newlocals;
    ierr = PetscMalloc1(bs*mat->rmap->n,&newlocals);CHKERRQ(ierr);
    for (i=0;i<mat->rmap->n;i++) {
      for (j=0;j<bs;j++) {
        newlocals[bs*i+j] = locals[i];
      }
    }
    ierr = PetscFree(locals);CHKERRQ(ierr);
    ierr = ISCreateGeneral(pcomm,bs*mat->rmap->n,newlocals,PETSC_OWN_POINTER,partitioning);CHKERRQ(ierr);
  } else {
    ierr = ISCreateGeneral(pcomm,mat->rmap->n,locals,PETSC_OWN_POINTER,partitioning);CHKERRQ(ierr);
  }
  if (useND) {
    IS ndis;

    if (bs > 1) {
      ierr = ISCreateBlock(pcomm,bs,mat->rmap->n,NDorder,PETSC_OWN_POINTER,&ndis);CHKERRQ(ierr);
    } else {
      ierr = ISCreateGeneral(pcomm,mat->rmap->n,NDorder,PETSC_OWN_POINTER,&ndis);CHKERRQ(ierr);
    }
    ierr = ISSetPermutation(ndis);CHKERRQ(ierr);
    ierr = PetscObjectCompose((PetscObject)(*partitioning),"_petsc_matpartitioning_ndorder",(PetscObject)ndis);CHKERRQ(ierr);
    ierr = ISDestroy(&ndis);CHKERRQ(ierr);
  }

  if (!flg) {
    ierr = MatDestroy(&mat);CHKERRQ(ierr);
  }
  PetscFunctionReturn(0);
}
コード例 #5
0
ファイル: gmap.c プロジェクト: Hartorn/AN304
int
main (
int                         argc,
char *                      argv[])
{
  SCOTCH_Graph        grafdat;                    /* Source graph              */
  SCOTCH_Num          grafflag;                   /* Source graph properties   */
  SCOTCH_Arch         archdat;                    /* Target architecture       */
  SCOTCH_Strat        stradat;                    /* Mapping strategy          */
  SCOTCH_Mapping      mapdat;                     /* Mapping data              */
  Clock               runtime[2];                 /* Timing variables          */
  double              kbalval;                    /* Imbalance tolerance value */
  int                 flagval;
  SCOTCH_Num          straval;
  char *              straptr;
  int                 i, j;

  flagval = C_FLAGNONE;                           /* Default behavior  */
  straval = 0;                                    /* No strategy flags */
  straptr = NULL;

  i = strlen (argv[0]);
  if ((i >= 5) && (strncmp (argv[0] + i - 5, "gpart", 5) == 0)) {
    flagval |= C_FLAGPART;
    C_paraNbr = 1;                                /* One more parameter       */
    C_fileNbr = 3;                                /* One less file to provide */
    errorProg ("gpart");
  }
  else
    errorProg ("gmap");

  intRandInit ();

  if ((argc >= 2) && (argv[1][0] == '?')) {       /* If need for help */
    usagePrint (stdout, C_usageList);
    return     (0);
  }

  grafflag = 0;                                   /* Use vertex and edge weights  */
  SCOTCH_stratInit (&stradat);                    /* Set default mapping strategy */

  kbalval = 0.01;                                 /* Set default load imbalance value */

  for (i = 0; i < C_FILENBR; i ++)                /* Set default stream pointers */
    C_fileTab[i].pntr = (C_fileTab[i].mode[0] == 'r') ? stdin : stdout;
  for (i = 1; i < argc; i ++) {                   /* Loop for all option codes                        */
    if ((argv[i][0] != '-') || (argv[i][1] == '\0') || (argv[i][1] == '.')) { /* If found a file name */
      if (C_paraNum < C_paraNbr) {                /* If number of parameters not reached              */
        if ((C_partNbr = atoi (argv[i])) < 1)     /* Get the number of parts                          */
          errorPrint ("main: invalid number of parts (\"%s\")", argv[i]);
        C_paraNum ++;
        continue;                                 /* Process the other parameters */
      }
      if (C_fileNum < C_fileNbr)                  /* A file name has been given */
        C_fileTab[C_fileNum ++].name = argv[i];
      else
        errorPrint ("main: too many file names given");
    }
    else {                                        /* If found an option name */
      switch (argv[i][1]) {
        case 'B' :
        case 'b' :
          flagval |= C_FLAGKBALVAL;
          kbalval = atof (&argv[i][2]);
          if ((kbalval < 0.0) ||
              (kbalval > 1.0) ||
              ((kbalval == 0.0) &&
               ((argv[i][2] != '0') && (argv[i][2] != '.')))) {
            errorPrint ("main: invalid load imbalance ratio");
          }
          break;
        case 'C' :
        case 'c' :                                /* Strategy selection parameters */
          for (j = 2; argv[i][j] != '\0'; j ++) {
            switch (argv[i][j]) {
              case 'B' :
              case 'b' :
                straval |= SCOTCH_STRATBALANCE;
                break;
              case 'Q' :
              case 'q' :
                straval |= SCOTCH_STRATQUALITY;
                break;
              case 'S' :
              case 's' :
                straval |= SCOTCH_STRATSPEED;
                break;
              case 'T' :
              case 't' :
                straval |= SCOTCH_STRATSAFETY;
                break;
              default :
                errorPrint ("main: invalid strategy selection option (\"%c\")", argv[i][j]);
            }
          }
          break;
        case 'H' :                                /* Give the usage message */
        case 'h' :
          usagePrint (stdout, C_usageList);
          return     (0);
        case 'M' :
        case 'm' :
          straptr = &argv[i][2];
          SCOTCH_stratExit (&stradat);
          SCOTCH_stratInit (&stradat);
          SCOTCH_stratGraphMap (&stradat, straptr);
          break;
        case 'S' :
        case 's' :                                /* Source graph parameters */
          for (j = 2; argv[i][j] != '\0'; j ++) {
            switch (argv[i][j]) {
              case 'E' :
              case 'e' :
                grafflag |= 2;                    /* Do not load edge weights */
                break;
              case 'V' :
              case 'v' :
                grafflag |= 1;                    /* Do not load vertex weights */
                break;
              default :
                errorPrint ("main: invalid source graph option (\"%c\")", argv[i][j]);
            }
          }
          break;
        case 'V' :
          fprintf (stderr, "gmap/gpart, version " SCOTCH_VERSION_STRING "\n");
          fprintf (stderr, "Copyright 2004,2007,2008,2010 ENSEIRB, INRIA & CNRS, France\n");
          fprintf (stderr, "This software is libre/free software under CeCILL-C -- see the user's manual for more information\n");
          return  (0);
        case 'v' :                                /* Output control info */
          for (j = 2; argv[i][j] != '\0'; j ++) {
            switch (argv[i][j]) {
              case 'M' :
              case 'm' :
                flagval |= C_FLAGVERBMAP;
                break;
              case 'S' :
              case 's' :
                flagval |= C_FLAGVERBSTR;
                break;
              case 'T' :
              case 't' :
                flagval |= C_FLAGVERBTIM;
                break;
              default :
                errorPrint ("main: unprocessed parameter \"%c\" in \"%s\"", argv[i][j], argv[i]);
            }
          }
          break;
        default :
          errorPrint ("main: unprocessed option (\"%s\")", argv[i]);
      }
    }
  }

  if ((flagval & C_FLAGPART) != 0) {              /* If program run as the partitioner            */
    C_fileTab[3].name = C_fileTab[2].name;        /* Put provided file names at their right place */
    C_fileTab[2].name = C_fileTab[1].name;
    C_fileTab[1].name = "-";
  }

  fileBlockOpen (C_fileTab, C_FILENBR);           /* Open all files */

  clockInit  (&runtime[0]);
  clockStart (&runtime[0]);

  SCOTCH_graphInit (&grafdat);                    /* Create graph structure         */
  SCOTCH_graphLoad (&grafdat, C_filepntrsrcinp, -1, grafflag); /* Read source graph */

  SCOTCH_archInit (&archdat);                     /* Create architecture structure          */
  if ((flagval & C_FLAGPART) != 0)                /* If program run as the partitioner      */
    SCOTCH_archCmplt (&archdat, C_partNbr);       /* Create a complete graph of proper size */
  else {
    SCOTCH_archLoad (&archdat, C_filepntrtgtinp); /* Read target architecture */
    C_partNbr = SCOTCH_archSize (&archdat);
  }

  if ((straval != 0) || ((flagval & C_FLAGKBALVAL) != 0)) {
    if (straptr != NULL)
      errorPrint ("main: options '-b' / '-c' and '-m' are exclusive");

    SCOTCH_stratGraphMapBuild (&stradat, straval, (SCOTCH_Num) C_partNbr, kbalval);
  }

  clockStop  (&runtime[0]);                       /* Get input time */
  clockInit  (&runtime[1]);
  clockStart (&runtime[1]);

  SCOTCH_graphMapInit    (&grafdat, &mapdat, &archdat, NULL);
  SCOTCH_graphMapCompute (&grafdat, &mapdat, &stradat); /* Perform mapping */

  clockStop  (&runtime[1]);                       /* Get computation time */
  clockStart (&runtime[0]);

  SCOTCH_graphMapSave (&grafdat, &mapdat, C_filepntrmapout); /* Write mapping */

  clockStop (&runtime[0]);                        /* Get output time */

  if (flagval & C_FLAGVERBSTR) {
    fprintf (C_filepntrlogout, "S\tStrat=");
    SCOTCH_stratSave (&stradat, C_filepntrlogout);
    putc ('\n', C_filepntrlogout);
  }
  if (flagval & C_FLAGVERBTIM) {
    fprintf (C_filepntrlogout, "T\tMapping\t\t%g\nT\tI/O\t\t%g\nT\tTotal\t\t%g\n",
             (double) clockVal (&runtime[1]),
             (double) clockVal (&runtime[0]),
             (double) clockVal (&runtime[0]) +
             (double) clockVal (&runtime[1]));
  }
  if (flagval & C_FLAGVERBMAP)
    SCOTCH_graphMapView (&grafdat, &mapdat, C_filepntrlogout);

  fileBlockClose (C_fileTab, C_FILENBR);          /* Always close explicitely to end eventual (un)compression tasks */

  SCOTCH_graphMapExit (&grafdat, &mapdat);
  SCOTCH_graphExit    (&grafdat);
  SCOTCH_stratExit    (&stradat);
  SCOTCH_archExit     (&archdat);

#ifdef COMMON_PTHREAD
  pthread_exit ((void *) 0);                      /* Allow potential (un)compression tasks to complete */
#endif /* COMMON_PTHREAD */
  return (0);
}
コード例 #6
0
// Call scotch with options from dictionary.
Foam::label Foam::ptscotchDecomp::decompose
(
    const fileName& meshPath,
    const List<int>& adjncy,
    const List<int>& xadj,
    const scalarField& cWeights,

    List<int>& finalDecomp
) const
{
    if (debug)
    {
        Pout<< "ptscotchDecomp : entering with xadj:" << xadj.size() << endl;
    }

    // Dump graph
    if (decompositionDict_.found("ptscotchCoeffs"))
    {
        const dictionary& scotchCoeffs =
            decompositionDict_.subDict("ptscotchCoeffs");

        if (scotchCoeffs.lookupOrDefault("writeGraph", false))
        {
            OFstream str
            (
               meshPath + "_" + Foam::name(Pstream::myProcNo()) + ".dgr"
            );

            Pout<< "Dumping Scotch graph file to " << str.name() << endl
                << "Use this in combination with dgpart." << endl;

            globalIndex globalCells(xadj.size()-1);

            // Distributed graph file (.grf)
            label version = 2;
            str << version << nl;
            // Number of files (procglbnbr)
            str << Pstream::nProcs();
            // My file number (procloc)
            str << ' ' << Pstream::myProcNo() << nl;

            // Total number of vertices (vertglbnbr)
            str << globalCells.size();
            // Total number of connections (edgeglbnbr)
            str << ' ' << returnReduce(xadj[xadj.size()-1], sumOp<label>())
                << nl;
            // Local number of vertices (vertlocnbr)
            str << xadj.size()-1;
            // Local number of connections (edgelocnbr)
            str << ' ' << xadj[xadj.size()-1] << nl;
            // Numbering starts from 0
            label baseval = 0;
            // 100*hasVertlabels+10*hasEdgeWeights+1*hasVertWeighs
            str << baseval << ' ' << "000" << nl;
            for (label cellI = 0; cellI < xadj.size()-1; cellI++)
            {
                label start = xadj[cellI];
                label end = xadj[cellI+1];
                str << end-start;

                for (label i = start; i < end; i++)
                {
                    str << ' ' << adjncy[i];
                }
                str << nl;
            }
        }
    }

    // Strategy
    // ~~~~~~~~

    // Default.
    SCOTCH_Strat stradat;
    check(SCOTCH_stratInit(&stradat), "SCOTCH_stratInit");

    if (decompositionDict_.found("scotchCoeffs"))
    {
        const dictionary& scotchCoeffs =
            decompositionDict_.subDict("scotchCoeffs");


        string strategy;
        if (scotchCoeffs.readIfPresent("strategy", strategy))
        {
            if (debug)
            {
                Info<< "ptscotchDecomp : Using strategy " << strategy << endl;
            }
            SCOTCH_stratDgraphMap(&stradat, strategy.c_str());
            //fprintf(stdout, "S\tStrat=");
            //SCOTCH_stratSave(&stradat, stdout);
            //fprintf(stdout, "\n");
        }
    }


    // Graph
    // ~~~~~

    List<int> velotab;


    // Check for externally provided cellweights and if so initialise weights
    scalar minWeights = gMin(cWeights);
    if (cWeights.size() > 0)
    {
        if (minWeights <= 0)
        {
            WarningIn
            (
                "ptscotchDecomp::decompose(..)"
            )   << "Illegal minimum weight " << minWeights
                << endl;
        }

        if (cWeights.size() != xadj.size()-1)
        {
            FatalErrorIn
            (
                "ptscotchDecomp::decompose(..)"
            )   << "Number of cell weights " << cWeights.size()
                << " does not equal number of cells " << xadj.size()-1
                << exit(FatalError);
        }

        // Convert to integers.
        velotab.setSize(cWeights.size());
        forAll(velotab, i)
        {
            velotab[i] = int(cWeights[i]/minWeights);
        }
    }



    if (debug)
    {
        Pout<< "SCOTCH_dgraphInit" << endl;
    }
    SCOTCH_Dgraph grafdat;
    check(SCOTCH_dgraphInit(&grafdat, MPI_COMM_WORLD), "SCOTCH_dgraphInit");


    if (debug)
    {
        Pout<< "SCOTCH_dgraphBuild with:" << nl
            << "xadj.size()-1   : " << xadj.size()-1 << nl
            << "xadj            : " << long(xadj.begin()) << nl
            << "velotab         : " << long(velotab.begin()) << nl
            << "adjncy.size()   : " << adjncy.size() << nl
            << "adjncy          : " << long(adjncy.begin()) << nl
            << endl;
    }

    check
    (
        SCOTCH_dgraphBuild
        (
            &grafdat,               // grafdat
            0,                      // baseval, c-style numbering
            xadj.size()-1,          // vertlocnbr, nCells
            xadj.size()-1,          // vertlocmax
            const_cast<SCOTCH_Num*>(xadj.begin()),
                                    // vertloctab, start index per cell into
                                    // adjncy
            const_cast<SCOTCH_Num*>(&xadj[1]),// vendloctab, end index  ,,

            const_cast<SCOTCH_Num*>(velotab.begin()),// veloloctab, vtx weights
            NULL,                   // vlblloctab

            adjncy.size(),          // edgelocnbr, number of arcs
            adjncy.size(),          // edgelocsiz
            const_cast<SCOTCH_Num*>(adjncy.begin()),         // edgeloctab
            NULL,                   // edgegsttab
            NULL                    // edlotab, edge weights
        ),
        "SCOTCH_dgraphBuild"
    );


    if (debug)
    {
        Pout<< "SCOTCH_dgraphCheck" << endl;
    }
    check(SCOTCH_dgraphCheck(&grafdat), "SCOTCH_dgraphCheck");


    // Architecture
    // ~~~~~~~~~~~~
    // (fully connected network topology since using switch)

    if (debug)
    {
        Pout<< "SCOTCH_archInit" << endl;
    }
    SCOTCH_Arch archdat;
    check(SCOTCH_archInit(&archdat), "SCOTCH_archInit");

    List<label> processorWeights;
    if (decompositionDict_.found("scotchCoeffs"))
    {
        const dictionary& scotchCoeffs =
            decompositionDict_.subDict("scotchCoeffs");

        scotchCoeffs.readIfPresent("processorWeights", processorWeights);
    }
    if (processorWeights.size())
    {
        if (debug)
        {
            Info<< "ptscotchDecomp : Using procesor weights "
                << processorWeights
                << endl;
        }
        check
        (
            SCOTCH_archCmpltw(&archdat, nProcessors_, processorWeights.begin()),
            "SCOTCH_archCmpltw"
        );
    }
    else
    {
        if (debug)
        {
            Pout<< "SCOTCH_archCmplt" << endl;
        }
        check
        (
            SCOTCH_archCmplt(&archdat, nProcessors_),
            "SCOTCH_archCmplt"
        );
    }


    //SCOTCH_Mapping mapdat;
    //SCOTCH_dgraphMapInit(&grafdat, &mapdat, &archdat, NULL);
    //SCOTCH_dgraphMapCompute(&grafdat, &mapdat, &stradat); /*Perform mapping*/
    //SCOTCHdgraphMapExit(&grafdat, &mapdat);


    // Hack:switch off fpu error trapping
#   ifdef LINUX_GNUC
    int oldExcepts = fedisableexcept
    (
        FE_DIVBYZERO
      | FE_INVALID
      | FE_OVERFLOW
    );
#   endif

    if (debug)
    {
        Pout<< "SCOTCH_dgraphMap" << endl;
    }
    finalDecomp.setSize(xadj.size()-1);
    finalDecomp = 0;
    check
    (
        SCOTCH_dgraphMap
        (
            &grafdat,
            &archdat,
            &stradat,           // const SCOTCH_Strat *
            finalDecomp.begin() // parttab
        ),
        "SCOTCH_graphMap"
    );

#   ifdef LINUX_GNUC
    feenableexcept(oldExcepts);
#   endif



    //finalDecomp.setSize(xadj.size()-1);
    //check
    //(
    //    SCOTCH_dgraphPart
    //    (
    //        &grafdat,
    //        nProcessors_,       // partnbr
    //        &stradat,           // const SCOTCH_Strat *
    //        finalDecomp.begin() // parttab
    //    ),
    //    "SCOTCH_graphPart"
    //);

    if (debug)
    {
        Pout<< "SCOTCH_dgraphExit" << endl;
    }
    // Release storage for graph
    SCOTCH_dgraphExit(&grafdat);
    // Release storage for strategy
    SCOTCH_stratExit(&stradat);
    // Release storage for network topology
    SCOTCH_archExit(&archdat);

    return 0;
}
コード例 #7
0
ファイル: partition_old.hpp プロジェクト: Lazin/RaftLib
      static
      void run_scotch( Container         &c,
                       MapContainer      &mapping,
                       const std::size_t cores,
                       weight_function_t weight_func,
                       void              *weight )
   {
#if 0
      static_assert( std::is_signed< 
         std::remove_reference< decltype( c.end() ) >::type >::value, 
            "Container must have signed types so that -1 may signify no mapping" );
#endif            
      raftgraph_t raft_graph;
      get_graph_info( c, 
                      raft_graph, 
                      weight_func, 
                      nullptr );
      SCOTCH_Graph graph;
      if( SCOTCH_graphInit( &graph ) != 0 )
      {
         /** TODO, add RaftLib Exception **/
         std::cerr << "Failed to initialize graph!!\n";
         exit( EXIT_FAILURE );
      }
      auto table( raft_graph.getScotchTables() );
      if( SCOTCH_graphBuild( 
            &graph                  /** graph ptr     **/,
            0                       /** base value    **/,
            table.num_vertices      /** vertex nmbr (zero indexed)   **/,
            table.vtable            /** vertex tab **/,
            &table.vtable[ 1 ]      /** vendtab **/,
            nullptr           /** velotab **/,
            nullptr           /** vlbltab **/,
            table.num_edges                 /** edge number **/,
            table.etable             /** edge tab **/,
            table.eweight         /** edlotab **/
          ) != 0 )
      {
         /** TODO, add RaftLib Exception **/
         std::cerr << "Failed to build graph\n";
         exit( EXIT_FAILURE );
      }
      if( SCOTCH_graphCheck( &graph ) != 0 )
      {
         /** TODO, add RaftLib Exception **/
         std::cerr << "Graph is inconsistent\n";
         std::remove_reference< decltype( table ) >::type::print( std::cerr, table );
         std::cerr << "\n";
         raft_graph.print( std::cerr );
         
         exit( EXIT_FAILURE );
      }
      /** TODO, we can do much more with this arch file **/
      SCOTCH_Arch archdat;
      if( SCOTCH_archInit( &archdat )  != 0 )
      {
         /** TODO, add RaftLib Exception **/
         std::cerr << "Architecture initialization failed\n";
         exit( EXIT_FAILURE );
      }
      /** core are equal **/
      if( SCOTCH_archCmplt( &archdat, cores /** num cores **/) != 0 )
      {
         /** TODO, add RaftLib Exception **/
         std::cerr << "Failed to create architecture file\n";
         exit( EXIT_FAILURE );
      }
      /** strategy **/
      SCOTCH_Strat stradat;
      if( SCOTCH_stratInit( &stradat ) != 0 )
      {
         /** TODO, add RaftLib Exception **/
         std::cerr << "Failed to init strategy!!\n";
         exit( EXIT_FAILURE );
      }
      /** build recursive strategy **/
      if( SCOTCH_stratGraphClusterBuild(
                                   &stradat,
                                   SCOTCH_STRATSPEED,
                                   cores,
                                   .75,
                                   .01) != 0 )
      {
         /** TODO, add RaftLib Exception **/
         std::cerr << "Failed to map strategy graph!!\n";
         exit( EXIT_FAILURE );
      }
      if( SCOTCH_graphMap( 
            &graph             /** graph ptr **/,
            &archdat,
            &stradat,
            table.partition    /** parttab **/
            ) != 0 )
      {
         /** TODO, add RaftLib Exception **/
         std::cerr << "Failed to map!!\n";
         exit( EXIT_FAILURE );
      }
      /**
       * first case is for if we've mapped all vertices, 
       * second is for when some of the kernels are innactive
       * in which case the number of vertices in the 
       * table will be less than the size of c in which case
       * we need to get which vertices (the actual number id
       * from the application) are mapped and to where, the 
       * returned table in mapping must include even the 
       * vertices that aren't active (indicated by a -1) so
       * that the returning loop can be as simple as possible
       */
      if( c.size() == table.num_vertices )
      {
         /** copy mapping **/ 
         for( auto i( 0 ); i < table.num_vertices; i++ )
         {
            mapping.emplace_back( table.partition[ i ] );
         }
      }
      else
      {
         const auto &vmapping( raft_graph.getVertexNumbersAtIndicies() );
         auto it_map_index( vmapping.cbegin() );
         auto table_index( 0 );
         const auto size( c.size() );
         for( auto i( 0 ); i < size; i++ )
         {
            if( i == (*it_map_index) &&  it_map_index != vmapping.cend() )
            {
               mapping.emplace_back( table.partition[ table_index++ ] );
               ++it_map_index;
            }
            else
            {
               mapping.emplace_back( -1 );
            }
         }
      }
      /** call exit graph **/
      SCOTCH_graphExit( &graph    );
      SCOTCH_stratExit( &stradat );
      SCOTCH_archExit ( &archdat );
      return;
   }