173e7a558SHong Zhang 2*26e093fcSHong Zhang #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 3*26e093fcSHong Zhang 473e7a558SHong Zhang #undef __FUNCT__ 573e7a558SHong Zhang #define __FUNCT__ "Mat_CheckCompressedRow" 673e7a558SHong Zhang /*@C 773e7a558SHong Zhang Mat_CheckCompressedRow - Determines whether the compressed row matrix format should be used. If 873e7a558SHong Zhang the format is to be used, this routine creates Mat_CompressedRow struct. 973e7a558SHong Zhang 1073e7a558SHong Zhang Compressed row format provides high performance routines by 1173e7a558SHong Zhang taking advantage of zero rows. 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 1973e7a558SHong Zhang - ratio - ratio of (num of zero rows)/m, used to determine if the compressed row format should be used 2073e7a558SHong Zhang 2173e7a558SHong Zhang Level: developer 2273e7a558SHong Zhang @*/ 2373e7a558SHong Zhang PetscErrorCode Mat_CheckCompressedRow(Mat A,Mat_CompressedRow *compressedrow,PetscInt *ai,PetscReal ratio) 2473e7a558SHong Zhang { 2573e7a558SHong Zhang PetscErrorCode ierr; 26*26e093fcSHong Zhang PetscInt nrows,*cpi=PETSC_NULL,*rindex=PETSC_NULL,nz,i,row,m=A->m/A->bs; 2773e7a558SHong Zhang 2873e7a558SHong Zhang PetscFunctionBegin; 2973e7a558SHong Zhang compressedrow->checked = PETSC_TRUE; 3073e7a558SHong Zhang 3173e7a558SHong Zhang /* compute number of zero rows */ 3273e7a558SHong Zhang nrows = 0; 3373e7a558SHong Zhang for (i=0; i<m; i++){ /* for each row */ 3473e7a558SHong Zhang nz = ai[i+1] - ai[i]; /* number of nonzeros */ 3573e7a558SHong Zhang if (nz == 0) nrows++; 3673e7a558SHong Zhang } 3773e7a558SHong Zhang /* if enough zero rows are found, use compressedrow data structure */ 3873e7a558SHong Zhang if (nrows < ratio*m) { 3973e7a558SHong Zhang compressedrow->use = PETSC_FALSE; 40*26e093fcSHong Zhang PetscLogInfo(A,"Mat_CheckCompressedRow: Found the ratio (num_zerorows %d)/(num_localrows %d) < %g. Do not use CompressedRow routines.\n",nrows,m,ratio); 4173e7a558SHong Zhang } else { 4273e7a558SHong Zhang compressedrow->use = PETSC_TRUE; 4373e7a558SHong Zhang PetscLogInfo(A,"Mat_CheckCompressedRow: Found the ratio (num_zerorows %d)/(num_localrows %d) > %g. Use CompressedRow routines.\n",nrows,m,ratio); 4473e7a558SHong Zhang 4573e7a558SHong Zhang /* set compressed row format */ 4673e7a558SHong Zhang nrows = m - nrows; /* num of non-zero rows */ 4773e7a558SHong Zhang ierr = PetscMalloc((2*nrows+1)*sizeof(PetscInt),&cpi);CHKERRQ(ierr); 4873e7a558SHong Zhang rindex = cpi + nrows + 1; 4973e7a558SHong Zhang row = 0; 5073e7a558SHong Zhang cpi[0] = 0; 5173e7a558SHong Zhang for (i=0; i<m; i++){ 5273e7a558SHong Zhang nz = ai[i+1] - ai[i]; 5373e7a558SHong Zhang if (nz == 0) continue; 5473e7a558SHong Zhang cpi[row+1] = ai[i+1]; /* compressed row pointer */ 5573e7a558SHong Zhang rindex[row++] = i; /* compressed row local index */ 5673e7a558SHong Zhang } 5773e7a558SHong Zhang compressedrow->nrows = nrows; 5873e7a558SHong Zhang compressedrow->i = cpi; 5973e7a558SHong Zhang compressedrow->rindex = rindex; 6073e7a558SHong Zhang } 6173e7a558SHong Zhang PetscFunctionReturn(0); 6273e7a558SHong Zhang } 63