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 */ 507529efd4SKris 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; 25063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af));CHKERRQ(ierr); 25163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af));CHKERRQ(ierr); 25263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af));CHKERRQ(ierr); 25363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:for best performance.\n"));CHKERRQ(ierr); 254431e710aSBarry Smith 2557529efd4SKris 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 } 35463ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 355464e8d44SSatish Balay if (ai[n] != 0) { 3561393bc97SHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 35763ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr); 35863ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af));CHKERRQ(ierr); 35963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af));CHKERRQ(ierr); 36063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr); 36105bf559cSSatish Balay } else { 36263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n"));CHKERRQ(ierr); 36305bf559cSSatish Balay } 36463ba0a88SBarry Smith #endif 3652fb3ed76SBarry Smith 366898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 367898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3681987afe7SBarry Smith 3691393bc97SHong Zhang /* destroy list of free space and other temporary array(s) */ 3701393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 3711393bc97SHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 3721393bc97SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 3731393bc97SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 374289bc588SBarry Smith 375289bc588SBarry Smith /* put together the new matrix */ 376f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr); 377f204ca49SKris Buschelman ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr); 378f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr); 37952e6d16bSBarry Smith ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr); 380416022c9SBarry Smith b = (Mat_SeqAIJ*)(*B)->data; 381606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 3827c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 383e8d4e0b9SBarry Smith /* the next line frees the default space generated by the Create() */ 384606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 385606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 3861393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 3871393bc97SHong Zhang b->j = bj; 3881393bc97SHong Zhang b->i = bi; 3891393bc97SHong Zhang b->diag = bdiag; 390416022c9SBarry Smith b->ilen = 0; 391416022c9SBarry Smith b->imax = 0; 392416022c9SBarry Smith b->row = isrow; 393416022c9SBarry Smith b->col = iscol; 394c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 395c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 39682bf6240SBarry Smith b->icol = isicol; 39787828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 3981393bc97SHong Zhang 3991393bc97SHong Zhang /* In b structure: Free imax, ilen, old a, old j. Allocate solve_work, new a, new j */ 4001393bc97SHong Zhang ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 4011393bc97SHong Zhang b->maxnz = b->nz = bi[n] ; 4027fd98bd6SLois Curfman McInnes 403329f5518SBarry Smith (*B)->factor = FACTOR_LU; 404418422e8SSatish Balay (*B)->info.factor_mallocs = reallocs; 405ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 406703deb49SSatish Balay 407e93ccc4dSBarry Smith if (ai[n] != 0) { 4081393bc97SHong Zhang (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 409eea2913aSSatish Balay } else { 410eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 411eea2913aSSatish Balay } 4124846f1f5SKris Buschelman ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr); 41371c2f376SKris Buschelman (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 4143a40ed3dSBarry Smith PetscFunctionReturn(0); 415289bc588SBarry Smith } 4161393bc97SHong Zhang 4170a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 4184a2ae208SSatish Balay #undef __FUNCT__ 4194a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 420af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 421289bc588SBarry Smith { 422416022c9SBarry Smith Mat C=*B; 423416022c9SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 42482bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 4256849ba73SBarry Smith PetscErrorCode ierr; 426b3bf805bSHong Zhang PetscInt *r,*ic,i,j,n=A->m,*bi=b->i,*bj=b->j; 427b3bf805bSHong Zhang PetscInt *ajtmp,*bjtmp,nz,row,*ics; 428d2276718SHong Zhang PetscInt *diag_offset = b->diag,diag,*pj; 42987828ca2SBarry Smith PetscScalar *rtmp,*v,*pc,multiplier,*pv,*rtmps; 430*6a7c8fc2SHong Zhang PetscScalar d; 431*6a7c8fc2SHong Zhang PetscReal rs; 432b3bf805bSHong Zhang LUShift_Ctx sctx; 433d2276718SHong Zhang PetscInt newshift; 434289bc588SBarry Smith 4353a40ed3dSBarry Smith PetscFunctionBegin; 43678b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 43778b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 43887828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 43987828ca2SBarry Smith ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 440010a20caSHong Zhang rtmps = rtmp; ics = ic; 441289bc588SBarry Smith 4426cc28720Svictorle if (!a->diag) { 4430cf777b8SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 4440cf777b8SBarry Smith } 445118f5789SBarry Smith /* if both shift schemes are chosen by user, only use info->shiftpd */ 446118f5789SBarry Smith if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0; 4471a890ab1SHong Zhang if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 44838baddfdSBarry Smith PetscInt *aai = a->i,*ddiag = a->diag; 449118f5789SBarry Smith sctx.shift_top = 0; 4506cc28720Svictorle for (i=0; i<n; i++) { 4519f95998fSHong Zhang /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 4529f95998fSHong Zhang d = (a->a)[ddiag[i]]; 4531a890ab1SHong Zhang rs = -PetscAbsScalar(d) - PetscRealPart(d); 454010a20caSHong Zhang v = a->a+aai[i]; 455e24cfe64SHong Zhang nz = aai[i+1] - aai[i]; 456e24cfe64SHong Zhang for (j=0; j<nz; j++) 4571a890ab1SHong Zhang rs += PetscAbsScalar(v[j]); 4581a890ab1SHong Zhang if (rs>sctx.shift_top) sctx.shift_top = rs; 4596cc28720Svictorle } 460118f5789SBarry Smith if (sctx.shift_top == 0.0) sctx.shift_top += 1.e-12; 461118f5789SBarry Smith sctx.shift_top *= 1.1; 462d2276718SHong Zhang sctx.nshift_max = 5; 463d2276718SHong Zhang sctx.shift_lo = 0.; 464d2276718SHong Zhang sctx.shift_hi = 1.; 465a0a356efSHong Zhang } 466a0a356efSHong Zhang 467a0a356efSHong Zhang sctx.shift_amount = 0; 468a0a356efSHong Zhang sctx.nshift = 0; 469085256b3SBarry Smith do { 470d2276718SHong Zhang sctx.lushift = PETSC_FALSE; 471289bc588SBarry Smith for (i=0; i<n; i++){ 472b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 473b3bf805bSHong Zhang bjtmp = bj + bi[i]; 474b3bf805bSHong Zhang for (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0; 475289bc588SBarry Smith 476289bc588SBarry Smith /* load in initial (unfactored row) */ 477416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 478b3bf805bSHong Zhang ajtmp = a->j + a->i[r[i]]; 479010a20caSHong Zhang v = a->a + a->i[r[i]]; 480085256b3SBarry Smith for (j=0; j<nz; j++) { 481b3bf805bSHong Zhang rtmp[ics[ajtmp[j]]] = v[j]; 482085256b3SBarry Smith } 483d2276718SHong Zhang rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 484289bc588SBarry Smith 485b3bf805bSHong Zhang row = *bjtmp++; 486f2582111SSatish Balay while (row < i) { 4878c37ef55SBarry Smith pc = rtmp + row; 4888c37ef55SBarry Smith if (*pc != 0.0) { 489010a20caSHong Zhang pv = b->a + diag_offset[row]; 490010a20caSHong Zhang pj = b->j + diag_offset[row] + 1; 49135aab85fSBarry Smith multiplier = *pc / *pv++; 4928c37ef55SBarry Smith *pc = multiplier; 493b3bf805bSHong Zhang nz = bi[row+1] - diag_offset[row] - 1; 494f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 495efee365bSSatish Balay ierr = PetscLogFlops(2*nz);CHKERRQ(ierr); 496289bc588SBarry Smith } 497b3bf805bSHong Zhang row = *bjtmp++; 498289bc588SBarry Smith } 499416022c9SBarry Smith /* finished row so stick it into b->a */ 500b3bf805bSHong Zhang pv = b->a + bi[i] ; 501b3bf805bSHong Zhang pj = b->j + bi[i] ; 502b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 503b3bf805bSHong Zhang diag = diag_offset[i] - bi[i]; 504d3d32019SBarry Smith rs = 0.0; 505d3d32019SBarry Smith for (j=0; j<nz; j++) { 506d3d32019SBarry Smith pv[j] = rtmps[pj[j]]; 507d3d32019SBarry Smith if (j != diag) rs += PetscAbsScalar(pv[j]); 508d3d32019SBarry Smith } 509d2276718SHong Zhang 510b3bf805bSHong Zhang /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 511d2276718SHong Zhang sctx.rs = rs; 512d2276718SHong Zhang sctx.pv = pv[diag]; 5133e8c821fSHong Zhang ierr = MatLUCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 514d2276718SHong Zhang if (newshift == 1){ 515b3bf805bSHong Zhang break; /* sctx.shift_amount is updated */ 516d2276718SHong Zhang } else if (newshift == -1){ 517d2276718SHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs); 518d708749eSBarry Smith } 51935aab85fSBarry Smith } 520d2276718SHong Zhang 521118f5789SBarry Smith if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 5226cc28720Svictorle /* 5236c037d1bSvictorle * if no shift in this attempt & shifting & started shifting & can refine, 5246cc28720Svictorle * then try lower shift 5256cc28720Svictorle */ 526d2276718SHong Zhang sctx.shift_hi = info->shift_fraction; 527d2276718SHong Zhang info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 528d2276718SHong Zhang sctx.shift_amount = info->shift_fraction * sctx.shift_top; 529d2276718SHong Zhang sctx.lushift = PETSC_TRUE; 530d2276718SHong Zhang sctx.nshift++; 5316cc28720Svictorle } 532d2276718SHong Zhang } while (sctx.lushift); 5330f11f9aeSSatish Balay 53435aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 53535aab85fSBarry Smith for (i=0; i<n; i++) { 536010a20caSHong Zhang b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]]; 53735aab85fSBarry Smith } 53835aab85fSBarry Smith 539606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 54078b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 54178b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 542416022c9SBarry Smith C->factor = FACTOR_LU; 5437dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 544c456f294SBarry Smith C->assembled = PETSC_TRUE; 545efee365bSSatish Balay ierr = PetscLogFlops(C->n);CHKERRQ(ierr); 546d2276718SHong Zhang if (sctx.nshift){ 547118f5789SBarry Smith if (info->shiftnz) { 54863ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatLUFactorNumeric_SeqAIJ: number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 549118f5789SBarry Smith } else if (info->shiftpd) { 55063ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatLUFactorNumeric_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,sctx.shift_top));CHKERRQ(ierr); 5516cc28720Svictorle } 5520697b6caSHong Zhang } 5533a40ed3dSBarry Smith PetscFunctionReturn(0); 554289bc588SBarry Smith } 5550889b5dcSvictorle 5560889b5dcSvictorle #undef __FUNCT__ 5570889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ" 558dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A) 5590889b5dcSvictorle { 5600889b5dcSvictorle PetscFunctionBegin; 5610889b5dcSvictorle A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ; 5620889b5dcSvictorle A->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 5630889b5dcSvictorle PetscFunctionReturn(0); 5640889b5dcSvictorle } 5650889b5dcSvictorle 5660889b5dcSvictorle 5670a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5684a2ae208SSatish Balay #undef __FUNCT__ 5694a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 570dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 571da3a660dSBarry Smith { 572dfbe8321SBarry Smith PetscErrorCode ierr; 573416022c9SBarry Smith Mat C; 574416022c9SBarry Smith 5753a40ed3dSBarry Smith PetscFunctionBegin; 5769e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 577af281ebdSHong Zhang ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 578273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 57952e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr); 5803a40ed3dSBarry Smith PetscFunctionReturn(0); 581da3a660dSBarry Smith } 5820a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5834a2ae208SSatish Balay #undef __FUNCT__ 5844a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 585dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 5868c37ef55SBarry Smith { 587416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 588416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 5896849ba73SBarry Smith PetscErrorCode ierr; 59038baddfdSBarry Smith PetscInt *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 59138baddfdSBarry Smith PetscInt nz,*rout,*cout; 59287828ca2SBarry Smith PetscScalar *x,*b,*tmp,*tmps,*aa = a->a,sum,*v; 5938c37ef55SBarry Smith 5943a40ed3dSBarry Smith PetscFunctionBegin; 5953a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 59691bf9adeSBarry Smith 5971ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 5981ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 599416022c9SBarry Smith tmp = a->solve_work; 6008c37ef55SBarry Smith 6013b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6023b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 6038c37ef55SBarry Smith 6048c37ef55SBarry Smith /* forward solve the lower triangular */ 6058c37ef55SBarry Smith tmp[0] = b[*r++]; 606010a20caSHong Zhang tmps = tmp; 6078c37ef55SBarry Smith for (i=1; i<n; i++) { 608010a20caSHong Zhang v = aa + ai[i] ; 609010a20caSHong Zhang vi = aj + ai[i] ; 61053e7425aSBarry Smith nz = a->diag[i] - ai[i]; 61153e7425aSBarry Smith sum = b[*r++]; 612ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 6138c37ef55SBarry Smith tmp[i] = sum; 6148c37ef55SBarry Smith } 6158c37ef55SBarry Smith 6168c37ef55SBarry Smith /* backward solve the upper triangular */ 6178c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 618010a20caSHong Zhang v = aa + a->diag[i] + 1; 619010a20caSHong Zhang vi = aj + a->diag[i] + 1; 620416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 6218c37ef55SBarry Smith sum = tmp[i]; 622ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 623010a20caSHong Zhang x[*c--] = tmp[i] = sum*aa[a->diag[i]]; 6248c37ef55SBarry Smith } 6258c37ef55SBarry Smith 6263b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6273b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 6281ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 6291ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 630efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr); 6313a40ed3dSBarry Smith PetscFunctionReturn(0); 6328c37ef55SBarry Smith } 633026e39d0SSatish Balay 634930ae53cSBarry Smith /* ----------------------------------------------------------- */ 6354a2ae208SSatish Balay #undef __FUNCT__ 6364a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 637dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 638930ae53cSBarry Smith { 639930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 6406849ba73SBarry Smith PetscErrorCode ierr; 64138baddfdSBarry Smith PetscInt n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag; 642362ced78SSatish Balay PetscScalar *x,*b,*aa = a->a; 643aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 64438baddfdSBarry Smith PetscInt adiag_i,i,*vi,nz,ai_i; 645362ced78SSatish Balay PetscScalar *v,sum; 646d85afc42SSatish Balay #endif 647930ae53cSBarry Smith 6483a40ed3dSBarry Smith PetscFunctionBegin; 6493a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 650930ae53cSBarry Smith 6511ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 6521ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 653930ae53cSBarry Smith 654aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 6551c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 6561c4feecaSSatish Balay #else 657930ae53cSBarry Smith /* forward solve the lower triangular */ 658930ae53cSBarry Smith x[0] = b[0]; 659930ae53cSBarry Smith for (i=1; i<n; i++) { 660930ae53cSBarry Smith ai_i = ai[i]; 661930ae53cSBarry Smith v = aa + ai_i; 662930ae53cSBarry Smith vi = aj + ai_i; 663930ae53cSBarry Smith nz = adiag[i] - ai_i; 664930ae53cSBarry Smith sum = b[i]; 665930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 666930ae53cSBarry Smith x[i] = sum; 667930ae53cSBarry Smith } 668930ae53cSBarry Smith 669930ae53cSBarry Smith /* backward solve the upper triangular */ 670930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 671930ae53cSBarry Smith adiag_i = adiag[i]; 672d708749eSBarry Smith v = aa + adiag_i + 1; 673d708749eSBarry Smith vi = aj + adiag_i + 1; 674930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 675930ae53cSBarry Smith sum = x[i]; 676930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 677930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 678930ae53cSBarry Smith } 6791c4feecaSSatish Balay #endif 680efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr); 6811ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 6821ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6833a40ed3dSBarry Smith PetscFunctionReturn(0); 684930ae53cSBarry Smith } 685930ae53cSBarry Smith 6864a2ae208SSatish Balay #undef __FUNCT__ 6874a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 688dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 689da3a660dSBarry Smith { 690416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 691416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 6926849ba73SBarry Smith PetscErrorCode ierr; 69338baddfdSBarry Smith PetscInt *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 69438baddfdSBarry Smith PetscInt nz,*rout,*cout; 69587828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,sum,*v; 696da3a660dSBarry Smith 6973a40ed3dSBarry Smith PetscFunctionBegin; 69878b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 699da3a660dSBarry Smith 7001ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 7011ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 702416022c9SBarry Smith tmp = a->solve_work; 703da3a660dSBarry Smith 7043b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7053b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 706da3a660dSBarry Smith 707da3a660dSBarry Smith /* forward solve the lower triangular */ 708da3a660dSBarry Smith tmp[0] = b[*r++]; 709da3a660dSBarry Smith for (i=1; i<n; i++) { 710010a20caSHong Zhang v = aa + ai[i] ; 711010a20caSHong Zhang vi = aj + ai[i] ; 712416022c9SBarry Smith nz = a->diag[i] - ai[i]; 713da3a660dSBarry Smith sum = b[*r++]; 714010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 715da3a660dSBarry Smith tmp[i] = sum; 716da3a660dSBarry Smith } 717da3a660dSBarry Smith 718da3a660dSBarry Smith /* backward solve the upper triangular */ 719da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 720010a20caSHong Zhang v = aa + a->diag[i] + 1; 721010a20caSHong Zhang vi = aj + a->diag[i] + 1; 722416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 723da3a660dSBarry Smith sum = tmp[i]; 724010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 725010a20caSHong Zhang tmp[i] = sum*aa[a->diag[i]]; 726da3a660dSBarry Smith x[*c--] += tmp[i]; 727da3a660dSBarry Smith } 728da3a660dSBarry Smith 7293b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7303b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7311ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 7321ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 733efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 734898183ecSLois Curfman McInnes 7353a40ed3dSBarry Smith PetscFunctionReturn(0); 736da3a660dSBarry Smith } 737da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 7384a2ae208SSatish Balay #undef __FUNCT__ 7394a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 740dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 741da3a660dSBarry Smith { 742416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7432235e667SBarry Smith IS iscol = a->col,isrow = a->row; 7446849ba73SBarry Smith PetscErrorCode ierr; 74538baddfdSBarry Smith PetscInt *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 74638baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 74787828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v,s1; 748da3a660dSBarry Smith 7493a40ed3dSBarry Smith PetscFunctionBegin; 7501ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 7511ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 752416022c9SBarry Smith tmp = a->solve_work; 753da3a660dSBarry Smith 7542235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7552235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 756da3a660dSBarry Smith 757da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 7582235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 759da3a660dSBarry Smith 760da3a660dSBarry Smith /* forward solve the U^T */ 761da3a660dSBarry Smith for (i=0; i<n; i++) { 762010a20caSHong Zhang v = aa + diag[i] ; 763010a20caSHong Zhang vi = aj + diag[i] + 1; 764f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 765c38d4ed2SBarry Smith s1 = tmp[i]; 766ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 767da3a660dSBarry Smith while (nz--) { 768010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*s1; 769da3a660dSBarry Smith } 770c38d4ed2SBarry Smith tmp[i] = s1; 771da3a660dSBarry Smith } 772da3a660dSBarry Smith 773da3a660dSBarry Smith /* backward solve the L^T */ 774da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 775010a20caSHong Zhang v = aa + diag[i] - 1 ; 776010a20caSHong Zhang vi = aj + diag[i] - 1 ; 777f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 778f1af5d2fSBarry Smith s1 = tmp[i]; 779da3a660dSBarry Smith while (nz--) { 780010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*s1; 781da3a660dSBarry Smith } 782da3a660dSBarry Smith } 783da3a660dSBarry Smith 784da3a660dSBarry Smith /* copy tmp into x according to permutation */ 785da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 786da3a660dSBarry Smith 7872235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7882235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7891ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 7901ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 791da3a660dSBarry Smith 792efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz-A->n);CHKERRQ(ierr); 7933a40ed3dSBarry Smith PetscFunctionReturn(0); 794da3a660dSBarry Smith } 795da3a660dSBarry Smith 7964a2ae208SSatish Balay #undef __FUNCT__ 7974a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 798dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 799da3a660dSBarry Smith { 800416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 8012235e667SBarry Smith IS iscol = a->col,isrow = a->row; 8026849ba73SBarry Smith PetscErrorCode ierr; 80338baddfdSBarry Smith PetscInt *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 80438baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 80587828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v; 8066abc6512SBarry Smith 8073a40ed3dSBarry Smith PetscFunctionBegin; 80823bc6035SMatthew Knepley if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);} 8096abc6512SBarry Smith 8101ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 8111ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 812416022c9SBarry Smith tmp = a->solve_work; 8136abc6512SBarry Smith 8142235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8152235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 8166abc6512SBarry Smith 8176abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 8182235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 8196abc6512SBarry Smith 8206abc6512SBarry Smith /* forward solve the U^T */ 8216abc6512SBarry Smith for (i=0; i<n; i++) { 822010a20caSHong Zhang v = aa + diag[i] ; 823010a20caSHong Zhang vi = aj + diag[i] + 1; 824f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 8256abc6512SBarry Smith tmp[i] *= *v++; 8266abc6512SBarry Smith while (nz--) { 827010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*tmp[i]; 8286abc6512SBarry Smith } 8296abc6512SBarry Smith } 8306abc6512SBarry Smith 8316abc6512SBarry Smith /* backward solve the L^T */ 8326abc6512SBarry Smith for (i=n-1; i>=0; i--){ 833010a20caSHong Zhang v = aa + diag[i] - 1 ; 834010a20caSHong Zhang vi = aj + diag[i] - 1 ; 835f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 8366abc6512SBarry Smith while (nz--) { 837010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*tmp[i]; 8386abc6512SBarry Smith } 8396abc6512SBarry Smith } 8406abc6512SBarry Smith 8416abc6512SBarry Smith /* copy tmp into x according to permutation */ 8426abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 8436abc6512SBarry Smith 8442235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8452235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 8461ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 8471ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8486abc6512SBarry Smith 849efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 8503a40ed3dSBarry Smith PetscFunctionReturn(0); 851da3a660dSBarry Smith } 8529e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 853dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat); 85408480c60SBarry Smith 8554a2ae208SSatish Balay #undef __FUNCT__ 8564a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 857dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 8589e25ed09SBarry Smith { 859416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 8609e25ed09SBarry Smith IS isicol; 8616849ba73SBarry Smith PetscErrorCode ierr; 8625a8e39fbSHong Zhang PetscInt *r,*ic,n=A->m,*ai=a->i,*aj=a->j; 8635a8e39fbSHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 8645a8e39fbSHong Zhang PetscInt *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0; 8655a8e39fbSHong Zhang PetscInt i,levels,diagonal_fill; 86677c4ece6SBarry Smith PetscTruth col_identity,row_identity; 867329f5518SBarry Smith PetscReal f; 8685a8e39fbSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 8695a8e39fbSHong Zhang PetscBT lnkbt; 8705a8e39fbSHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 8715a8e39fbSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 8725a8e39fbSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 8739e25ed09SBarry Smith 8743a40ed3dSBarry Smith PetscFunctionBegin; 8756cf997b0SBarry Smith f = info->fill; 87638baddfdSBarry Smith levels = (PetscInt)info->levels; 87738baddfdSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 8784c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 87982bf6240SBarry Smith 88008480c60SBarry Smith /* special case that simply copies fill pattern */ 881be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 882be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 88386bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 8842e8a6d31SBarry Smith ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 88508480c60SBarry Smith (*fact)->factor = FACTOR_LU; 88608480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 88708480c60SBarry Smith if (!b->diag) { 8887c922b88SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 88908480c60SBarry Smith } 8907c922b88SBarry Smith ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 89108480c60SBarry Smith b->row = isrow; 89208480c60SBarry Smith b->col = iscol; 89382bf6240SBarry Smith b->icol = isicol; 89487828ca2SBarry Smith ierr = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 895f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 896c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 897c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 8983a40ed3dSBarry Smith PetscFunctionReturn(0); 89908480c60SBarry Smith } 90008480c60SBarry Smith 901898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 902898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 9039e25ed09SBarry Smith 9049e25ed09SBarry Smith /* get new row pointers */ 9055a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 9065a8e39fbSHong Zhang bi[0] = 0; 9075a8e39fbSHong Zhang /* bdiag is location of diagonal in factor */ 9085a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 9095a8e39fbSHong Zhang bdiag[0] = 0; 9106cf997b0SBarry Smith 9115a8e39fbSHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr); 9125a8e39fbSHong Zhang bjlvl_ptr = (PetscInt**)(bj_ptr + n); 9135a8e39fbSHong Zhang 9145a8e39fbSHong Zhang /* create a linked list for storing column indices of the active row */ 9155a8e39fbSHong Zhang nlnk = n + 1; 9165a8e39fbSHong Zhang ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 9175a8e39fbSHong Zhang 9185a8e39fbSHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 9195a8e39fbSHong Zhang ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 9205a8e39fbSHong Zhang current_space = free_space; 9215a8e39fbSHong Zhang ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 9225a8e39fbSHong Zhang current_space_lvl = free_space_lvl; 9235a8e39fbSHong Zhang 9245a8e39fbSHong Zhang for (i=0; i<n; i++) { 9255a8e39fbSHong Zhang nzi = 0; 9266cf997b0SBarry Smith /* copy current row into linked list */ 9275a8e39fbSHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 9285a8e39fbSHong Zhang if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 9295a8e39fbSHong Zhang cols = aj + ai[r[i]]; 9305a8e39fbSHong Zhang lnk[i] = -1; /* marker to indicate if diagonal exists */ 9315a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 9325a8e39fbSHong Zhang nzi += nlnk; 9336cf997b0SBarry Smith 9346cf997b0SBarry Smith /* make sure diagonal entry is included */ 9355a8e39fbSHong Zhang if (diagonal_fill && lnk[i] == -1) { 9366cf997b0SBarry Smith fm = n; 9375a8e39fbSHong Zhang while (lnk[fm] < i) fm = lnk[fm]; 9385a8e39fbSHong Zhang lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 9395a8e39fbSHong Zhang lnk[fm] = i; 9405a8e39fbSHong Zhang lnk_lvl[i] = 0; 9415a8e39fbSHong Zhang nzi++; dcount++; 9426cf997b0SBarry Smith } 9436cf997b0SBarry Smith 9445a8e39fbSHong Zhang /* add pivot rows into the active row */ 9455a8e39fbSHong Zhang nzbd = 0; 9465a8e39fbSHong Zhang prow = lnk[n]; 9475a8e39fbSHong Zhang while (prow < i) { 9485a8e39fbSHong Zhang nnz = bdiag[prow]; 9495a8e39fbSHong Zhang cols = bj_ptr[prow] + nnz + 1; 9505a8e39fbSHong Zhang cols_lvl = bjlvl_ptr[prow] + nnz + 1; 9515a8e39fbSHong Zhang nnz = bi[prow+1] - bi[prow] - nnz - 1; 9525a8e39fbSHong Zhang ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 9535a8e39fbSHong Zhang nzi += nlnk; 9545a8e39fbSHong Zhang prow = lnk[prow]; 9555a8e39fbSHong Zhang nzbd++; 95656beaf04SBarry Smith } 9575a8e39fbSHong Zhang bdiag[i] = nzbd; 9585a8e39fbSHong Zhang bi[i+1] = bi[i] + nzi; 959ecf371e4SBarry Smith 9605a8e39fbSHong Zhang /* if free space is not available, make more free space */ 9615a8e39fbSHong Zhang if (current_space->local_remaining<nzi) { 9625a8e39fbSHong Zhang nnz = nzi*(n - i); /* estimated and max additional space needed */ 9635a8e39fbSHong Zhang ierr = GetMoreSpace(nnz,¤t_space);CHKERRQ(ierr); 9645a8e39fbSHong Zhang ierr = GetMoreSpace(nnz,¤t_space_lvl);CHKERRQ(ierr); 9655a8e39fbSHong Zhang reallocs++; 9665a8e39fbSHong Zhang } 967ecf371e4SBarry Smith 9685a8e39fbSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 9695a8e39fbSHong Zhang ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 9705a8e39fbSHong Zhang bj_ptr[i] = current_space->array; 9715a8e39fbSHong Zhang bjlvl_ptr[i] = current_space_lvl->array; 9725a8e39fbSHong Zhang 9735a8e39fbSHong Zhang /* make sure the active row i has diagonal entry */ 9745a8e39fbSHong Zhang if (*(bj_ptr[i]+bdiag[i]) != i) { 97577431f27SBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 9765a8e39fbSHong Zhang try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i); 9776cf997b0SBarry Smith } 9785a8e39fbSHong Zhang 9795a8e39fbSHong Zhang current_space->array += nzi; 9805a8e39fbSHong Zhang current_space->local_used += nzi; 9815a8e39fbSHong Zhang current_space->local_remaining -= nzi; 9825a8e39fbSHong Zhang current_space_lvl->array += nzi; 9835a8e39fbSHong Zhang current_space_lvl->local_used += nzi; 9845a8e39fbSHong Zhang current_space_lvl->local_remaining -= nzi; 9859e25ed09SBarry Smith } 9865a8e39fbSHong Zhang 987898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 988898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 9895a8e39fbSHong Zhang 9905a8e39fbSHong Zhang /* destroy list of free space and other temporary arrays */ 9915a8e39fbSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 9925a8e39fbSHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 9935a8e39fbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 9945a8e39fbSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 9955a8e39fbSHong Zhang ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 9969e25ed09SBarry Smith 99763ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 998f64065bbSBarry Smith { 9995a8e39fbSHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 100063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr); 100163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af));CHKERRQ(ierr); 100263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af));CHKERRQ(ierr); 100363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr); 1004335d9088SBarry Smith if (diagonal_fill) { 100563ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount));CHKERRQ(ierr); 1006335d9088SBarry Smith } 1007f64065bbSBarry Smith } 100863ba0a88SBarry Smith #endif 1009bbb6d6a8SBarry Smith 10109e25ed09SBarry Smith /* put together the new matrix */ 1011f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr); 1012f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 1013f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr); 101452e6d16bSBarry Smith ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr); 1015416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 1016606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 10177c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 10185a8e39fbSHong Zhang len = (bi[n] )*sizeof(PetscScalar); 10199e25ed09SBarry Smith /* the next line frees the default space generated by the Create() */ 1020606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 1021606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 1022b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 10235a8e39fbSHong Zhang b->j = bj; 10245a8e39fbSHong Zhang b->i = bi; 10255a8e39fbSHong Zhang for (i=0; i<n; i++) bdiag[i] += bi[i]; 10265a8e39fbSHong Zhang b->diag = bdiag; 1027416022c9SBarry Smith b->ilen = 0; 1028416022c9SBarry Smith b->imax = 0; 1029416022c9SBarry Smith b->row = isrow; 1030416022c9SBarry Smith b->col = iscol; 1031c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1032c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 103382bf6240SBarry Smith b->icol = isicol; 103487828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1035416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 10365a8e39fbSHong Zhang Allocate bdiag, solve_work, new a, new j */ 10375a8e39fbSHong Zhang ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 10385a8e39fbSHong Zhang b->maxnz = b->nz = bi[n] ; 10399e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 1040418422e8SSatish Balay (*fact)->info.factor_mallocs = reallocs; 1041ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 10425a8e39fbSHong Zhang (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 104371c2f376SKris Buschelman 104454e71613SHong Zhang ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr); 104571c2f376SKris Buschelman (*fact)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 10464846f1f5SKris Buschelman 10473a40ed3dSBarry Smith PetscFunctionReturn(0); 10489e25ed09SBarry Smith } 1049d5d45c9bSBarry Smith 10503c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1051a6175056SHong Zhang #undef __FUNCT__ 1052a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1053af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 1054a6175056SHong Zhang { 10552ed133dbSHong Zhang Mat C = *B; 10560968510aSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 10572ed133dbSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 10582ed133dbSHong Zhang IS ip=b->row; 10592ed133dbSHong Zhang PetscErrorCode ierr; 10602ed133dbSHong Zhang PetscInt *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol; 10612ed133dbSHong Zhang PetscInt *ai=a->i,*aj=a->j; 10622ed133dbSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1063622e2028SHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1064540703acSHong Zhang PetscReal zeropivot,rs,shiftnz; 1065540703acSHong Zhang PetscTruth shiftpd; 1066540703acSHong Zhang ChShift_Ctx sctx; 1067540703acSHong Zhang PetscInt newshift; 1068d5d45c9bSBarry Smith 1069a6175056SHong Zhang PetscFunctionBegin; 1070540703acSHong Zhang shiftnz = info->shiftnz; 1071540703acSHong Zhang shiftpd = info->shiftpd; 1072ee45ca4aSHong Zhang zeropivot = info->zeropivot; 1073ee45ca4aSHong Zhang 10742ed133dbSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 10752ed133dbSHong Zhang 10762ed133dbSHong Zhang /* initialization */ 10772ed133dbSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 10782ed133dbSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 10792ed133dbSHong Zhang jl = il + mbs; 10802ed133dbSHong Zhang rtmp = (MatScalar*)(jl + mbs); 10812ed133dbSHong Zhang 1082540703acSHong Zhang sctx.shift_amount = 0; 1083540703acSHong Zhang sctx.nshift = 0; 10842ed133dbSHong Zhang do { 1085540703acSHong Zhang sctx.chshift = PETSC_FALSE; 10862ed133dbSHong Zhang for (i=0; i<mbs; i++) { 10872ed133dbSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1088a6175056SHong Zhang } 10892ed133dbSHong Zhang 10902ed133dbSHong Zhang for (k = 0; k<mbs; k++){ 1091c05c3958SHong Zhang bval = ba + bi[k]; 10922ed133dbSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 10932ed133dbSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 10942ed133dbSHong Zhang for (j = jmin; j < jmax; j++){ 10952ed133dbSHong Zhang col = rip[aj[j]]; 10962ed133dbSHong Zhang if (col >= k){ /* only take upper triangular entry */ 10972ed133dbSHong Zhang rtmp[col] = aa[j]; 1098c05c3958SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 10992ed133dbSHong Zhang } 11002ed133dbSHong Zhang } 110139e3d630SHong Zhang /* shift the diagonal of the matrix */ 1102540703acSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 11032ed133dbSHong Zhang 11042ed133dbSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 11052ed133dbSHong Zhang dk = rtmp[k]; 11062ed133dbSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 11072ed133dbSHong Zhang 11082ed133dbSHong Zhang while (i < k){ 11092ed133dbSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 11102ed133dbSHong Zhang 11112ed133dbSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 11122ed133dbSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 11132ed133dbSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 11142ed133dbSHong Zhang dk += uikdi*ba[ili]; 11152ed133dbSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 11162ed133dbSHong Zhang 11172ed133dbSHong Zhang /* add multiple of row i to k-th row */ 11182ed133dbSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 11192ed133dbSHong Zhang if (jmin < jmax){ 11202ed133dbSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 11212ed133dbSHong Zhang /* update il and jl for row i */ 11222ed133dbSHong Zhang il[i] = jmin; 11232ed133dbSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 11242ed133dbSHong Zhang } 11252ed133dbSHong Zhang i = nexti; 11262ed133dbSHong Zhang } 11272ed133dbSHong Zhang 11280697b6caSHong Zhang /* shift the diagonals when zero pivot is detected */ 11290697b6caSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 11300697b6caSHong Zhang rs = 0.0; 11313c9e8598SHong Zhang jmin = bi[k]+1; 11323c9e8598SHong Zhang nz = bi[k+1] - jmin; 11333c9e8598SHong Zhang if (nz){ 11343c9e8598SHong Zhang bcol = bj + jmin; 11353c9e8598SHong Zhang while (nz--){ 113689f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 113789f845c8SHong Zhang bcol++; 11383c9e8598SHong Zhang } 11393c9e8598SHong Zhang } 1140540703acSHong Zhang 1141540703acSHong Zhang sctx.rs = rs; 1142540703acSHong Zhang sctx.pv = dk; 11433e8c821fSHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 1144540703acSHong Zhang if (newshift == 1){ 1145540703acSHong Zhang break; /* sctx.shift_amount is updated */ 1146540703acSHong Zhang } else if (newshift == -1){ 11470697b6caSHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 11483c9e8598SHong Zhang } 11493c9e8598SHong Zhang 11503c9e8598SHong Zhang /* copy data into U(k,:) */ 115139e3d630SHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 115239e3d630SHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 115339e3d630SHong Zhang if (jmin < jmax) { 115439e3d630SHong Zhang for (j=jmin; j<jmax; j++){ 115539e3d630SHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 11563c9e8598SHong Zhang } 115739e3d630SHong Zhang /* add the k-th row into il and jl */ 11583c9e8598SHong Zhang il[k] = jmin; 11593c9e8598SHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 11603c9e8598SHong Zhang } 11613c9e8598SHong Zhang } 1162540703acSHong Zhang } while (sctx.chshift); 11633c9e8598SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 11643c9e8598SHong Zhang 116539e3d630SHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 11663c9e8598SHong Zhang C->factor = FACTOR_CHOLESKY; 11673c9e8598SHong Zhang C->assembled = PETSC_TRUE; 11683c9e8598SHong Zhang C->preallocated = PETSC_TRUE; 1169efee365bSSatish Balay ierr = PetscLogFlops(C->m);CHKERRQ(ierr); 1170540703acSHong Zhang if (sctx.nshift){ 1171540703acSHong Zhang if (shiftnz) { 117263ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 1173540703acSHong Zhang } else if (shiftpd) { 117463ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 11750697b6caSHong Zhang } 11763c9e8598SHong Zhang } 11773c9e8598SHong Zhang PetscFunctionReturn(0); 11783c9e8598SHong Zhang } 1179a6175056SHong Zhang 1180a6175056SHong Zhang #undef __FUNCT__ 1181a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 1182dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1183a6175056SHong Zhang { 11840968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1185ed59401aSHong Zhang Mat_SeqSBAIJ *b; 1186ed59401aSHong Zhang Mat B; 1187dfbe8321SBarry Smith PetscErrorCode ierr; 1188653a6975SHong Zhang PetscTruth perm_identity; 1189622e2028SHong Zhang PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui; 1190ed59401aSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 1191391eac55SHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 11925a8e39fbSHong Zhang PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 1193ed59401aSHong Zhang PetscReal fill=info->fill,levels=info->levels; 11947a48dd6fSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 11957a48dd6fSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 11967a48dd6fSHong Zhang PetscBT lnkbt; 1197a6175056SHong Zhang 1198a6175056SHong Zhang PetscFunctionBegin; 1199653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 12007a48dd6fSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 12017a48dd6fSHong Zhang 120239e3d630SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 120339e3d630SHong Zhang ui[0] = 0; 120439e3d630SHong Zhang 1205622e2028SHong Zhang /* special case that simply copies fill pattern */ 1206622e2028SHong Zhang if (!levels && perm_identity) { 1207622e2028SHong Zhang ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 1208ed59401aSHong Zhang for (i=0; i<am; i++) { 120939e3d630SHong Zhang ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 1210ed59401aSHong Zhang } 121139e3d630SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 121239e3d630SHong Zhang cols = uj; 1213ed59401aSHong Zhang for (i=0; i<am; i++) { 1214ed59401aSHong Zhang aj = a->j + a->diag[i]; 121539e3d630SHong Zhang ncols = ui[i+1] - ui[i]; 121639e3d630SHong Zhang for (j=0; j<ncols; j++) *cols++ = *aj++; 1217ed59401aSHong Zhang } 121839e3d630SHong Zhang } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 12197a48dd6fSHong Zhang /* initialization */ 12205a8e39fbSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 12217a48dd6fSHong Zhang 12227a48dd6fSHong Zhang /* jl: linked list for storing indices of the pivot rows 12237a48dd6fSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 12241393bc97SHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 12257a48dd6fSHong Zhang il = jl + am; 12267a48dd6fSHong Zhang uj_ptr = (PetscInt**)(il + am); 12277a48dd6fSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 12287a48dd6fSHong Zhang for (i=0; i<am; i++){ 12297a48dd6fSHong Zhang jl[i] = am; il[i] = 0; 12307a48dd6fSHong Zhang } 12317a48dd6fSHong Zhang 12327a48dd6fSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 12337a48dd6fSHong Zhang nlnk = am + 1; 12342ed133dbSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12357a48dd6fSHong Zhang 12367a48dd6fSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 12377a48dd6fSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 12387a48dd6fSHong Zhang current_space = free_space; 12397a48dd6fSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 12407a48dd6fSHong Zhang current_space_lvl = free_space_lvl; 12417a48dd6fSHong Zhang 12427a48dd6fSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 12437a48dd6fSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 12447a48dd6fSHong Zhang nzk = 0; 1245622e2028SHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 1246622e2028SHong Zhang ncols_upper = 0; 1247622e2028SHong Zhang for (j=0; j<ncols; j++){ 12485a8e39fbSHong Zhang i = *(aj + ai[rip[k]] + j); 12495a8e39fbSHong Zhang if (rip[i] >= k){ /* only take upper triangular entry */ 12505a8e39fbSHong Zhang ajtmp[ncols_upper] = i; 1251622e2028SHong Zhang ncols_upper++; 1252622e2028SHong Zhang } 1253622e2028SHong Zhang } 12545a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12557a48dd6fSHong Zhang nzk += nlnk; 12567a48dd6fSHong Zhang 12577a48dd6fSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 12587a48dd6fSHong Zhang prow = jl[k]; /* 1st pivot row */ 12597a48dd6fSHong Zhang 12607a48dd6fSHong Zhang while (prow < k){ 12617a48dd6fSHong Zhang nextprow = jl[prow]; 12627a48dd6fSHong Zhang 12637a48dd6fSHong Zhang /* merge prow into k-th row */ 12647a48dd6fSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 12657a48dd6fSHong Zhang jmax = ui[prow+1]; 12667a48dd6fSHong Zhang ncols = jmax-jmin; 1267ed59401aSHong Zhang i = jmin - ui[prow]; 1268ed59401aSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 12695a8e39fbSHong Zhang uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 12705a8e39fbSHong Zhang j = *(uj - 1); 12715a8e39fbSHong Zhang ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 12727a48dd6fSHong Zhang nzk += nlnk; 12737a48dd6fSHong Zhang 12747a48dd6fSHong Zhang /* update il and jl for prow */ 12757a48dd6fSHong Zhang if (jmin < jmax){ 12767a48dd6fSHong Zhang il[prow] = jmin; 12777a48dd6fSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 12787a48dd6fSHong Zhang } 12797a48dd6fSHong Zhang prow = nextprow; 12807a48dd6fSHong Zhang } 12817a48dd6fSHong Zhang 12827a48dd6fSHong Zhang /* if free space is not available, make more free space */ 12837a48dd6fSHong Zhang if (current_space->local_remaining<nzk) { 12847a48dd6fSHong Zhang i = am - k + 1; /* num of unfactored rows */ 12857a48dd6fSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 12867a48dd6fSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 12877a48dd6fSHong Zhang ierr = GetMoreSpace(i,¤t_space_lvl);CHKERRQ(ierr); 12887a48dd6fSHong Zhang reallocs++; 12897a48dd6fSHong Zhang } 12907a48dd6fSHong Zhang 12917a48dd6fSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 12922ed133dbSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 12937a48dd6fSHong Zhang 12947a48dd6fSHong Zhang /* add the k-th row into il and jl */ 129539e3d630SHong Zhang if (nzk > 1){ 12967a48dd6fSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 12977a48dd6fSHong Zhang jl[k] = jl[i]; jl[i] = k; 12987a48dd6fSHong Zhang il[k] = ui[k] + 1; 12997a48dd6fSHong Zhang } 13007a48dd6fSHong Zhang uj_ptr[k] = current_space->array; 13017a48dd6fSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 13027a48dd6fSHong Zhang 13037a48dd6fSHong Zhang current_space->array += nzk; 13047a48dd6fSHong Zhang current_space->local_used += nzk; 13057a48dd6fSHong Zhang current_space->local_remaining -= nzk; 13067a48dd6fSHong Zhang 13077a48dd6fSHong Zhang current_space_lvl->array += nzk; 13087a48dd6fSHong Zhang current_space_lvl->local_used += nzk; 13097a48dd6fSHong Zhang current_space_lvl->local_remaining -= nzk; 13107a48dd6fSHong Zhang 13117a48dd6fSHong Zhang ui[k+1] = ui[k] + nzk; 13127a48dd6fSHong Zhang } 13137a48dd6fSHong Zhang 131463ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 13157a48dd6fSHong Zhang if (ai[am] != 0) { 131639e3d630SHong Zhang PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 131763ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr); 131863ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr); 131963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr); 13207a48dd6fSHong Zhang } else { 132163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr); 13227a48dd6fSHong Zhang } 132363ba0a88SBarry Smith #endif 13247a48dd6fSHong Zhang 13257a48dd6fSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 13267a48dd6fSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 13275a8e39fbSHong Zhang ierr = PetscFree(ajtmp);CHKERRQ(ierr); 13287a48dd6fSHong Zhang 13297a48dd6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 13307a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 13317a48dd6fSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 13322ed133dbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 13337a48dd6fSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 13347a48dd6fSHong Zhang 133539e3d630SHong Zhang } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 133639e3d630SHong Zhang 13377a48dd6fSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 13387a48dd6fSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 1339ed59401aSHong Zhang B = *fact; 1340ed59401aSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 1341ed59401aSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 13427a48dd6fSHong Zhang 1343ed59401aSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 13447a48dd6fSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 13457a48dd6fSHong Zhang b->singlemalloc = PETSC_FALSE; 13467a48dd6fSHong Zhang /* the next line frees the default space generated by the Create() */ 13477a48dd6fSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 13487a48dd6fSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 13497a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 13507a48dd6fSHong Zhang b->j = uj; 13517a48dd6fSHong Zhang b->i = ui; 13527a48dd6fSHong Zhang b->diag = 0; 13537a48dd6fSHong Zhang b->ilen = 0; 13547a48dd6fSHong Zhang b->imax = 0; 13557a48dd6fSHong Zhang b->row = perm; 13567a48dd6fSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 13577a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 13587a48dd6fSHong Zhang b->icol = perm; 13597a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 13607a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 136152e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 13627a48dd6fSHong Zhang b->maxnz = b->nz = ui[am]; 13637a48dd6fSHong Zhang 1364ed59401aSHong Zhang B->factor = FACTOR_CHOLESKY; 1365ed59401aSHong Zhang B->info.factor_mallocs = reallocs; 1366ed59401aSHong Zhang B->info.fill_ratio_given = fill; 13677a48dd6fSHong Zhang if (ai[am] != 0) { 1368ed59401aSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 13697a48dd6fSHong Zhang } else { 1370ed59401aSHong Zhang B->info.fill_ratio_needed = 0.0; 13717a48dd6fSHong Zhang } 137239e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1373622e2028SHong Zhang if (perm_identity){ 1374ed59401aSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1375ed59401aSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1376622e2028SHong Zhang } 1377a6175056SHong Zhang PetscFunctionReturn(0); 1378a6175056SHong Zhang } 1379d5d45c9bSBarry Smith 1380f76d2b81SHong Zhang #undef __FUNCT__ 1381f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1382dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1383f76d2b81SHong Zhang { 1384f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 138510c27e3dSHong Zhang Mat_SeqSBAIJ *b; 13862ed133dbSHong Zhang Mat B; 1387dfbe8321SBarry Smith PetscErrorCode ierr; 1388f76d2b81SHong Zhang PetscTruth perm_identity; 138910c27e3dSHong Zhang PetscReal fill = info->fill; 13901393bc97SHong Zhang PetscInt *rip,*riip,i,am=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow; 139110c27e3dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 13922ed133dbSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 139310c27e3dSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 139410c27e3dSHong Zhang PetscBT lnkbt; 13952ed133dbSHong Zhang IS iperm; 1396f76d2b81SHong Zhang 1397f76d2b81SHong Zhang PetscFunctionBegin; 13982ed133dbSHong Zhang /* check whether perm is the identity mapping */ 1399f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 14002ed133dbSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 14012ed133dbSHong Zhang 1402f76d2b81SHong Zhang if (!perm_identity){ 14032ed133dbSHong Zhang /* check if perm is symmetric! */ 14042ed133dbSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 14052ed133dbSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 14061393bc97SHong Zhang for (i=0; i<am; i++) { 14072ed133dbSHong Zhang if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation"); 14082ed133dbSHong Zhang } 14092ed133dbSHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 14102ed133dbSHong Zhang ierr = ISDestroy(iperm);CHKERRQ(ierr); 1411f76d2b81SHong Zhang } 141210c27e3dSHong Zhang 141310c27e3dSHong Zhang /* initialization */ 14141393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 141510c27e3dSHong Zhang ui[0] = 0; 141610c27e3dSHong Zhang 141710c27e3dSHong Zhang /* jl: linked list for storing indices of the pivot rows 14181393bc97SHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 14191393bc97SHong Zhang ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 14201393bc97SHong Zhang il = jl + am; 14211393bc97SHong Zhang cols = il + am; 14221393bc97SHong Zhang ui_ptr = (PetscInt**)(cols + am); 14231393bc97SHong Zhang for (i=0; i<am; i++){ 14241393bc97SHong Zhang jl[i] = am; il[i] = 0; 1425f76d2b81SHong Zhang } 142610c27e3dSHong Zhang 142710c27e3dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 14281393bc97SHong Zhang nlnk = am + 1; 14291393bc97SHong Zhang ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 143010c27e3dSHong Zhang 14311393bc97SHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 14321393bc97SHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 143310c27e3dSHong Zhang current_space = free_space; 143410c27e3dSHong Zhang 14351393bc97SHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 143610c27e3dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 143710c27e3dSHong Zhang nzk = 0; 14382ed133dbSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 14392ed133dbSHong Zhang ncols_upper = 0; 14402ed133dbSHong Zhang for (j=0; j<ncols; j++){ 1441622e2028SHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 14422ed133dbSHong Zhang if (i >= k){ /* only take upper triangular entry */ 14432ed133dbSHong Zhang cols[ncols_upper] = i; 14442ed133dbSHong Zhang ncols_upper++; 14452ed133dbSHong Zhang } 14462ed133dbSHong Zhang } 14471393bc97SHong Zhang ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 144810c27e3dSHong Zhang nzk += nlnk; 144910c27e3dSHong Zhang 145010c27e3dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 145110c27e3dSHong Zhang prow = jl[k]; /* 1st pivot row */ 145210c27e3dSHong Zhang 145310c27e3dSHong Zhang while (prow < k){ 145410c27e3dSHong Zhang nextprow = jl[prow]; 145510c27e3dSHong Zhang /* merge prow into k-th row */ 14561393bc97SHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 145710c27e3dSHong Zhang jmax = ui[prow+1]; 145810c27e3dSHong Zhang ncols = jmax-jmin; 14591393bc97SHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 14605a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 146110c27e3dSHong Zhang nzk += nlnk; 146210c27e3dSHong Zhang 146310c27e3dSHong Zhang /* update il and jl for prow */ 146410c27e3dSHong Zhang if (jmin < jmax){ 146510c27e3dSHong Zhang il[prow] = jmin; 14662ed133dbSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 146710c27e3dSHong Zhang } 146810c27e3dSHong Zhang prow = nextprow; 146910c27e3dSHong Zhang } 147010c27e3dSHong Zhang 147110c27e3dSHong Zhang /* if free space is not available, make more free space */ 147210c27e3dSHong Zhang if (current_space->local_remaining<nzk) { 14731393bc97SHong Zhang i = am - k + 1; /* num of unfactored rows */ 147410c27e3dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 147510c27e3dSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 147610c27e3dSHong Zhang reallocs++; 147710c27e3dSHong Zhang } 147810c27e3dSHong Zhang 147910c27e3dSHong Zhang /* copy data into free space, then initialize lnk */ 14801393bc97SHong Zhang ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 148110c27e3dSHong Zhang 148210c27e3dSHong Zhang /* add the k-th row into il and jl */ 148310c27e3dSHong Zhang if (nzk-1 > 0){ 14841393bc97SHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 148510c27e3dSHong Zhang jl[k] = jl[i]; jl[i] = k; 148610c27e3dSHong Zhang il[k] = ui[k] + 1; 148710c27e3dSHong Zhang } 14882ed133dbSHong Zhang ui_ptr[k] = current_space->array; 148910c27e3dSHong Zhang current_space->array += nzk; 149010c27e3dSHong Zhang current_space->local_used += nzk; 149110c27e3dSHong Zhang current_space->local_remaining -= nzk; 149210c27e3dSHong Zhang 149310c27e3dSHong Zhang ui[k+1] = ui[k] + nzk; 149410c27e3dSHong Zhang } 149510c27e3dSHong Zhang 149663ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 14971393bc97SHong Zhang if (ai[am] != 0) { 14981393bc97SHong Zhang PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]); 149963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr); 150063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr); 150163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr); 150210c27e3dSHong Zhang } else { 150363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr); 150410c27e3dSHong Zhang } 150563ba0a88SBarry Smith #endif 150610c27e3dSHong Zhang 150710c27e3dSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 150810c27e3dSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 150910c27e3dSHong Zhang 151010c27e3dSHong Zhang /* destroy list of free space and other temporary array(s) */ 15111393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 151210c27e3dSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 151310c27e3dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 151410c27e3dSHong Zhang 151510c27e3dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 15161393bc97SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 15172ed133dbSHong Zhang B = *fact; 15182ed133dbSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 15192ed133dbSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 152010c27e3dSHong Zhang 15212ed133dbSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 152210c27e3dSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 152310c27e3dSHong Zhang b->singlemalloc = PETSC_FALSE; 152410c27e3dSHong Zhang /* the next line frees the default space generated by the Create() */ 152510c27e3dSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 152610c27e3dSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 15271393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 152810c27e3dSHong Zhang b->j = uj; 152910c27e3dSHong Zhang b->i = ui; 153010c27e3dSHong Zhang b->diag = 0; 153110c27e3dSHong Zhang b->ilen = 0; 153210c27e3dSHong Zhang b->imax = 0; 153310c27e3dSHong Zhang b->row = perm; 153410c27e3dSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 153510c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 153610c27e3dSHong Zhang b->icol = perm; 153710c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 15381393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 15391393bc97SHong Zhang ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 15401393bc97SHong Zhang b->maxnz = b->nz = ui[am]; 154110c27e3dSHong Zhang 15422ed133dbSHong Zhang B->factor = FACTOR_CHOLESKY; 15432ed133dbSHong Zhang B->info.factor_mallocs = reallocs; 15442ed133dbSHong Zhang B->info.fill_ratio_given = fill; 15451393bc97SHong Zhang if (ai[am] != 0) { 15461393bc97SHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 154710c27e3dSHong Zhang } else { 15482ed133dbSHong Zhang B->info.fill_ratio_needed = 0.0; 154910c27e3dSHong Zhang } 155039e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 15512ed133dbSHong Zhang if (perm_identity){ 155210c27e3dSHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 155310c27e3dSHong Zhang (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 15542ed133dbSHong Zhang } 1555f76d2b81SHong Zhang PetscFunctionReturn(0); 1556f76d2b81SHong Zhang } 1557