Actual source code: ij.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2:  #include <../src/mat/impls/aij/seq/aij.h>

  4: /*
  5:   MatToSymmetricIJ_SeqAIJ - Convert a (generally nonsymmetric) sparse AIJ matrix
  6:            to IJ format (ignore the "A" part) Allocates the space needed. Uses only
  7:            the lower triangular part of the matrix.

  9:     Description:
 10:     Take the data in the row-oriented sparse storage and build the
 11:     IJ data for the Matrix.  Return 0 on success,row + 1 on failure
 12:     at that row. Produces the ij for a symmetric matrix by using
 13:     the lower triangular part of the matrix if lower_triangular is PETSC_TRUE;
 14:     it uses the upper triangular otherwise.

 16:     Input Parameters:
 17: .   Matrix - matrix to convert
 18: .   lower_triangular - symmetrize the lower triangular part
 19: .   shiftin - the shift for the original matrix (0 or 1)
 20: .   shiftout - the shift required for the ordering routine (0 or 1)

 22:     Output Parameters:
 23: .   ia     - ia part of IJ representation (row information)
 24: .   ja     - ja part (column indices)

 26:     Notes:
 27:     Both ia and ja may be freed with PetscFree();
 28:     This routine is provided for ordering routines that require a
 29:     symmetric structure.  It is required since those routines call
 30:     SparsePak routines that expect a symmetric  matrix.
 31: */
 32: PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt m,PetscInt *ai,PetscInt *aj,PetscBool lower_triangular,PetscInt shiftin,PetscInt shiftout,PetscInt **iia,PetscInt **jja)
 33: {
 35:   PetscInt       *work,*ia,*ja,*j,i,nz,row,col;

 38:   /* allocate space for row pointers */
 39:   PetscMalloc1(m+1,&ia);
 40:   *iia = ia;
 41:   PetscMemzero(ia,(m+1)*sizeof(PetscInt));
 42:   PetscMalloc1(m+1,&work);

 44:   /* determine the number of columns in each row */
 45:   ia[0] = shiftout;
 46:   for (row = 0; row < m; row++) {
 47:     nz = ai[row+1] - ai[row];
 48:     j  = aj + ai[row] + shiftin;
 49:     while (nz--) {
 50:       col = *j++ + shiftin;
 51:       if (lower_triangular) {
 52:         if (col > row) break;
 53:       } else {
 54:         if (col < row) break;
 55:       }
 56:       if (col != row) ia[row+1]++;
 57:       ia[col+1]++;
 58:     }
 59:   }

 61:   /* shiftin ia[i] to point to next row */
 62:   for (i=1; i<m+1; i++) {
 63:     row       = ia[i-1];
 64:     ia[i]    += row;
 65:     work[i-1] = row - shiftout;
 66:   }

 68:   /* allocate space for column pointers */
 69:   nz   = ia[m] + (!shiftin);
 70:   PetscMalloc1(nz,&ja);
 71:   *jja = ja;

 73:   /* loop over lower triangular part putting into ja */
 74:   for (row = 0; row < m; row++) {
 75:     nz = ai[row+1] - ai[row];
 76:     j  = aj + ai[row] + shiftin;
 77:     while (nz--) {
 78:       col = *j++ + shiftin;
 79:       if (lower_triangular) {
 80:         if (col > row) break;
 81:       } else {
 82:         if (col < row) break;
 83:       }
 84:       if (col != row) ja[work[col]++] = row + shiftout;
 85:       ja[work[row]++] = col + shiftout;
 86:     }
 87:   }
 88:   PetscFree(work);
 89:   return(0);
 90: }