2: #include <petsc-private/matimpl.h> /*I "petscmat.h" I*/
6: /*@C
7: MatCheckCompressedRow - Determines whether the compressed row matrix format should be used.
8: If the format is to be used, this routine creates Mat_CompressedRow struct.
9: Compressed row format provides high performance routines by taking advantage of zero rows.
10: Supported types are MATAIJ, MATBAIJ and MATSBAIJ.
12: Collective
14: Input Parameters:
15: + A - the matrix
16: . compressedrow - pointer to the struct Mat_CompressedRow
17: . ai - row pointer used by seqaij and seqbaij
18: . mbs - number of (block) rows represented by ai
19: - ratio - ratio of (num of zero rows)/m, used to determine if the compressed row format should be used
21: Notes: By default PETSc will not check for compressed rows on sequential matrices. Call MatSetOption(Mat,MAT_CHECK_COMPRESSED_ROW,PETSC_TRUE); before
22: MatAssemblyBegin() to have it check.
24: Developer Note: The reason this takes the compressedrow, ai and mbs arguments is because it is called by both the SeqAIJ and SEQBAIJ matrices and
25: the values are not therefore obtained by directly taking the values from the matrix object.
26: Level: developer
27: @*/
28: PetscErrorCodeMatCheckCompressedRow(Mat A,Mat_CompressedRow *compressedrow,PetscInt *ai,PetscInt mbs,PetscReal ratio) 29: {
31: PetscInt nrows,*cpi=PETSC_NULL,*ridx=PETSC_NULL,nz,i,row;
34: if (!compressedrow->check) return(0);
36: /* in case this is being reused, delete old space */
37: PetscFree2(compressedrow->i,compressedrow->rindex);
38: compressedrow->i = PETSC_NULL;
39: compressedrow->rindex = PETSC_NULL;
42: /* compute number of zero rows */
43: nrows = 0;
44: for (i=0; i<mbs; i++){ /* for each row */
45: nz = ai[i+1] - ai[i]; /* number of nonzeros */
46: if (nz == 0) nrows++;
47: }
49: /* if a large number of zero rows is found, use compressedrow data structure */
50: if (nrows < ratio*mbs) {
51: compressedrow->use = PETSC_FALSE;
52: PetscInfo3(A,"Found the ratio (num_zerorows %d)/(num_localrows %d) < %G. Do not use CompressedRow routines.\n",nrows,mbs,ratio);
53: } else {
54: compressedrow->use = PETSC_TRUE;
55: PetscInfo3(A,"Found the ratio (num_zerorows %d)/(num_localrows %d) > %G. Use CompressedRow routines.\n",nrows,mbs,ratio);
57: /* set compressed row format */
58: nrows = mbs - nrows; /* num of non-zero rows */
59: PetscMalloc2(nrows+1,PetscInt,&cpi,nrows,PetscInt,&ridx);
60: row = 0;
61: cpi[0] = 0;
62: for (i=0; i<mbs; i++){
63: nz = ai[i+1] - ai[i];
64: if (nz == 0) continue;
65: cpi[row+1] = ai[i+1]; /* compressed row pointer */
66: ridx[row++] = i; /* compressed row local index */
67: }
68: compressedrow->nrows = nrows;
69: compressedrow->i = cpi;
70: compressedrow->rindex = ridx;
71: }
72: 73: return(0);
74: }