xref: /petsc/src/mat/utils/compressedrow.c (revision 8865f1ea4ea560cd84ab8db62e98b7095cdff96f)
173e7a558SHong Zhang 
2b45d2f2cSJed Brown #include <petsc-private/matimpl.h>  /*I   "petscmat.h"  I*/
326e093fcSHong Zhang 
473e7a558SHong Zhang #undef __FUNCT__
5cd6b891eSBarry Smith #define __FUNCT__ "MatCheckCompressedRow"
673e7a558SHong Zhang /*@C
7cd6b891eSBarry Smith    MatCheckCompressedRow - Determines whether the compressed row matrix format should be used.
8f38b99b6SHong Zhang       If the format is to be used, this routine creates Mat_CompressedRow struct.
9f38b99b6SHong Zhang       Compressed row format provides high performance routines by taking advantage of zero rows.
10f38b99b6SHong Zhang       Supported types are MATAIJ, MATBAIJ and MATSBAIJ.
1173e7a558SHong Zhang 
1273e7a558SHong Zhang    Collective
1373e7a558SHong Zhang 
1473e7a558SHong Zhang    Input Parameters:
1573e7a558SHong Zhang +  A             - the matrix
1673e7a558SHong Zhang .  compressedrow - pointer to the struct Mat_CompressedRow
1773e7a558SHong Zhang .  ai            - row pointer used by seqaij and seqbaij
18317fbc4cSHong Zhang .  mbs           - number of (block) rows represented by ai
1973e7a558SHong Zhang -  ratio         - ratio of (num of zero rows)/m, used to determine if the compressed row format should be used
2073e7a558SHong Zhang 
21cd6b891eSBarry Smith    Notes: By default PETSc will not check for compressed rows on sequential matrices. Call MatSetOption(Mat,MAT_CHECK_COMPRESSED_ROW,PETSC_TRUE); before
22cd6b891eSBarry Smith           MatAssemblyBegin() to have it check.
23cd6b891eSBarry Smith 
24340b11eeSBarry Smith    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
25340b11eeSBarry Smith                    the values are not therefore obtained by directly taking the values from the matrix object.
2673e7a558SHong Zhang    Level: developer
2773e7a558SHong Zhang @*/
28cd6b891eSBarry Smith PetscErrorCode MatCheckCompressedRow(Mat A,Mat_CompressedRow *compressedrow,PetscInt *ai,PetscInt mbs,PetscReal ratio)
2973e7a558SHong Zhang {
3073e7a558SHong Zhang   PetscErrorCode ierr;
31317fbc4cSHong Zhang   PetscInt       nrows,*cpi=PETSC_NULL,*ridx=PETSC_NULL,nz,i,row;
3273e7a558SHong Zhang 
3373e7a558SHong Zhang   PetscFunctionBegin;
34cd6b891eSBarry Smith   if (!compressedrow->check) PetscFunctionReturn(0);
35cd6b891eSBarry Smith 
36cd6b891eSBarry Smith   /* in case this is being reused, delete old space */
370e83c824SBarry Smith   ierr = PetscFree2(compressedrow->i,compressedrow->rindex);CHKERRQ(ierr);
38*8865f1eaSKarl Rupp 
390e83c824SBarry Smith   compressedrow->i      = PETSC_NULL;
4088e51ccdSHong Zhang   compressedrow->rindex = PETSC_NULL;
41cd6b891eSBarry Smith 
4273e7a558SHong Zhang 
4373e7a558SHong Zhang   /* compute number of zero rows */
4473e7a558SHong Zhang   nrows = 0;
45317fbc4cSHong Zhang   for (i=0; i<mbs; i++) {        /* for each row */
4673e7a558SHong Zhang     nz = ai[i+1] - ai[i];       /* number of nonzeros */
4773e7a558SHong Zhang     if (nz == 0) nrows++;
4873e7a558SHong Zhang   }
49baa4e9faSMark F. Adams 
50317fbc4cSHong Zhang   /* if a large number of zero rows is found, use compressedrow data structure */
51317fbc4cSHong Zhang   if (nrows < ratio*mbs) {
5273e7a558SHong Zhang     compressedrow->use = PETSC_FALSE;
53*8865f1eaSKarl Rupp 
54ae15b995SBarry Smith     ierr = PetscInfo3(A,"Found the ratio (num_zerorows %d)/(num_localrows %d) < %G. Do not use CompressedRow routines.\n",nrows,mbs,ratio);CHKERRQ(ierr);
5573e7a558SHong Zhang   } else {
5673e7a558SHong Zhang     compressedrow->use = PETSC_TRUE;
57*8865f1eaSKarl Rupp 
58ae15b995SBarry Smith     ierr = PetscInfo3(A,"Found the ratio (num_zerorows %d)/(num_localrows %d) > %G. Use CompressedRow routines.\n",nrows,mbs,ratio);CHKERRQ(ierr);
5973e7a558SHong Zhang 
6073e7a558SHong Zhang     /* set compressed row format */
61317fbc4cSHong Zhang     nrows  = mbs - nrows; /* num of non-zero rows */
620e83c824SBarry Smith     ierr   = PetscMalloc2(nrows+1,PetscInt,&cpi,nrows,PetscInt,&ridx);CHKERRQ(ierr);
6373e7a558SHong Zhang     row    = 0;
6473e7a558SHong Zhang     cpi[0] = 0;
65317fbc4cSHong Zhang     for (i=0; i<mbs; i++) {
6673e7a558SHong Zhang       nz = ai[i+1] - ai[i];
6773e7a558SHong Zhang       if (nz == 0) continue;
6873e7a558SHong Zhang       cpi[row+1]  = ai[i+1];    /* compressed row pointer */
697b2bb3b9SHong Zhang       ridx[row++] = i;          /* compressed row local index */
7073e7a558SHong Zhang     }
7173e7a558SHong Zhang     compressedrow->nrows  = nrows;
7273e7a558SHong Zhang     compressedrow->i      = cpi;
737b2bb3b9SHong Zhang     compressedrow->rindex = ridx;
7473e7a558SHong Zhang   }
7573e7a558SHong Zhang   PetscFunctionReturn(0);
7673e7a558SHong Zhang }
77