fixes for MatILUDTFactor_SeqAIJ
Lisandro Dalcin
dalcinl at gmail.com
Tue Sep 15 18:27:06 CDT 2009
Can someone review the patch below?
Perhaps the dt & dtcol initialization should also go to PCCreate_ILU ?
$ hg diff
diff -r 3af956418390 src/mat/impls/aij/seq/aijfact.c
--- a/src/mat/impls/aij/seq/aijfact.c Mon Sep 14 22:53:57 2009 -0500
+++ b/src/mat/impls/aij/seq/aijfact.c Tue Sep 15 20:18:12 2009 -0300
@@ -2245,18 +2245,22 @@
const PetscInt *r,*ic;
PetscInt i,n=A->rmap->n,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
PetscInt *bi,*bj,*bdiag,*bdiag_rev;
- PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
+ PetscInt row,nzi,nzi_bl,nzi_bu,*im,nzi_al,nzi_au;
PetscInt nlnk,*lnk;
PetscBT lnkbt;
PetscTruth row_identity,icol_identity,both_identity;
MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,multiplier,*vtmp,diag_tmp;
const PetscInt *ics;
PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp;
- PetscReal dt=info->dt,shift=info->shiftinblocks;
- PetscInt nnz_max;
+ PetscReal dt=info->dt,dtcol=info->dtcol,shift=info->shiftinblocks;
+ PetscInt dtcount=(PetscInt)info->dtcount,nnz_max;
PetscTruth missing;
PetscFunctionBegin;
+
+ if (dt == PETSC_DEFAULT) dt = 0.005;
+ if (dtcol == PETSC_DEFAULT) dtcol = 0.01; /* XXX unused! */
+
/* ------- symbolic factorization, can be reused ---------*/
ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing
diagonal entry %D",i);
@@ -2272,10 +2276,10 @@
ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
/* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
- dtcount = (PetscInt)info->dtcount;
+ if (dtcount == PETSC_DEFAULT) dtcount = (PetscInt)(1.5*a->rmax);
if (dtcount > n-1) dtcount = n-1; /* diagonal is excluded */
nnz_max = ai[n]+2*n*dtcount+2;
-
+
ierr = PetscMalloc((nnz_max+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
ierr = PetscMalloc((nnz_max+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
--
Lisandro Dalcín
---------------
Centro Internacional de Métodos Computacionales en Ingeniería (CIMEC)
Instituto de Desarrollo Tecnológico para la Industria Química (INTEC)
Consejo Nacional de Investigaciones Científicas y Técnicas (CONICET)
PTLC - Güemes 3450, (3000) Santa Fe, Argentina
Tel/Fax: +54-(0)342-451.1594
More information about the petsc-dev
mailing list