xref: /petsc/src/mat/utils/compressedrow.c (revision c6db04a5321582041def2b1e244c75985478b3ef)
173e7a558SHong Zhang 
2*c6db04a5SJed Brown #include <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 
2473e7a558SHong Zhang    Level: developer
2573e7a558SHong Zhang @*/
26cd6b891eSBarry Smith PetscErrorCode MatCheckCompressedRow(Mat A,Mat_CompressedRow *compressedrow,PetscInt *ai,PetscInt mbs,PetscReal ratio)
2773e7a558SHong Zhang {
2873e7a558SHong Zhang   PetscErrorCode ierr;
29317fbc4cSHong Zhang   PetscInt       nrows,*cpi=PETSC_NULL,*ridx=PETSC_NULL,nz,i,row;
3073e7a558SHong Zhang 
3173e7a558SHong Zhang   PetscFunctionBegin;
32cd6b891eSBarry Smith   if (!compressedrow->check) PetscFunctionReturn(0);
33cd6b891eSBarry Smith 
34cd6b891eSBarry Smith   /* in case this is being reused, delete old space */
350e83c824SBarry Smith   ierr = PetscFree2(compressedrow->i,compressedrow->rindex);CHKERRQ(ierr);
360e83c824SBarry Smith   compressedrow->i      = PETSC_NULL;
3788e51ccdSHong Zhang   compressedrow->rindex = PETSC_NULL;
38cd6b891eSBarry Smith 
3973e7a558SHong Zhang 
4073e7a558SHong Zhang   /* compute number of zero rows */
4173e7a558SHong Zhang   nrows = 0;
42317fbc4cSHong Zhang   for (i=0; i<mbs; i++){        /* for each row */
4373e7a558SHong Zhang     nz = ai[i+1] - ai[i];       /* number of nonzeros */
4473e7a558SHong Zhang     if (nz == 0) nrows++;
4573e7a558SHong Zhang   }
46317fbc4cSHong Zhang   /* if a large number of zero rows is found, use compressedrow data structure */
47317fbc4cSHong Zhang   if (nrows < ratio*mbs) {
4873e7a558SHong Zhang     compressedrow->use = PETSC_FALSE;
49ae15b995SBarry 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);
5073e7a558SHong Zhang   } else {
5173e7a558SHong Zhang     compressedrow->use = PETSC_TRUE;
52ae15b995SBarry Smith     ierr = PetscInfo3(A,"Found the ratio (num_zerorows %d)/(num_localrows %d) > %G. Use CompressedRow routines.\n",nrows,mbs,ratio);CHKERRQ(ierr);
5373e7a558SHong Zhang 
5473e7a558SHong Zhang     /* set compressed row format */
55317fbc4cSHong Zhang     nrows = mbs - nrows; /* num of non-zero rows */
560e83c824SBarry Smith     ierr = PetscMalloc2(nrows+1,PetscInt,&cpi,nrows,PetscInt,&ridx);CHKERRQ(ierr);
5773e7a558SHong Zhang     row    = 0;
5873e7a558SHong Zhang     cpi[0] = 0;
59317fbc4cSHong Zhang     for (i=0; i<mbs; i++){
6073e7a558SHong Zhang       nz = ai[i+1] - ai[i];
6173e7a558SHong Zhang       if (nz == 0) continue;
6273e7a558SHong Zhang       cpi[row+1]  = ai[i+1];    /* compressed row pointer */
637b2bb3b9SHong Zhang       ridx[row++] = i;          /* compressed row local index */
6473e7a558SHong Zhang     }
6573e7a558SHong Zhang     compressedrow->nrows  = nrows;
6673e7a558SHong Zhang     compressedrow->i      = cpi;
677b2bb3b9SHong Zhang     compressedrow->rindex = ridx;
6873e7a558SHong Zhang   }
69f38b99b6SHong Zhang 
7073e7a558SHong Zhang   PetscFunctionReturn(0);
7173e7a558SHong Zhang }
72