Actual source code: gennd.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2: /* gennd.f -- translated by f2c (version 19931217).*/

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

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

 12:   /* Local variables */
 13:   PetscInt swap,i,m,in;

 16:   /* Parameter adjustments */
 17:   --perm;

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


 32: /*****************************************************************/
 33: /*********     GENND ..... GENERAL NESTED DISSECTION     *********/
 34: /*****************************************************************/

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

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

 43: /*    OUTPUT PARAMETERS -*/
 44: /*       PERM - THE NESTED DISSECTION ORDERING.*/

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

 52: /*    PROGRAM SUBROUTINES -*/
 53: /*       FNDSEP, REVRSE.*/
 54: /*****************************************************************/

 56: PetscErrorCode SPARSEPACKgennd(const PetscInt *neqns,const PetscInt *xadj,const PetscInt *adjncy,PetscInt *mask,PetscInt *perm,PetscInt *xls,PetscInt *ls)
 57: {
 58:   /* System generated locals */
 59:   PetscInt i__1;

 61:   /* Local variables */
 62:   PetscInt nsep,root,i;
 63:   PetscInt num;

 66:   /* Parameter adjustments */
 67:   --ls;
 68:   --xls;
 69:   --perm;
 70:   --mask;
 71:   --adjncy;
 72:   --xadj;

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