Actual source code: gennd.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: /* gennd.f -- translated by f2c (version 19931217).*/

  4: #include <petscsys.h>
  5: #include <petsc/private/matorderimpl.h>

  9: PetscErrorCode SPARSEPACKrevrse(const PetscInt *n,PetscInt *perm)
 10: {
 11:   /* System generated locals */
 12:   PetscInt i__1;

 14:   /* Local variables */
 15:   PetscInt swap,i,m,in;

 18:   /* Parameter adjustments */
 19:   --perm;

 21:   in   = *n;
 22:   m    = *n / 2;
 23:   i__1 = m;
 24:   for (i = 1; i <= i__1; ++i) {
 25:     swap     = perm[i];
 26:     perm[i]  = perm[in];
 27:     perm[in] = swap;
 28:     --in;
 29:   }
 30:   return(0);
 31: }


 34: /*****************************************************************/
 35: /*********     GENND ..... GENERAL NESTED DISSECTION     *********/
 36: /*****************************************************************/

 38: /*    PURPOSE - SUBROUTINE GENND FINDS A NESTED DISSECTION*/
 39: /*       ORDERING FOR A GENERAL GRAPH.*/

 41: /*    INPUT PARAMETERS -*/
 42: /*       NEQNS - NUMBER OF EQUATIONS.*/
 43: /*       (XADJ, ADJNCY) - ADJACENCY STRUCTURE PAIR.*/

 45: /*    OUTPUT PARAMETERS -*/
 46: /*       PERM - THE NESTED DISSECTION ORDERING.*/

 48: /*    WORKING PARAMETERS -*/
 49: /*       MASK - IS USED TO MASK OFF VARIABLES THAT HAVE*/
 50: /*              BEEN NUMBERED DURING THE ORDERNG PROCESS.*/
 51: /*       (XLS, LS) - THIS LEVEL STRUCTURE PAIR IS USED AS*/
 52: /*              TEMPORARY STORAGE BY FN../../...*/

 54: /*    PROGRAM SUBROUTINES -*/
 55: /*       FNDSEP, REVRSE.*/
 56: /*****************************************************************/

 60: PetscErrorCode SPARSEPACKgennd(const PetscInt *neqns,const PetscInt *xadj,const PetscInt *adjncy,PetscInt *mask,PetscInt *perm,PetscInt *xls,PetscInt *ls)
 61: {
 62:   /* System generated locals */
 63:   PetscInt i__1;

 65:   /* Local variables */
 66:   PetscInt nsep,root,i;
 67:   PetscInt num;

 70:   /* Parameter adjustments */
 71:   --ls;
 72:   --xls;
 73:   --perm;
 74:   --mask;
 75:   --adjncy;
 76:   --xadj;

 78:   i__1 = *neqns;
 79:   for (i = 1; i <= i__1; ++i) mask[i] = 1;
 80:   num  = 0;
 81:   i__1 = *neqns;
 82:   for (i = 1; i <= i__1; ++i) {
 83: /*           FOR EACH MASKED COMPONENT ...*/
 84: L200:
 85:     if (!mask[i]) goto L300;
 86:     root = i;
 87: /*              FIND A SEPARATOR AND NUMBER THE NODES NEXT.*/
 88:     SPARSEPACKfndsep(&root,&xadj[1],&adjncy[1],&mask[1],&nsep,&perm[num + 1],&xls[1],&ls[1]);
 89:     num += nsep;
 90:     if (num >= *neqns) goto L400;
 91:     goto L200;
 92: L300:
 93:     ;
 94:   }
 95: /*        SINCE SEPARATORS FOUND FIRST SHOULD BE ORDERED*/
 96: /*        LAST, ROUTINE REVRSE IS CALLED TO ADJUST THE*/
 97: /*        ORDERING VECTOR.*/
 98: L400:
 99:   SPARSEPACKrevrse(neqns,&perm[1]);
100:   return(0);
101: }