1289bc588SBarry Smith 270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 31c4feecaSSatish Balay #include "src/inline/dot.h" 4ed480e8bSBarry Smith #include "src/inline/spops.h" 51393bc97SHong Zhang #include "petscbt.h" 61393bc97SHong Zhang #include "src/mat/utils/freespace.h" 73b2fbd54SBarry Smith 84a2ae208SSatish Balay #undef __FUNCT__ 94a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ" 10dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol) 113b2fbd54SBarry Smith { 123a40ed3dSBarry Smith PetscFunctionBegin; 133a40ed3dSBarry Smith 1429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not written"); 15aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG) 16d64ed03dSBarry Smith PetscFunctionReturn(0); 17d64ed03dSBarry Smith #endif 183b2fbd54SBarry Smith } 193b2fbd54SBarry Smith 2086bacbd2SBarry Smith 21dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 2286bacbd2SBarry Smith 2338baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 249cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*); 2538baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*); 2686bacbd2SBarry Smith 274a2ae208SSatish Balay #undef __FUNCT__ 284a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 2986bacbd2SBarry Smith /* ------------------------------------------------------------ 3086bacbd2SBarry Smith 3186bacbd2SBarry Smith This interface was contribed by Tony Caola 3286bacbd2SBarry Smith 3386bacbd2SBarry Smith This routine is an interface to the pivoting drop-tolerance 345bd2ddc7SBarry Smith ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of 355bd2ddc7SBarry Smith SPARSEKIT2. 365bd2ddc7SBarry Smith 375bd2ddc7SBarry Smith The SPARSEKIT2 routines used here are covered by the GNU 385bd2ddc7SBarry Smith copyright; see the file gnu in this directory. 3986bacbd2SBarry Smith 4086bacbd2SBarry Smith Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their 4186bacbd2SBarry Smith help in getting this routine ironed out. 4286bacbd2SBarry Smith 435bd2ddc7SBarry Smith The major drawback to this routine is that if info->fill is 445bd2ddc7SBarry Smith not large enough it fails rather than allocating more space; 455bd2ddc7SBarry Smith this can be fixed by hacking/improving the f2c version of 465bd2ddc7SBarry Smith Yousef Saad's code. 4786bacbd2SBarry Smith 4886bacbd2SBarry Smith ------------------------------------------------------------ 4986bacbd2SBarry Smith */ 50*7529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 5186bacbd2SBarry Smith { 5260ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE) 5360ec86bdSBarry Smith PetscFunctionBegin; 54e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\ 5560ec86bdSBarry Smith You can obtain the drop tolerance routines by installing PETSc from\n\ 5660ec86bdSBarry Smith www.mcs.anl.gov/petsc\n"); 5760ec86bdSBarry Smith #else 5886bacbd2SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 5907d2056aSBarry Smith IS iscolf,isicol,isirow; 6086bacbd2SBarry Smith PetscTruth reorder; 619cc1f4e3SBarry Smith PetscErrorCode ierr,sierr; 629cc1f4e3SBarry Smith PetscInt *c,*r,*ic,i,n = A->m; 6338baddfdSBarry Smith PetscInt *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j; 6438baddfdSBarry Smith PetscInt *ordcol,*iwk,*iperm,*jw; 6538baddfdSBarry Smith PetscInt jmax,lfill,job,*o_i,*o_j; 6687828ca2SBarry Smith PetscScalar *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a; 6754a8161fSBarry Smith PetscReal af; 6886bacbd2SBarry Smith 6986bacbd2SBarry Smith PetscFunctionBegin; 7086bacbd2SBarry Smith 7186bacbd2SBarry Smith if (info->dt == PETSC_DEFAULT) info->dt = .005; 7238baddfdSBarry Smith if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax); 7386bacbd2SBarry Smith if (info->dtcol == PETSC_DEFAULT) info->dtcol = .01; 7433259db3SBarry Smith if (info->fill == PETSC_DEFAULT) info->fill = ((double)(n*(info->dtcount+1)))/a->nz; 7538baddfdSBarry Smith lfill = (PetscInt)(info->dtcount/2.0); 7638baddfdSBarry Smith jmax = (PetscInt)(info->fill*a->nz); 7786bacbd2SBarry Smith 7886bacbd2SBarry Smith 7986bacbd2SBarry Smith /* ------------------------------------------------------------ 8086bacbd2SBarry Smith If reorder=.TRUE., then the original matrix has to be 8186bacbd2SBarry Smith reordered to reflect the user selected ordering scheme, and 8286bacbd2SBarry Smith then de-reordered so it is in it's original format. 8386bacbd2SBarry Smith Because Saad's dperm() is NOT in place, we have to copy 8486bacbd2SBarry Smith the original matrix and allocate more storage. . . 8586bacbd2SBarry Smith ------------------------------------------------------------ 8686bacbd2SBarry Smith */ 8786bacbd2SBarry Smith 8886bacbd2SBarry Smith /* set reorder to true if either isrow or iscol is not identity */ 8986bacbd2SBarry Smith ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr); 9086bacbd2SBarry Smith if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);} 9186bacbd2SBarry Smith reorder = PetscNot(reorder); 9286bacbd2SBarry Smith 9386bacbd2SBarry Smith 9486bacbd2SBarry Smith /* storage for ilu factor */ 9538baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr); 9638baddfdSBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr); 9787828ca2SBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr); 9838baddfdSBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr); 9986bacbd2SBarry Smith 10086bacbd2SBarry Smith /* ------------------------------------------------------------ 10186bacbd2SBarry Smith Make sure that everything is Fortran formatted (1-Based) 10286bacbd2SBarry Smith ------------------------------------------------------------ 10386bacbd2SBarry Smith */ 104b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 105b19713cbSBarry Smith old_j[i]++; 10686bacbd2SBarry Smith } 107b19713cbSBarry Smith for(i=0;i<n+1;i++) { 108b19713cbSBarry Smith old_i[i]++; 109b19713cbSBarry Smith }; 110010a20caSHong Zhang 11186bacbd2SBarry Smith 11286bacbd2SBarry Smith if (reorder) { 113c0c2c81eSBarry Smith ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr); 114c0c2c81eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 115c0c2c81eSBarry Smith for(i=0;i<n;i++) { 116c0c2c81eSBarry Smith r[i] = r[i]+1; 117c0c2c81eSBarry Smith c[i] = c[i]+1; 118c0c2c81eSBarry Smith } 11938baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr); 12038baddfdSBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr); 12187828ca2SBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr); 1225bd2ddc7SBarry Smith job = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job); 123c0c2c81eSBarry Smith for (i=0;i<n;i++) { 124c0c2c81eSBarry Smith r[i] = r[i]-1; 125c0c2c81eSBarry Smith c[i] = c[i]-1; 126c0c2c81eSBarry Smith } 127c0c2c81eSBarry Smith ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 128c0c2c81eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 129b19713cbSBarry Smith o_a = old_a2; 130b19713cbSBarry Smith o_j = old_j2; 131b19713cbSBarry Smith o_i = old_i2; 132b19713cbSBarry Smith } else { 133b19713cbSBarry Smith o_a = old_a; 134b19713cbSBarry Smith o_j = old_j; 135b19713cbSBarry Smith o_i = old_i; 13686bacbd2SBarry Smith } 13786bacbd2SBarry Smith 13886bacbd2SBarry Smith /* ------------------------------------------------------------ 13986bacbd2SBarry Smith Call Saad's ilutp() routine to generate the factorization 14086bacbd2SBarry Smith ------------------------------------------------------------ 14186bacbd2SBarry Smith */ 14286bacbd2SBarry Smith 14338baddfdSBarry Smith ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr); 14438baddfdSBarry Smith ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr); 14587828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr); 14686bacbd2SBarry Smith 14754a8161fSBarry 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); 1486849ba73SBarry Smith if (sierr) { 1496849ba73SBarry Smith switch (sierr) { 15077431f27SBarry Smith case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax); 15177431f27SBarry Smith case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax); 152e005ede5SBarry Smith case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered"); 153e005ede5SBarry Smith case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong"); 15477431f27SBarry Smith case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax); 15577431f27SBarry Smith default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr); 15686bacbd2SBarry Smith } 15786bacbd2SBarry Smith } 15886bacbd2SBarry Smith 15986bacbd2SBarry Smith ierr = PetscFree(w);CHKERRQ(ierr); 16086bacbd2SBarry Smith ierr = PetscFree(jw);CHKERRQ(ierr); 16186bacbd2SBarry Smith 16286bacbd2SBarry Smith /* ------------------------------------------------------------ 16386bacbd2SBarry Smith Saad's routine gives the result in Modified Sparse Row (msr) 16486bacbd2SBarry Smith Convert to Compressed Sparse Row format (csr) 16586bacbd2SBarry Smith ------------------------------------------------------------ 16686bacbd2SBarry Smith */ 16786bacbd2SBarry Smith 16887828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr); 16938baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr); 17086bacbd2SBarry Smith 1715bd2ddc7SBarry Smith SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk); 17286bacbd2SBarry Smith 17386bacbd2SBarry Smith ierr = PetscFree(iwk);CHKERRQ(ierr); 17486bacbd2SBarry Smith ierr = PetscFree(wk);CHKERRQ(ierr); 17586bacbd2SBarry Smith 17686bacbd2SBarry Smith if (reorder) { 17786bacbd2SBarry Smith ierr = PetscFree(old_a2);CHKERRQ(ierr); 17886bacbd2SBarry Smith ierr = PetscFree(old_j2);CHKERRQ(ierr); 17986bacbd2SBarry Smith ierr = PetscFree(old_i2);CHKERRQ(ierr); 180313ae024SBarry Smith } else { 181b19713cbSBarry Smith /* fix permutation of old_j that the factorization introduced */ 182f170005cSBarry Smith for (i=old_i[0]; i<old_i[n]; i++) { 183b19713cbSBarry Smith old_j[i-1] = iperm[old_j[i-1]-1]; 184b19713cbSBarry Smith } 185b19713cbSBarry Smith } 186b19713cbSBarry Smith 187b801d0f9SBarry Smith /* get rid of the shift to indices starting at 1 */ 18886bacbd2SBarry Smith for (i=0; i<n+1; i++) { 189b19713cbSBarry Smith old_i[i]--; 190b19713cbSBarry Smith } 191b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 192b19713cbSBarry Smith old_j[i]--; 193b19713cbSBarry Smith } 19486bacbd2SBarry Smith 195b19713cbSBarry Smith /* Make the factored matrix 0-based */ 19686bacbd2SBarry Smith for (i=0; i<n+1; i++) { 197b19713cbSBarry Smith new_i[i]--; 198b19713cbSBarry Smith } 199b19713cbSBarry Smith for (i=new_i[0];i<new_i[n];i++) { 200b19713cbSBarry Smith new_j[i]--; 201b19713cbSBarry Smith } 20286bacbd2SBarry Smith 20386bacbd2SBarry Smith /*-- due to the pivoting, we need to reorder iscol to correctly --*/ 20486bacbd2SBarry Smith /*-- permute the right-hand-side and solution vectors --*/ 2054c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2064c49b128SBarry Smith ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr); 207c0c2c81eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 20886bacbd2SBarry Smith for(i=0; i<n; i++) { 20986bacbd2SBarry Smith ordcol[i] = ic[iperm[i]-1]; 21086bacbd2SBarry Smith }; 21186bacbd2SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 21207d2056aSBarry Smith ierr = ISDestroy(isicol);CHKERRQ(ierr); 21386bacbd2SBarry Smith 214c0c2c81eSBarry Smith ierr = PetscFree(iperm);CHKERRQ(ierr); 215c0c2c81eSBarry Smith 216329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr); 21707d2056aSBarry Smith ierr = PetscFree(ordcol);CHKERRQ(ierr); 21886bacbd2SBarry Smith 21986bacbd2SBarry Smith /*----- put together the new matrix -----*/ 22086bacbd2SBarry Smith 221f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr); 222f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 223f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr); 22486bacbd2SBarry Smith (*fact)->factor = FACTOR_LU; 22586bacbd2SBarry Smith (*fact)->assembled = PETSC_TRUE; 22686bacbd2SBarry Smith 22786bacbd2SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 22886bacbd2SBarry Smith ierr = PetscFree(b->imax);CHKERRQ(ierr); 22986bacbd2SBarry Smith b->sorted = PETSC_FALSE; 23007d2056aSBarry Smith b->singlemalloc = PETSC_FALSE; 231b19713cbSBarry Smith /* the next line frees the default space generated by the MatCreate() */ 23286bacbd2SBarry Smith ierr = PetscFree(b->a);CHKERRQ(ierr); 23386bacbd2SBarry Smith ierr = PetscFree(b->ilen);CHKERRQ(ierr); 23486bacbd2SBarry Smith b->a = new_a; 23586bacbd2SBarry Smith b->j = new_j; 23686bacbd2SBarry Smith b->i = new_i; 23786bacbd2SBarry Smith b->ilen = 0; 23886bacbd2SBarry Smith b->imax = 0; 2391f9e874aSBarry Smith /* I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */ 240313ae024SBarry Smith b->row = isirow; 24186bacbd2SBarry Smith b->col = iscolf; 24287828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 24386bacbd2SBarry Smith b->maxnz = b->nz = new_i[n]; 24486bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 24586bacbd2SBarry Smith (*fact)->info.factor_mallocs = 0; 24686bacbd2SBarry Smith 24786bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 24886bacbd2SBarry Smith 249431e710aSBarry Smith af = ((double)b->nz)/((double)a->nz) + .001; 250b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af); 251b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af); 252b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af); 253b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n"); 254431e710aSBarry Smith 255*7529efd4SKris Buschelman ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr); 25671c2f376SKris Buschelman 25786bacbd2SBarry Smith PetscFunctionReturn(0); 25860ec86bdSBarry Smith #endif 25986bacbd2SBarry Smith } 26086bacbd2SBarry Smith 2614a2ae208SSatish Balay #undef __FUNCT__ 262b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ" 263dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B) 264289bc588SBarry Smith { 265416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 266289bc588SBarry Smith IS isicol; 2676849ba73SBarry Smith PetscErrorCode ierr; 26838baddfdSBarry Smith PetscInt *r,*ic,i,n=A->m,*ai=a->i,*aj=a->j; 2691393bc97SHong Zhang PetscInt *bi,*bj,*ajtmp; 2701393bc97SHong Zhang PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 2719e046878SBarry Smith PetscReal f; 2721393bc97SHong Zhang PetscInt nlnk,*lnk,k,*cols,**bi_ptr; 2731393bc97SHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2741393bc97SHong Zhang PetscBT lnkbt; 275289bc588SBarry Smith 2763a40ed3dSBarry Smith PetscFunctionBegin; 27729bbc08cSBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 2784c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 279ac121b3dSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 280ac121b3dSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 281289bc588SBarry Smith 282289bc588SBarry Smith /* get new row pointers */ 2831393bc97SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 2841393bc97SHong Zhang bi[0] = 0; 2851393bc97SHong Zhang 2861393bc97SHong Zhang /* bdiag is location of diagonal in factor */ 2871393bc97SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 2881393bc97SHong Zhang bdiag[0] = 0; 2891393bc97SHong Zhang 2901393bc97SHong Zhang /* linked list for storing column indices of the active row */ 2911393bc97SHong Zhang nlnk = n + 1; 2921393bc97SHong Zhang ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2931393bc97SHong Zhang 2941393bc97SHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt)+n*sizeof(PetscInt**),&cols);CHKERRQ(ierr); 2951393bc97SHong Zhang im = cols + n; 2961393bc97SHong Zhang bi_ptr = (PetscInt**)(im + n); 2971393bc97SHong Zhang 2981393bc97SHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 299d3d32019SBarry Smith f = info->fill; 3001393bc97SHong Zhang ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 3011393bc97SHong Zhang current_space = free_space; 302289bc588SBarry Smith 303289bc588SBarry Smith for (i=0; i<n; i++) { 3041393bc97SHong Zhang /* copy previous fill into linked list */ 305289bc588SBarry Smith nzi = 0; 3061393bc97SHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 3071393bc97SHong Zhang if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 3081393bc97SHong Zhang ajtmp = aj + ai[r[i]]; 3091393bc97SHong Zhang for (k=0; k<nnz; k++) cols[k] = ic[*(ajtmp+k)]; /* note: cols is not sorted when iscol!=indentity */ 3101393bc97SHong Zhang ierr = PetscLLAdd(nnz,cols,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 3111393bc97SHong Zhang nzi += nlnk; 3121393bc97SHong Zhang 3131393bc97SHong Zhang /* add pivot rows into linked list */ 3141393bc97SHong Zhang row = lnk[n]; 3151393bc97SHong Zhang while (row < i) { 3161393bc97SHong Zhang nzbd = bdiag[row] - bi[row] + 1; 3171393bc97SHong Zhang ajtmp = bi_ptr[row] + nzbd; 3181393bc97SHong Zhang nnz = im[row] - nzbd; /* num of columns with row<indices<=i */ 3191393bc97SHong Zhang im[row] = nzbd; 3201393bc97SHong Zhang ierr = PetscLLAddSortedLU(nnz,ajtmp,row,nlnk,lnk,lnkbt,i,nzbd);CHKERRQ(ierr); 3211393bc97SHong Zhang nzi += nlnk; 3221393bc97SHong Zhang im[row] += nzbd; /* update im[row]: num of cols with index<=i */ 3231393bc97SHong Zhang 3241393bc97SHong Zhang row = lnk[row]; 325289bc588SBarry Smith } 3261393bc97SHong Zhang 3271393bc97SHong Zhang bi[i+1] = bi[i] + nzi; 3281393bc97SHong Zhang im[i] = nzi; 3291393bc97SHong Zhang 3301393bc97SHong Zhang /* mark bdiag */ 3311393bc97SHong Zhang nzbd = 0; 3321393bc97SHong Zhang nnz = nzi; 3331393bc97SHong Zhang k = lnk[n]; 3341393bc97SHong Zhang while (nnz-- && k < i){ 3351393bc97SHong Zhang nzbd++; 3361393bc97SHong Zhang k = lnk[k]; 3371393bc97SHong Zhang } 3381393bc97SHong Zhang bdiag[i] = bi[i] + nzbd; 3391393bc97SHong Zhang 3401393bc97SHong Zhang /* if free space is not available, make more free space */ 3411393bc97SHong Zhang if (current_space->local_remaining<nzi) { 3421393bc97SHong Zhang nnz = (n - i)*nzi; /* estimated and max additional space needed */ 3431393bc97SHong Zhang ierr = GetMoreSpace(nnz,¤t_space);CHKERRQ(ierr); 3441393bc97SHong Zhang reallocs++; 3451393bc97SHong Zhang } 3461393bc97SHong Zhang 3471393bc97SHong Zhang /* copy data into free space, then initialize lnk */ 3481393bc97SHong Zhang ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 3491393bc97SHong Zhang bi_ptr[i] = current_space->array; 3501393bc97SHong Zhang current_space->array += nzi; 3511393bc97SHong Zhang current_space->local_used += nzi; 3521393bc97SHong Zhang current_space->local_remaining -= nzi; 353289bc588SBarry Smith } 354464e8d44SSatish Balay if (ai[n] != 0) { 3551393bc97SHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 356418422e8SSatish Balay PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af); 357b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af); 358b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af); 359b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"); 36005bf559cSSatish Balay } else { 361b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n"); 36205bf559cSSatish Balay } 3632fb3ed76SBarry Smith 364898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 365898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3661987afe7SBarry Smith 3671393bc97SHong Zhang /* destroy list of free space and other temporary array(s) */ 3681393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 3691393bc97SHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 3701393bc97SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 3711393bc97SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 372289bc588SBarry Smith 373289bc588SBarry Smith /* put together the new matrix */ 374f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr); 375f204ca49SKris Buschelman ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr); 376f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr); 37752e6d16bSBarry Smith ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr); 378416022c9SBarry Smith b = (Mat_SeqAIJ*)(*B)->data; 379606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 3807c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 381e8d4e0b9SBarry Smith /* the next line frees the default space generated by the Create() */ 382606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 383606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 3841393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 3851393bc97SHong Zhang b->j = bj; 3861393bc97SHong Zhang b->i = bi; 3871393bc97SHong Zhang b->diag = bdiag; 388416022c9SBarry Smith b->ilen = 0; 389416022c9SBarry Smith b->imax = 0; 390416022c9SBarry Smith b->row = isrow; 391416022c9SBarry Smith b->col = iscol; 392c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 393c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 39482bf6240SBarry Smith b->icol = isicol; 39587828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 3961393bc97SHong Zhang 3971393bc97SHong Zhang /* In b structure: Free imax, ilen, old a, old j. Allocate solve_work, new a, new j */ 3981393bc97SHong Zhang ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 3991393bc97SHong Zhang b->maxnz = b->nz = bi[n] ; 4007fd98bd6SLois Curfman McInnes 401329f5518SBarry Smith (*B)->factor = FACTOR_LU; 402418422e8SSatish Balay (*B)->info.factor_mallocs = reallocs; 403ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 404703deb49SSatish Balay 405e93ccc4dSBarry Smith if (ai[n] != 0) { 4061393bc97SHong Zhang (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 407eea2913aSSatish Balay } else { 408eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 409eea2913aSSatish Balay } 4104846f1f5SKris Buschelman ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr); 41171c2f376SKris Buschelman (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 4123a40ed3dSBarry Smith PetscFunctionReturn(0); 413289bc588SBarry Smith } 4141393bc97SHong Zhang 4150a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 4164a2ae208SSatish Balay #undef __FUNCT__ 4174a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 418af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 419289bc588SBarry Smith { 420416022c9SBarry Smith Mat C=*B; 421416022c9SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 42282bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 4236849ba73SBarry Smith PetscErrorCode ierr; 424b3bf805bSHong Zhang PetscInt *r,*ic,i,j,n=A->m,*bi=b->i,*bj=b->j; 425b3bf805bSHong Zhang PetscInt *ajtmp,*bjtmp,nz,row,*ics; 426d2276718SHong Zhang PetscInt *diag_offset = b->diag,diag,*pj; 42787828ca2SBarry Smith PetscScalar *rtmp,*v,*pc,multiplier,*pv,*rtmps; 428540703acSHong Zhang PetscReal zeropivot,rs,d,shift_nz; 429d2276718SHong Zhang PetscReal row_shift,shift_top=0.; 430d2276718SHong Zhang PetscTruth shift_pd; 431b3bf805bSHong Zhang LUShift_Ctx sctx; 432d2276718SHong Zhang PetscInt newshift; 433289bc588SBarry Smith 4343a40ed3dSBarry Smith PetscFunctionBegin; 435540703acSHong Zhang shift_nz = info->shiftnz; 4360a29876aSHong Zhang shift_pd = info->shiftpd; 4370a29876aSHong Zhang zeropivot = info->zeropivot; 4380a29876aSHong Zhang 43978b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 44078b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 44187828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 44287828ca2SBarry Smith ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 443010a20caSHong Zhang rtmps = rtmp; ics = ic; 444289bc588SBarry Smith 4456cc28720Svictorle if (!a->diag) { 4460cf777b8SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 4470cf777b8SBarry Smith } 448e24cfe64SHong Zhang /* if both shift schemes are chosen by user, only use shift_pd */ 449540703acSHong Zhang if (shift_pd && shift_nz) shift_nz = 0.0; 450e24cfe64SHong Zhang if (shift_pd) { /* set shift_top=max{row_shift} */ 45138baddfdSBarry Smith PetscInt *aai = a->i,*ddiag = a->diag; 4520cf777b8SBarry Smith shift_top = 0; 4536cc28720Svictorle for (i=0; i<n; i++) { 454a0a356efSHong Zhang d = PetscRealPart((a->a)[ddiag[i]]); 4556cc28720Svictorle /* calculate amt of shift needed for this row */ 4566c037d1bSvictorle if (d<=0) { 4573aeef889SHong Zhang row_shift = 0; 4583aeef889SHong Zhang } else { 4593aeef889SHong Zhang row_shift = -2*d; 4603aeef889SHong Zhang } 461010a20caSHong Zhang v = a->a+aai[i]; 462e24cfe64SHong Zhang nz = aai[i+1] - aai[i]; 463e24cfe64SHong Zhang for (j=0; j<nz; j++) 4646cc28720Svictorle row_shift += PetscAbsScalar(v[j]); 4656cc28720Svictorle if (row_shift>shift_top) shift_top = row_shift; 4666cc28720Svictorle } 467a0a356efSHong Zhang if (shift_top == 0.0) shift_top += 1.e-12; 468d2276718SHong Zhang sctx.shift_top = shift_top; 469d2276718SHong Zhang sctx.nshift_max = 5; 470d2276718SHong Zhang sctx.shift_lo = 0.; 471d2276718SHong Zhang sctx.shift_hi = 1.; 472a0a356efSHong Zhang } 473a0a356efSHong Zhang 474a0a356efSHong Zhang sctx.shift_amount = 0; 475a0a356efSHong Zhang sctx.nshift = 0; 476085256b3SBarry Smith do { 477d2276718SHong Zhang sctx.lushift = PETSC_FALSE; 478289bc588SBarry Smith for (i=0; i<n; i++){ 479b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 480b3bf805bSHong Zhang bjtmp = bj + bi[i]; 481b3bf805bSHong Zhang for (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0; 482289bc588SBarry Smith 483289bc588SBarry Smith /* load in initial (unfactored row) */ 484416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 485b3bf805bSHong Zhang ajtmp = a->j + a->i[r[i]]; 486010a20caSHong Zhang v = a->a + a->i[r[i]]; 487085256b3SBarry Smith for (j=0; j<nz; j++) { 488b3bf805bSHong Zhang rtmp[ics[ajtmp[j]]] = v[j]; 489085256b3SBarry Smith } 490d2276718SHong Zhang rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 491289bc588SBarry Smith 492b3bf805bSHong Zhang row = *bjtmp++; 493f2582111SSatish Balay while (row < i) { 4948c37ef55SBarry Smith pc = rtmp + row; 4958c37ef55SBarry Smith if (*pc != 0.0) { 496010a20caSHong Zhang pv = b->a + diag_offset[row]; 497010a20caSHong Zhang pj = b->j + diag_offset[row] + 1; 49835aab85fSBarry Smith multiplier = *pc / *pv++; 4998c37ef55SBarry Smith *pc = multiplier; 500b3bf805bSHong Zhang nz = bi[row+1] - diag_offset[row] - 1; 501f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 502efee365bSSatish Balay ierr = PetscLogFlops(2*nz);CHKERRQ(ierr); 503289bc588SBarry Smith } 504b3bf805bSHong Zhang row = *bjtmp++; 505289bc588SBarry Smith } 506416022c9SBarry Smith /* finished row so stick it into b->a */ 507b3bf805bSHong Zhang pv = b->a + bi[i] ; 508b3bf805bSHong Zhang pj = b->j + bi[i] ; 509b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 510b3bf805bSHong Zhang diag = diag_offset[i] - bi[i]; 511d3d32019SBarry Smith rs = 0.0; 512d3d32019SBarry Smith for (j=0; j<nz; j++) { 513d3d32019SBarry Smith pv[j] = rtmps[pj[j]]; 514d3d32019SBarry Smith if (j != diag) rs += PetscAbsScalar(pv[j]); 515d3d32019SBarry Smith } 516d2276718SHong Zhang 517b3bf805bSHong Zhang /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 518d2276718SHong Zhang sctx.rs = rs; 519d2276718SHong Zhang sctx.pv = pv[diag]; 5203e8c821fSHong Zhang ierr = MatLUCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 521d2276718SHong Zhang if (newshift == 1){ 522b3bf805bSHong Zhang break; /* sctx.shift_amount is updated */ 523d2276718SHong Zhang } else if (newshift == -1){ 524d2276718SHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs); 525d708749eSBarry Smith } 52635aab85fSBarry Smith } 527d2276718SHong Zhang 528d2276718SHong Zhang if (shift_pd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 5296cc28720Svictorle /* 5306c037d1bSvictorle * if no shift in this attempt & shifting & started shifting & can refine, 5316cc28720Svictorle * then try lower shift 5326cc28720Svictorle */ 533d2276718SHong Zhang sctx.shift_hi = info->shift_fraction; 534d2276718SHong Zhang info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 535d2276718SHong Zhang sctx.shift_amount = info->shift_fraction * sctx.shift_top; 536d2276718SHong Zhang sctx.lushift = PETSC_TRUE; 537d2276718SHong Zhang sctx.nshift++; 5386cc28720Svictorle } 539d2276718SHong Zhang } while (sctx.lushift); 5400f11f9aeSSatish Balay 54135aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 54235aab85fSBarry Smith for (i=0; i<n; i++) { 543010a20caSHong Zhang b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]]; 54435aab85fSBarry Smith } 54535aab85fSBarry Smith 546606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 54778b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 54878b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 549416022c9SBarry Smith C->factor = FACTOR_LU; 5507dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 551c456f294SBarry Smith C->assembled = PETSC_TRUE; 552efee365bSSatish Balay ierr = PetscLogFlops(C->n);CHKERRQ(ierr); 553d2276718SHong Zhang if (sctx.nshift){ 554540703acSHong Zhang if (shift_nz) { 555540703acSHong Zhang PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 5560697b6caSHong Zhang } else if (shift_pd) { 557d2276718SHong Zhang b->lu_shift_fraction = info->shift_fraction; 558a0a356efSHong Zhang PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: 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,shift_top); 5596cc28720Svictorle } 5600697b6caSHong Zhang } 5613a40ed3dSBarry Smith PetscFunctionReturn(0); 562289bc588SBarry Smith } 5630889b5dcSvictorle 5640889b5dcSvictorle #undef __FUNCT__ 5650889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ" 566dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A) 5670889b5dcSvictorle { 5680889b5dcSvictorle PetscFunctionBegin; 5690889b5dcSvictorle A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ; 5700889b5dcSvictorle A->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 5710889b5dcSvictorle PetscFunctionReturn(0); 5720889b5dcSvictorle } 5730889b5dcSvictorle 5740889b5dcSvictorle 5750a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5764a2ae208SSatish Balay #undef __FUNCT__ 5774a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 578dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 579da3a660dSBarry Smith { 580dfbe8321SBarry Smith PetscErrorCode ierr; 581416022c9SBarry Smith Mat C; 582416022c9SBarry Smith 5833a40ed3dSBarry Smith PetscFunctionBegin; 5849e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 585af281ebdSHong Zhang ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 586273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 58752e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr); 5883a40ed3dSBarry Smith PetscFunctionReturn(0); 589da3a660dSBarry Smith } 5900a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5914a2ae208SSatish Balay #undef __FUNCT__ 5924a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 593dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 5948c37ef55SBarry Smith { 595416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 596416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 5976849ba73SBarry Smith PetscErrorCode ierr; 59838baddfdSBarry Smith PetscInt *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 59938baddfdSBarry Smith PetscInt nz,*rout,*cout; 60087828ca2SBarry Smith PetscScalar *x,*b,*tmp,*tmps,*aa = a->a,sum,*v; 6018c37ef55SBarry Smith 6023a40ed3dSBarry Smith PetscFunctionBegin; 6033a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 60491bf9adeSBarry Smith 6051ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 6061ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 607416022c9SBarry Smith tmp = a->solve_work; 6088c37ef55SBarry Smith 6093b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6103b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 6118c37ef55SBarry Smith 6128c37ef55SBarry Smith /* forward solve the lower triangular */ 6138c37ef55SBarry Smith tmp[0] = b[*r++]; 614010a20caSHong Zhang tmps = tmp; 6158c37ef55SBarry Smith for (i=1; i<n; i++) { 616010a20caSHong Zhang v = aa + ai[i] ; 617010a20caSHong Zhang vi = aj + ai[i] ; 61853e7425aSBarry Smith nz = a->diag[i] - ai[i]; 61953e7425aSBarry Smith sum = b[*r++]; 620ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 6218c37ef55SBarry Smith tmp[i] = sum; 6228c37ef55SBarry Smith } 6238c37ef55SBarry Smith 6248c37ef55SBarry Smith /* backward solve the upper triangular */ 6258c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 626010a20caSHong Zhang v = aa + a->diag[i] + 1; 627010a20caSHong Zhang vi = aj + a->diag[i] + 1; 628416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 6298c37ef55SBarry Smith sum = tmp[i]; 630ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 631010a20caSHong Zhang x[*c--] = tmp[i] = sum*aa[a->diag[i]]; 6328c37ef55SBarry Smith } 6338c37ef55SBarry Smith 6343b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6353b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 6361ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 6371ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 638efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr); 6393a40ed3dSBarry Smith PetscFunctionReturn(0); 6408c37ef55SBarry Smith } 641026e39d0SSatish Balay 642930ae53cSBarry Smith /* ----------------------------------------------------------- */ 6434a2ae208SSatish Balay #undef __FUNCT__ 6444a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 645dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 646930ae53cSBarry Smith { 647930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 6486849ba73SBarry Smith PetscErrorCode ierr; 64938baddfdSBarry Smith PetscInt n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag; 650362ced78SSatish Balay PetscScalar *x,*b,*aa = a->a; 651aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 65238baddfdSBarry Smith PetscInt adiag_i,i,*vi,nz,ai_i; 653362ced78SSatish Balay PetscScalar *v,sum; 654d85afc42SSatish Balay #endif 655930ae53cSBarry Smith 6563a40ed3dSBarry Smith PetscFunctionBegin; 6573a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 658930ae53cSBarry Smith 6591ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 6601ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 661930ae53cSBarry Smith 662aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 6631c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 6641c4feecaSSatish Balay #else 665930ae53cSBarry Smith /* forward solve the lower triangular */ 666930ae53cSBarry Smith x[0] = b[0]; 667930ae53cSBarry Smith for (i=1; i<n; i++) { 668930ae53cSBarry Smith ai_i = ai[i]; 669930ae53cSBarry Smith v = aa + ai_i; 670930ae53cSBarry Smith vi = aj + ai_i; 671930ae53cSBarry Smith nz = adiag[i] - ai_i; 672930ae53cSBarry Smith sum = b[i]; 673930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 674930ae53cSBarry Smith x[i] = sum; 675930ae53cSBarry Smith } 676930ae53cSBarry Smith 677930ae53cSBarry Smith /* backward solve the upper triangular */ 678930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 679930ae53cSBarry Smith adiag_i = adiag[i]; 680d708749eSBarry Smith v = aa + adiag_i + 1; 681d708749eSBarry Smith vi = aj + adiag_i + 1; 682930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 683930ae53cSBarry Smith sum = x[i]; 684930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 685930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 686930ae53cSBarry Smith } 6871c4feecaSSatish Balay #endif 688efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr); 6891ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 6901ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6913a40ed3dSBarry Smith PetscFunctionReturn(0); 692930ae53cSBarry Smith } 693930ae53cSBarry Smith 6944a2ae208SSatish Balay #undef __FUNCT__ 6954a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 696dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 697da3a660dSBarry Smith { 698416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 699416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 7006849ba73SBarry Smith PetscErrorCode ierr; 70138baddfdSBarry Smith PetscInt *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 70238baddfdSBarry Smith PetscInt nz,*rout,*cout; 70387828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,sum,*v; 704da3a660dSBarry Smith 7053a40ed3dSBarry Smith PetscFunctionBegin; 70678b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 707da3a660dSBarry Smith 7081ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 7091ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 710416022c9SBarry Smith tmp = a->solve_work; 711da3a660dSBarry Smith 7123b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7133b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 714da3a660dSBarry Smith 715da3a660dSBarry Smith /* forward solve the lower triangular */ 716da3a660dSBarry Smith tmp[0] = b[*r++]; 717da3a660dSBarry Smith for (i=1; i<n; i++) { 718010a20caSHong Zhang v = aa + ai[i] ; 719010a20caSHong Zhang vi = aj + ai[i] ; 720416022c9SBarry Smith nz = a->diag[i] - ai[i]; 721da3a660dSBarry Smith sum = b[*r++]; 722010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 723da3a660dSBarry Smith tmp[i] = sum; 724da3a660dSBarry Smith } 725da3a660dSBarry Smith 726da3a660dSBarry Smith /* backward solve the upper triangular */ 727da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 728010a20caSHong Zhang v = aa + a->diag[i] + 1; 729010a20caSHong Zhang vi = aj + a->diag[i] + 1; 730416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 731da3a660dSBarry Smith sum = tmp[i]; 732010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 733010a20caSHong Zhang tmp[i] = sum*aa[a->diag[i]]; 734da3a660dSBarry Smith x[*c--] += tmp[i]; 735da3a660dSBarry Smith } 736da3a660dSBarry Smith 7373b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7383b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7391ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 7401ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 741efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 742898183ecSLois Curfman McInnes 7433a40ed3dSBarry Smith PetscFunctionReturn(0); 744da3a660dSBarry Smith } 745da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 7464a2ae208SSatish Balay #undef __FUNCT__ 7474a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 748dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 749da3a660dSBarry Smith { 750416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7512235e667SBarry Smith IS iscol = a->col,isrow = a->row; 7526849ba73SBarry Smith PetscErrorCode ierr; 75338baddfdSBarry Smith PetscInt *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 75438baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 75587828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v,s1; 756da3a660dSBarry Smith 7573a40ed3dSBarry Smith PetscFunctionBegin; 7581ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 7591ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 760416022c9SBarry Smith tmp = a->solve_work; 761da3a660dSBarry Smith 7622235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7632235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 764da3a660dSBarry Smith 765da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 7662235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 767da3a660dSBarry Smith 768da3a660dSBarry Smith /* forward solve the U^T */ 769da3a660dSBarry Smith for (i=0; i<n; i++) { 770010a20caSHong Zhang v = aa + diag[i] ; 771010a20caSHong Zhang vi = aj + diag[i] + 1; 772f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 773c38d4ed2SBarry Smith s1 = tmp[i]; 774ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 775da3a660dSBarry Smith while (nz--) { 776010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*s1; 777da3a660dSBarry Smith } 778c38d4ed2SBarry Smith tmp[i] = s1; 779da3a660dSBarry Smith } 780da3a660dSBarry Smith 781da3a660dSBarry Smith /* backward solve the L^T */ 782da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 783010a20caSHong Zhang v = aa + diag[i] - 1 ; 784010a20caSHong Zhang vi = aj + diag[i] - 1 ; 785f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 786f1af5d2fSBarry Smith s1 = tmp[i]; 787da3a660dSBarry Smith while (nz--) { 788010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*s1; 789da3a660dSBarry Smith } 790da3a660dSBarry Smith } 791da3a660dSBarry Smith 792da3a660dSBarry Smith /* copy tmp into x according to permutation */ 793da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 794da3a660dSBarry Smith 7952235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7962235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7971ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 7981ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 799da3a660dSBarry Smith 800efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz-A->n);CHKERRQ(ierr); 8013a40ed3dSBarry Smith PetscFunctionReturn(0); 802da3a660dSBarry Smith } 803da3a660dSBarry Smith 8044a2ae208SSatish Balay #undef __FUNCT__ 8054a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 806dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 807da3a660dSBarry Smith { 808416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 8092235e667SBarry Smith IS iscol = a->col,isrow = a->row; 8106849ba73SBarry Smith PetscErrorCode ierr; 81138baddfdSBarry Smith PetscInt *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 81238baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 81387828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v; 8146abc6512SBarry Smith 8153a40ed3dSBarry Smith PetscFunctionBegin; 81623bc6035SMatthew Knepley if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);} 8176abc6512SBarry Smith 8181ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 8191ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 820416022c9SBarry Smith tmp = a->solve_work; 8216abc6512SBarry Smith 8222235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8232235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 8246abc6512SBarry Smith 8256abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 8262235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 8276abc6512SBarry Smith 8286abc6512SBarry Smith /* forward solve the U^T */ 8296abc6512SBarry Smith for (i=0; i<n; i++) { 830010a20caSHong Zhang v = aa + diag[i] ; 831010a20caSHong Zhang vi = aj + diag[i] + 1; 832f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 8336abc6512SBarry Smith tmp[i] *= *v++; 8346abc6512SBarry Smith while (nz--) { 835010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*tmp[i]; 8366abc6512SBarry Smith } 8376abc6512SBarry Smith } 8386abc6512SBarry Smith 8396abc6512SBarry Smith /* backward solve the L^T */ 8406abc6512SBarry Smith for (i=n-1; i>=0; i--){ 841010a20caSHong Zhang v = aa + diag[i] - 1 ; 842010a20caSHong Zhang vi = aj + diag[i] - 1 ; 843f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 8446abc6512SBarry Smith while (nz--) { 845010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*tmp[i]; 8466abc6512SBarry Smith } 8476abc6512SBarry Smith } 8486abc6512SBarry Smith 8496abc6512SBarry Smith /* copy tmp into x according to permutation */ 8506abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 8516abc6512SBarry Smith 8522235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8532235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 8541ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 8551ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8566abc6512SBarry Smith 857efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 8583a40ed3dSBarry Smith PetscFunctionReturn(0); 859da3a660dSBarry Smith } 8609e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 861dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat); 86208480c60SBarry Smith 8634a2ae208SSatish Balay #undef __FUNCT__ 8644a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 865dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 8669e25ed09SBarry Smith { 867416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 8689e25ed09SBarry Smith IS isicol; 8696849ba73SBarry Smith PetscErrorCode ierr; 8705a8e39fbSHong Zhang PetscInt *r,*ic,n=A->m,*ai=a->i,*aj=a->j; 8715a8e39fbSHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 8725a8e39fbSHong Zhang PetscInt *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0; 8735a8e39fbSHong Zhang PetscInt i,levels,diagonal_fill; 87477c4ece6SBarry Smith PetscTruth col_identity,row_identity; 875329f5518SBarry Smith PetscReal f; 8765a8e39fbSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 8775a8e39fbSHong Zhang PetscBT lnkbt; 8785a8e39fbSHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 8795a8e39fbSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 8805a8e39fbSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 8819e25ed09SBarry Smith 8823a40ed3dSBarry Smith PetscFunctionBegin; 8836cf997b0SBarry Smith f = info->fill; 88438baddfdSBarry Smith levels = (PetscInt)info->levels; 88538baddfdSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 8864c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 88782bf6240SBarry Smith 88808480c60SBarry Smith /* special case that simply copies fill pattern */ 889be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 890be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 89186bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 8922e8a6d31SBarry Smith ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 89308480c60SBarry Smith (*fact)->factor = FACTOR_LU; 89408480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 89508480c60SBarry Smith if (!b->diag) { 8967c922b88SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 89708480c60SBarry Smith } 8987c922b88SBarry Smith ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 89908480c60SBarry Smith b->row = isrow; 90008480c60SBarry Smith b->col = iscol; 90182bf6240SBarry Smith b->icol = isicol; 90287828ca2SBarry Smith ierr = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 903f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 904c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 905c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 9063a40ed3dSBarry Smith PetscFunctionReturn(0); 90708480c60SBarry Smith } 90808480c60SBarry Smith 909898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 910898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 9119e25ed09SBarry Smith 9129e25ed09SBarry Smith /* get new row pointers */ 9135a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 9145a8e39fbSHong Zhang bi[0] = 0; 9155a8e39fbSHong Zhang /* bdiag is location of diagonal in factor */ 9165a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 9175a8e39fbSHong Zhang bdiag[0] = 0; 9186cf997b0SBarry Smith 9195a8e39fbSHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr); 9205a8e39fbSHong Zhang bjlvl_ptr = (PetscInt**)(bj_ptr + n); 9215a8e39fbSHong Zhang 9225a8e39fbSHong Zhang /* create a linked list for storing column indices of the active row */ 9235a8e39fbSHong Zhang nlnk = n + 1; 9245a8e39fbSHong Zhang ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 9255a8e39fbSHong Zhang 9265a8e39fbSHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 9275a8e39fbSHong Zhang ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 9285a8e39fbSHong Zhang current_space = free_space; 9295a8e39fbSHong Zhang ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 9305a8e39fbSHong Zhang current_space_lvl = free_space_lvl; 9315a8e39fbSHong Zhang 9325a8e39fbSHong Zhang for (i=0; i<n; i++) { 9335a8e39fbSHong Zhang nzi = 0; 9346cf997b0SBarry Smith /* copy current row into linked list */ 9355a8e39fbSHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 9365a8e39fbSHong Zhang if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 9375a8e39fbSHong Zhang cols = aj + ai[r[i]]; 9385a8e39fbSHong Zhang lnk[i] = -1; /* marker to indicate if diagonal exists */ 9395a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 9405a8e39fbSHong Zhang nzi += nlnk; 9416cf997b0SBarry Smith 9426cf997b0SBarry Smith /* make sure diagonal entry is included */ 9435a8e39fbSHong Zhang if (diagonal_fill && lnk[i] == -1) { 9446cf997b0SBarry Smith fm = n; 9455a8e39fbSHong Zhang while (lnk[fm] < i) fm = lnk[fm]; 9465a8e39fbSHong Zhang lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 9475a8e39fbSHong Zhang lnk[fm] = i; 9485a8e39fbSHong Zhang lnk_lvl[i] = 0; 9495a8e39fbSHong Zhang nzi++; dcount++; 9506cf997b0SBarry Smith } 9516cf997b0SBarry Smith 9525a8e39fbSHong Zhang /* add pivot rows into the active row */ 9535a8e39fbSHong Zhang nzbd = 0; 9545a8e39fbSHong Zhang prow = lnk[n]; 9555a8e39fbSHong Zhang while (prow < i) { 9565a8e39fbSHong Zhang nnz = bdiag[prow]; 9575a8e39fbSHong Zhang cols = bj_ptr[prow] + nnz + 1; 9585a8e39fbSHong Zhang cols_lvl = bjlvl_ptr[prow] + nnz + 1; 9595a8e39fbSHong Zhang nnz = bi[prow+1] - bi[prow] - nnz - 1; 9605a8e39fbSHong Zhang ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 9615a8e39fbSHong Zhang nzi += nlnk; 9625a8e39fbSHong Zhang prow = lnk[prow]; 9635a8e39fbSHong Zhang nzbd++; 96456beaf04SBarry Smith } 9655a8e39fbSHong Zhang bdiag[i] = nzbd; 9665a8e39fbSHong Zhang bi[i+1] = bi[i] + nzi; 967ecf371e4SBarry Smith 9685a8e39fbSHong Zhang /* if free space is not available, make more free space */ 9695a8e39fbSHong Zhang if (current_space->local_remaining<nzi) { 9705a8e39fbSHong Zhang nnz = nzi*(n - i); /* estimated and max additional space needed */ 9715a8e39fbSHong Zhang ierr = GetMoreSpace(nnz,¤t_space);CHKERRQ(ierr); 9725a8e39fbSHong Zhang ierr = GetMoreSpace(nnz,¤t_space_lvl);CHKERRQ(ierr); 9735a8e39fbSHong Zhang reallocs++; 9745a8e39fbSHong Zhang } 975ecf371e4SBarry Smith 9765a8e39fbSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 9775a8e39fbSHong Zhang ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 9785a8e39fbSHong Zhang bj_ptr[i] = current_space->array; 9795a8e39fbSHong Zhang bjlvl_ptr[i] = current_space_lvl->array; 9805a8e39fbSHong Zhang 9815a8e39fbSHong Zhang /* make sure the active row i has diagonal entry */ 9825a8e39fbSHong Zhang if (*(bj_ptr[i]+bdiag[i]) != i) { 98377431f27SBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 9845a8e39fbSHong Zhang try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i); 9856cf997b0SBarry Smith } 9865a8e39fbSHong Zhang 9875a8e39fbSHong Zhang current_space->array += nzi; 9885a8e39fbSHong Zhang current_space->local_used += nzi; 9895a8e39fbSHong Zhang current_space->local_remaining -= nzi; 9905a8e39fbSHong Zhang current_space_lvl->array += nzi; 9915a8e39fbSHong Zhang current_space_lvl->local_used += nzi; 9925a8e39fbSHong Zhang current_space_lvl->local_remaining -= nzi; 9939e25ed09SBarry Smith } 9945a8e39fbSHong Zhang 995898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 996898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 9975a8e39fbSHong Zhang 9985a8e39fbSHong Zhang /* destroy list of free space and other temporary arrays */ 9995a8e39fbSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 10005a8e39fbSHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 10015a8e39fbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 10025a8e39fbSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 10035a8e39fbSHong Zhang ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 10049e25ed09SBarry Smith 1005f64065bbSBarry Smith { 10065a8e39fbSHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1007418422e8SSatish Balay PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af); 1008c0d6ac57SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af); 1009c0d6ac57SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af); 1010b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"); 1011335d9088SBarry Smith if (diagonal_fill) { 101277431f27SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount); 1013335d9088SBarry Smith } 1014f64065bbSBarry Smith } 1015bbb6d6a8SBarry Smith 10169e25ed09SBarry Smith /* put together the new matrix */ 1017f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr); 1018f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 1019f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr); 102052e6d16bSBarry Smith ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr); 1021416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 1022606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 10237c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 10245a8e39fbSHong Zhang len = (bi[n] )*sizeof(PetscScalar); 10259e25ed09SBarry Smith /* the next line frees the default space generated by the Create() */ 1026606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 1027606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 1028b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 10295a8e39fbSHong Zhang b->j = bj; 10305a8e39fbSHong Zhang b->i = bi; 10315a8e39fbSHong Zhang for (i=0; i<n; i++) bdiag[i] += bi[i]; 10325a8e39fbSHong Zhang b->diag = bdiag; 1033416022c9SBarry Smith b->ilen = 0; 1034416022c9SBarry Smith b->imax = 0; 1035416022c9SBarry Smith b->row = isrow; 1036416022c9SBarry Smith b->col = iscol; 1037c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1038c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 103982bf6240SBarry Smith b->icol = isicol; 104087828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1041416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 10425a8e39fbSHong Zhang Allocate bdiag, solve_work, new a, new j */ 10435a8e39fbSHong Zhang ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 10445a8e39fbSHong Zhang b->maxnz = b->nz = bi[n] ; 10459e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 1046418422e8SSatish Balay (*fact)->info.factor_mallocs = reallocs; 1047ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 10485a8e39fbSHong Zhang (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 104971c2f376SKris Buschelman 1050fd62cfe2SHong Zhang ierr = MatILUFactorSymbolic_Inode(A,info,isrow,iscol,fact);CHKERRQ(ierr); 105171c2f376SKris Buschelman (*fact)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 10524846f1f5SKris Buschelman 10533a40ed3dSBarry Smith PetscFunctionReturn(0); 10549e25ed09SBarry Smith } 1055d5d45c9bSBarry Smith 10563c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1057a6175056SHong Zhang #undef __FUNCT__ 1058a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1059af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 1060a6175056SHong Zhang { 10612ed133dbSHong Zhang Mat C = *B; 10620968510aSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 10632ed133dbSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 10642ed133dbSHong Zhang IS ip=b->row; 10652ed133dbSHong Zhang PetscErrorCode ierr; 10662ed133dbSHong Zhang PetscInt *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol; 10672ed133dbSHong Zhang PetscInt *ai=a->i,*aj=a->j; 10682ed133dbSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1069622e2028SHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1070540703acSHong Zhang PetscReal zeropivot,rs,shiftnz; 1071540703acSHong Zhang PetscTruth shiftpd; 1072540703acSHong Zhang ChShift_Ctx sctx; 1073540703acSHong Zhang PetscInt newshift; 1074d5d45c9bSBarry Smith 1075a6175056SHong Zhang PetscFunctionBegin; 1076540703acSHong Zhang shiftnz = info->shiftnz; 1077540703acSHong Zhang shiftpd = info->shiftpd; 1078ee45ca4aSHong Zhang zeropivot = info->zeropivot; 1079ee45ca4aSHong Zhang 10802ed133dbSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 10812ed133dbSHong Zhang 10822ed133dbSHong Zhang /* initialization */ 10832ed133dbSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 10842ed133dbSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 10852ed133dbSHong Zhang jl = il + mbs; 10862ed133dbSHong Zhang rtmp = (MatScalar*)(jl + mbs); 10872ed133dbSHong Zhang 1088540703acSHong Zhang sctx.shift_amount = 0; 1089540703acSHong Zhang sctx.nshift = 0; 10902ed133dbSHong Zhang do { 1091540703acSHong Zhang sctx.chshift = PETSC_FALSE; 10922ed133dbSHong Zhang for (i=0; i<mbs; i++) { 10932ed133dbSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1094a6175056SHong Zhang } 10952ed133dbSHong Zhang 10962ed133dbSHong Zhang for (k = 0; k<mbs; k++){ 1097c05c3958SHong Zhang bval = ba + bi[k]; 10982ed133dbSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 10992ed133dbSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 11002ed133dbSHong Zhang for (j = jmin; j < jmax; j++){ 11012ed133dbSHong Zhang col = rip[aj[j]]; 11022ed133dbSHong Zhang if (col >= k){ /* only take upper triangular entry */ 11032ed133dbSHong Zhang rtmp[col] = aa[j]; 1104c05c3958SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 11052ed133dbSHong Zhang } 11062ed133dbSHong Zhang } 110739e3d630SHong Zhang /* shift the diagonal of the matrix */ 1108540703acSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 11092ed133dbSHong Zhang 11102ed133dbSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 11112ed133dbSHong Zhang dk = rtmp[k]; 11122ed133dbSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 11132ed133dbSHong Zhang 11142ed133dbSHong Zhang while (i < k){ 11152ed133dbSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 11162ed133dbSHong Zhang 11172ed133dbSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 11182ed133dbSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 11192ed133dbSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 11202ed133dbSHong Zhang dk += uikdi*ba[ili]; 11212ed133dbSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 11222ed133dbSHong Zhang 11232ed133dbSHong Zhang /* add multiple of row i to k-th row */ 11242ed133dbSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 11252ed133dbSHong Zhang if (jmin < jmax){ 11262ed133dbSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 11272ed133dbSHong Zhang /* update il and jl for row i */ 11282ed133dbSHong Zhang il[i] = jmin; 11292ed133dbSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 11302ed133dbSHong Zhang } 11312ed133dbSHong Zhang i = nexti; 11322ed133dbSHong Zhang } 11332ed133dbSHong Zhang 11340697b6caSHong Zhang /* shift the diagonals when zero pivot is detected */ 11350697b6caSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 11360697b6caSHong Zhang rs = 0.0; 11373c9e8598SHong Zhang jmin = bi[k]+1; 11383c9e8598SHong Zhang nz = bi[k+1] - jmin; 11393c9e8598SHong Zhang if (nz){ 11403c9e8598SHong Zhang bcol = bj + jmin; 11413c9e8598SHong Zhang while (nz--){ 114289f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 114389f845c8SHong Zhang bcol++; 11443c9e8598SHong Zhang } 11453c9e8598SHong Zhang } 1146540703acSHong Zhang 1147540703acSHong Zhang sctx.rs = rs; 1148540703acSHong Zhang sctx.pv = dk; 11493e8c821fSHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 1150540703acSHong Zhang if (newshift == 1){ 1151540703acSHong Zhang break; /* sctx.shift_amount is updated */ 1152540703acSHong Zhang } else if (newshift == -1){ 11530697b6caSHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 11543c9e8598SHong Zhang } 11553c9e8598SHong Zhang 11563c9e8598SHong Zhang /* copy data into U(k,:) */ 115739e3d630SHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 115839e3d630SHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 115939e3d630SHong Zhang if (jmin < jmax) { 116039e3d630SHong Zhang for (j=jmin; j<jmax; j++){ 116139e3d630SHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 11623c9e8598SHong Zhang } 116339e3d630SHong Zhang /* add the k-th row into il and jl */ 11643c9e8598SHong Zhang il[k] = jmin; 11653c9e8598SHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 11663c9e8598SHong Zhang } 11673c9e8598SHong Zhang } 1168540703acSHong Zhang } while (sctx.chshift); 11693c9e8598SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 11703c9e8598SHong Zhang 117139e3d630SHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 11723c9e8598SHong Zhang C->factor = FACTOR_CHOLESKY; 11733c9e8598SHong Zhang C->assembled = PETSC_TRUE; 11743c9e8598SHong Zhang C->preallocated = PETSC_TRUE; 1175efee365bSSatish Balay ierr = PetscLogFlops(C->m);CHKERRQ(ierr); 1176540703acSHong Zhang if (sctx.nshift){ 1177540703acSHong Zhang if (shiftnz) { 1178540703acSHong Zhang PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 1179540703acSHong Zhang } else if (shiftpd) { 1180540703acSHong Zhang PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 11810697b6caSHong Zhang } 11823c9e8598SHong Zhang } 11833c9e8598SHong Zhang PetscFunctionReturn(0); 11843c9e8598SHong Zhang } 1185a6175056SHong Zhang 1186a6175056SHong Zhang #undef __FUNCT__ 1187a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 1188dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1189a6175056SHong Zhang { 11900968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1191ed59401aSHong Zhang Mat_SeqSBAIJ *b; 1192ed59401aSHong Zhang Mat B; 1193dfbe8321SBarry Smith PetscErrorCode ierr; 1194653a6975SHong Zhang PetscTruth perm_identity; 1195622e2028SHong Zhang PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui; 1196ed59401aSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 1197391eac55SHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 11985a8e39fbSHong Zhang PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 1199ed59401aSHong Zhang PetscReal fill=info->fill,levels=info->levels; 12007a48dd6fSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 12017a48dd6fSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 12027a48dd6fSHong Zhang PetscBT lnkbt; 1203a6175056SHong Zhang 1204a6175056SHong Zhang PetscFunctionBegin; 1205653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 12067a48dd6fSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 12077a48dd6fSHong Zhang 120839e3d630SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 120939e3d630SHong Zhang ui[0] = 0; 121039e3d630SHong Zhang 1211622e2028SHong Zhang /* special case that simply copies fill pattern */ 1212622e2028SHong Zhang if (!levels && perm_identity) { 1213622e2028SHong Zhang ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 1214ed59401aSHong Zhang for (i=0; i<am; i++) { 121539e3d630SHong Zhang ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 1216ed59401aSHong Zhang } 121739e3d630SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 121839e3d630SHong Zhang cols = uj; 1219ed59401aSHong Zhang for (i=0; i<am; i++) { 1220ed59401aSHong Zhang aj = a->j + a->diag[i]; 122139e3d630SHong Zhang ncols = ui[i+1] - ui[i]; 122239e3d630SHong Zhang for (j=0; j<ncols; j++) *cols++ = *aj++; 1223ed59401aSHong Zhang } 122439e3d630SHong Zhang } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 12257a48dd6fSHong Zhang /* initialization */ 12265a8e39fbSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 12277a48dd6fSHong Zhang 12287a48dd6fSHong Zhang /* jl: linked list for storing indices of the pivot rows 12297a48dd6fSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 12301393bc97SHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 12317a48dd6fSHong Zhang il = jl + am; 12327a48dd6fSHong Zhang uj_ptr = (PetscInt**)(il + am); 12337a48dd6fSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 12347a48dd6fSHong Zhang for (i=0; i<am; i++){ 12357a48dd6fSHong Zhang jl[i] = am; il[i] = 0; 12367a48dd6fSHong Zhang } 12377a48dd6fSHong Zhang 12387a48dd6fSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 12397a48dd6fSHong Zhang nlnk = am + 1; 12402ed133dbSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12417a48dd6fSHong Zhang 12427a48dd6fSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 12437a48dd6fSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 12447a48dd6fSHong Zhang current_space = free_space; 12457a48dd6fSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 12467a48dd6fSHong Zhang current_space_lvl = free_space_lvl; 12477a48dd6fSHong Zhang 12487a48dd6fSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 12497a48dd6fSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 12507a48dd6fSHong Zhang nzk = 0; 1251622e2028SHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 1252622e2028SHong Zhang ncols_upper = 0; 1253622e2028SHong Zhang for (j=0; j<ncols; j++){ 12545a8e39fbSHong Zhang i = *(aj + ai[rip[k]] + j); 12555a8e39fbSHong Zhang if (rip[i] >= k){ /* only take upper triangular entry */ 12565a8e39fbSHong Zhang ajtmp[ncols_upper] = i; 1257622e2028SHong Zhang ncols_upper++; 1258622e2028SHong Zhang } 1259622e2028SHong Zhang } 12605a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12617a48dd6fSHong Zhang nzk += nlnk; 12627a48dd6fSHong Zhang 12637a48dd6fSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 12647a48dd6fSHong Zhang prow = jl[k]; /* 1st pivot row */ 12657a48dd6fSHong Zhang 12667a48dd6fSHong Zhang while (prow < k){ 12677a48dd6fSHong Zhang nextprow = jl[prow]; 12687a48dd6fSHong Zhang 12697a48dd6fSHong Zhang /* merge prow into k-th row */ 12707a48dd6fSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 12717a48dd6fSHong Zhang jmax = ui[prow+1]; 12727a48dd6fSHong Zhang ncols = jmax-jmin; 1273ed59401aSHong Zhang i = jmin - ui[prow]; 1274ed59401aSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 12755a8e39fbSHong Zhang uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 12765a8e39fbSHong Zhang j = *(uj - 1); 12775a8e39fbSHong Zhang ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 12787a48dd6fSHong Zhang nzk += nlnk; 12797a48dd6fSHong Zhang 12807a48dd6fSHong Zhang /* update il and jl for prow */ 12817a48dd6fSHong Zhang if (jmin < jmax){ 12827a48dd6fSHong Zhang il[prow] = jmin; 12837a48dd6fSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 12847a48dd6fSHong Zhang } 12857a48dd6fSHong Zhang prow = nextprow; 12867a48dd6fSHong Zhang } 12877a48dd6fSHong Zhang 12887a48dd6fSHong Zhang /* if free space is not available, make more free space */ 12897a48dd6fSHong Zhang if (current_space->local_remaining<nzk) { 12907a48dd6fSHong Zhang i = am - k + 1; /* num of unfactored rows */ 12917a48dd6fSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 12927a48dd6fSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 12937a48dd6fSHong Zhang ierr = GetMoreSpace(i,¤t_space_lvl);CHKERRQ(ierr); 12947a48dd6fSHong Zhang reallocs++; 12957a48dd6fSHong Zhang } 12967a48dd6fSHong Zhang 12977a48dd6fSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 12982ed133dbSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 12997a48dd6fSHong Zhang 13007a48dd6fSHong Zhang /* add the k-th row into il and jl */ 130139e3d630SHong Zhang if (nzk > 1){ 13027a48dd6fSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 13037a48dd6fSHong Zhang jl[k] = jl[i]; jl[i] = k; 13047a48dd6fSHong Zhang il[k] = ui[k] + 1; 13057a48dd6fSHong Zhang } 13067a48dd6fSHong Zhang uj_ptr[k] = current_space->array; 13077a48dd6fSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 13087a48dd6fSHong Zhang 13097a48dd6fSHong Zhang current_space->array += nzk; 13107a48dd6fSHong Zhang current_space->local_used += nzk; 13117a48dd6fSHong Zhang current_space->local_remaining -= nzk; 13127a48dd6fSHong Zhang 13137a48dd6fSHong Zhang current_space_lvl->array += nzk; 13147a48dd6fSHong Zhang current_space_lvl->local_used += nzk; 13157a48dd6fSHong Zhang current_space_lvl->local_remaining -= nzk; 13167a48dd6fSHong Zhang 13177a48dd6fSHong Zhang ui[k+1] = ui[k] + nzk; 13187a48dd6fSHong Zhang } 13197a48dd6fSHong Zhang 13207a48dd6fSHong Zhang if (ai[am] != 0) { 132139e3d630SHong Zhang PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 13227a48dd6fSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af); 13237a48dd6fSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af); 13247a48dd6fSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af); 13257a48dd6fSHong Zhang } else { 1326ed59401aSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n"); 13277a48dd6fSHong Zhang } 13287a48dd6fSHong Zhang 13297a48dd6fSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 13307a48dd6fSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 13315a8e39fbSHong Zhang ierr = PetscFree(ajtmp);CHKERRQ(ierr); 13327a48dd6fSHong Zhang 13337a48dd6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 13347a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 13357a48dd6fSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 13362ed133dbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 13377a48dd6fSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 13387a48dd6fSHong Zhang 133939e3d630SHong Zhang } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 134039e3d630SHong Zhang 13417a48dd6fSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 13427a48dd6fSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 1343ed59401aSHong Zhang B = *fact; 1344ed59401aSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 1345ed59401aSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 13467a48dd6fSHong Zhang 1347ed59401aSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 13487a48dd6fSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 13497a48dd6fSHong Zhang b->singlemalloc = PETSC_FALSE; 13507a48dd6fSHong Zhang /* the next line frees the default space generated by the Create() */ 13517a48dd6fSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 13527a48dd6fSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 13537a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 13547a48dd6fSHong Zhang b->j = uj; 13557a48dd6fSHong Zhang b->i = ui; 13567a48dd6fSHong Zhang b->diag = 0; 13577a48dd6fSHong Zhang b->ilen = 0; 13587a48dd6fSHong Zhang b->imax = 0; 13597a48dd6fSHong Zhang b->row = perm; 13607a48dd6fSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 13617a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 13627a48dd6fSHong Zhang b->icol = perm; 13637a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 13647a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 136552e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 13667a48dd6fSHong Zhang b->maxnz = b->nz = ui[am]; 13677a48dd6fSHong Zhang 1368ed59401aSHong Zhang B->factor = FACTOR_CHOLESKY; 1369ed59401aSHong Zhang B->info.factor_mallocs = reallocs; 1370ed59401aSHong Zhang B->info.fill_ratio_given = fill; 13717a48dd6fSHong Zhang if (ai[am] != 0) { 1372ed59401aSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 13737a48dd6fSHong Zhang } else { 1374ed59401aSHong Zhang B->info.fill_ratio_needed = 0.0; 13757a48dd6fSHong Zhang } 137639e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1377622e2028SHong Zhang if (perm_identity){ 1378ed59401aSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1379ed59401aSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1380622e2028SHong Zhang } 1381a6175056SHong Zhang PetscFunctionReturn(0); 1382a6175056SHong Zhang } 1383d5d45c9bSBarry Smith 1384f76d2b81SHong Zhang #undef __FUNCT__ 1385f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1386dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1387f76d2b81SHong Zhang { 1388f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 138910c27e3dSHong Zhang Mat_SeqSBAIJ *b; 13902ed133dbSHong Zhang Mat B; 1391dfbe8321SBarry Smith PetscErrorCode ierr; 1392f76d2b81SHong Zhang PetscTruth perm_identity; 139310c27e3dSHong Zhang PetscReal fill = info->fill; 13941393bc97SHong Zhang PetscInt *rip,*riip,i,am=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow; 139510c27e3dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 13962ed133dbSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 139710c27e3dSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 139810c27e3dSHong Zhang PetscBT lnkbt; 13992ed133dbSHong Zhang IS iperm; 1400f76d2b81SHong Zhang 1401f76d2b81SHong Zhang PetscFunctionBegin; 14022ed133dbSHong Zhang /* check whether perm is the identity mapping */ 1403f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 14042ed133dbSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 14052ed133dbSHong Zhang 1406f76d2b81SHong Zhang if (!perm_identity){ 14072ed133dbSHong Zhang /* check if perm is symmetric! */ 14082ed133dbSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 14092ed133dbSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 14101393bc97SHong Zhang for (i=0; i<am; i++) { 14112ed133dbSHong Zhang if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation"); 14122ed133dbSHong Zhang } 14132ed133dbSHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 14142ed133dbSHong Zhang ierr = ISDestroy(iperm);CHKERRQ(ierr); 1415f76d2b81SHong Zhang } 141610c27e3dSHong Zhang 141710c27e3dSHong Zhang /* initialization */ 14181393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 141910c27e3dSHong Zhang ui[0] = 0; 142010c27e3dSHong Zhang 142110c27e3dSHong Zhang /* jl: linked list for storing indices of the pivot rows 14221393bc97SHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 14231393bc97SHong Zhang ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 14241393bc97SHong Zhang il = jl + am; 14251393bc97SHong Zhang cols = il + am; 14261393bc97SHong Zhang ui_ptr = (PetscInt**)(cols + am); 14271393bc97SHong Zhang for (i=0; i<am; i++){ 14281393bc97SHong Zhang jl[i] = am; il[i] = 0; 1429f76d2b81SHong Zhang } 143010c27e3dSHong Zhang 143110c27e3dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 14321393bc97SHong Zhang nlnk = am + 1; 14331393bc97SHong Zhang ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 143410c27e3dSHong Zhang 14351393bc97SHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 14361393bc97SHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 143710c27e3dSHong Zhang current_space = free_space; 143810c27e3dSHong Zhang 14391393bc97SHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 144010c27e3dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 144110c27e3dSHong Zhang nzk = 0; 14422ed133dbSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 14432ed133dbSHong Zhang ncols_upper = 0; 14442ed133dbSHong Zhang for (j=0; j<ncols; j++){ 1445622e2028SHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 14462ed133dbSHong Zhang if (i >= k){ /* only take upper triangular entry */ 14472ed133dbSHong Zhang cols[ncols_upper] = i; 14482ed133dbSHong Zhang ncols_upper++; 14492ed133dbSHong Zhang } 14502ed133dbSHong Zhang } 14511393bc97SHong Zhang ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 145210c27e3dSHong Zhang nzk += nlnk; 145310c27e3dSHong Zhang 145410c27e3dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 145510c27e3dSHong Zhang prow = jl[k]; /* 1st pivot row */ 145610c27e3dSHong Zhang 145710c27e3dSHong Zhang while (prow < k){ 145810c27e3dSHong Zhang nextprow = jl[prow]; 145910c27e3dSHong Zhang /* merge prow into k-th row */ 14601393bc97SHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 146110c27e3dSHong Zhang jmax = ui[prow+1]; 146210c27e3dSHong Zhang ncols = jmax-jmin; 14631393bc97SHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 14645a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 146510c27e3dSHong Zhang nzk += nlnk; 146610c27e3dSHong Zhang 146710c27e3dSHong Zhang /* update il and jl for prow */ 146810c27e3dSHong Zhang if (jmin < jmax){ 146910c27e3dSHong Zhang il[prow] = jmin; 14702ed133dbSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 147110c27e3dSHong Zhang } 147210c27e3dSHong Zhang prow = nextprow; 147310c27e3dSHong Zhang } 147410c27e3dSHong Zhang 147510c27e3dSHong Zhang /* if free space is not available, make more free space */ 147610c27e3dSHong Zhang if (current_space->local_remaining<nzk) { 14771393bc97SHong Zhang i = am - k + 1; /* num of unfactored rows */ 147810c27e3dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 147910c27e3dSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 148010c27e3dSHong Zhang reallocs++; 148110c27e3dSHong Zhang } 148210c27e3dSHong Zhang 148310c27e3dSHong Zhang /* copy data into free space, then initialize lnk */ 14841393bc97SHong Zhang ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 148510c27e3dSHong Zhang 148610c27e3dSHong Zhang /* add the k-th row into il and jl */ 148710c27e3dSHong Zhang if (nzk-1 > 0){ 14881393bc97SHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 148910c27e3dSHong Zhang jl[k] = jl[i]; jl[i] = k; 149010c27e3dSHong Zhang il[k] = ui[k] + 1; 149110c27e3dSHong Zhang } 14922ed133dbSHong Zhang ui_ptr[k] = current_space->array; 149310c27e3dSHong Zhang current_space->array += nzk; 149410c27e3dSHong Zhang current_space->local_used += nzk; 149510c27e3dSHong Zhang current_space->local_remaining -= nzk; 149610c27e3dSHong Zhang 149710c27e3dSHong Zhang ui[k+1] = ui[k] + nzk; 149810c27e3dSHong Zhang } 149910c27e3dSHong Zhang 15001393bc97SHong Zhang if (ai[am] != 0) { 15011393bc97SHong Zhang PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]); 150278910aadSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af); 150378910aadSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af); 150478910aadSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af); 150510c27e3dSHong Zhang } else { 150678910aadSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n"); 150710c27e3dSHong Zhang } 150810c27e3dSHong Zhang 150910c27e3dSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 151010c27e3dSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 151110c27e3dSHong Zhang 151210c27e3dSHong Zhang /* destroy list of free space and other temporary array(s) */ 15131393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 151410c27e3dSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 151510c27e3dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 151610c27e3dSHong Zhang 151710c27e3dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 15181393bc97SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 15192ed133dbSHong Zhang B = *fact; 15202ed133dbSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 15212ed133dbSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 152210c27e3dSHong Zhang 15232ed133dbSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 152410c27e3dSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 152510c27e3dSHong Zhang b->singlemalloc = PETSC_FALSE; 152610c27e3dSHong Zhang /* the next line frees the default space generated by the Create() */ 152710c27e3dSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 152810c27e3dSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 15291393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 153010c27e3dSHong Zhang b->j = uj; 153110c27e3dSHong Zhang b->i = ui; 153210c27e3dSHong Zhang b->diag = 0; 153310c27e3dSHong Zhang b->ilen = 0; 153410c27e3dSHong Zhang b->imax = 0; 153510c27e3dSHong Zhang b->row = perm; 153610c27e3dSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 153710c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 153810c27e3dSHong Zhang b->icol = perm; 153910c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 15401393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 15411393bc97SHong Zhang ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 15421393bc97SHong Zhang b->maxnz = b->nz = ui[am]; 154310c27e3dSHong Zhang 15442ed133dbSHong Zhang B->factor = FACTOR_CHOLESKY; 15452ed133dbSHong Zhang B->info.factor_mallocs = reallocs; 15462ed133dbSHong Zhang B->info.fill_ratio_given = fill; 15471393bc97SHong Zhang if (ai[am] != 0) { 15481393bc97SHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 154910c27e3dSHong Zhang } else { 15502ed133dbSHong Zhang B->info.fill_ratio_needed = 0.0; 155110c27e3dSHong Zhang } 155239e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 15532ed133dbSHong Zhang if (perm_identity){ 155410c27e3dSHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 155510c27e3dSHong Zhang (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 15562ed133dbSHong Zhang } 1557f76d2b81SHong Zhang PetscFunctionReturn(0); 1558f76d2b81SHong Zhang } 1559