#define PETSCMAT_DLL #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ #undef __FUNCT__ #define __FUNCT__ "Mat_CheckCompressedRow" /*@C Mat_CheckCompressedRow - Determines whether the compressed row matrix format should be used. If the format is to be used, this routine creates Mat_CompressedRow struct. Compressed row format provides high performance routines by taking advantage of zero rows. Supported types are MATAIJ, MATBAIJ and MATSBAIJ. Collective Input Parameters: + A - the matrix . compressedrow - pointer to the struct Mat_CompressedRow . ai - row pointer used by seqaij and seqbaij . mbs - number of (block) rows represented by ai - ratio - ratio of (num of zero rows)/m, used to determine if the compressed row format should be used Level: developer @*/ PetscErrorCode Mat_CheckCompressedRow(Mat A,Mat_CompressedRow *compressedrow,PetscInt *ai,PetscInt mbs,PetscReal ratio) { PetscErrorCode ierr; PetscInt nrows,*cpi=PETSC_NULL,*ridx=PETSC_NULL,nz,i,row; PetscFunctionBegin; if (!compressedrow->use) PetscFunctionReturn(0); if (compressedrow->checked){ if (!A->same_nonzero){ ierr = PetscFree(compressedrow->i);CHKERRQ(ierr); compressedrow->rindex = PETSC_NULL; ierr = PetscLogInfo((A,"Mat_CheckCompressedRow: Mat structure might be changed. Free memory and recheck.\n"));CHKERRQ(ierr); } else if (compressedrow->i == PETSC_NULL) { /* Don't know why this occures. For safe, recheck. */ ierr = PetscLogInfo((A,"Mat_CheckCompressedRow: compressedrow.checked, but compressedrow.i==null. Recheck.\n"));CHKERRQ(ierr); } else { /* use compressedrow, checked, A->same_nonzero = PETSC_TRUE. Skip check */ ierr = PetscLogInfo((A,"Mat_CheckCompressedRow: Skip check. m: %d, n: %d,M: %d, N: %d,nrows: %d, ii: %p, type: %s\n",A->m,A->n,A->M,A->N,compressedrow->nrows,compressedrow->i,A->type_name));CHKERRQ(ierr); PetscFunctionReturn(0); } } compressedrow->checked = PETSC_TRUE; /* compute number of zero rows */ nrows = 0; for (i=0; iuse = PETSC_FALSE; ierr = PetscLogInfo((A,"Mat_CheckCompressedRow: Found the ratio (num_zerorows %d)/(num_localrows %d) < %g. Do not use CompressedRow routines.\n",nrows,mbs,ratio));CHKERRQ(ierr); } else { compressedrow->use = PETSC_TRUE; ierr = PetscLogInfo((A,"Mat_CheckCompressedRow: Found the ratio (num_zerorows %d)/(num_localrows %d) > %g. Use CompressedRow routines.\n",nrows,mbs,ratio));CHKERRQ(ierr); /* set compressed row format */ nrows = mbs - nrows; /* num of non-zero rows */ ierr = PetscMalloc((2*nrows+1)*sizeof(PetscInt),&cpi);CHKERRQ(ierr); ridx = cpi + nrows + 1; row = 0; cpi[0] = 0; for (i=0; inrows = nrows; compressedrow->i = cpi; compressedrow->rindex = ridx; } PetscFunctionReturn(0); }