1be1d678aSKris Buschelman #define PETSCMAT_DLL 2289bc588SBarry Smith 370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 41c4feecaSSatish Balay #include "src/inline/dot.h" 5ed480e8bSBarry Smith #include "src/inline/spops.h" 61393bc97SHong Zhang #include "petscbt.h" 71393bc97SHong Zhang #include "src/mat/utils/freespace.h" 83b2fbd54SBarry Smith 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ" 11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol) 123b2fbd54SBarry Smith { 133a40ed3dSBarry Smith PetscFunctionBegin; 143a40ed3dSBarry Smith 1529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not written"); 16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG) 17d64ed03dSBarry Smith PetscFunctionReturn(0); 18d64ed03dSBarry Smith #endif 193b2fbd54SBarry Smith } 203b2fbd54SBarry Smith 2186bacbd2SBarry Smith 22e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE) 23a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 24a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscErrorCode*); 25a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,MatScalar*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*); 26e34fafa9SBarry Smith #endif 2786bacbd2SBarry Smith 284a2ae208SSatish Balay #undef __FUNCT__ 294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 3086bacbd2SBarry Smith /* ------------------------------------------------------------ 3186bacbd2SBarry Smith 3286bacbd2SBarry Smith This interface was contribed by Tony Caola 3386bacbd2SBarry Smith 3486bacbd2SBarry Smith This routine is an interface to the pivoting drop-tolerance 355bd2ddc7SBarry Smith ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of 365bd2ddc7SBarry Smith SPARSEKIT2. 375bd2ddc7SBarry Smith 385bd2ddc7SBarry Smith The SPARSEKIT2 routines used here are covered by the GNU 395bd2ddc7SBarry Smith copyright; see the file gnu in this directory. 4086bacbd2SBarry Smith 4186bacbd2SBarry Smith Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their 4286bacbd2SBarry Smith help in getting this routine ironed out. 4386bacbd2SBarry Smith 445bd2ddc7SBarry Smith The major drawback to this routine is that if info->fill is 455bd2ddc7SBarry Smith not large enough it fails rather than allocating more space; 465bd2ddc7SBarry Smith this can be fixed by hacking/improving the f2c version of 475bd2ddc7SBarry Smith Yousef Saad's code. 4886bacbd2SBarry Smith 4986bacbd2SBarry Smith ------------------------------------------------------------ 5086bacbd2SBarry Smith */ 517529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 5286bacbd2SBarry Smith { 5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE) 5460ec86bdSBarry Smith PetscFunctionBegin; 55e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\ 5660ec86bdSBarry Smith You can obtain the drop tolerance routines by installing PETSc from\n\ 5760ec86bdSBarry Smith www.mcs.anl.gov/petsc\n"); 5860ec86bdSBarry Smith #else 5986bacbd2SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 6007d2056aSBarry Smith IS iscolf,isicol,isirow; 6186bacbd2SBarry Smith PetscTruth reorder; 629cc1f4e3SBarry Smith PetscErrorCode ierr,sierr; 63899cda47SBarry Smith PetscInt *c,*r,*ic,i,n = A->rmap.n; 6438baddfdSBarry Smith PetscInt *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j; 6538baddfdSBarry Smith PetscInt *ordcol,*iwk,*iperm,*jw; 6638baddfdSBarry Smith PetscInt jmax,lfill,job,*o_i,*o_j; 6754f21887SBarry Smith MatScalar *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a; 6854a8161fSBarry Smith PetscReal af; 6986bacbd2SBarry Smith 7086bacbd2SBarry Smith PetscFunctionBegin; 7186bacbd2SBarry Smith 7286bacbd2SBarry Smith if (info->dt == PETSC_DEFAULT) info->dt = .005; 7338baddfdSBarry Smith if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax); 7486bacbd2SBarry Smith if (info->dtcol == PETSC_DEFAULT) info->dtcol = .01; 7533259db3SBarry Smith if (info->fill == PETSC_DEFAULT) info->fill = ((double)(n*(info->dtcount+1)))/a->nz; 7638baddfdSBarry Smith lfill = (PetscInt)(info->dtcount/2.0); 7738baddfdSBarry Smith jmax = (PetscInt)(info->fill*a->nz); 7886bacbd2SBarry Smith 7986bacbd2SBarry Smith 8086bacbd2SBarry Smith /* ------------------------------------------------------------ 8186bacbd2SBarry Smith If reorder=.TRUE., then the original matrix has to be 8286bacbd2SBarry Smith reordered to reflect the user selected ordering scheme, and 8386bacbd2SBarry Smith then de-reordered so it is in it's original format. 8486bacbd2SBarry Smith Because Saad's dperm() is NOT in place, we have to copy 8586bacbd2SBarry Smith the original matrix and allocate more storage. . . 8686bacbd2SBarry Smith ------------------------------------------------------------ 8786bacbd2SBarry Smith */ 8886bacbd2SBarry Smith 8986bacbd2SBarry Smith /* set reorder to true if either isrow or iscol is not identity */ 9086bacbd2SBarry Smith ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr); 9186bacbd2SBarry Smith if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);} 9286bacbd2SBarry Smith reorder = PetscNot(reorder); 9386bacbd2SBarry Smith 9486bacbd2SBarry Smith 9586bacbd2SBarry Smith /* storage for ilu factor */ 9638baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr); 9738baddfdSBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr); 98a77337e4SBarry Smith ierr = PetscMalloc(jmax*sizeof(MatScalar),&new_a);CHKERRQ(ierr); 9938baddfdSBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr); 10086bacbd2SBarry Smith 10186bacbd2SBarry Smith /* ------------------------------------------------------------ 10286bacbd2SBarry Smith Make sure that everything is Fortran formatted (1-Based) 10386bacbd2SBarry Smith ------------------------------------------------------------ 10486bacbd2SBarry Smith */ 105b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 106b19713cbSBarry Smith old_j[i]++; 10786bacbd2SBarry Smith } 108b19713cbSBarry Smith for(i=0;i<n+1;i++) { 109b19713cbSBarry Smith old_i[i]++; 110b19713cbSBarry Smith }; 111010a20caSHong Zhang 11286bacbd2SBarry Smith 11386bacbd2SBarry Smith if (reorder) { 114c0c2c81eSBarry Smith ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr); 115c0c2c81eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 116c0c2c81eSBarry Smith for(i=0;i<n;i++) { 117c0c2c81eSBarry Smith r[i] = r[i]+1; 118c0c2c81eSBarry Smith c[i] = c[i]+1; 119c0c2c81eSBarry Smith } 12038baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr); 12138baddfdSBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr); 122a77337e4SBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(MatScalar),&old_a2);CHKERRQ(ierr); 1235bd2ddc7SBarry Smith job = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job); 124c0c2c81eSBarry Smith for (i=0;i<n;i++) { 125c0c2c81eSBarry Smith r[i] = r[i]-1; 126c0c2c81eSBarry Smith c[i] = c[i]-1; 127c0c2c81eSBarry Smith } 128c0c2c81eSBarry Smith ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 129c0c2c81eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 130b19713cbSBarry Smith o_a = old_a2; 131b19713cbSBarry Smith o_j = old_j2; 132b19713cbSBarry Smith o_i = old_i2; 133b19713cbSBarry Smith } else { 134b19713cbSBarry Smith o_a = old_a; 135b19713cbSBarry Smith o_j = old_j; 136b19713cbSBarry Smith o_i = old_i; 13786bacbd2SBarry Smith } 13886bacbd2SBarry Smith 13986bacbd2SBarry Smith /* ------------------------------------------------------------ 14086bacbd2SBarry Smith Call Saad's ilutp() routine to generate the factorization 14186bacbd2SBarry Smith ------------------------------------------------------------ 14286bacbd2SBarry Smith */ 14386bacbd2SBarry Smith 14438baddfdSBarry Smith ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr); 14538baddfdSBarry Smith ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr); 14687828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr); 14786bacbd2SBarry Smith 14854a8161fSBarry Smith SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&info->dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr); 1496849ba73SBarry Smith if (sierr) { 1506849ba73SBarry Smith switch (sierr) { 151a83599f4SBarry Smith case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax); 152a83599f4SBarry Smith case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax); 153e005ede5SBarry Smith case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered"); 154e005ede5SBarry Smith case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong"); 15577431f27SBarry Smith case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax); 15677431f27SBarry Smith default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr); 15786bacbd2SBarry Smith } 15886bacbd2SBarry Smith } 15986bacbd2SBarry Smith 16086bacbd2SBarry Smith ierr = PetscFree(w);CHKERRQ(ierr); 16186bacbd2SBarry Smith ierr = PetscFree(jw);CHKERRQ(ierr); 16286bacbd2SBarry Smith 16386bacbd2SBarry Smith /* ------------------------------------------------------------ 16486bacbd2SBarry Smith Saad's routine gives the result in Modified Sparse Row (msr) 16586bacbd2SBarry Smith Convert to Compressed Sparse Row format (csr) 16686bacbd2SBarry Smith ------------------------------------------------------------ 16786bacbd2SBarry Smith */ 16886bacbd2SBarry Smith 16987828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr); 17038baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr); 17186bacbd2SBarry Smith 1725bd2ddc7SBarry Smith SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk); 17386bacbd2SBarry Smith 17486bacbd2SBarry Smith ierr = PetscFree(iwk);CHKERRQ(ierr); 17586bacbd2SBarry Smith ierr = PetscFree(wk);CHKERRQ(ierr); 17686bacbd2SBarry Smith 17786bacbd2SBarry Smith if (reorder) { 17886bacbd2SBarry Smith ierr = PetscFree(old_a2);CHKERRQ(ierr); 17986bacbd2SBarry Smith ierr = PetscFree(old_j2);CHKERRQ(ierr); 18086bacbd2SBarry Smith ierr = PetscFree(old_i2);CHKERRQ(ierr); 181313ae024SBarry Smith } else { 182b19713cbSBarry Smith /* fix permutation of old_j that the factorization introduced */ 183f170005cSBarry Smith for (i=old_i[0]; i<old_i[n]; i++) { 184b19713cbSBarry Smith old_j[i-1] = iperm[old_j[i-1]-1]; 185b19713cbSBarry Smith } 186b19713cbSBarry Smith } 187b19713cbSBarry Smith 188b801d0f9SBarry Smith /* get rid of the shift to indices starting at 1 */ 18986bacbd2SBarry Smith for (i=0; i<n+1; i++) { 190b19713cbSBarry Smith old_i[i]--; 191b19713cbSBarry Smith } 192b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 193b19713cbSBarry Smith old_j[i]--; 194b19713cbSBarry Smith } 19586bacbd2SBarry Smith 196b19713cbSBarry Smith /* Make the factored matrix 0-based */ 19786bacbd2SBarry Smith for (i=0; i<n+1; i++) { 198b19713cbSBarry Smith new_i[i]--; 199b19713cbSBarry Smith } 200b19713cbSBarry Smith for (i=new_i[0];i<new_i[n];i++) { 201b19713cbSBarry Smith new_j[i]--; 202b19713cbSBarry Smith } 20386bacbd2SBarry Smith 20486bacbd2SBarry Smith /*-- due to the pivoting, we need to reorder iscol to correctly --*/ 20586bacbd2SBarry Smith /*-- permute the right-hand-side and solution vectors --*/ 2064c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2074c49b128SBarry Smith ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr); 208c0c2c81eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 20986bacbd2SBarry Smith for(i=0; i<n; i++) { 21086bacbd2SBarry Smith ordcol[i] = ic[iperm[i]-1]; 21186bacbd2SBarry Smith }; 21286bacbd2SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 21307d2056aSBarry Smith ierr = ISDestroy(isicol);CHKERRQ(ierr); 21486bacbd2SBarry Smith 215c0c2c81eSBarry Smith ierr = PetscFree(iperm);CHKERRQ(ierr); 216c0c2c81eSBarry Smith 217329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr); 21807d2056aSBarry Smith ierr = PetscFree(ordcol);CHKERRQ(ierr); 21986bacbd2SBarry Smith 22086bacbd2SBarry Smith /*----- put together the new matrix -----*/ 22186bacbd2SBarry Smith 2227adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 223f69a0ea3SMatthew Knepley ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr); 2247adad957SLisandro Dalcin ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 225ab93d7beSBarry Smith ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 22686bacbd2SBarry Smith (*fact)->factor = FACTOR_LU; 22786bacbd2SBarry Smith (*fact)->assembled = PETSC_TRUE; 22886bacbd2SBarry Smith 22986bacbd2SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 230e6b907acSBarry Smith b->free_a = PETSC_TRUE; 231e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 23207d2056aSBarry Smith b->singlemalloc = PETSC_FALSE; 23386bacbd2SBarry Smith b->a = new_a; 23486bacbd2SBarry Smith b->j = new_j; 23586bacbd2SBarry Smith b->i = new_i; 23686bacbd2SBarry Smith b->ilen = 0; 23786bacbd2SBarry Smith b->imax = 0; 2381f9e874aSBarry Smith /* I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */ 239313ae024SBarry Smith b->row = isirow; 24086bacbd2SBarry Smith b->col = iscolf; 24187828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 24286bacbd2SBarry Smith b->maxnz = b->nz = new_i[n]; 24386bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 24486bacbd2SBarry Smith (*fact)->info.factor_mallocs = 0; 24586bacbd2SBarry Smith 246431e710aSBarry Smith af = ((double)b->nz)/((double)a->nz) + .001; 247ae15b995SBarry Smith ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",info->fill,af);CHKERRQ(ierr); 248ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 249ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 250ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 251431e710aSBarry Smith 2527529efd4SKris Buschelman ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr); 25371c2f376SKris Buschelman 25486bacbd2SBarry Smith PetscFunctionReturn(0); 25560ec86bdSBarry Smith #endif 25686bacbd2SBarry Smith } 25786bacbd2SBarry Smith 2584a2ae208SSatish Balay #undef __FUNCT__ 259*b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_petsc" 260*b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B) 261*b24902e0SBarry Smith { 262*b24902e0SBarry Smith PetscInt n = A->rmap.n; 263*b24902e0SBarry Smith PetscErrorCode ierr; 264*b24902e0SBarry Smith 265*b24902e0SBarry Smith PetscFunctionBegin; 266*b24902e0SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 267*b24902e0SBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 268*b24902e0SBarry Smith if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 269*b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr); 270*b24902e0SBarry Smith } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 271*b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 272*b24902e0SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 273*b24902e0SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 274*b24902e0SBarry Smith PetscFunctionReturn(0); 275*b24902e0SBarry Smith } 276*b24902e0SBarry Smith 277*b24902e0SBarry Smith #undef __FUNCT__ 278b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ" 279dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B) 280289bc588SBarry Smith { 281416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 282289bc588SBarry Smith IS isicol; 2836849ba73SBarry Smith PetscErrorCode ierr; 284899cda47SBarry Smith PetscInt *r,*ic,i,n=A->rmap.n,*ai=a->i,*aj=a->j; 2851393bc97SHong Zhang PetscInt *bi,*bj,*ajtmp; 2861393bc97SHong Zhang PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 2879e046878SBarry Smith PetscReal f; 2884875ba38SHong Zhang PetscInt nlnk,*lnk,k,**bi_ptr; 289a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2901393bc97SHong Zhang PetscBT lnkbt; 291289bc588SBarry Smith 2923a40ed3dSBarry Smith PetscFunctionBegin; 293899cda47SBarry Smith if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 2944c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 295ac121b3dSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 296ac121b3dSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 297289bc588SBarry Smith 298289bc588SBarry Smith /* get new row pointers */ 2991393bc97SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 3001393bc97SHong Zhang bi[0] = 0; 3011393bc97SHong Zhang 3021393bc97SHong Zhang /* bdiag is location of diagonal in factor */ 3031393bc97SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 3041393bc97SHong Zhang bdiag[0] = 0; 3051393bc97SHong Zhang 3061393bc97SHong Zhang /* linked list for storing column indices of the active row */ 3071393bc97SHong Zhang nlnk = n + 1; 3081393bc97SHong Zhang ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 3091393bc97SHong Zhang 31035baf27eSSatish Balay ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 3111393bc97SHong Zhang 3121393bc97SHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 313d3d32019SBarry Smith f = info->fill; 314a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 3151393bc97SHong Zhang current_space = free_space; 316289bc588SBarry Smith 317289bc588SBarry Smith for (i=0; i<n; i++) { 3181393bc97SHong Zhang /* copy previous fill into linked list */ 319289bc588SBarry Smith nzi = 0; 3201393bc97SHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 3211393bc97SHong Zhang if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 3221393bc97SHong Zhang ajtmp = aj + ai[r[i]]; 323afcc9904SHong Zhang ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 3241393bc97SHong Zhang nzi += nlnk; 3251393bc97SHong Zhang 3261393bc97SHong Zhang /* add pivot rows into linked list */ 3271393bc97SHong Zhang row = lnk[n]; 3281393bc97SHong Zhang while (row < i) { 329d1170560SHong Zhang nzbd = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */ 330d1170560SHong Zhang ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 331d1170560SHong Zhang ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 3321393bc97SHong Zhang nzi += nlnk; 3331393bc97SHong Zhang row = lnk[row]; 334289bc588SBarry Smith } 3351393bc97SHong Zhang bi[i+1] = bi[i] + nzi; 3361393bc97SHong Zhang im[i] = nzi; 3371393bc97SHong Zhang 3381393bc97SHong Zhang /* mark bdiag */ 3391393bc97SHong Zhang nzbd = 0; 3401393bc97SHong Zhang nnz = nzi; 3411393bc97SHong Zhang k = lnk[n]; 3421393bc97SHong Zhang while (nnz-- && k < i){ 3431393bc97SHong Zhang nzbd++; 3441393bc97SHong Zhang k = lnk[k]; 3451393bc97SHong Zhang } 3461393bc97SHong Zhang bdiag[i] = bi[i] + nzbd; 3471393bc97SHong Zhang 3481393bc97SHong Zhang /* if free space is not available, make more free space */ 3491393bc97SHong Zhang if (current_space->local_remaining<nzi) { 3501393bc97SHong Zhang nnz = (n - i)*nzi; /* estimated and max additional space needed */ 351a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 3521393bc97SHong Zhang reallocs++; 3531393bc97SHong Zhang } 3541393bc97SHong Zhang 3551393bc97SHong Zhang /* copy data into free space, then initialize lnk */ 3561393bc97SHong Zhang ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 3571393bc97SHong Zhang bi_ptr[i] = current_space->array; 3581393bc97SHong Zhang current_space->array += nzi; 3591393bc97SHong Zhang current_space->local_used += nzi; 3601393bc97SHong Zhang current_space->local_remaining -= nzi; 361289bc588SBarry Smith } 3626cf91177SBarry Smith #if defined(PETSC_USE_INFO) 363464e8d44SSatish Balay if (ai[n] != 0) { 3641393bc97SHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 365ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 366ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 367ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 368ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 36905bf559cSSatish Balay } else { 370ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 37105bf559cSSatish Balay } 37263ba0a88SBarry Smith #endif 3732fb3ed76SBarry Smith 374898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 375898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3761987afe7SBarry Smith 3771393bc97SHong Zhang /* destroy list of free space and other temporary array(s) */ 3781393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 379a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 3801393bc97SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 38135baf27eSSatish Balay ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 382289bc588SBarry Smith 383289bc588SBarry Smith /* put together the new matrix */ 38452e6d16bSBarry Smith ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr); 385416022c9SBarry Smith b = (Mat_SeqAIJ*)(*B)->data; 386e6b907acSBarry Smith b->free_a = PETSC_TRUE; 387e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 3887c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 3891393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 3901393bc97SHong Zhang b->j = bj; 3911393bc97SHong Zhang b->i = bi; 3921393bc97SHong Zhang b->diag = bdiag; 393416022c9SBarry Smith b->ilen = 0; 394416022c9SBarry Smith b->imax = 0; 395416022c9SBarry Smith b->row = isrow; 396416022c9SBarry Smith b->col = iscol; 397c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 398c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 39982bf6240SBarry Smith b->icol = isicol; 40087828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 4011393bc97SHong Zhang 4021393bc97SHong Zhang /* In b structure: Free imax, ilen, old a, old j. Allocate solve_work, new a, new j */ 4031393bc97SHong Zhang ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 4041393bc97SHong Zhang b->maxnz = b->nz = bi[n] ; 4057fd98bd6SLois Curfman McInnes 406329f5518SBarry Smith (*B)->factor = FACTOR_LU; 407418422e8SSatish Balay (*B)->info.factor_mallocs = reallocs; 408ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 409703deb49SSatish Balay 410e93ccc4dSBarry Smith if (ai[n] != 0) { 4111393bc97SHong Zhang (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 412eea2913aSSatish Balay } else { 413eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 414eea2913aSSatish Balay } 4154846f1f5SKris Buschelman ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr); 41671c2f376SKris Buschelman (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 4173a40ed3dSBarry Smith PetscFunctionReturn(0); 418289bc588SBarry Smith } 4191393bc97SHong Zhang 4205b5bc046SBarry Smith /* 4215b5bc046SBarry Smith Trouble in factorization, should we dump the original matrix? 4225b5bc046SBarry Smith */ 4235b5bc046SBarry Smith #undef __FUNCT__ 4245b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix" 4255b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A) 4265b5bc046SBarry Smith { 4275b5bc046SBarry Smith PetscErrorCode ierr; 4285b5bc046SBarry Smith PetscTruth flg; 4295b5bc046SBarry Smith 4305b5bc046SBarry Smith PetscFunctionBegin; 4315b5bc046SBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_factor_dump_on_error",&flg);CHKERRQ(ierr); 4325b5bc046SBarry Smith if (flg) { 4335b5bc046SBarry Smith PetscViewer viewer; 4345b5bc046SBarry Smith char filename[PETSC_MAX_PATH_LEN]; 4355b5bc046SBarry Smith 4365b5bc046SBarry Smith ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr); 4377adad957SLisandro Dalcin ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 4385b5bc046SBarry Smith ierr = MatView(A,viewer);CHKERRQ(ierr); 4395b5bc046SBarry Smith ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 4405b5bc046SBarry Smith } 4415b5bc046SBarry Smith PetscFunctionReturn(0); 4425b5bc046SBarry Smith } 4435b5bc046SBarry Smith 4440a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 4454a2ae208SSatish Balay #undef __FUNCT__ 4464a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 447af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 448289bc588SBarry Smith { 449416022c9SBarry Smith Mat C=*B; 450416022c9SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 45182bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 4526849ba73SBarry Smith PetscErrorCode ierr; 453899cda47SBarry Smith PetscInt *r,*ic,i,j,n=A->rmap.n,*bi=b->i,*bj=b->j; 454b3bf805bSHong Zhang PetscInt *ajtmp,*bjtmp,nz,row,*ics; 455d2276718SHong Zhang PetscInt *diag_offset = b->diag,diag,*pj; 45654f21887SBarry Smith PetscScalar *rtmp,*pc,multiplier,*rtmps; 45754f21887SBarry Smith MatScalar *v,*pv; 4586a7c8fc2SHong Zhang PetscScalar d; 4596a7c8fc2SHong Zhang PetscReal rs; 460b3bf805bSHong Zhang LUShift_Ctx sctx; 46142898933SHong Zhang PetscInt newshift,*ddiag; 462289bc588SBarry Smith 4633a40ed3dSBarry Smith PetscFunctionBegin; 46478b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 46578b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 46687828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 46787828ca2SBarry Smith ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 468010a20caSHong Zhang rtmps = rtmp; ics = ic; 469289bc588SBarry Smith 470ada7143aSSatish Balay sctx.shift_top = 0; 471ada7143aSSatish Balay sctx.nshift_max = 0; 472ada7143aSSatish Balay sctx.shift_lo = 0; 473ada7143aSSatish Balay sctx.shift_hi = 0; 474ada7143aSSatish Balay 475118f5789SBarry Smith /* if both shift schemes are chosen by user, only use info->shiftpd */ 476118f5789SBarry Smith if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0; 4771a890ab1SHong Zhang if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 4789deb78dcSHong Zhang PetscInt *aai = a->i; 47942898933SHong Zhang ddiag = a->diag; 480118f5789SBarry Smith sctx.shift_top = 0; 4816cc28720Svictorle for (i=0; i<n; i++) { 4829f95998fSHong Zhang /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 4839f95998fSHong Zhang d = (a->a)[ddiag[i]]; 4841a890ab1SHong Zhang rs = -PetscAbsScalar(d) - PetscRealPart(d); 485010a20caSHong Zhang v = a->a+aai[i]; 486e24cfe64SHong Zhang nz = aai[i+1] - aai[i]; 487e24cfe64SHong Zhang for (j=0; j<nz; j++) 4881a890ab1SHong Zhang rs += PetscAbsScalar(v[j]); 4891a890ab1SHong Zhang if (rs>sctx.shift_top) sctx.shift_top = rs; 4906cc28720Svictorle } 491c3af1dc1SHong Zhang if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot; 492118f5789SBarry Smith sctx.shift_top *= 1.1; 493d2276718SHong Zhang sctx.nshift_max = 5; 494d2276718SHong Zhang sctx.shift_lo = 0.; 495d2276718SHong Zhang sctx.shift_hi = 1.; 496a0a356efSHong Zhang } 497a0a356efSHong Zhang 498a0a356efSHong Zhang sctx.shift_amount = 0; 499a0a356efSHong Zhang sctx.nshift = 0; 500085256b3SBarry Smith do { 501d2276718SHong Zhang sctx.lushift = PETSC_FALSE; 502289bc588SBarry Smith for (i=0; i<n; i++){ 503b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 504b3bf805bSHong Zhang bjtmp = bj + bi[i]; 505b3bf805bSHong Zhang for (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0; 506289bc588SBarry Smith 507289bc588SBarry Smith /* load in initial (unfactored row) */ 508416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 509b3bf805bSHong Zhang ajtmp = a->j + a->i[r[i]]; 510010a20caSHong Zhang v = a->a + a->i[r[i]]; 511085256b3SBarry Smith for (j=0; j<nz; j++) { 512b3bf805bSHong Zhang rtmp[ics[ajtmp[j]]] = v[j]; 513085256b3SBarry Smith } 514d2276718SHong Zhang rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 515289bc588SBarry Smith 516b3bf805bSHong Zhang row = *bjtmp++; 517f2582111SSatish Balay while (row < i) { 5188c37ef55SBarry Smith pc = rtmp + row; 5198c37ef55SBarry Smith if (*pc != 0.0) { 520010a20caSHong Zhang pv = b->a + diag_offset[row]; 521010a20caSHong Zhang pj = b->j + diag_offset[row] + 1; 52235aab85fSBarry Smith multiplier = *pc / *pv++; 5238c37ef55SBarry Smith *pc = multiplier; 524b3bf805bSHong Zhang nz = bi[row+1] - diag_offset[row] - 1; 525f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 526efee365bSSatish Balay ierr = PetscLogFlops(2*nz);CHKERRQ(ierr); 527289bc588SBarry Smith } 528b3bf805bSHong Zhang row = *bjtmp++; 529289bc588SBarry Smith } 530416022c9SBarry Smith /* finished row so stick it into b->a */ 531b3bf805bSHong Zhang pv = b->a + bi[i] ; 532b3bf805bSHong Zhang pj = b->j + bi[i] ; 533b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 534b3bf805bSHong Zhang diag = diag_offset[i] - bi[i]; 535d3d32019SBarry Smith rs = 0.0; 536d3d32019SBarry Smith for (j=0; j<nz; j++) { 537d3d32019SBarry Smith pv[j] = rtmps[pj[j]]; 538d3d32019SBarry Smith if (j != diag) rs += PetscAbsScalar(pv[j]); 539d3d32019SBarry Smith } 540d2276718SHong Zhang 541b3bf805bSHong Zhang /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 542d2276718SHong Zhang sctx.rs = rs; 543d2276718SHong Zhang sctx.pv = pv[diag]; 544c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr); 5455b5bc046SBarry Smith if (newshift == 1) break; 54635aab85fSBarry Smith } 547d2276718SHong Zhang 548118f5789SBarry Smith if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 5496cc28720Svictorle /* 5506c037d1bSvictorle * if no shift in this attempt & shifting & started shifting & can refine, 5516cc28720Svictorle * then try lower shift 5526cc28720Svictorle */ 553d2276718SHong Zhang sctx.shift_hi = info->shift_fraction; 554d2276718SHong Zhang info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 555d2276718SHong Zhang sctx.shift_amount = info->shift_fraction * sctx.shift_top; 556d2276718SHong Zhang sctx.lushift = PETSC_TRUE; 557d2276718SHong Zhang sctx.nshift++; 5586cc28720Svictorle } 559d2276718SHong Zhang } while (sctx.lushift); 5600f11f9aeSSatish Balay 56135aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 56235aab85fSBarry Smith for (i=0; i<n; i++) { 563010a20caSHong Zhang b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]]; 56435aab85fSBarry Smith } 56535aab85fSBarry Smith 566606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 56778b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 56878b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 569416022c9SBarry Smith C->factor = FACTOR_LU; 5707dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 571c456f294SBarry Smith C->assembled = PETSC_TRUE; 572899cda47SBarry Smith ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr); 573d2276718SHong Zhang if (sctx.nshift){ 574118f5789SBarry Smith if (info->shiftnz) { 5751e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 576118f5789SBarry Smith } else if (info->shiftpd) { 5771e2582c4SBarry Smith ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr); 5786cc28720Svictorle } 5790697b6caSHong Zhang } 5803a40ed3dSBarry Smith PetscFunctionReturn(0); 581289bc588SBarry Smith } 5820889b5dcSvictorle 583137fb511SHong Zhang /* 584137fb511SHong Zhang This routine implements inplace ILU(0) with row or/and column permutations. 585137fb511SHong Zhang Input: 586137fb511SHong Zhang A - original matrix 587137fb511SHong Zhang Output; 588137fb511SHong Zhang A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i] 589137fb511SHong Zhang a->j (col index) is permuted by the inverse of colperm, then sorted 590137fb511SHong Zhang a->a reordered accordingly with a->j 591137fb511SHong Zhang a->diag (ptr to diagonal elements) is updated. 592137fb511SHong Zhang */ 593137fb511SHong Zhang #undef __FUNCT__ 594137fb511SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm" 595137fb511SHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat A,MatFactorInfo *info,Mat *B) 596137fb511SHong Zhang { 597b051339dSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 598b051339dSHong Zhang IS isrow = a->row,isicol = a->icol; 599137fb511SHong Zhang PetscErrorCode ierr; 600b051339dSHong Zhang PetscInt *r,*ic,i,j,n=A->rmap.n,*ai=a->i,*aj=a->j; 601b051339dSHong Zhang PetscInt *ajtmp,nz,row,*ics; 602b051339dSHong Zhang PetscInt *diag = a->diag,nbdiag,*pj; 603a77337e4SBarry Smith PetscScalar *rtmp,*pc,multiplier,d; 604a77337e4SBarry Smith MatScalar *v,*pv; 605137fb511SHong Zhang PetscReal rs; 606137fb511SHong Zhang LUShift_Ctx sctx; 607b051339dSHong Zhang PetscInt newshift; 608137fb511SHong Zhang 609137fb511SHong Zhang PetscFunctionBegin; 610b051339dSHong Zhang if (A != *B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address"); 611137fb511SHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 612137fb511SHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 613137fb511SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 614137fb511SHong Zhang ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 615b051339dSHong Zhang ics = ic; 616137fb511SHong Zhang 617137fb511SHong Zhang sctx.shift_top = 0; 618137fb511SHong Zhang sctx.nshift_max = 0; 619137fb511SHong Zhang sctx.shift_lo = 0; 620137fb511SHong Zhang sctx.shift_hi = 0; 621137fb511SHong Zhang 622137fb511SHong Zhang /* if both shift schemes are chosen by user, only use info->shiftpd */ 623137fb511SHong Zhang if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0; 624137fb511SHong Zhang if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 625137fb511SHong Zhang sctx.shift_top = 0; 626137fb511SHong Zhang for (i=0; i<n; i++) { 627137fb511SHong Zhang /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 628b051339dSHong Zhang d = (a->a)[diag[i]]; 629137fb511SHong Zhang rs = -PetscAbsScalar(d) - PetscRealPart(d); 630b051339dSHong Zhang v = a->a+ai[i]; 631b051339dSHong Zhang nz = ai[i+1] - ai[i]; 632137fb511SHong Zhang for (j=0; j<nz; j++) 633137fb511SHong Zhang rs += PetscAbsScalar(v[j]); 634137fb511SHong Zhang if (rs>sctx.shift_top) sctx.shift_top = rs; 635137fb511SHong Zhang } 636137fb511SHong Zhang if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot; 637137fb511SHong Zhang sctx.shift_top *= 1.1; 638137fb511SHong Zhang sctx.nshift_max = 5; 639137fb511SHong Zhang sctx.shift_lo = 0.; 640137fb511SHong Zhang sctx.shift_hi = 1.; 641137fb511SHong Zhang } 642137fb511SHong Zhang 643137fb511SHong Zhang sctx.shift_amount = 0; 644137fb511SHong Zhang sctx.nshift = 0; 645137fb511SHong Zhang do { 646137fb511SHong Zhang sctx.lushift = PETSC_FALSE; 647137fb511SHong Zhang for (i=0; i<n; i++){ 648b051339dSHong Zhang /* load in initial unfactored row */ 649b051339dSHong Zhang nz = ai[r[i]+1] - ai[r[i]]; 650b051339dSHong Zhang ajtmp = aj + ai[r[i]]; 651b051339dSHong Zhang v = a->a + ai[r[i]]; 652137fb511SHong Zhang /* sort permuted ajtmp and values v accordingly */ 653b051339dSHong Zhang for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]]; 654137fb511SHong Zhang ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr); 655137fb511SHong Zhang 656b051339dSHong Zhang diag[r[i]] = ai[r[i]]; 657137fb511SHong Zhang for (j=0; j<nz; j++) { 658137fb511SHong Zhang rtmp[ajtmp[j]] = v[j]; 659b051339dSHong Zhang if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */ 660137fb511SHong Zhang } 661137fb511SHong Zhang rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 662137fb511SHong Zhang 663b051339dSHong Zhang row = *ajtmp++; 664137fb511SHong Zhang while (row < i) { 665137fb511SHong Zhang pc = rtmp + row; 666137fb511SHong Zhang if (*pc != 0.0) { 667b051339dSHong Zhang pv = a->a + diag[r[row]]; 668b051339dSHong Zhang pj = aj + diag[r[row]] + 1; 669137fb511SHong Zhang 670137fb511SHong Zhang multiplier = *pc / *pv++; 671137fb511SHong Zhang *pc = multiplier; 672b051339dSHong Zhang nz = ai[r[row]+1] - diag[r[row]] - 1; 673b051339dSHong Zhang for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 674137fb511SHong Zhang ierr = PetscLogFlops(2*nz);CHKERRQ(ierr); 675137fb511SHong Zhang } 676b051339dSHong Zhang row = *ajtmp++; 677137fb511SHong Zhang } 678b051339dSHong Zhang /* finished row so overwrite it onto a->a */ 679b051339dSHong Zhang pv = a->a + ai[r[i]] ; 680b051339dSHong Zhang pj = aj + ai[r[i]] ; 681b051339dSHong Zhang nz = ai[r[i]+1] - ai[r[i]]; 682b051339dSHong Zhang nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */ 683137fb511SHong Zhang 684137fb511SHong Zhang rs = 0.0; 685137fb511SHong Zhang for (j=0; j<nz; j++) { 686b051339dSHong Zhang pv[j] = rtmp[pj[j]]; 687b051339dSHong Zhang if (j != nbdiag) rs += PetscAbsScalar(pv[j]); 688137fb511SHong Zhang } 689137fb511SHong Zhang 690137fb511SHong Zhang /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 691137fb511SHong Zhang sctx.rs = rs; 692b051339dSHong Zhang sctx.pv = pv[nbdiag]; 693c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr); 694137fb511SHong Zhang if (newshift == 1) break; 695137fb511SHong Zhang } 696137fb511SHong Zhang 697137fb511SHong Zhang if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 698137fb511SHong Zhang /* 699137fb511SHong Zhang * if no shift in this attempt & shifting & started shifting & can refine, 700137fb511SHong Zhang * then try lower shift 701137fb511SHong Zhang */ 702137fb511SHong Zhang sctx.shift_hi = info->shift_fraction; 703137fb511SHong Zhang info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 704137fb511SHong Zhang sctx.shift_amount = info->shift_fraction * sctx.shift_top; 705137fb511SHong Zhang sctx.lushift = PETSC_TRUE; 706137fb511SHong Zhang sctx.nshift++; 707137fb511SHong Zhang } 708137fb511SHong Zhang } while (sctx.lushift); 709137fb511SHong Zhang 710137fb511SHong Zhang /* invert diagonal entries for simplier triangular solves */ 711137fb511SHong Zhang for (i=0; i<n; i++) { 712b051339dSHong Zhang a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]]; 713137fb511SHong Zhang } 714137fb511SHong Zhang 715137fb511SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 716137fb511SHong Zhang ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 717137fb511SHong Zhang ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 718b051339dSHong Zhang A->factor = FACTOR_LU; 719b051339dSHong Zhang A->ops->solve = MatSolve_SeqAIJ_InplaceWithPerm; 720b051339dSHong Zhang A->assembled = PETSC_TRUE; 721b051339dSHong Zhang ierr = PetscLogFlops(A->cmap.n);CHKERRQ(ierr); 722137fb511SHong Zhang if (sctx.nshift){ 723137fb511SHong Zhang if (info->shiftnz) { 7241e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 725137fb511SHong Zhang } else if (info->shiftpd) { 7261e2582c4SBarry Smith ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr); 727137fb511SHong Zhang } 728137fb511SHong Zhang } 729137fb511SHong Zhang PetscFunctionReturn(0); 730137fb511SHong Zhang } 731137fb511SHong Zhang 7320a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 7334a2ae208SSatish Balay #undef __FUNCT__ 7344a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 735dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 736da3a660dSBarry Smith { 737dfbe8321SBarry Smith PetscErrorCode ierr; 738416022c9SBarry Smith Mat C; 739416022c9SBarry Smith 7403a40ed3dSBarry Smith PetscFunctionBegin; 741*b24902e0SBarry Smith ierr = MatGetFactor(A,"petsc",MAT_FACTOR_LU,&C);CHKERRQ(ierr); 7429e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 743af281ebdSHong Zhang ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 744273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 74552e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr); 7463a40ed3dSBarry Smith PetscFunctionReturn(0); 747da3a660dSBarry Smith } 7480a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 7494a2ae208SSatish Balay #undef __FUNCT__ 7504a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 751dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 7528c37ef55SBarry Smith { 753416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 754416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 7556849ba73SBarry Smith PetscErrorCode ierr; 756899cda47SBarry Smith PetscInt *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j; 75738baddfdSBarry Smith PetscInt nz,*rout,*cout; 758dd6ea824SBarry Smith PetscScalar *x,*tmp,*tmps,sum; 759d9fead3dSBarry Smith const PetscScalar *b; 760dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 7618c37ef55SBarry Smith 7623a40ed3dSBarry Smith PetscFunctionBegin; 7633a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 76491bf9adeSBarry Smith 765d9fead3dSBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 7661ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 767416022c9SBarry Smith tmp = a->solve_work; 7688c37ef55SBarry Smith 7693b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7703b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 7718c37ef55SBarry Smith 7728c37ef55SBarry Smith /* forward solve the lower triangular */ 7738c37ef55SBarry Smith tmp[0] = b[*r++]; 774010a20caSHong Zhang tmps = tmp; 7758c37ef55SBarry Smith for (i=1; i<n; i++) { 776010a20caSHong Zhang v = aa + ai[i] ; 777010a20caSHong Zhang vi = aj + ai[i] ; 77853e7425aSBarry Smith nz = a->diag[i] - ai[i]; 77953e7425aSBarry Smith sum = b[*r++]; 780ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 7818c37ef55SBarry Smith tmp[i] = sum; 7828c37ef55SBarry Smith } 7838c37ef55SBarry Smith 7848c37ef55SBarry Smith /* backward solve the upper triangular */ 7858c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 786010a20caSHong Zhang v = aa + a->diag[i] + 1; 787010a20caSHong Zhang vi = aj + a->diag[i] + 1; 788416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 7898c37ef55SBarry Smith sum = tmp[i]; 790ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 791010a20caSHong Zhang x[*c--] = tmp[i] = sum*aa[a->diag[i]]; 7928c37ef55SBarry Smith } 7938c37ef55SBarry Smith 7943b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7953b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 796d9fead3dSBarry Smith ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 7971ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 798899cda47SBarry Smith ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr); 7993a40ed3dSBarry Smith PetscFunctionReturn(0); 8008c37ef55SBarry Smith } 801026e39d0SSatish Balay 8022bebee5dSHong Zhang #undef __FUNCT__ 8032bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ" 8042bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X) 8052bebee5dSHong Zhang { 8062bebee5dSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 8072bebee5dSHong Zhang IS iscol = a->col,isrow = a->row; 8082bebee5dSHong Zhang PetscErrorCode ierr; 809899cda47SBarry Smith PetscInt *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j; 8102bebee5dSHong Zhang PetscInt nz,*rout,*cout,neq; 811dd6ea824SBarry Smith PetscScalar *x,*b,*tmp,*tmps,sum; 812dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 8132bebee5dSHong Zhang 8142bebee5dSHong Zhang PetscFunctionBegin; 8152bebee5dSHong Zhang if (!n) PetscFunctionReturn(0); 8162bebee5dSHong Zhang 8172bebee5dSHong Zhang ierr = MatGetArray(B,&b);CHKERRQ(ierr); 8182bebee5dSHong Zhang ierr = MatGetArray(X,&x);CHKERRQ(ierr); 8192bebee5dSHong Zhang 8202bebee5dSHong Zhang tmp = a->solve_work; 8212bebee5dSHong Zhang ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8222bebee5dSHong Zhang ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 8232bebee5dSHong Zhang 8242bebee5dSHong Zhang for (neq=0; neq<n; neq++){ 8252bebee5dSHong Zhang /* forward solve the lower triangular */ 8262bebee5dSHong Zhang tmp[0] = b[r[0]]; 8272bebee5dSHong Zhang tmps = tmp; 8282bebee5dSHong Zhang for (i=1; i<n; i++) { 8292bebee5dSHong Zhang v = aa + ai[i] ; 8302bebee5dSHong Zhang vi = aj + ai[i] ; 8312bebee5dSHong Zhang nz = a->diag[i] - ai[i]; 8322bebee5dSHong Zhang sum = b[r[i]]; 8332bebee5dSHong Zhang SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 8342bebee5dSHong Zhang tmp[i] = sum; 8352bebee5dSHong Zhang } 8362bebee5dSHong Zhang /* backward solve the upper triangular */ 8372bebee5dSHong Zhang for (i=n-1; i>=0; i--){ 8382bebee5dSHong Zhang v = aa + a->diag[i] + 1; 8392bebee5dSHong Zhang vi = aj + a->diag[i] + 1; 8402bebee5dSHong Zhang nz = ai[i+1] - a->diag[i] - 1; 8412bebee5dSHong Zhang sum = tmp[i]; 8422bebee5dSHong Zhang SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 8432bebee5dSHong Zhang x[c[i]] = tmp[i] = sum*aa[a->diag[i]]; 8442bebee5dSHong Zhang } 8452bebee5dSHong Zhang 8462bebee5dSHong Zhang b += n; 8472bebee5dSHong Zhang x += n; 8482bebee5dSHong Zhang } 8492bebee5dSHong Zhang ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8502bebee5dSHong Zhang ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 8512bebee5dSHong Zhang ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 8522bebee5dSHong Zhang ierr = MatRestoreArray(X,&x);CHKERRQ(ierr); 8532bebee5dSHong Zhang ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr); 8542bebee5dSHong Zhang PetscFunctionReturn(0); 8552bebee5dSHong Zhang } 8562bebee5dSHong Zhang 857137fb511SHong Zhang #undef __FUNCT__ 858137fb511SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm" 859137fb511SHong Zhang PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx) 860137fb511SHong Zhang { 861137fb511SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 862137fb511SHong Zhang IS iscol = a->col,isrow = a->row; 863137fb511SHong Zhang PetscErrorCode ierr; 864137fb511SHong Zhang PetscInt *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j; 865137fb511SHong Zhang PetscInt nz,*rout,*cout,row; 866dd6ea824SBarry Smith PetscScalar *x,*b,*tmp,*tmps,sum; 867dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 868137fb511SHong Zhang 869137fb511SHong Zhang PetscFunctionBegin; 870137fb511SHong Zhang if (!n) PetscFunctionReturn(0); 871137fb511SHong Zhang 872137fb511SHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 873137fb511SHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 874137fb511SHong Zhang tmp = a->solve_work; 875137fb511SHong Zhang 876137fb511SHong Zhang ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 877137fb511SHong Zhang ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 878137fb511SHong Zhang 879137fb511SHong Zhang /* forward solve the lower triangular */ 880137fb511SHong Zhang tmp[0] = b[*r++]; 881137fb511SHong Zhang tmps = tmp; 882137fb511SHong Zhang for (row=1; row<n; row++) { 883137fb511SHong Zhang i = rout[row]; /* permuted row */ 884137fb511SHong Zhang v = aa + ai[i] ; 885137fb511SHong Zhang vi = aj + ai[i] ; 886137fb511SHong Zhang nz = a->diag[i] - ai[i]; 887137fb511SHong Zhang sum = b[*r++]; 888137fb511SHong Zhang SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 889137fb511SHong Zhang tmp[row] = sum; 890137fb511SHong Zhang } 891137fb511SHong Zhang 892137fb511SHong Zhang /* backward solve the upper triangular */ 893137fb511SHong Zhang for (row=n-1; row>=0; row--){ 894137fb511SHong Zhang i = rout[row]; /* permuted row */ 895137fb511SHong Zhang v = aa + a->diag[i] + 1; 896137fb511SHong Zhang vi = aj + a->diag[i] + 1; 897137fb511SHong Zhang nz = ai[i+1] - a->diag[i] - 1; 898137fb511SHong Zhang sum = tmp[row]; 899137fb511SHong Zhang SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 900137fb511SHong Zhang x[*c--] = tmp[row] = sum*aa[a->diag[i]]; 901137fb511SHong Zhang } 902137fb511SHong Zhang 903137fb511SHong Zhang ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 904137fb511SHong Zhang ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 905137fb511SHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 906137fb511SHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 907137fb511SHong Zhang ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr); 908137fb511SHong Zhang PetscFunctionReturn(0); 909137fb511SHong Zhang } 910137fb511SHong Zhang 911930ae53cSBarry Smith /* ----------------------------------------------------------- */ 9124a2ae208SSatish Balay #undef __FUNCT__ 9134a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 914dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 915930ae53cSBarry Smith { 916930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 9176849ba73SBarry Smith PetscErrorCode ierr; 918899cda47SBarry Smith PetscInt n = A->rmap.n,*ai = a->i,*aj = a->j,*adiag = a->diag; 919dd6ea824SBarry Smith PetscScalar *x,*b; 920dd6ea824SBarry Smith const MatScalar *aa = a->a; 921aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 92238baddfdSBarry Smith PetscInt adiag_i,i,*vi,nz,ai_i; 923dd6ea824SBarry Smith const MatScalar *v; 924dd6ea824SBarry Smith PetscScalar sum; 925d85afc42SSatish Balay #endif 926930ae53cSBarry Smith 9273a40ed3dSBarry Smith PetscFunctionBegin; 9283a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 929930ae53cSBarry Smith 9301ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 9311ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 932930ae53cSBarry Smith 933aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 9341c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 9351c4feecaSSatish Balay #else 936930ae53cSBarry Smith /* forward solve the lower triangular */ 937930ae53cSBarry Smith x[0] = b[0]; 938930ae53cSBarry Smith for (i=1; i<n; i++) { 939930ae53cSBarry Smith ai_i = ai[i]; 940930ae53cSBarry Smith v = aa + ai_i; 941930ae53cSBarry Smith vi = aj + ai_i; 942930ae53cSBarry Smith nz = adiag[i] - ai_i; 943930ae53cSBarry Smith sum = b[i]; 944930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 945930ae53cSBarry Smith x[i] = sum; 946930ae53cSBarry Smith } 947930ae53cSBarry Smith 948930ae53cSBarry Smith /* backward solve the upper triangular */ 949930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 950930ae53cSBarry Smith adiag_i = adiag[i]; 951d708749eSBarry Smith v = aa + adiag_i + 1; 952d708749eSBarry Smith vi = aj + adiag_i + 1; 953930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 954930ae53cSBarry Smith sum = x[i]; 955930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 956930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 957930ae53cSBarry Smith } 9581c4feecaSSatish Balay #endif 959899cda47SBarry Smith ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr); 9601ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 9611ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 9623a40ed3dSBarry Smith PetscFunctionReturn(0); 963930ae53cSBarry Smith } 964930ae53cSBarry Smith 9654a2ae208SSatish Balay #undef __FUNCT__ 9664a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 967dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 968da3a660dSBarry Smith { 969416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 970416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 9716849ba73SBarry Smith PetscErrorCode ierr; 972899cda47SBarry Smith PetscInt *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j; 97338baddfdSBarry Smith PetscInt nz,*rout,*cout; 974dd6ea824SBarry Smith PetscScalar *x,*b,*tmp,sum; 975dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 976da3a660dSBarry Smith 9773a40ed3dSBarry Smith PetscFunctionBegin; 97878b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 979da3a660dSBarry Smith 9801ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 9811ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 982416022c9SBarry Smith tmp = a->solve_work; 983da3a660dSBarry Smith 9843b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 9853b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 986da3a660dSBarry Smith 987da3a660dSBarry Smith /* forward solve the lower triangular */ 988da3a660dSBarry Smith tmp[0] = b[*r++]; 989da3a660dSBarry Smith for (i=1; i<n; i++) { 990010a20caSHong Zhang v = aa + ai[i] ; 991010a20caSHong Zhang vi = aj + ai[i] ; 992416022c9SBarry Smith nz = a->diag[i] - ai[i]; 993da3a660dSBarry Smith sum = b[*r++]; 994010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 995da3a660dSBarry Smith tmp[i] = sum; 996da3a660dSBarry Smith } 997da3a660dSBarry Smith 998da3a660dSBarry Smith /* backward solve the upper triangular */ 999da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 1000010a20caSHong Zhang v = aa + a->diag[i] + 1; 1001010a20caSHong Zhang vi = aj + a->diag[i] + 1; 1002416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 1003da3a660dSBarry Smith sum = tmp[i]; 1004010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 1005010a20caSHong Zhang tmp[i] = sum*aa[a->diag[i]]; 1006da3a660dSBarry Smith x[*c--] += tmp[i]; 1007da3a660dSBarry Smith } 1008da3a660dSBarry Smith 10093b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 10103b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 10111ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 10121ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1013efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 1014898183ecSLois Curfman McInnes 10153a40ed3dSBarry Smith PetscFunctionReturn(0); 1016da3a660dSBarry Smith } 1017da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 10184a2ae208SSatish Balay #undef __FUNCT__ 10194a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 1020dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 1021da3a660dSBarry Smith { 1022416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 10232235e667SBarry Smith IS iscol = a->col,isrow = a->row; 10246849ba73SBarry Smith PetscErrorCode ierr; 1025899cda47SBarry Smith PetscInt *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j; 102638baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 1027dd6ea824SBarry Smith PetscScalar *x,*b,*tmp,s1; 1028dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 1029da3a660dSBarry Smith 10303a40ed3dSBarry Smith PetscFunctionBegin; 10311ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 10321ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1033416022c9SBarry Smith tmp = a->solve_work; 1034da3a660dSBarry Smith 10352235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 10362235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 1037da3a660dSBarry Smith 1038da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 10392235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 1040da3a660dSBarry Smith 1041da3a660dSBarry Smith /* forward solve the U^T */ 1042da3a660dSBarry Smith for (i=0; i<n; i++) { 1043010a20caSHong Zhang v = aa + diag[i] ; 1044010a20caSHong Zhang vi = aj + diag[i] + 1; 1045f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 1046c38d4ed2SBarry Smith s1 = tmp[i]; 1047ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 1048da3a660dSBarry Smith while (nz--) { 1049010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*s1; 1050da3a660dSBarry Smith } 1051c38d4ed2SBarry Smith tmp[i] = s1; 1052da3a660dSBarry Smith } 1053da3a660dSBarry Smith 1054da3a660dSBarry Smith /* backward solve the L^T */ 1055da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 1056010a20caSHong Zhang v = aa + diag[i] - 1 ; 1057010a20caSHong Zhang vi = aj + diag[i] - 1 ; 1058f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 1059f1af5d2fSBarry Smith s1 = tmp[i]; 1060da3a660dSBarry Smith while (nz--) { 1061010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*s1; 1062da3a660dSBarry Smith } 1063da3a660dSBarry Smith } 1064da3a660dSBarry Smith 1065da3a660dSBarry Smith /* copy tmp into x according to permutation */ 1066da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 1067da3a660dSBarry Smith 10682235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 10692235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 10701ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 10711ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1072da3a660dSBarry Smith 1073899cda47SBarry Smith ierr = PetscLogFlops(2*a->nz-A->cmap.n);CHKERRQ(ierr); 10743a40ed3dSBarry Smith PetscFunctionReturn(0); 1075da3a660dSBarry Smith } 1076da3a660dSBarry Smith 10774a2ae208SSatish Balay #undef __FUNCT__ 10784a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 1079dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 1080da3a660dSBarry Smith { 1081416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 10822235e667SBarry Smith IS iscol = a->col,isrow = a->row; 10836849ba73SBarry Smith PetscErrorCode ierr; 1084899cda47SBarry Smith PetscInt *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j; 108538baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 1086dd6ea824SBarry Smith PetscScalar *x,*b,*tmp; 1087dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 10886abc6512SBarry Smith 10893a40ed3dSBarry Smith PetscFunctionBegin; 109023bc6035SMatthew Knepley if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);} 10916abc6512SBarry Smith 10921ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 10931ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1094416022c9SBarry Smith tmp = a->solve_work; 10956abc6512SBarry Smith 10962235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 10972235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 10986abc6512SBarry Smith 10996abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 11002235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 11016abc6512SBarry Smith 11026abc6512SBarry Smith /* forward solve the U^T */ 11036abc6512SBarry Smith for (i=0; i<n; i++) { 1104010a20caSHong Zhang v = aa + diag[i] ; 1105010a20caSHong Zhang vi = aj + diag[i] + 1; 1106f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 11076abc6512SBarry Smith tmp[i] *= *v++; 11086abc6512SBarry Smith while (nz--) { 1109010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*tmp[i]; 11106abc6512SBarry Smith } 11116abc6512SBarry Smith } 11126abc6512SBarry Smith 11136abc6512SBarry Smith /* backward solve the L^T */ 11146abc6512SBarry Smith for (i=n-1; i>=0; i--){ 1115010a20caSHong Zhang v = aa + diag[i] - 1 ; 1116010a20caSHong Zhang vi = aj + diag[i] - 1 ; 1117f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 11186abc6512SBarry Smith while (nz--) { 1119010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*tmp[i]; 11206abc6512SBarry Smith } 11216abc6512SBarry Smith } 11226abc6512SBarry Smith 11236abc6512SBarry Smith /* copy tmp into x according to permutation */ 11246abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 11256abc6512SBarry Smith 11262235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 11272235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 11281ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 11291ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 11306abc6512SBarry Smith 1131efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 11323a40ed3dSBarry Smith PetscFunctionReturn(0); 1133da3a660dSBarry Smith } 11349e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 11353c5fc038SBarry Smith EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth); 1136*b24902e0SBarry Smith EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 113708480c60SBarry Smith 11384a2ae208SSatish Balay #undef __FUNCT__ 11394a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 1140dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 11419e25ed09SBarry Smith { 1142416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 11439e25ed09SBarry Smith IS isicol; 11446849ba73SBarry Smith PetscErrorCode ierr; 114509f38230SBarry Smith PetscInt *r,*ic,n=A->rmap.n,*ai=a->i,*aj=a->j,d; 11465a8e39fbSHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 11475a8e39fbSHong Zhang PetscInt *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0; 11485a8e39fbSHong Zhang PetscInt i,levels,diagonal_fill; 114977c4ece6SBarry Smith PetscTruth col_identity,row_identity; 1150329f5518SBarry Smith PetscReal f; 11515a8e39fbSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 11525a8e39fbSHong Zhang PetscBT lnkbt; 11535a8e39fbSHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 1154a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1155a1a86e44SBarry Smith PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 115609f38230SBarry Smith PetscTruth missing; 11579e25ed09SBarry Smith 11583a40ed3dSBarry Smith PetscFunctionBegin; 11596cf997b0SBarry Smith f = info->fill; 116038baddfdSBarry Smith levels = (PetscInt)info->levels; 116138baddfdSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 11624c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 116382bf6240SBarry Smith 116408480c60SBarry Smith /* special case that simply copies fill pattern */ 1165be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 1166be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 116786bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 1168*b24902e0SBarry Smith ierr = MatDuplicateNoCreate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 116908480c60SBarry Smith (*fact)->factor = FACTOR_LU; 1170f3a39becSBarry Smith (*fact)->info.factor_mallocs = 0; 1171f3a39becSBarry Smith (*fact)->info.fill_ratio_given = info->fill; 1172f3a39becSBarry Smith (*fact)->info.fill_ratio_needed = 1.0; 117308480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 11748ef3462fSBarry Smith ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 117509f38230SBarry Smith if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 117608480c60SBarry Smith b->row = isrow; 117708480c60SBarry Smith b->col = iscol; 117882bf6240SBarry Smith b->icol = isicol; 1179899cda47SBarry Smith ierr = PetscMalloc(((*fact)->rmap.n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1180f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 1181c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1182c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 11833c5fc038SBarry Smith ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 11843a40ed3dSBarry Smith PetscFunctionReturn(0); 118508480c60SBarry Smith } 118608480c60SBarry Smith 1187898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1188898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 11899e25ed09SBarry Smith 11909e25ed09SBarry Smith /* get new row pointers */ 11915a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 11925a8e39fbSHong Zhang bi[0] = 0; 11935a8e39fbSHong Zhang /* bdiag is location of diagonal in factor */ 11945a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 11955a8e39fbSHong Zhang bdiag[0] = 0; 11966cf997b0SBarry Smith 11975a8e39fbSHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr); 11985a8e39fbSHong Zhang bjlvl_ptr = (PetscInt**)(bj_ptr + n); 11995a8e39fbSHong Zhang 12005a8e39fbSHong Zhang /* create a linked list for storing column indices of the active row */ 12015a8e39fbSHong Zhang nlnk = n + 1; 12025a8e39fbSHong Zhang ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12035a8e39fbSHong Zhang 12045a8e39fbSHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 1205a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 12065a8e39fbSHong Zhang current_space = free_space; 1207a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 12085a8e39fbSHong Zhang current_space_lvl = free_space_lvl; 12095a8e39fbSHong Zhang 12105a8e39fbSHong Zhang for (i=0; i<n; i++) { 12115a8e39fbSHong Zhang nzi = 0; 12126cf997b0SBarry Smith /* copy current row into linked list */ 12135a8e39fbSHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 12145a8e39fbSHong Zhang if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 12155a8e39fbSHong Zhang cols = aj + ai[r[i]]; 12165a8e39fbSHong Zhang lnk[i] = -1; /* marker to indicate if diagonal exists */ 12175a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12185a8e39fbSHong Zhang nzi += nlnk; 12196cf997b0SBarry Smith 12206cf997b0SBarry Smith /* make sure diagonal entry is included */ 12215a8e39fbSHong Zhang if (diagonal_fill && lnk[i] == -1) { 12226cf997b0SBarry Smith fm = n; 12235a8e39fbSHong Zhang while (lnk[fm] < i) fm = lnk[fm]; 12245a8e39fbSHong Zhang lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 12255a8e39fbSHong Zhang lnk[fm] = i; 12265a8e39fbSHong Zhang lnk_lvl[i] = 0; 12275a8e39fbSHong Zhang nzi++; dcount++; 12286cf997b0SBarry Smith } 12296cf997b0SBarry Smith 12305a8e39fbSHong Zhang /* add pivot rows into the active row */ 12315a8e39fbSHong Zhang nzbd = 0; 12325a8e39fbSHong Zhang prow = lnk[n]; 12335a8e39fbSHong Zhang while (prow < i) { 12345a8e39fbSHong Zhang nnz = bdiag[prow]; 12355a8e39fbSHong Zhang cols = bj_ptr[prow] + nnz + 1; 12365a8e39fbSHong Zhang cols_lvl = bjlvl_ptr[prow] + nnz + 1; 12375a8e39fbSHong Zhang nnz = bi[prow+1] - bi[prow] - nnz - 1; 12385a8e39fbSHong Zhang ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 12395a8e39fbSHong Zhang nzi += nlnk; 12405a8e39fbSHong Zhang prow = lnk[prow]; 12415a8e39fbSHong Zhang nzbd++; 124256beaf04SBarry Smith } 12435a8e39fbSHong Zhang bdiag[i] = nzbd; 12445a8e39fbSHong Zhang bi[i+1] = bi[i] + nzi; 1245ecf371e4SBarry Smith 12465a8e39fbSHong Zhang /* if free space is not available, make more free space */ 12475a8e39fbSHong Zhang if (current_space->local_remaining<nzi) { 12485a8e39fbSHong Zhang nnz = nzi*(n - i); /* estimated and max additional space needed */ 1249a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 1250a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 12515a8e39fbSHong Zhang reallocs++; 12525a8e39fbSHong Zhang } 1253ecf371e4SBarry Smith 12545a8e39fbSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 12555a8e39fbSHong Zhang ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 12565a8e39fbSHong Zhang bj_ptr[i] = current_space->array; 12575a8e39fbSHong Zhang bjlvl_ptr[i] = current_space_lvl->array; 12585a8e39fbSHong Zhang 12595a8e39fbSHong Zhang /* make sure the active row i has diagonal entry */ 12605a8e39fbSHong Zhang if (*(bj_ptr[i]+bdiag[i]) != i) { 126177431f27SBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 1262d7ee0231SBarry Smith try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i); 12636cf997b0SBarry Smith } 12645a8e39fbSHong Zhang 12655a8e39fbSHong Zhang current_space->array += nzi; 12665a8e39fbSHong Zhang current_space->local_used += nzi; 12675a8e39fbSHong Zhang current_space->local_remaining -= nzi; 12685a8e39fbSHong Zhang current_space_lvl->array += nzi; 12695a8e39fbSHong Zhang current_space_lvl->local_used += nzi; 12705a8e39fbSHong Zhang current_space_lvl->local_remaining -= nzi; 12719e25ed09SBarry Smith } 12725a8e39fbSHong Zhang 1273898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1274898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 12755a8e39fbSHong Zhang 12765a8e39fbSHong Zhang /* destroy list of free space and other temporary arrays */ 12775a8e39fbSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1278a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 12795a8e39fbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1280a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 12815a8e39fbSHong Zhang ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 12829e25ed09SBarry Smith 12836cf91177SBarry Smith #if defined(PETSC_USE_INFO) 1284f64065bbSBarry Smith { 12855a8e39fbSHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1286ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 1287ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1288ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr); 1289ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 1290335d9088SBarry Smith if (diagonal_fill) { 1291ae15b995SBarry Smith ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr); 1292335d9088SBarry Smith } 1293f64065bbSBarry Smith } 129463ba0a88SBarry Smith #endif 1295bbb6d6a8SBarry Smith 12969e25ed09SBarry Smith /* put together the new matrix */ 129752e6d16bSBarry Smith ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr); 1298416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 1299e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1300e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 13017c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 13025a8e39fbSHong Zhang len = (bi[n] )*sizeof(PetscScalar); 1303b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 13045a8e39fbSHong Zhang b->j = bj; 13055a8e39fbSHong Zhang b->i = bi; 13065a8e39fbSHong Zhang for (i=0; i<n; i++) bdiag[i] += bi[i]; 13075a8e39fbSHong Zhang b->diag = bdiag; 1308416022c9SBarry Smith b->ilen = 0; 1309416022c9SBarry Smith b->imax = 0; 1310416022c9SBarry Smith b->row = isrow; 1311416022c9SBarry Smith b->col = iscol; 1312c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1313c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 131482bf6240SBarry Smith b->icol = isicol; 131587828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1316416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 13175a8e39fbSHong Zhang Allocate bdiag, solve_work, new a, new j */ 13185a8e39fbSHong Zhang ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 13195a8e39fbSHong Zhang b->maxnz = b->nz = bi[n] ; 13209e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 1321418422e8SSatish Balay (*fact)->info.factor_mallocs = reallocs; 1322ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 13235a8e39fbSHong Zhang (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 132471c2f376SKris Buschelman 132554e71613SHong Zhang ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr); 132671c2f376SKris Buschelman (*fact)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 13274846f1f5SKris Buschelman 13283a40ed3dSBarry Smith PetscFunctionReturn(0); 13299e25ed09SBarry Smith } 1330d5d45c9bSBarry Smith 13313c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1332a6175056SHong Zhang #undef __FUNCT__ 1333a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1334af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 1335a6175056SHong Zhang { 13362ed133dbSHong Zhang Mat C = *B; 13370968510aSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 13382ed133dbSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 13399bfd6278SHong Zhang IS ip=b->row,iip = b->icol; 13402ed133dbSHong Zhang PetscErrorCode ierr; 13419bfd6278SHong Zhang PetscInt *rip,*riip,i,j,mbs=A->rmap.n,*bi=b->i,*bj=b->j,*bcol; 13422ed133dbSHong Zhang PetscInt *ai=a->i,*aj=a->j; 13432ed133dbSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1344622e2028SHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1345540703acSHong Zhang PetscReal zeropivot,rs,shiftnz; 1346fbf22428SSatish Balay PetscReal shiftpd; 1347540703acSHong Zhang ChShift_Ctx sctx; 1348540703acSHong Zhang PetscInt newshift; 1349d5d45c9bSBarry Smith 1350a6175056SHong Zhang PetscFunctionBegin; 1351117f1891Stmunson 1352540703acSHong Zhang shiftnz = info->shiftnz; 1353540703acSHong Zhang shiftpd = info->shiftpd; 1354ee45ca4aSHong Zhang zeropivot = info->zeropivot; 1355ee45ca4aSHong Zhang 13562ed133dbSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 13579bfd6278SHong Zhang ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr); 13582ed133dbSHong Zhang 13592ed133dbSHong Zhang /* initialization */ 13602ed133dbSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 13612ed133dbSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 13622ed133dbSHong Zhang jl = il + mbs; 13632ed133dbSHong Zhang rtmp = (MatScalar*)(jl + mbs); 13642ed133dbSHong Zhang 1365540703acSHong Zhang sctx.shift_amount = 0; 1366540703acSHong Zhang sctx.nshift = 0; 13672ed133dbSHong Zhang do { 1368540703acSHong Zhang sctx.chshift = PETSC_FALSE; 13692ed133dbSHong Zhang for (i=0; i<mbs; i++) { 13702ed133dbSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1371a6175056SHong Zhang } 13722ed133dbSHong Zhang 13732ed133dbSHong Zhang for (k = 0; k<mbs; k++){ 1374c05c3958SHong Zhang bval = ba + bi[k]; 13752ed133dbSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 13762ed133dbSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 13772ed133dbSHong Zhang for (j = jmin; j < jmax; j++){ 13789bfd6278SHong Zhang col = riip[aj[j]]; 13792ed133dbSHong Zhang if (col >= k){ /* only take upper triangular entry */ 13802ed133dbSHong Zhang rtmp[col] = aa[j]; 1381c05c3958SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 13822ed133dbSHong Zhang } 13832ed133dbSHong Zhang } 138439e3d630SHong Zhang /* shift the diagonal of the matrix */ 1385540703acSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 13862ed133dbSHong Zhang 13872ed133dbSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 13882ed133dbSHong Zhang dk = rtmp[k]; 13892ed133dbSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 13902ed133dbSHong Zhang 13912ed133dbSHong Zhang while (i < k){ 13922ed133dbSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 13932ed133dbSHong Zhang 13942ed133dbSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 13952ed133dbSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 13962ed133dbSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 13972ed133dbSHong Zhang dk += uikdi*ba[ili]; 13982ed133dbSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 13992ed133dbSHong Zhang 14002ed133dbSHong Zhang /* add multiple of row i to k-th row */ 14012ed133dbSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 14022ed133dbSHong Zhang if (jmin < jmax){ 14032ed133dbSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 14042ed133dbSHong Zhang /* update il and jl for row i */ 14052ed133dbSHong Zhang il[i] = jmin; 14062ed133dbSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 14072ed133dbSHong Zhang } 14082ed133dbSHong Zhang i = nexti; 14092ed133dbSHong Zhang } 14102ed133dbSHong Zhang 14110697b6caSHong Zhang /* shift the diagonals when zero pivot is detected */ 14120697b6caSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 14130697b6caSHong Zhang rs = 0.0; 14143c9e8598SHong Zhang jmin = bi[k]+1; 14153c9e8598SHong Zhang nz = bi[k+1] - jmin; 14163c9e8598SHong Zhang bcol = bj + jmin; 14173c9e8598SHong Zhang while (nz--){ 141889f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 141989f845c8SHong Zhang bcol++; 14203c9e8598SHong Zhang } 1421540703acSHong Zhang 1422540703acSHong Zhang sctx.rs = rs; 1423540703acSHong Zhang sctx.pv = dk; 14245b5bc046SBarry Smith ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 1425117f1891Stmunson 1426117f1891Stmunson if (newshift == 1) { 1427117f1891Stmunson if (!sctx.shift_amount) { 1428117f1891Stmunson sctx.shift_amount = 1e-5; 1429117f1891Stmunson } 1430117f1891Stmunson break; 1431117f1891Stmunson } 14323c9e8598SHong Zhang 14333c9e8598SHong Zhang /* copy data into U(k,:) */ 143439e3d630SHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 143539e3d630SHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 143639e3d630SHong Zhang if (jmin < jmax) { 143739e3d630SHong Zhang for (j=jmin; j<jmax; j++){ 143839e3d630SHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 14393c9e8598SHong Zhang } 144039e3d630SHong Zhang /* add the k-th row into il and jl */ 14413c9e8598SHong Zhang il[k] = jmin; 14423c9e8598SHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 14433c9e8598SHong Zhang } 14443c9e8598SHong Zhang } 1445540703acSHong Zhang } while (sctx.chshift); 14463c9e8598SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 14473c9e8598SHong Zhang 144839e3d630SHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 14499bfd6278SHong Zhang ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr); 14503c9e8598SHong Zhang C->factor = FACTOR_CHOLESKY; 14513c9e8598SHong Zhang C->assembled = PETSC_TRUE; 14523c9e8598SHong Zhang C->preallocated = PETSC_TRUE; 1453899cda47SBarry Smith ierr = PetscLogFlops(C->rmap.n);CHKERRQ(ierr); 1454540703acSHong Zhang if (sctx.nshift){ 1455540703acSHong Zhang if (shiftnz) { 14561e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1457540703acSHong Zhang } else if (shiftpd) { 14581e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 14590697b6caSHong Zhang } 14603c9e8598SHong Zhang } 14613c9e8598SHong Zhang PetscFunctionReturn(0); 14623c9e8598SHong Zhang } 1463a6175056SHong Zhang 1464a6175056SHong Zhang #undef __FUNCT__ 1465a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 1466dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1467a6175056SHong Zhang { 14680968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1469ed59401aSHong Zhang Mat_SeqSBAIJ *b; 1470dfbe8321SBarry Smith PetscErrorCode ierr; 147158ebbce7SBarry Smith PetscTruth perm_identity,missing; 1472b635d36bSHong Zhang PetscInt reallocs=0,*rip,*riip,i,*ai=a->i,*aj=a->j,am=A->rmap.n,*ui; 1473ed59401aSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 147458ebbce7SBarry Smith PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,d; 14755a8e39fbSHong Zhang PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 1476ed59401aSHong Zhang PetscReal fill=info->fill,levels=info->levels; 1477a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1478a1a86e44SBarry Smith PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 14797a48dd6fSHong Zhang PetscBT lnkbt; 1480b635d36bSHong Zhang IS iperm; 1481a6175056SHong Zhang 1482a6175056SHong Zhang PetscFunctionBegin; 148358ebbce7SBarry Smith ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 148458ebbce7SBarry Smith if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 1485653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1486b635d36bSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 14877a48dd6fSHong Zhang 148839e3d630SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 148939e3d630SHong Zhang ui[0] = 0; 149039e3d630SHong Zhang 1491b635d36bSHong Zhang /* ICC(0) without matrix ordering: simply copies fill pattern */ 1492622e2028SHong Zhang if (!levels && perm_identity) { 149358ebbce7SBarry Smith 1494ed59401aSHong Zhang for (i=0; i<am; i++) { 149539e3d630SHong Zhang ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 1496ed59401aSHong Zhang } 149739e3d630SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 149839e3d630SHong Zhang cols = uj; 1499ed59401aSHong Zhang for (i=0; i<am; i++) { 1500ed59401aSHong Zhang aj = a->j + a->diag[i]; 150139e3d630SHong Zhang ncols = ui[i+1] - ui[i]; 150239e3d630SHong Zhang for (j=0; j<ncols; j++) *cols++ = *aj++; 1503ed59401aSHong Zhang } 150439e3d630SHong Zhang } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 1505b635d36bSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 1506b635d36bSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 1507b635d36bSHong Zhang 15087a48dd6fSHong Zhang /* initialization */ 15095a8e39fbSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 15107a48dd6fSHong Zhang 15117a48dd6fSHong Zhang /* jl: linked list for storing indices of the pivot rows 15127a48dd6fSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 15131393bc97SHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 15147a48dd6fSHong Zhang il = jl + am; 15157a48dd6fSHong Zhang uj_ptr = (PetscInt**)(il + am); 15167a48dd6fSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 15177a48dd6fSHong Zhang for (i=0; i<am; i++){ 15187a48dd6fSHong Zhang jl[i] = am; il[i] = 0; 15197a48dd6fSHong Zhang } 15207a48dd6fSHong Zhang 15217a48dd6fSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 15227a48dd6fSHong Zhang nlnk = am + 1; 15232ed133dbSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 15247a48dd6fSHong Zhang 15257a48dd6fSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 1526a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 15277a48dd6fSHong Zhang current_space = free_space; 1528a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 15297a48dd6fSHong Zhang current_space_lvl = free_space_lvl; 15307a48dd6fSHong Zhang 15317a48dd6fSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 15327a48dd6fSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 15337a48dd6fSHong Zhang nzk = 0; 1534622e2028SHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 1535b635d36bSHong Zhang if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix"); 1536622e2028SHong Zhang ncols_upper = 0; 1537622e2028SHong Zhang for (j=0; j<ncols; j++){ 1538b635d36bSHong Zhang i = *(aj + ai[rip[k]] + j); /* unpermuted column index */ 1539b635d36bSHong Zhang if (riip[i] >= k){ /* only take upper triangular entry */ 15405a8e39fbSHong Zhang ajtmp[ncols_upper] = i; 1541622e2028SHong Zhang ncols_upper++; 1542622e2028SHong Zhang } 1543622e2028SHong Zhang } 1544b635d36bSHong Zhang ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 15457a48dd6fSHong Zhang nzk += nlnk; 15467a48dd6fSHong Zhang 15477a48dd6fSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 15487a48dd6fSHong Zhang prow = jl[k]; /* 1st pivot row */ 15497a48dd6fSHong Zhang 15507a48dd6fSHong Zhang while (prow < k){ 15517a48dd6fSHong Zhang nextprow = jl[prow]; 15527a48dd6fSHong Zhang 15537a48dd6fSHong Zhang /* merge prow into k-th row */ 15547a48dd6fSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 15557a48dd6fSHong Zhang jmax = ui[prow+1]; 15567a48dd6fSHong Zhang ncols = jmax-jmin; 1557ed59401aSHong Zhang i = jmin - ui[prow]; 1558ed59401aSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 15595a8e39fbSHong Zhang uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 15605a8e39fbSHong Zhang j = *(uj - 1); 15615a8e39fbSHong Zhang ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 15627a48dd6fSHong Zhang nzk += nlnk; 15637a48dd6fSHong Zhang 15647a48dd6fSHong Zhang /* update il and jl for prow */ 15657a48dd6fSHong Zhang if (jmin < jmax){ 15667a48dd6fSHong Zhang il[prow] = jmin; 15677a48dd6fSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 15687a48dd6fSHong Zhang } 15697a48dd6fSHong Zhang prow = nextprow; 15707a48dd6fSHong Zhang } 15717a48dd6fSHong Zhang 15727a48dd6fSHong Zhang /* if free space is not available, make more free space */ 15737a48dd6fSHong Zhang if (current_space->local_remaining<nzk) { 15747a48dd6fSHong Zhang i = am - k + 1; /* num of unfactored rows */ 15757a48dd6fSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1576a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1577a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 15787a48dd6fSHong Zhang reallocs++; 15797a48dd6fSHong Zhang } 15807a48dd6fSHong Zhang 15817a48dd6fSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 1582b635d36bSHong Zhang if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k); 15832ed133dbSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 15847a48dd6fSHong Zhang 15857a48dd6fSHong Zhang /* add the k-th row into il and jl */ 158639e3d630SHong Zhang if (nzk > 1){ 15877a48dd6fSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 15887a48dd6fSHong Zhang jl[k] = jl[i]; jl[i] = k; 15897a48dd6fSHong Zhang il[k] = ui[k] + 1; 15907a48dd6fSHong Zhang } 15917a48dd6fSHong Zhang uj_ptr[k] = current_space->array; 15927a48dd6fSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 15937a48dd6fSHong Zhang 15947a48dd6fSHong Zhang current_space->array += nzk; 15957a48dd6fSHong Zhang current_space->local_used += nzk; 15967a48dd6fSHong Zhang current_space->local_remaining -= nzk; 15977a48dd6fSHong Zhang 15987a48dd6fSHong Zhang current_space_lvl->array += nzk; 15997a48dd6fSHong Zhang current_space_lvl->local_used += nzk; 16007a48dd6fSHong Zhang current_space_lvl->local_remaining -= nzk; 16017a48dd6fSHong Zhang 16027a48dd6fSHong Zhang ui[k+1] = ui[k] + nzk; 16037a48dd6fSHong Zhang } 16047a48dd6fSHong Zhang 16056cf91177SBarry Smith #if defined(PETSC_USE_INFO) 16067a48dd6fSHong Zhang if (ai[am] != 0) { 160739e3d630SHong Zhang PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 1608ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1609ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1610ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 16117a48dd6fSHong Zhang } else { 1612ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 16137a48dd6fSHong Zhang } 161463ba0a88SBarry Smith #endif 16157a48dd6fSHong Zhang 16167a48dd6fSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1617b635d36bSHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 16187a48dd6fSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 16195a8e39fbSHong Zhang ierr = PetscFree(ajtmp);CHKERRQ(ierr); 16207a48dd6fSHong Zhang 16217a48dd6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 16227a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1623a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 16242ed133dbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1625a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 16267a48dd6fSHong Zhang 162739e3d630SHong Zhang } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 162839e3d630SHong Zhang 16297a48dd6fSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 16307a48dd6fSHong Zhang 1631*b24902e0SBarry Smith b = (Mat_SeqSBAIJ*)(*fact)->data; 16327a48dd6fSHong Zhang b->singlemalloc = PETSC_FALSE; 16337a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 16347a48dd6fSHong Zhang b->j = uj; 16357a48dd6fSHong Zhang b->i = ui; 16367a48dd6fSHong Zhang b->diag = 0; 16377a48dd6fSHong Zhang b->ilen = 0; 16387a48dd6fSHong Zhang b->imax = 0; 16397a48dd6fSHong Zhang b->row = perm; 1640b635d36bSHong Zhang b->col = perm; 1641b635d36bSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1642b635d36bSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1643b635d36bSHong Zhang b->icol = iperm; 16447a48dd6fSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 16457a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1646*b24902e0SBarry Smith ierr = PetscLogObjectMemory((*fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 16477a48dd6fSHong Zhang b->maxnz = b->nz = ui[am]; 1648e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1649e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 16507a48dd6fSHong Zhang 1651*b24902e0SBarry Smith (*fact)->factor = FACTOR_CHOLESKY; 1652*b24902e0SBarry Smith (*fact)->info.factor_mallocs = reallocs; 1653*b24902e0SBarry Smith (*fact)->info.fill_ratio_given = fill; 16547a48dd6fSHong Zhang if (ai[am] != 0) { 1655*b24902e0SBarry Smith (*fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 16567a48dd6fSHong Zhang } else { 1657*b24902e0SBarry Smith (*fact)->info.fill_ratio_needed = 0.0; 16587a48dd6fSHong Zhang } 165939e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1660622e2028SHong Zhang if (perm_identity){ 1661*b24902e0SBarry Smith (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1662*b24902e0SBarry Smith (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1663622e2028SHong Zhang } 1664a6175056SHong Zhang PetscFunctionReturn(0); 1665a6175056SHong Zhang } 1666d5d45c9bSBarry Smith 1667f76d2b81SHong Zhang #undef __FUNCT__ 1668f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1669dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1670f76d2b81SHong Zhang { 1671f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 167210c27e3dSHong Zhang Mat_SeqSBAIJ *b; 1673dfbe8321SBarry Smith PetscErrorCode ierr; 1674f76d2b81SHong Zhang PetscTruth perm_identity; 167510c27e3dSHong Zhang PetscReal fill = info->fill; 1676899cda47SBarry Smith PetscInt *rip,*riip,i,am=A->rmap.n,*ai=a->i,*aj=a->j,reallocs=0,prow; 167710c27e3dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 16782ed133dbSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 1679a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 168010c27e3dSHong Zhang PetscBT lnkbt; 16812ed133dbSHong Zhang IS iperm; 1682f76d2b81SHong Zhang 1683f76d2b81SHong Zhang PetscFunctionBegin; 16842ed133dbSHong Zhang /* check whether perm is the identity mapping */ 1685f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 16862ed133dbSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 16872ed133dbSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 16889bfd6278SHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 168910c27e3dSHong Zhang 169010c27e3dSHong Zhang /* initialization */ 16911393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 169210c27e3dSHong Zhang ui[0] = 0; 169310c27e3dSHong Zhang 169410c27e3dSHong Zhang /* jl: linked list for storing indices of the pivot rows 16951393bc97SHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 16961393bc97SHong Zhang ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 16971393bc97SHong Zhang il = jl + am; 16981393bc97SHong Zhang cols = il + am; 16991393bc97SHong Zhang ui_ptr = (PetscInt**)(cols + am); 17001393bc97SHong Zhang for (i=0; i<am; i++){ 17011393bc97SHong Zhang jl[i] = am; il[i] = 0; 1702f76d2b81SHong Zhang } 170310c27e3dSHong Zhang 170410c27e3dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 17051393bc97SHong Zhang nlnk = am + 1; 17061393bc97SHong Zhang ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 170710c27e3dSHong Zhang 17081393bc97SHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 1709a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 171010c27e3dSHong Zhang current_space = free_space; 171110c27e3dSHong Zhang 17121393bc97SHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 171310c27e3dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 171410c27e3dSHong Zhang nzk = 0; 17152ed133dbSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 17169bfd6278SHong Zhang if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix"); 17172ed133dbSHong Zhang ncols_upper = 0; 17182ed133dbSHong Zhang for (j=0; j<ncols; j++){ 17199bfd6278SHong Zhang i = riip[*(aj + ai[rip[k]] + j)]; 17202ed133dbSHong Zhang if (i >= k){ /* only take upper triangular entry */ 17212ed133dbSHong Zhang cols[ncols_upper] = i; 17222ed133dbSHong Zhang ncols_upper++; 17232ed133dbSHong Zhang } 17242ed133dbSHong Zhang } 17251393bc97SHong Zhang ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 172610c27e3dSHong Zhang nzk += nlnk; 172710c27e3dSHong Zhang 172810c27e3dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 172910c27e3dSHong Zhang prow = jl[k]; /* 1st pivot row */ 173010c27e3dSHong Zhang 173110c27e3dSHong Zhang while (prow < k){ 173210c27e3dSHong Zhang nextprow = jl[prow]; 173310c27e3dSHong Zhang /* merge prow into k-th row */ 17341393bc97SHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 173510c27e3dSHong Zhang jmax = ui[prow+1]; 173610c27e3dSHong Zhang ncols = jmax-jmin; 17371393bc97SHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 17385a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 173910c27e3dSHong Zhang nzk += nlnk; 174010c27e3dSHong Zhang 174110c27e3dSHong Zhang /* update il and jl for prow */ 174210c27e3dSHong Zhang if (jmin < jmax){ 174310c27e3dSHong Zhang il[prow] = jmin; 17442ed133dbSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 174510c27e3dSHong Zhang } 174610c27e3dSHong Zhang prow = nextprow; 174710c27e3dSHong Zhang } 174810c27e3dSHong Zhang 174910c27e3dSHong Zhang /* if free space is not available, make more free space */ 175010c27e3dSHong Zhang if (current_space->local_remaining<nzk) { 17511393bc97SHong Zhang i = am - k + 1; /* num of unfactored rows */ 175210c27e3dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1753a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 175410c27e3dSHong Zhang reallocs++; 175510c27e3dSHong Zhang } 175610c27e3dSHong Zhang 175710c27e3dSHong Zhang /* copy data into free space, then initialize lnk */ 17581393bc97SHong Zhang ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 175910c27e3dSHong Zhang 176010c27e3dSHong Zhang /* add the k-th row into il and jl */ 176110c27e3dSHong Zhang if (nzk-1 > 0){ 17621393bc97SHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 176310c27e3dSHong Zhang jl[k] = jl[i]; jl[i] = k; 176410c27e3dSHong Zhang il[k] = ui[k] + 1; 176510c27e3dSHong Zhang } 17662ed133dbSHong Zhang ui_ptr[k] = current_space->array; 176710c27e3dSHong Zhang current_space->array += nzk; 176810c27e3dSHong Zhang current_space->local_used += nzk; 176910c27e3dSHong Zhang current_space->local_remaining -= nzk; 177010c27e3dSHong Zhang 177110c27e3dSHong Zhang ui[k+1] = ui[k] + nzk; 177210c27e3dSHong Zhang } 177310c27e3dSHong Zhang 17746cf91177SBarry Smith #if defined(PETSC_USE_INFO) 17751393bc97SHong Zhang if (ai[am] != 0) { 17761393bc97SHong Zhang PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]); 1777ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1778ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1779ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 178010c27e3dSHong Zhang } else { 1781ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 178210c27e3dSHong Zhang } 178363ba0a88SBarry Smith #endif 178410c27e3dSHong Zhang 178510c27e3dSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 17869bfd6278SHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 178710c27e3dSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 178810c27e3dSHong Zhang 178910c27e3dSHong Zhang /* destroy list of free space and other temporary array(s) */ 17901393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1791a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 179210c27e3dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 179310c27e3dSHong Zhang 179410c27e3dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 179510c27e3dSHong Zhang 1796*b24902e0SBarry Smith b = (Mat_SeqSBAIJ*)(*fact)->data; 179710c27e3dSHong Zhang b->singlemalloc = PETSC_FALSE; 1798e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1799e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 18001393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 180110c27e3dSHong Zhang b->j = uj; 180210c27e3dSHong Zhang b->i = ui; 180310c27e3dSHong Zhang b->diag = 0; 180410c27e3dSHong Zhang b->ilen = 0; 180510c27e3dSHong Zhang b->imax = 0; 180610c27e3dSHong Zhang b->row = perm; 18079bfd6278SHong Zhang b->col = perm; 18089bfd6278SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 18099bfd6278SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 18109bfd6278SHong Zhang b->icol = iperm; 181110c27e3dSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 18121393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1813*b24902e0SBarry Smith ierr = PetscLogObjectMemory(*fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 18141393bc97SHong Zhang b->maxnz = b->nz = ui[am]; 181510c27e3dSHong Zhang 1816*b24902e0SBarry Smith (*fact)->factor = FACTOR_CHOLESKY; 1817*b24902e0SBarry Smith (*fact)->info.factor_mallocs = reallocs; 1818*b24902e0SBarry Smith (*fact)->info.fill_ratio_given = fill; 18191393bc97SHong Zhang if (ai[am] != 0) { 1820*b24902e0SBarry Smith (*fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 182110c27e3dSHong Zhang } else { 1822*b24902e0SBarry Smith (*fact)->info.fill_ratio_needed = 0.0; 182310c27e3dSHong Zhang } 182439e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 18252ed133dbSHong Zhang if (perm_identity){ 182610c27e3dSHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 182710c27e3dSHong Zhang (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1828fd531716SHong Zhang (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 1829fd531716SHong Zhang (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 1830fd531716SHong Zhang } else { 1831fd531716SHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1; 1832fd531716SHong Zhang (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 1833fd531716SHong Zhang (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1; 1834fd531716SHong Zhang (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1; 18352ed133dbSHong Zhang } 1836f76d2b81SHong Zhang PetscFunctionReturn(0); 1837f76d2b81SHong Zhang } 1838