xref: /petsc/src/mat/utils/compressedrow.c (revision cd6b891e1504900ecbee365dbcdfc13910f27c49)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
273e7a558SHong Zhang 
37c4f633dSBarry Smith #include "private/matimpl.h"  /*I   "petscmat.h"  I*/
426e093fcSHong Zhang 
573e7a558SHong Zhang #undef __FUNCT__
6*cd6b891eSBarry Smith #define __FUNCT__ "MatCheckCompressedRow"
773e7a558SHong Zhang /*@C
8*cd6b891eSBarry Smith    MatCheckCompressedRow - Determines whether the compressed row matrix format should be used.
9f38b99b6SHong Zhang       If the format is to be used, this routine creates Mat_CompressedRow struct.
10f38b99b6SHong Zhang       Compressed row format provides high performance routines by taking advantage of zero rows.
11f38b99b6SHong Zhang       Supported types are MATAIJ, MATBAIJ and MATSBAIJ.
1273e7a558SHong Zhang 
1373e7a558SHong Zhang    Collective
1473e7a558SHong Zhang 
1573e7a558SHong Zhang    Input Parameters:
1673e7a558SHong Zhang +  A             - the matrix
1773e7a558SHong Zhang .  compressedrow - pointer to the struct Mat_CompressedRow
1873e7a558SHong Zhang .  ai            - row pointer used by seqaij and seqbaij
19317fbc4cSHong Zhang .  mbs           - number of (block) rows represented by ai
2073e7a558SHong Zhang -  ratio         - ratio of (num of zero rows)/m, used to determine if the compressed row format should be used
2173e7a558SHong Zhang 
22*cd6b891eSBarry Smith    Notes: By default PETSc will not check for compressed rows on sequential matrices. Call MatSetOption(Mat,MAT_CHECK_COMPRESSED_ROW,PETSC_TRUE); before
23*cd6b891eSBarry Smith           MatAssemblyBegin() to have it check.
24*cd6b891eSBarry Smith 
2573e7a558SHong Zhang    Level: developer
2673e7a558SHong Zhang @*/
27*cd6b891eSBarry Smith PetscErrorCode MatCheckCompressedRow(Mat A,Mat_CompressedRow *compressedrow,PetscInt *ai,PetscInt mbs,PetscReal ratio)
2873e7a558SHong Zhang {
2973e7a558SHong Zhang   PetscErrorCode ierr;
30317fbc4cSHong Zhang   PetscInt       nrows,*cpi=PETSC_NULL,*ridx=PETSC_NULL,nz,i,row;
3173e7a558SHong Zhang 
3273e7a558SHong Zhang   PetscFunctionBegin;
33*cd6b891eSBarry Smith   if (!compressedrow->check) PetscFunctionReturn(0);
34*cd6b891eSBarry Smith 
35*cd6b891eSBarry Smith   /* in case this is being reused, delete old space */
360e83c824SBarry Smith   ierr = PetscFree2(compressedrow->i,compressedrow->rindex);CHKERRQ(ierr);
370e83c824SBarry Smith   compressedrow->i      = PETSC_NULL;
3888e51ccdSHong Zhang   compressedrow->rindex = PETSC_NULL;
39*cd6b891eSBarry Smith 
4073e7a558SHong Zhang 
4173e7a558SHong Zhang   /* compute number of zero rows */
4273e7a558SHong Zhang   nrows = 0;
43317fbc4cSHong Zhang   for (i=0; i<mbs; i++){        /* for each row */
4473e7a558SHong Zhang     nz = ai[i+1] - ai[i];       /* number of nonzeros */
4573e7a558SHong Zhang     if (nz == 0) nrows++;
4673e7a558SHong Zhang   }
47317fbc4cSHong Zhang   /* if a large number of zero rows is found, use compressedrow data structure */
48317fbc4cSHong Zhang   if (nrows < ratio*mbs) {
4973e7a558SHong Zhang     compressedrow->use = PETSC_FALSE;
50ae15b995SBarry 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);
5173e7a558SHong Zhang   } else {
5273e7a558SHong Zhang     compressedrow->use = PETSC_TRUE;
53ae15b995SBarry Smith     ierr = PetscInfo3(A,"Found the ratio (num_zerorows %d)/(num_localrows %d) > %G. Use CompressedRow routines.\n",nrows,mbs,ratio);CHKERRQ(ierr);
5473e7a558SHong Zhang 
5573e7a558SHong Zhang     /* set compressed row format */
56317fbc4cSHong Zhang     nrows = mbs - nrows; /* num of non-zero rows */
570e83c824SBarry Smith     ierr = PetscMalloc2(nrows+1,PetscInt,&cpi,nrows,PetscInt,&ridx);CHKERRQ(ierr);
5873e7a558SHong Zhang     row    = 0;
5973e7a558SHong Zhang     cpi[0] = 0;
60317fbc4cSHong Zhang     for (i=0; i<mbs; i++){
6173e7a558SHong Zhang       nz = ai[i+1] - ai[i];
6273e7a558SHong Zhang       if (nz == 0) continue;
6373e7a558SHong Zhang       cpi[row+1]  = ai[i+1];    /* compressed row pointer */
647b2bb3b9SHong Zhang       ridx[row++] = i;          /* compressed row local index */
6573e7a558SHong Zhang     }
6673e7a558SHong Zhang     compressedrow->nrows  = nrows;
6773e7a558SHong Zhang     compressedrow->i      = cpi;
687b2bb3b9SHong Zhang     compressedrow->rindex = ridx;
6973e7a558SHong Zhang   }
70f38b99b6SHong Zhang 
7173e7a558SHong Zhang   PetscFunctionReturn(0);
7273e7a558SHong Zhang }
73