1*be1d678aSKris Buschelman #define PETSCMAT_DLL 2289bc588SBarry Smith 370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 41c4feecaSSatish Balay #include "src/inline/dot.h" 5ed480e8bSBarry Smith #include "src/inline/spops.h" 61393bc97SHong Zhang #include "petscbt.h" 71393bc97SHong Zhang #include "src/mat/utils/freespace.h" 83b2fbd54SBarry Smith 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ" 11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol) 123b2fbd54SBarry Smith { 133a40ed3dSBarry Smith PetscFunctionBegin; 143a40ed3dSBarry Smith 1529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not written"); 16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG) 17d64ed03dSBarry Smith PetscFunctionReturn(0); 18d64ed03dSBarry Smith #endif 193b2fbd54SBarry Smith } 203b2fbd54SBarry Smith 2186bacbd2SBarry Smith 22dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 2386bacbd2SBarry Smith 2438baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 259cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*); 2638baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*); 2786bacbd2SBarry Smith 284a2ae208SSatish Balay #undef __FUNCT__ 294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 3086bacbd2SBarry Smith /* ------------------------------------------------------------ 3186bacbd2SBarry Smith 3286bacbd2SBarry Smith This interface was contribed by Tony Caola 3386bacbd2SBarry Smith 3486bacbd2SBarry Smith This routine is an interface to the pivoting drop-tolerance 355bd2ddc7SBarry Smith ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of 365bd2ddc7SBarry Smith SPARSEKIT2. 375bd2ddc7SBarry Smith 385bd2ddc7SBarry Smith The SPARSEKIT2 routines used here are covered by the GNU 395bd2ddc7SBarry Smith copyright; see the file gnu in this directory. 4086bacbd2SBarry Smith 4186bacbd2SBarry Smith Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their 4286bacbd2SBarry Smith help in getting this routine ironed out. 4386bacbd2SBarry Smith 445bd2ddc7SBarry Smith The major drawback to this routine is that if info->fill is 455bd2ddc7SBarry Smith not large enough it fails rather than allocating more space; 465bd2ddc7SBarry Smith this can be fixed by hacking/improving the f2c version of 475bd2ddc7SBarry Smith Yousef Saad's code. 4886bacbd2SBarry Smith 4986bacbd2SBarry Smith ------------------------------------------------------------ 5086bacbd2SBarry Smith */ 517529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 5286bacbd2SBarry Smith { 5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE) 5460ec86bdSBarry Smith PetscFunctionBegin; 55e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\ 5660ec86bdSBarry Smith You can obtain the drop tolerance routines by installing PETSc from\n\ 5760ec86bdSBarry Smith www.mcs.anl.gov/petsc\n"); 5860ec86bdSBarry Smith #else 5986bacbd2SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 6007d2056aSBarry Smith IS iscolf,isicol,isirow; 6186bacbd2SBarry Smith PetscTruth reorder; 629cc1f4e3SBarry Smith PetscErrorCode ierr,sierr; 639cc1f4e3SBarry Smith PetscInt *c,*r,*ic,i,n = A->m; 6438baddfdSBarry Smith PetscInt *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j; 6538baddfdSBarry Smith PetscInt *ordcol,*iwk,*iperm,*jw; 6638baddfdSBarry Smith PetscInt jmax,lfill,job,*o_i,*o_j; 6787828ca2SBarry Smith PetscScalar *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a; 6854a8161fSBarry Smith PetscReal af; 6986bacbd2SBarry Smith 7086bacbd2SBarry Smith PetscFunctionBegin; 7186bacbd2SBarry Smith 7286bacbd2SBarry Smith if (info->dt == PETSC_DEFAULT) info->dt = .005; 7338baddfdSBarry Smith if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax); 7486bacbd2SBarry Smith if (info->dtcol == PETSC_DEFAULT) info->dtcol = .01; 7533259db3SBarry Smith if (info->fill == PETSC_DEFAULT) info->fill = ((double)(n*(info->dtcount+1)))/a->nz; 7638baddfdSBarry Smith lfill = (PetscInt)(info->dtcount/2.0); 7738baddfdSBarry Smith jmax = (PetscInt)(info->fill*a->nz); 7886bacbd2SBarry Smith 7986bacbd2SBarry Smith 8086bacbd2SBarry Smith /* ------------------------------------------------------------ 8186bacbd2SBarry Smith If reorder=.TRUE., then the original matrix has to be 8286bacbd2SBarry Smith reordered to reflect the user selected ordering scheme, and 8386bacbd2SBarry Smith then de-reordered so it is in it's original format. 8486bacbd2SBarry Smith Because Saad's dperm() is NOT in place, we have to copy 8586bacbd2SBarry Smith the original matrix and allocate more storage. . . 8686bacbd2SBarry Smith ------------------------------------------------------------ 8786bacbd2SBarry Smith */ 8886bacbd2SBarry Smith 8986bacbd2SBarry Smith /* set reorder to true if either isrow or iscol is not identity */ 9086bacbd2SBarry Smith ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr); 9186bacbd2SBarry Smith if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);} 9286bacbd2SBarry Smith reorder = PetscNot(reorder); 9386bacbd2SBarry Smith 9486bacbd2SBarry Smith 9586bacbd2SBarry Smith /* storage for ilu factor */ 9638baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr); 9738baddfdSBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr); 9887828ca2SBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr); 9938baddfdSBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr); 10086bacbd2SBarry Smith 10186bacbd2SBarry Smith /* ------------------------------------------------------------ 10286bacbd2SBarry Smith Make sure that everything is Fortran formatted (1-Based) 10386bacbd2SBarry Smith ------------------------------------------------------------ 10486bacbd2SBarry Smith */ 105b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 106b19713cbSBarry Smith old_j[i]++; 10786bacbd2SBarry Smith } 108b19713cbSBarry Smith for(i=0;i<n+1;i++) { 109b19713cbSBarry Smith old_i[i]++; 110b19713cbSBarry Smith }; 111010a20caSHong Zhang 11286bacbd2SBarry Smith 11386bacbd2SBarry Smith if (reorder) { 114c0c2c81eSBarry Smith ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr); 115c0c2c81eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 116c0c2c81eSBarry Smith for(i=0;i<n;i++) { 117c0c2c81eSBarry Smith r[i] = r[i]+1; 118c0c2c81eSBarry Smith c[i] = c[i]+1; 119c0c2c81eSBarry Smith } 12038baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr); 12138baddfdSBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr); 12287828ca2SBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr); 1235bd2ddc7SBarry Smith job = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job); 124c0c2c81eSBarry Smith for (i=0;i<n;i++) { 125c0c2c81eSBarry Smith r[i] = r[i]-1; 126c0c2c81eSBarry Smith c[i] = c[i]-1; 127c0c2c81eSBarry Smith } 128c0c2c81eSBarry Smith ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 129c0c2c81eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 130b19713cbSBarry Smith o_a = old_a2; 131b19713cbSBarry Smith o_j = old_j2; 132b19713cbSBarry Smith o_i = old_i2; 133b19713cbSBarry Smith } else { 134b19713cbSBarry Smith o_a = old_a; 135b19713cbSBarry Smith o_j = old_j; 136b19713cbSBarry Smith o_i = old_i; 13786bacbd2SBarry Smith } 13886bacbd2SBarry Smith 13986bacbd2SBarry Smith /* ------------------------------------------------------------ 14086bacbd2SBarry Smith Call Saad's ilutp() routine to generate the factorization 14186bacbd2SBarry Smith ------------------------------------------------------------ 14286bacbd2SBarry Smith */ 14386bacbd2SBarry Smith 14438baddfdSBarry Smith ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr); 14538baddfdSBarry Smith ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr); 14687828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr); 14786bacbd2SBarry Smith 14854a8161fSBarry Smith SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&info->dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr); 1496849ba73SBarry Smith if (sierr) { 1506849ba73SBarry Smith switch (sierr) { 15177431f27SBarry Smith case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax); 15277431f27SBarry Smith case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax); 153e005ede5SBarry Smith case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered"); 154e005ede5SBarry Smith case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong"); 15577431f27SBarry Smith case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax); 15677431f27SBarry Smith default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr); 15786bacbd2SBarry Smith } 15886bacbd2SBarry Smith } 15986bacbd2SBarry Smith 16086bacbd2SBarry Smith ierr = PetscFree(w);CHKERRQ(ierr); 16186bacbd2SBarry Smith ierr = PetscFree(jw);CHKERRQ(ierr); 16286bacbd2SBarry Smith 16386bacbd2SBarry Smith /* ------------------------------------------------------------ 16486bacbd2SBarry Smith Saad's routine gives the result in Modified Sparse Row (msr) 16586bacbd2SBarry Smith Convert to Compressed Sparse Row format (csr) 16686bacbd2SBarry Smith ------------------------------------------------------------ 16786bacbd2SBarry Smith */ 16886bacbd2SBarry Smith 16987828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr); 17038baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr); 17186bacbd2SBarry Smith 1725bd2ddc7SBarry Smith SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk); 17386bacbd2SBarry Smith 17486bacbd2SBarry Smith ierr = PetscFree(iwk);CHKERRQ(ierr); 17586bacbd2SBarry Smith ierr = PetscFree(wk);CHKERRQ(ierr); 17686bacbd2SBarry Smith 17786bacbd2SBarry Smith if (reorder) { 17886bacbd2SBarry Smith ierr = PetscFree(old_a2);CHKERRQ(ierr); 17986bacbd2SBarry Smith ierr = PetscFree(old_j2);CHKERRQ(ierr); 18086bacbd2SBarry Smith ierr = PetscFree(old_i2);CHKERRQ(ierr); 181313ae024SBarry Smith } else { 182b19713cbSBarry Smith /* fix permutation of old_j that the factorization introduced */ 183f170005cSBarry Smith for (i=old_i[0]; i<old_i[n]; i++) { 184b19713cbSBarry Smith old_j[i-1] = iperm[old_j[i-1]-1]; 185b19713cbSBarry Smith } 186b19713cbSBarry Smith } 187b19713cbSBarry Smith 188b801d0f9SBarry Smith /* get rid of the shift to indices starting at 1 */ 18986bacbd2SBarry Smith for (i=0; i<n+1; i++) { 190b19713cbSBarry Smith old_i[i]--; 191b19713cbSBarry Smith } 192b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 193b19713cbSBarry Smith old_j[i]--; 194b19713cbSBarry Smith } 19586bacbd2SBarry Smith 196b19713cbSBarry Smith /* Make the factored matrix 0-based */ 19786bacbd2SBarry Smith for (i=0; i<n+1; i++) { 198b19713cbSBarry Smith new_i[i]--; 199b19713cbSBarry Smith } 200b19713cbSBarry Smith for (i=new_i[0];i<new_i[n];i++) { 201b19713cbSBarry Smith new_j[i]--; 202b19713cbSBarry Smith } 20386bacbd2SBarry Smith 20486bacbd2SBarry Smith /*-- due to the pivoting, we need to reorder iscol to correctly --*/ 20586bacbd2SBarry Smith /*-- permute the right-hand-side and solution vectors --*/ 2064c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2074c49b128SBarry Smith ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr); 208c0c2c81eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 20986bacbd2SBarry Smith for(i=0; i<n; i++) { 21086bacbd2SBarry Smith ordcol[i] = ic[iperm[i]-1]; 21186bacbd2SBarry Smith }; 21286bacbd2SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 21307d2056aSBarry Smith ierr = ISDestroy(isicol);CHKERRQ(ierr); 21486bacbd2SBarry Smith 215c0c2c81eSBarry Smith ierr = PetscFree(iperm);CHKERRQ(ierr); 216c0c2c81eSBarry Smith 217329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr); 21807d2056aSBarry Smith ierr = PetscFree(ordcol);CHKERRQ(ierr); 21986bacbd2SBarry Smith 22086bacbd2SBarry Smith /*----- put together the new matrix -----*/ 22186bacbd2SBarry Smith 222f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr); 223f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 224f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr); 22586bacbd2SBarry Smith (*fact)->factor = FACTOR_LU; 22686bacbd2SBarry Smith (*fact)->assembled = PETSC_TRUE; 22786bacbd2SBarry Smith 22886bacbd2SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 22986bacbd2SBarry Smith ierr = PetscFree(b->imax);CHKERRQ(ierr); 23086bacbd2SBarry Smith b->sorted = PETSC_FALSE; 23107d2056aSBarry Smith b->singlemalloc = PETSC_FALSE; 232b19713cbSBarry Smith /* the next line frees the default space generated by the MatCreate() */ 23386bacbd2SBarry Smith ierr = PetscFree(b->a);CHKERRQ(ierr); 23486bacbd2SBarry Smith ierr = PetscFree(b->ilen);CHKERRQ(ierr); 23586bacbd2SBarry Smith b->a = new_a; 23686bacbd2SBarry Smith b->j = new_j; 23786bacbd2SBarry Smith b->i = new_i; 23886bacbd2SBarry Smith b->ilen = 0; 23986bacbd2SBarry Smith b->imax = 0; 2401f9e874aSBarry Smith /* I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */ 241313ae024SBarry Smith b->row = isirow; 24286bacbd2SBarry Smith b->col = iscolf; 24387828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 24486bacbd2SBarry Smith b->maxnz = b->nz = new_i[n]; 24586bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 24686bacbd2SBarry Smith (*fact)->info.factor_mallocs = 0; 24786bacbd2SBarry Smith 24886bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 24986bacbd2SBarry Smith 250431e710aSBarry Smith af = ((double)b->nz)/((double)a->nz) + .001; 25163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af));CHKERRQ(ierr); 25263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af));CHKERRQ(ierr); 25363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af));CHKERRQ(ierr); 25463ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:for best performance.\n"));CHKERRQ(ierr); 255431e710aSBarry Smith 2567529efd4SKris Buschelman ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr); 25771c2f376SKris Buschelman 25886bacbd2SBarry Smith PetscFunctionReturn(0); 25960ec86bdSBarry Smith #endif 26086bacbd2SBarry Smith } 26186bacbd2SBarry Smith 2624a2ae208SSatish Balay #undef __FUNCT__ 263b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ" 264dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B) 265289bc588SBarry Smith { 266416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 267289bc588SBarry Smith IS isicol; 2686849ba73SBarry Smith PetscErrorCode ierr; 26938baddfdSBarry Smith PetscInt *r,*ic,i,n=A->m,*ai=a->i,*aj=a->j; 2701393bc97SHong Zhang PetscInt *bi,*bj,*ajtmp; 2711393bc97SHong Zhang PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 2729e046878SBarry Smith PetscReal f; 2731393bc97SHong Zhang PetscInt nlnk,*lnk,k,*cols,**bi_ptr; 2741393bc97SHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2751393bc97SHong Zhang PetscBT lnkbt; 276289bc588SBarry Smith 2773a40ed3dSBarry Smith PetscFunctionBegin; 27829bbc08cSBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 2794c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 280ac121b3dSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 281ac121b3dSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 282289bc588SBarry Smith 283289bc588SBarry Smith /* get new row pointers */ 2841393bc97SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 2851393bc97SHong Zhang bi[0] = 0; 2861393bc97SHong Zhang 2871393bc97SHong Zhang /* bdiag is location of diagonal in factor */ 2881393bc97SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 2891393bc97SHong Zhang bdiag[0] = 0; 2901393bc97SHong Zhang 2911393bc97SHong Zhang /* linked list for storing column indices of the active row */ 2921393bc97SHong Zhang nlnk = n + 1; 2931393bc97SHong Zhang ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2941393bc97SHong Zhang 2951393bc97SHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt)+n*sizeof(PetscInt**),&cols);CHKERRQ(ierr); 2961393bc97SHong Zhang im = cols + n; 2971393bc97SHong Zhang bi_ptr = (PetscInt**)(im + n); 2981393bc97SHong Zhang 2991393bc97SHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 300d3d32019SBarry Smith f = info->fill; 3011393bc97SHong Zhang ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 3021393bc97SHong Zhang current_space = free_space; 303289bc588SBarry Smith 304289bc588SBarry Smith for (i=0; i<n; i++) { 3051393bc97SHong Zhang /* copy previous fill into linked list */ 306289bc588SBarry Smith nzi = 0; 3071393bc97SHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 3081393bc97SHong Zhang if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 3091393bc97SHong Zhang ajtmp = aj + ai[r[i]]; 3101393bc97SHong Zhang for (k=0; k<nnz; k++) cols[k] = ic[*(ajtmp+k)]; /* note: cols is not sorted when iscol!=indentity */ 3111393bc97SHong Zhang ierr = PetscLLAdd(nnz,cols,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 3121393bc97SHong Zhang nzi += nlnk; 3131393bc97SHong Zhang 3141393bc97SHong Zhang /* add pivot rows into linked list */ 3151393bc97SHong Zhang row = lnk[n]; 3161393bc97SHong Zhang while (row < i) { 3171393bc97SHong Zhang nzbd = bdiag[row] - bi[row] + 1; 3181393bc97SHong Zhang ajtmp = bi_ptr[row] + nzbd; 3191393bc97SHong Zhang nnz = im[row] - nzbd; /* num of columns with row<indices<=i */ 3201393bc97SHong Zhang im[row] = nzbd; 3211393bc97SHong Zhang ierr = PetscLLAddSortedLU(nnz,ajtmp,row,nlnk,lnk,lnkbt,i,nzbd);CHKERRQ(ierr); 3221393bc97SHong Zhang nzi += nlnk; 3231393bc97SHong Zhang im[row] += nzbd; /* update im[row]: num of cols with index<=i */ 3241393bc97SHong Zhang 3251393bc97SHong Zhang row = lnk[row]; 326289bc588SBarry Smith } 3271393bc97SHong Zhang 3281393bc97SHong Zhang bi[i+1] = bi[i] + nzi; 3291393bc97SHong Zhang im[i] = nzi; 3301393bc97SHong Zhang 3311393bc97SHong Zhang /* mark bdiag */ 3321393bc97SHong Zhang nzbd = 0; 3331393bc97SHong Zhang nnz = nzi; 3341393bc97SHong Zhang k = lnk[n]; 3351393bc97SHong Zhang while (nnz-- && k < i){ 3361393bc97SHong Zhang nzbd++; 3371393bc97SHong Zhang k = lnk[k]; 3381393bc97SHong Zhang } 3391393bc97SHong Zhang bdiag[i] = bi[i] + nzbd; 3401393bc97SHong Zhang 3411393bc97SHong Zhang /* if free space is not available, make more free space */ 3421393bc97SHong Zhang if (current_space->local_remaining<nzi) { 3431393bc97SHong Zhang nnz = (n - i)*nzi; /* estimated and max additional space needed */ 3441393bc97SHong Zhang ierr = GetMoreSpace(nnz,¤t_space);CHKERRQ(ierr); 3451393bc97SHong Zhang reallocs++; 3461393bc97SHong Zhang } 3471393bc97SHong Zhang 3481393bc97SHong Zhang /* copy data into free space, then initialize lnk */ 3491393bc97SHong Zhang ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 3501393bc97SHong Zhang bi_ptr[i] = current_space->array; 3511393bc97SHong Zhang current_space->array += nzi; 3521393bc97SHong Zhang current_space->local_used += nzi; 3531393bc97SHong Zhang current_space->local_remaining -= nzi; 354289bc588SBarry Smith } 35563ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 356464e8d44SSatish Balay if (ai[n] != 0) { 3571393bc97SHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 35863ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr); 35963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af));CHKERRQ(ierr); 36063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af));CHKERRQ(ierr); 36163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr); 36205bf559cSSatish Balay } else { 36363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n"));CHKERRQ(ierr); 36405bf559cSSatish Balay } 36563ba0a88SBarry Smith #endif 3662fb3ed76SBarry Smith 367898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 368898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3691987afe7SBarry Smith 3701393bc97SHong Zhang /* destroy list of free space and other temporary array(s) */ 3711393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 3721393bc97SHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 3731393bc97SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 3741393bc97SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 375289bc588SBarry Smith 376289bc588SBarry Smith /* put together the new matrix */ 377f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr); 378f204ca49SKris Buschelman ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr); 379f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr); 38052e6d16bSBarry Smith ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr); 381416022c9SBarry Smith b = (Mat_SeqAIJ*)(*B)->data; 382606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 3837c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 384e8d4e0b9SBarry Smith /* the next line frees the default space generated by the Create() */ 385606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 386606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 3871393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 3881393bc97SHong Zhang b->j = bj; 3891393bc97SHong Zhang b->i = bi; 3901393bc97SHong Zhang b->diag = bdiag; 391416022c9SBarry Smith b->ilen = 0; 392416022c9SBarry Smith b->imax = 0; 393416022c9SBarry Smith b->row = isrow; 394416022c9SBarry Smith b->col = iscol; 395c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 396c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 39782bf6240SBarry Smith b->icol = isicol; 39887828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 3991393bc97SHong Zhang 4001393bc97SHong Zhang /* In b structure: Free imax, ilen, old a, old j. Allocate solve_work, new a, new j */ 4011393bc97SHong Zhang ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 4021393bc97SHong Zhang b->maxnz = b->nz = bi[n] ; 4037fd98bd6SLois Curfman McInnes 404329f5518SBarry Smith (*B)->factor = FACTOR_LU; 405418422e8SSatish Balay (*B)->info.factor_mallocs = reallocs; 406ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 407703deb49SSatish Balay 408e93ccc4dSBarry Smith if (ai[n] != 0) { 4091393bc97SHong Zhang (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 410eea2913aSSatish Balay } else { 411eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 412eea2913aSSatish Balay } 4134846f1f5SKris Buschelman ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr); 41471c2f376SKris Buschelman (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 4153a40ed3dSBarry Smith PetscFunctionReturn(0); 416289bc588SBarry Smith } 4171393bc97SHong Zhang 4180a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 4194a2ae208SSatish Balay #undef __FUNCT__ 4204a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 421af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 422289bc588SBarry Smith { 423416022c9SBarry Smith Mat C=*B; 424416022c9SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 42582bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 4266849ba73SBarry Smith PetscErrorCode ierr; 427b3bf805bSHong Zhang PetscInt *r,*ic,i,j,n=A->m,*bi=b->i,*bj=b->j; 428b3bf805bSHong Zhang PetscInt *ajtmp,*bjtmp,nz,row,*ics; 429d2276718SHong Zhang PetscInt *diag_offset = b->diag,diag,*pj; 43087828ca2SBarry Smith PetscScalar *rtmp,*v,*pc,multiplier,*pv,*rtmps; 4316a7c8fc2SHong Zhang PetscScalar d; 4326a7c8fc2SHong Zhang PetscReal rs; 433b3bf805bSHong Zhang LUShift_Ctx sctx; 434d2276718SHong Zhang PetscInt newshift; 435289bc588SBarry Smith 4363a40ed3dSBarry Smith PetscFunctionBegin; 43778b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 43878b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 43987828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 44087828ca2SBarry Smith ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 441010a20caSHong Zhang rtmps = rtmp; ics = ic; 442289bc588SBarry Smith 4436cc28720Svictorle if (!a->diag) { 4440cf777b8SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 4450cf777b8SBarry Smith } 446118f5789SBarry Smith /* if both shift schemes are chosen by user, only use info->shiftpd */ 447118f5789SBarry Smith if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0; 4481a890ab1SHong Zhang if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 44938baddfdSBarry Smith PetscInt *aai = a->i,*ddiag = a->diag; 450118f5789SBarry Smith sctx.shift_top = 0; 4516cc28720Svictorle for (i=0; i<n; i++) { 4529f95998fSHong Zhang /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 4539f95998fSHong Zhang d = (a->a)[ddiag[i]]; 4541a890ab1SHong Zhang rs = -PetscAbsScalar(d) - PetscRealPart(d); 455010a20caSHong Zhang v = a->a+aai[i]; 456e24cfe64SHong Zhang nz = aai[i+1] - aai[i]; 457e24cfe64SHong Zhang for (j=0; j<nz; j++) 4581a890ab1SHong Zhang rs += PetscAbsScalar(v[j]); 4591a890ab1SHong Zhang if (rs>sctx.shift_top) sctx.shift_top = rs; 4606cc28720Svictorle } 461118f5789SBarry Smith if (sctx.shift_top == 0.0) sctx.shift_top += 1.e-12; 462118f5789SBarry Smith sctx.shift_top *= 1.1; 463d2276718SHong Zhang sctx.nshift_max = 5; 464d2276718SHong Zhang sctx.shift_lo = 0.; 465d2276718SHong Zhang sctx.shift_hi = 1.; 466a0a356efSHong Zhang } 467a0a356efSHong Zhang 468a0a356efSHong Zhang sctx.shift_amount = 0; 469a0a356efSHong Zhang sctx.nshift = 0; 470085256b3SBarry Smith do { 471d2276718SHong Zhang sctx.lushift = PETSC_FALSE; 472289bc588SBarry Smith for (i=0; i<n; i++){ 473b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 474b3bf805bSHong Zhang bjtmp = bj + bi[i]; 475b3bf805bSHong Zhang for (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0; 476289bc588SBarry Smith 477289bc588SBarry Smith /* load in initial (unfactored row) */ 478416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 479b3bf805bSHong Zhang ajtmp = a->j + a->i[r[i]]; 480010a20caSHong Zhang v = a->a + a->i[r[i]]; 481085256b3SBarry Smith for (j=0; j<nz; j++) { 482b3bf805bSHong Zhang rtmp[ics[ajtmp[j]]] = v[j]; 483085256b3SBarry Smith } 484d2276718SHong Zhang rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 485289bc588SBarry Smith 486b3bf805bSHong Zhang row = *bjtmp++; 487f2582111SSatish Balay while (row < i) { 4888c37ef55SBarry Smith pc = rtmp + row; 4898c37ef55SBarry Smith if (*pc != 0.0) { 490010a20caSHong Zhang pv = b->a + diag_offset[row]; 491010a20caSHong Zhang pj = b->j + diag_offset[row] + 1; 49235aab85fSBarry Smith multiplier = *pc / *pv++; 4938c37ef55SBarry Smith *pc = multiplier; 494b3bf805bSHong Zhang nz = bi[row+1] - diag_offset[row] - 1; 495f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 496efee365bSSatish Balay ierr = PetscLogFlops(2*nz);CHKERRQ(ierr); 497289bc588SBarry Smith } 498b3bf805bSHong Zhang row = *bjtmp++; 499289bc588SBarry Smith } 500416022c9SBarry Smith /* finished row so stick it into b->a */ 501b3bf805bSHong Zhang pv = b->a + bi[i] ; 502b3bf805bSHong Zhang pj = b->j + bi[i] ; 503b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 504b3bf805bSHong Zhang diag = diag_offset[i] - bi[i]; 505d3d32019SBarry Smith rs = 0.0; 506d3d32019SBarry Smith for (j=0; j<nz; j++) { 507d3d32019SBarry Smith pv[j] = rtmps[pj[j]]; 508d3d32019SBarry Smith if (j != diag) rs += PetscAbsScalar(pv[j]); 509d3d32019SBarry Smith } 510d2276718SHong Zhang 511b3bf805bSHong Zhang /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 512d2276718SHong Zhang sctx.rs = rs; 513d2276718SHong Zhang sctx.pv = pv[diag]; 5143e8c821fSHong Zhang ierr = MatLUCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 515d2276718SHong Zhang if (newshift == 1){ 516b3bf805bSHong Zhang break; /* sctx.shift_amount is updated */ 517d2276718SHong Zhang } else if (newshift == -1){ 518d2276718SHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs); 519d708749eSBarry Smith } 52035aab85fSBarry Smith } 521d2276718SHong Zhang 522118f5789SBarry Smith if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 5236cc28720Svictorle /* 5246c037d1bSvictorle * if no shift in this attempt & shifting & started shifting & can refine, 5256cc28720Svictorle * then try lower shift 5266cc28720Svictorle */ 527d2276718SHong Zhang sctx.shift_hi = info->shift_fraction; 528d2276718SHong Zhang info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 529d2276718SHong Zhang sctx.shift_amount = info->shift_fraction * sctx.shift_top; 530d2276718SHong Zhang sctx.lushift = PETSC_TRUE; 531d2276718SHong Zhang sctx.nshift++; 5326cc28720Svictorle } 533d2276718SHong Zhang } while (sctx.lushift); 5340f11f9aeSSatish Balay 53535aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 53635aab85fSBarry Smith for (i=0; i<n; i++) { 537010a20caSHong Zhang b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]]; 53835aab85fSBarry Smith } 53935aab85fSBarry Smith 540606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 54178b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 54278b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 543416022c9SBarry Smith C->factor = FACTOR_LU; 5447dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 545c456f294SBarry Smith C->assembled = PETSC_TRUE; 546efee365bSSatish Balay ierr = PetscLogFlops(C->n);CHKERRQ(ierr); 547d2276718SHong Zhang if (sctx.nshift){ 548118f5789SBarry Smith if (info->shiftnz) { 54963ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatLUFactorNumeric_SeqAIJ: number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 550118f5789SBarry Smith } else if (info->shiftpd) { 55163ba0a88SBarry 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); 5526cc28720Svictorle } 5530697b6caSHong Zhang } 5543a40ed3dSBarry Smith PetscFunctionReturn(0); 555289bc588SBarry Smith } 5560889b5dcSvictorle 5570889b5dcSvictorle #undef __FUNCT__ 5580889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ" 559dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A) 5600889b5dcSvictorle { 5610889b5dcSvictorle PetscFunctionBegin; 5620889b5dcSvictorle A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ; 5630889b5dcSvictorle A->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 5640889b5dcSvictorle PetscFunctionReturn(0); 5650889b5dcSvictorle } 5660889b5dcSvictorle 5670889b5dcSvictorle 5680a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5694a2ae208SSatish Balay #undef __FUNCT__ 5704a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 571dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 572da3a660dSBarry Smith { 573dfbe8321SBarry Smith PetscErrorCode ierr; 574416022c9SBarry Smith Mat C; 575416022c9SBarry Smith 5763a40ed3dSBarry Smith PetscFunctionBegin; 5779e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 578af281ebdSHong Zhang ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 579273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 58052e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr); 5813a40ed3dSBarry Smith PetscFunctionReturn(0); 582da3a660dSBarry Smith } 5830a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5844a2ae208SSatish Balay #undef __FUNCT__ 5854a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 586dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 5878c37ef55SBarry Smith { 588416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 589416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 5906849ba73SBarry Smith PetscErrorCode ierr; 59138baddfdSBarry Smith PetscInt *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 59238baddfdSBarry Smith PetscInt nz,*rout,*cout; 59387828ca2SBarry Smith PetscScalar *x,*b,*tmp,*tmps,*aa = a->a,sum,*v; 5948c37ef55SBarry Smith 5953a40ed3dSBarry Smith PetscFunctionBegin; 5963a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 59791bf9adeSBarry Smith 5981ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 5991ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 600416022c9SBarry Smith tmp = a->solve_work; 6018c37ef55SBarry Smith 6023b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6033b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 6048c37ef55SBarry Smith 6058c37ef55SBarry Smith /* forward solve the lower triangular */ 6068c37ef55SBarry Smith tmp[0] = b[*r++]; 607010a20caSHong Zhang tmps = tmp; 6088c37ef55SBarry Smith for (i=1; i<n; i++) { 609010a20caSHong Zhang v = aa + ai[i] ; 610010a20caSHong Zhang vi = aj + ai[i] ; 61153e7425aSBarry Smith nz = a->diag[i] - ai[i]; 61253e7425aSBarry Smith sum = b[*r++]; 613ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 6148c37ef55SBarry Smith tmp[i] = sum; 6158c37ef55SBarry Smith } 6168c37ef55SBarry Smith 6178c37ef55SBarry Smith /* backward solve the upper triangular */ 6188c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 619010a20caSHong Zhang v = aa + a->diag[i] + 1; 620010a20caSHong Zhang vi = aj + a->diag[i] + 1; 621416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 6228c37ef55SBarry Smith sum = tmp[i]; 623ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 624010a20caSHong Zhang x[*c--] = tmp[i] = sum*aa[a->diag[i]]; 6258c37ef55SBarry Smith } 6268c37ef55SBarry Smith 6273b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6283b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 6291ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 6301ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 631efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr); 6323a40ed3dSBarry Smith PetscFunctionReturn(0); 6338c37ef55SBarry Smith } 634026e39d0SSatish Balay 635930ae53cSBarry Smith /* ----------------------------------------------------------- */ 6364a2ae208SSatish Balay #undef __FUNCT__ 6374a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 638dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 639930ae53cSBarry Smith { 640930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 6416849ba73SBarry Smith PetscErrorCode ierr; 64238baddfdSBarry Smith PetscInt n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag; 643362ced78SSatish Balay PetscScalar *x,*b,*aa = a->a; 644aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 64538baddfdSBarry Smith PetscInt adiag_i,i,*vi,nz,ai_i; 646362ced78SSatish Balay PetscScalar *v,sum; 647d85afc42SSatish Balay #endif 648930ae53cSBarry Smith 6493a40ed3dSBarry Smith PetscFunctionBegin; 6503a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 651930ae53cSBarry Smith 6521ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 6531ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 654930ae53cSBarry Smith 655aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 6561c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 6571c4feecaSSatish Balay #else 658930ae53cSBarry Smith /* forward solve the lower triangular */ 659930ae53cSBarry Smith x[0] = b[0]; 660930ae53cSBarry Smith for (i=1; i<n; i++) { 661930ae53cSBarry Smith ai_i = ai[i]; 662930ae53cSBarry Smith v = aa + ai_i; 663930ae53cSBarry Smith vi = aj + ai_i; 664930ae53cSBarry Smith nz = adiag[i] - ai_i; 665930ae53cSBarry Smith sum = b[i]; 666930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 667930ae53cSBarry Smith x[i] = sum; 668930ae53cSBarry Smith } 669930ae53cSBarry Smith 670930ae53cSBarry Smith /* backward solve the upper triangular */ 671930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 672930ae53cSBarry Smith adiag_i = adiag[i]; 673d708749eSBarry Smith v = aa + adiag_i + 1; 674d708749eSBarry Smith vi = aj + adiag_i + 1; 675930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 676930ae53cSBarry Smith sum = x[i]; 677930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 678930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 679930ae53cSBarry Smith } 6801c4feecaSSatish Balay #endif 681efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr); 6821ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 6831ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6843a40ed3dSBarry Smith PetscFunctionReturn(0); 685930ae53cSBarry Smith } 686930ae53cSBarry Smith 6874a2ae208SSatish Balay #undef __FUNCT__ 6884a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 689dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 690da3a660dSBarry Smith { 691416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 692416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 6936849ba73SBarry Smith PetscErrorCode ierr; 69438baddfdSBarry Smith PetscInt *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 69538baddfdSBarry Smith PetscInt nz,*rout,*cout; 69687828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,sum,*v; 697da3a660dSBarry Smith 6983a40ed3dSBarry Smith PetscFunctionBegin; 69978b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 700da3a660dSBarry Smith 7011ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 7021ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 703416022c9SBarry Smith tmp = a->solve_work; 704da3a660dSBarry Smith 7053b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7063b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 707da3a660dSBarry Smith 708da3a660dSBarry Smith /* forward solve the lower triangular */ 709da3a660dSBarry Smith tmp[0] = b[*r++]; 710da3a660dSBarry Smith for (i=1; i<n; i++) { 711010a20caSHong Zhang v = aa + ai[i] ; 712010a20caSHong Zhang vi = aj + ai[i] ; 713416022c9SBarry Smith nz = a->diag[i] - ai[i]; 714da3a660dSBarry Smith sum = b[*r++]; 715010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 716da3a660dSBarry Smith tmp[i] = sum; 717da3a660dSBarry Smith } 718da3a660dSBarry Smith 719da3a660dSBarry Smith /* backward solve the upper triangular */ 720da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 721010a20caSHong Zhang v = aa + a->diag[i] + 1; 722010a20caSHong Zhang vi = aj + a->diag[i] + 1; 723416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 724da3a660dSBarry Smith sum = tmp[i]; 725010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 726010a20caSHong Zhang tmp[i] = sum*aa[a->diag[i]]; 727da3a660dSBarry Smith x[*c--] += tmp[i]; 728da3a660dSBarry Smith } 729da3a660dSBarry Smith 7303b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7313b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7321ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 7331ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 734efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 735898183ecSLois Curfman McInnes 7363a40ed3dSBarry Smith PetscFunctionReturn(0); 737da3a660dSBarry Smith } 738da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 7394a2ae208SSatish Balay #undef __FUNCT__ 7404a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 741dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 742da3a660dSBarry Smith { 743416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7442235e667SBarry Smith IS iscol = a->col,isrow = a->row; 7456849ba73SBarry Smith PetscErrorCode ierr; 74638baddfdSBarry Smith PetscInt *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 74738baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 74887828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v,s1; 749da3a660dSBarry Smith 7503a40ed3dSBarry Smith PetscFunctionBegin; 7511ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 7521ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 753416022c9SBarry Smith tmp = a->solve_work; 754da3a660dSBarry Smith 7552235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7562235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 757da3a660dSBarry Smith 758da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 7592235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 760da3a660dSBarry Smith 761da3a660dSBarry Smith /* forward solve the U^T */ 762da3a660dSBarry Smith for (i=0; i<n; i++) { 763010a20caSHong Zhang v = aa + diag[i] ; 764010a20caSHong Zhang vi = aj + diag[i] + 1; 765f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 766c38d4ed2SBarry Smith s1 = tmp[i]; 767ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 768da3a660dSBarry Smith while (nz--) { 769010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*s1; 770da3a660dSBarry Smith } 771c38d4ed2SBarry Smith tmp[i] = s1; 772da3a660dSBarry Smith } 773da3a660dSBarry Smith 774da3a660dSBarry Smith /* backward solve the L^T */ 775da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 776010a20caSHong Zhang v = aa + diag[i] - 1 ; 777010a20caSHong Zhang vi = aj + diag[i] - 1 ; 778f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 779f1af5d2fSBarry Smith s1 = tmp[i]; 780da3a660dSBarry Smith while (nz--) { 781010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*s1; 782da3a660dSBarry Smith } 783da3a660dSBarry Smith } 784da3a660dSBarry Smith 785da3a660dSBarry Smith /* copy tmp into x according to permutation */ 786da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 787da3a660dSBarry Smith 7882235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7892235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7901ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 7911ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 792da3a660dSBarry Smith 793efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz-A->n);CHKERRQ(ierr); 7943a40ed3dSBarry Smith PetscFunctionReturn(0); 795da3a660dSBarry Smith } 796da3a660dSBarry Smith 7974a2ae208SSatish Balay #undef __FUNCT__ 7984a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 799dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 800da3a660dSBarry Smith { 801416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 8022235e667SBarry Smith IS iscol = a->col,isrow = a->row; 8036849ba73SBarry Smith PetscErrorCode ierr; 80438baddfdSBarry Smith PetscInt *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 80538baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 80687828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v; 8076abc6512SBarry Smith 8083a40ed3dSBarry Smith PetscFunctionBegin; 80923bc6035SMatthew Knepley if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);} 8106abc6512SBarry Smith 8111ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 8121ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 813416022c9SBarry Smith tmp = a->solve_work; 8146abc6512SBarry Smith 8152235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8162235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 8176abc6512SBarry Smith 8186abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 8192235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 8206abc6512SBarry Smith 8216abc6512SBarry Smith /* forward solve the U^T */ 8226abc6512SBarry Smith for (i=0; i<n; i++) { 823010a20caSHong Zhang v = aa + diag[i] ; 824010a20caSHong Zhang vi = aj + diag[i] + 1; 825f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 8266abc6512SBarry Smith tmp[i] *= *v++; 8276abc6512SBarry Smith while (nz--) { 828010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*tmp[i]; 8296abc6512SBarry Smith } 8306abc6512SBarry Smith } 8316abc6512SBarry Smith 8326abc6512SBarry Smith /* backward solve the L^T */ 8336abc6512SBarry Smith for (i=n-1; i>=0; i--){ 834010a20caSHong Zhang v = aa + diag[i] - 1 ; 835010a20caSHong Zhang vi = aj + diag[i] - 1 ; 836f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 8376abc6512SBarry Smith while (nz--) { 838010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*tmp[i]; 8396abc6512SBarry Smith } 8406abc6512SBarry Smith } 8416abc6512SBarry Smith 8426abc6512SBarry Smith /* copy tmp into x according to permutation */ 8436abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 8446abc6512SBarry Smith 8452235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8462235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 8471ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 8481ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8496abc6512SBarry Smith 850efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 8513a40ed3dSBarry Smith PetscFunctionReturn(0); 852da3a660dSBarry Smith } 8539e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 854dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat); 85508480c60SBarry Smith 8564a2ae208SSatish Balay #undef __FUNCT__ 8574a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 858dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 8599e25ed09SBarry Smith { 860416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 8619e25ed09SBarry Smith IS isicol; 8626849ba73SBarry Smith PetscErrorCode ierr; 8635a8e39fbSHong Zhang PetscInt *r,*ic,n=A->m,*ai=a->i,*aj=a->j; 8645a8e39fbSHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 8655a8e39fbSHong Zhang PetscInt *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0; 8665a8e39fbSHong Zhang PetscInt i,levels,diagonal_fill; 86777c4ece6SBarry Smith PetscTruth col_identity,row_identity; 868329f5518SBarry Smith PetscReal f; 8695a8e39fbSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 8705a8e39fbSHong Zhang PetscBT lnkbt; 8715a8e39fbSHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 8725a8e39fbSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 8735a8e39fbSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 8749e25ed09SBarry Smith 8753a40ed3dSBarry Smith PetscFunctionBegin; 8766cf997b0SBarry Smith f = info->fill; 87738baddfdSBarry Smith levels = (PetscInt)info->levels; 87838baddfdSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 8794c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 88082bf6240SBarry Smith 88108480c60SBarry Smith /* special case that simply copies fill pattern */ 882be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 883be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 88486bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 8852e8a6d31SBarry Smith ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 88608480c60SBarry Smith (*fact)->factor = FACTOR_LU; 88708480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 88808480c60SBarry Smith if (!b->diag) { 8897c922b88SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 89008480c60SBarry Smith } 8917c922b88SBarry Smith ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 89208480c60SBarry Smith b->row = isrow; 89308480c60SBarry Smith b->col = iscol; 89482bf6240SBarry Smith b->icol = isicol; 89587828ca2SBarry Smith ierr = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 896f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 897c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 898c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 8993a40ed3dSBarry Smith PetscFunctionReturn(0); 90008480c60SBarry Smith } 90108480c60SBarry Smith 902898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 903898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 9049e25ed09SBarry Smith 9059e25ed09SBarry Smith /* get new row pointers */ 9065a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 9075a8e39fbSHong Zhang bi[0] = 0; 9085a8e39fbSHong Zhang /* bdiag is location of diagonal in factor */ 9095a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 9105a8e39fbSHong Zhang bdiag[0] = 0; 9116cf997b0SBarry Smith 9125a8e39fbSHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr); 9135a8e39fbSHong Zhang bjlvl_ptr = (PetscInt**)(bj_ptr + n); 9145a8e39fbSHong Zhang 9155a8e39fbSHong Zhang /* create a linked list for storing column indices of the active row */ 9165a8e39fbSHong Zhang nlnk = n + 1; 9175a8e39fbSHong Zhang ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 9185a8e39fbSHong Zhang 9195a8e39fbSHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 9205a8e39fbSHong Zhang ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 9215a8e39fbSHong Zhang current_space = free_space; 9225a8e39fbSHong Zhang ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 9235a8e39fbSHong Zhang current_space_lvl = free_space_lvl; 9245a8e39fbSHong Zhang 9255a8e39fbSHong Zhang for (i=0; i<n; i++) { 9265a8e39fbSHong Zhang nzi = 0; 9276cf997b0SBarry Smith /* copy current row into linked list */ 9285a8e39fbSHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 9295a8e39fbSHong Zhang if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 9305a8e39fbSHong Zhang cols = aj + ai[r[i]]; 9315a8e39fbSHong Zhang lnk[i] = -1; /* marker to indicate if diagonal exists */ 9325a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 9335a8e39fbSHong Zhang nzi += nlnk; 9346cf997b0SBarry Smith 9356cf997b0SBarry Smith /* make sure diagonal entry is included */ 9365a8e39fbSHong Zhang if (diagonal_fill && lnk[i] == -1) { 9376cf997b0SBarry Smith fm = n; 9385a8e39fbSHong Zhang while (lnk[fm] < i) fm = lnk[fm]; 9395a8e39fbSHong Zhang lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 9405a8e39fbSHong Zhang lnk[fm] = i; 9415a8e39fbSHong Zhang lnk_lvl[i] = 0; 9425a8e39fbSHong Zhang nzi++; dcount++; 9436cf997b0SBarry Smith } 9446cf997b0SBarry Smith 9455a8e39fbSHong Zhang /* add pivot rows into the active row */ 9465a8e39fbSHong Zhang nzbd = 0; 9475a8e39fbSHong Zhang prow = lnk[n]; 9485a8e39fbSHong Zhang while (prow < i) { 9495a8e39fbSHong Zhang nnz = bdiag[prow]; 9505a8e39fbSHong Zhang cols = bj_ptr[prow] + nnz + 1; 9515a8e39fbSHong Zhang cols_lvl = bjlvl_ptr[prow] + nnz + 1; 9525a8e39fbSHong Zhang nnz = bi[prow+1] - bi[prow] - nnz - 1; 9535a8e39fbSHong Zhang ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 9545a8e39fbSHong Zhang nzi += nlnk; 9555a8e39fbSHong Zhang prow = lnk[prow]; 9565a8e39fbSHong Zhang nzbd++; 95756beaf04SBarry Smith } 9585a8e39fbSHong Zhang bdiag[i] = nzbd; 9595a8e39fbSHong Zhang bi[i+1] = bi[i] + nzi; 960ecf371e4SBarry Smith 9615a8e39fbSHong Zhang /* if free space is not available, make more free space */ 9625a8e39fbSHong Zhang if (current_space->local_remaining<nzi) { 9635a8e39fbSHong Zhang nnz = nzi*(n - i); /* estimated and max additional space needed */ 9645a8e39fbSHong Zhang ierr = GetMoreSpace(nnz,¤t_space);CHKERRQ(ierr); 9655a8e39fbSHong Zhang ierr = GetMoreSpace(nnz,¤t_space_lvl);CHKERRQ(ierr); 9665a8e39fbSHong Zhang reallocs++; 9675a8e39fbSHong Zhang } 968ecf371e4SBarry Smith 9695a8e39fbSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 9705a8e39fbSHong Zhang ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 9715a8e39fbSHong Zhang bj_ptr[i] = current_space->array; 9725a8e39fbSHong Zhang bjlvl_ptr[i] = current_space_lvl->array; 9735a8e39fbSHong Zhang 9745a8e39fbSHong Zhang /* make sure the active row i has diagonal entry */ 9755a8e39fbSHong Zhang if (*(bj_ptr[i]+bdiag[i]) != i) { 97677431f27SBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 9775a8e39fbSHong Zhang try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i); 9786cf997b0SBarry Smith } 9795a8e39fbSHong Zhang 9805a8e39fbSHong Zhang current_space->array += nzi; 9815a8e39fbSHong Zhang current_space->local_used += nzi; 9825a8e39fbSHong Zhang current_space->local_remaining -= nzi; 9835a8e39fbSHong Zhang current_space_lvl->array += nzi; 9845a8e39fbSHong Zhang current_space_lvl->local_used += nzi; 9855a8e39fbSHong Zhang current_space_lvl->local_remaining -= nzi; 9869e25ed09SBarry Smith } 9875a8e39fbSHong Zhang 988898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 989898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 9905a8e39fbSHong Zhang 9915a8e39fbSHong Zhang /* destroy list of free space and other temporary arrays */ 9925a8e39fbSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 9935a8e39fbSHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 9945a8e39fbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 9955a8e39fbSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 9965a8e39fbSHong Zhang ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 9979e25ed09SBarry Smith 99863ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 999f64065bbSBarry Smith { 10005a8e39fbSHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 100163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr); 100263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af));CHKERRQ(ierr); 100363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af));CHKERRQ(ierr); 100463ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr); 1005335d9088SBarry Smith if (diagonal_fill) { 100663ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount));CHKERRQ(ierr); 1007335d9088SBarry Smith } 1008f64065bbSBarry Smith } 100963ba0a88SBarry Smith #endif 1010bbb6d6a8SBarry Smith 10119e25ed09SBarry Smith /* put together the new matrix */ 1012f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr); 1013f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 1014f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr); 101552e6d16bSBarry Smith ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr); 1016416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 1017606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 10187c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 10195a8e39fbSHong Zhang len = (bi[n] )*sizeof(PetscScalar); 10209e25ed09SBarry Smith /* the next line frees the default space generated by the Create() */ 1021606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 1022606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 1023b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 10245a8e39fbSHong Zhang b->j = bj; 10255a8e39fbSHong Zhang b->i = bi; 10265a8e39fbSHong Zhang for (i=0; i<n; i++) bdiag[i] += bi[i]; 10275a8e39fbSHong Zhang b->diag = bdiag; 1028416022c9SBarry Smith b->ilen = 0; 1029416022c9SBarry Smith b->imax = 0; 1030416022c9SBarry Smith b->row = isrow; 1031416022c9SBarry Smith b->col = iscol; 1032c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1033c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 103482bf6240SBarry Smith b->icol = isicol; 103587828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1036416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 10375a8e39fbSHong Zhang Allocate bdiag, solve_work, new a, new j */ 10385a8e39fbSHong Zhang ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 10395a8e39fbSHong Zhang b->maxnz = b->nz = bi[n] ; 10409e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 1041418422e8SSatish Balay (*fact)->info.factor_mallocs = reallocs; 1042ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 10435a8e39fbSHong Zhang (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 104471c2f376SKris Buschelman 104554e71613SHong Zhang ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr); 104671c2f376SKris Buschelman (*fact)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 10474846f1f5SKris Buschelman 10483a40ed3dSBarry Smith PetscFunctionReturn(0); 10499e25ed09SBarry Smith } 1050d5d45c9bSBarry Smith 10513c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1052a6175056SHong Zhang #undef __FUNCT__ 1053a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1054af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 1055a6175056SHong Zhang { 10562ed133dbSHong Zhang Mat C = *B; 10570968510aSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 10582ed133dbSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 10592ed133dbSHong Zhang IS ip=b->row; 10602ed133dbSHong Zhang PetscErrorCode ierr; 10612ed133dbSHong Zhang PetscInt *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol; 10622ed133dbSHong Zhang PetscInt *ai=a->i,*aj=a->j; 10632ed133dbSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1064622e2028SHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1065540703acSHong Zhang PetscReal zeropivot,rs,shiftnz; 1066540703acSHong Zhang PetscTruth shiftpd; 1067540703acSHong Zhang ChShift_Ctx sctx; 1068540703acSHong Zhang PetscInt newshift; 1069d5d45c9bSBarry Smith 1070a6175056SHong Zhang PetscFunctionBegin; 1071540703acSHong Zhang shiftnz = info->shiftnz; 1072540703acSHong Zhang shiftpd = info->shiftpd; 1073ee45ca4aSHong Zhang zeropivot = info->zeropivot; 1074ee45ca4aSHong Zhang 10752ed133dbSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 10762ed133dbSHong Zhang 10772ed133dbSHong Zhang /* initialization */ 10782ed133dbSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 10792ed133dbSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 10802ed133dbSHong Zhang jl = il + mbs; 10812ed133dbSHong Zhang rtmp = (MatScalar*)(jl + mbs); 10822ed133dbSHong Zhang 1083540703acSHong Zhang sctx.shift_amount = 0; 1084540703acSHong Zhang sctx.nshift = 0; 10852ed133dbSHong Zhang do { 1086540703acSHong Zhang sctx.chshift = PETSC_FALSE; 10872ed133dbSHong Zhang for (i=0; i<mbs; i++) { 10882ed133dbSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1089a6175056SHong Zhang } 10902ed133dbSHong Zhang 10912ed133dbSHong Zhang for (k = 0; k<mbs; k++){ 1092c05c3958SHong Zhang bval = ba + bi[k]; 10932ed133dbSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 10942ed133dbSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 10952ed133dbSHong Zhang for (j = jmin; j < jmax; j++){ 10962ed133dbSHong Zhang col = rip[aj[j]]; 10972ed133dbSHong Zhang if (col >= k){ /* only take upper triangular entry */ 10982ed133dbSHong Zhang rtmp[col] = aa[j]; 1099c05c3958SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 11002ed133dbSHong Zhang } 11012ed133dbSHong Zhang } 110239e3d630SHong Zhang /* shift the diagonal of the matrix */ 1103540703acSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 11042ed133dbSHong Zhang 11052ed133dbSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 11062ed133dbSHong Zhang dk = rtmp[k]; 11072ed133dbSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 11082ed133dbSHong Zhang 11092ed133dbSHong Zhang while (i < k){ 11102ed133dbSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 11112ed133dbSHong Zhang 11122ed133dbSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 11132ed133dbSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 11142ed133dbSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 11152ed133dbSHong Zhang dk += uikdi*ba[ili]; 11162ed133dbSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 11172ed133dbSHong Zhang 11182ed133dbSHong Zhang /* add multiple of row i to k-th row */ 11192ed133dbSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 11202ed133dbSHong Zhang if (jmin < jmax){ 11212ed133dbSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 11222ed133dbSHong Zhang /* update il and jl for row i */ 11232ed133dbSHong Zhang il[i] = jmin; 11242ed133dbSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 11252ed133dbSHong Zhang } 11262ed133dbSHong Zhang i = nexti; 11272ed133dbSHong Zhang } 11282ed133dbSHong Zhang 11290697b6caSHong Zhang /* shift the diagonals when zero pivot is detected */ 11300697b6caSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 11310697b6caSHong Zhang rs = 0.0; 11323c9e8598SHong Zhang jmin = bi[k]+1; 11333c9e8598SHong Zhang nz = bi[k+1] - jmin; 11343c9e8598SHong Zhang if (nz){ 11353c9e8598SHong Zhang bcol = bj + jmin; 11363c9e8598SHong Zhang while (nz--){ 113789f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 113889f845c8SHong Zhang bcol++; 11393c9e8598SHong Zhang } 11403c9e8598SHong Zhang } 1141540703acSHong Zhang 1142540703acSHong Zhang sctx.rs = rs; 1143540703acSHong Zhang sctx.pv = dk; 11443e8c821fSHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 1145540703acSHong Zhang if (newshift == 1){ 1146540703acSHong Zhang break; /* sctx.shift_amount is updated */ 1147540703acSHong Zhang } else if (newshift == -1){ 11480697b6caSHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 11493c9e8598SHong Zhang } 11503c9e8598SHong Zhang 11513c9e8598SHong Zhang /* copy data into U(k,:) */ 115239e3d630SHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 115339e3d630SHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 115439e3d630SHong Zhang if (jmin < jmax) { 115539e3d630SHong Zhang for (j=jmin; j<jmax; j++){ 115639e3d630SHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 11573c9e8598SHong Zhang } 115839e3d630SHong Zhang /* add the k-th row into il and jl */ 11593c9e8598SHong Zhang il[k] = jmin; 11603c9e8598SHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 11613c9e8598SHong Zhang } 11623c9e8598SHong Zhang } 1163540703acSHong Zhang } while (sctx.chshift); 11643c9e8598SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 11653c9e8598SHong Zhang 116639e3d630SHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 11673c9e8598SHong Zhang C->factor = FACTOR_CHOLESKY; 11683c9e8598SHong Zhang C->assembled = PETSC_TRUE; 11693c9e8598SHong Zhang C->preallocated = PETSC_TRUE; 1170efee365bSSatish Balay ierr = PetscLogFlops(C->m);CHKERRQ(ierr); 1171540703acSHong Zhang if (sctx.nshift){ 1172540703acSHong Zhang if (shiftnz) { 117363ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 1174540703acSHong Zhang } else if (shiftpd) { 117563ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 11760697b6caSHong Zhang } 11773c9e8598SHong Zhang } 11783c9e8598SHong Zhang PetscFunctionReturn(0); 11793c9e8598SHong Zhang } 1180a6175056SHong Zhang 1181a6175056SHong Zhang #undef __FUNCT__ 1182a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 1183dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1184a6175056SHong Zhang { 11850968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1186ed59401aSHong Zhang Mat_SeqSBAIJ *b; 1187ed59401aSHong Zhang Mat B; 1188dfbe8321SBarry Smith PetscErrorCode ierr; 1189653a6975SHong Zhang PetscTruth perm_identity; 1190622e2028SHong Zhang PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui; 1191ed59401aSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 1192391eac55SHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 11935a8e39fbSHong Zhang PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 1194ed59401aSHong Zhang PetscReal fill=info->fill,levels=info->levels; 11957a48dd6fSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 11967a48dd6fSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 11977a48dd6fSHong Zhang PetscBT lnkbt; 1198a6175056SHong Zhang 1199a6175056SHong Zhang PetscFunctionBegin; 1200653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 12017a48dd6fSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 12027a48dd6fSHong Zhang 120339e3d630SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 120439e3d630SHong Zhang ui[0] = 0; 120539e3d630SHong Zhang 1206622e2028SHong Zhang /* special case that simply copies fill pattern */ 1207622e2028SHong Zhang if (!levels && perm_identity) { 1208622e2028SHong Zhang ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 1209ed59401aSHong Zhang for (i=0; i<am; i++) { 121039e3d630SHong Zhang ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 1211ed59401aSHong Zhang } 121239e3d630SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 121339e3d630SHong Zhang cols = uj; 1214ed59401aSHong Zhang for (i=0; i<am; i++) { 1215ed59401aSHong Zhang aj = a->j + a->diag[i]; 121639e3d630SHong Zhang ncols = ui[i+1] - ui[i]; 121739e3d630SHong Zhang for (j=0; j<ncols; j++) *cols++ = *aj++; 1218ed59401aSHong Zhang } 121939e3d630SHong Zhang } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 12207a48dd6fSHong Zhang /* initialization */ 12215a8e39fbSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 12227a48dd6fSHong Zhang 12237a48dd6fSHong Zhang /* jl: linked list for storing indices of the pivot rows 12247a48dd6fSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 12251393bc97SHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 12267a48dd6fSHong Zhang il = jl + am; 12277a48dd6fSHong Zhang uj_ptr = (PetscInt**)(il + am); 12287a48dd6fSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 12297a48dd6fSHong Zhang for (i=0; i<am; i++){ 12307a48dd6fSHong Zhang jl[i] = am; il[i] = 0; 12317a48dd6fSHong Zhang } 12327a48dd6fSHong Zhang 12337a48dd6fSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 12347a48dd6fSHong Zhang nlnk = am + 1; 12352ed133dbSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12367a48dd6fSHong Zhang 12377a48dd6fSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 12387a48dd6fSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 12397a48dd6fSHong Zhang current_space = free_space; 12407a48dd6fSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 12417a48dd6fSHong Zhang current_space_lvl = free_space_lvl; 12427a48dd6fSHong Zhang 12437a48dd6fSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 12447a48dd6fSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 12457a48dd6fSHong Zhang nzk = 0; 1246622e2028SHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 1247622e2028SHong Zhang ncols_upper = 0; 1248622e2028SHong Zhang for (j=0; j<ncols; j++){ 12495a8e39fbSHong Zhang i = *(aj + ai[rip[k]] + j); 12505a8e39fbSHong Zhang if (rip[i] >= k){ /* only take upper triangular entry */ 12515a8e39fbSHong Zhang ajtmp[ncols_upper] = i; 1252622e2028SHong Zhang ncols_upper++; 1253622e2028SHong Zhang } 1254622e2028SHong Zhang } 12555a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12567a48dd6fSHong Zhang nzk += nlnk; 12577a48dd6fSHong Zhang 12587a48dd6fSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 12597a48dd6fSHong Zhang prow = jl[k]; /* 1st pivot row */ 12607a48dd6fSHong Zhang 12617a48dd6fSHong Zhang while (prow < k){ 12627a48dd6fSHong Zhang nextprow = jl[prow]; 12637a48dd6fSHong Zhang 12647a48dd6fSHong Zhang /* merge prow into k-th row */ 12657a48dd6fSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 12667a48dd6fSHong Zhang jmax = ui[prow+1]; 12677a48dd6fSHong Zhang ncols = jmax-jmin; 1268ed59401aSHong Zhang i = jmin - ui[prow]; 1269ed59401aSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 12705a8e39fbSHong Zhang uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 12715a8e39fbSHong Zhang j = *(uj - 1); 12725a8e39fbSHong Zhang ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 12737a48dd6fSHong Zhang nzk += nlnk; 12747a48dd6fSHong Zhang 12757a48dd6fSHong Zhang /* update il and jl for prow */ 12767a48dd6fSHong Zhang if (jmin < jmax){ 12777a48dd6fSHong Zhang il[prow] = jmin; 12787a48dd6fSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 12797a48dd6fSHong Zhang } 12807a48dd6fSHong Zhang prow = nextprow; 12817a48dd6fSHong Zhang } 12827a48dd6fSHong Zhang 12837a48dd6fSHong Zhang /* if free space is not available, make more free space */ 12847a48dd6fSHong Zhang if (current_space->local_remaining<nzk) { 12857a48dd6fSHong Zhang i = am - k + 1; /* num of unfactored rows */ 12867a48dd6fSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 12877a48dd6fSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 12887a48dd6fSHong Zhang ierr = GetMoreSpace(i,¤t_space_lvl);CHKERRQ(ierr); 12897a48dd6fSHong Zhang reallocs++; 12907a48dd6fSHong Zhang } 12917a48dd6fSHong Zhang 12927a48dd6fSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 12932ed133dbSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 12947a48dd6fSHong Zhang 12957a48dd6fSHong Zhang /* add the k-th row into il and jl */ 129639e3d630SHong Zhang if (nzk > 1){ 12977a48dd6fSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 12987a48dd6fSHong Zhang jl[k] = jl[i]; jl[i] = k; 12997a48dd6fSHong Zhang il[k] = ui[k] + 1; 13007a48dd6fSHong Zhang } 13017a48dd6fSHong Zhang uj_ptr[k] = current_space->array; 13027a48dd6fSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 13037a48dd6fSHong Zhang 13047a48dd6fSHong Zhang current_space->array += nzk; 13057a48dd6fSHong Zhang current_space->local_used += nzk; 13067a48dd6fSHong Zhang current_space->local_remaining -= nzk; 13077a48dd6fSHong Zhang 13087a48dd6fSHong Zhang current_space_lvl->array += nzk; 13097a48dd6fSHong Zhang current_space_lvl->local_used += nzk; 13107a48dd6fSHong Zhang current_space_lvl->local_remaining -= nzk; 13117a48dd6fSHong Zhang 13127a48dd6fSHong Zhang ui[k+1] = ui[k] + nzk; 13137a48dd6fSHong Zhang } 13147a48dd6fSHong Zhang 131563ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 13167a48dd6fSHong Zhang if (ai[am] != 0) { 131739e3d630SHong Zhang PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 131863ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr); 131963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr); 132063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr); 13217a48dd6fSHong Zhang } else { 132263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr); 13237a48dd6fSHong Zhang } 132463ba0a88SBarry Smith #endif 13257a48dd6fSHong Zhang 13267a48dd6fSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 13277a48dd6fSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 13285a8e39fbSHong Zhang ierr = PetscFree(ajtmp);CHKERRQ(ierr); 13297a48dd6fSHong Zhang 13307a48dd6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 13317a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 13327a48dd6fSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 13332ed133dbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 13347a48dd6fSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 13357a48dd6fSHong Zhang 133639e3d630SHong Zhang } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 133739e3d630SHong Zhang 13387a48dd6fSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 13397a48dd6fSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 1340ed59401aSHong Zhang B = *fact; 1341ed59401aSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 1342ed59401aSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 13437a48dd6fSHong Zhang 1344ed59401aSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 13457a48dd6fSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 13467a48dd6fSHong Zhang b->singlemalloc = PETSC_FALSE; 13477a48dd6fSHong Zhang /* the next line frees the default space generated by the Create() */ 13487a48dd6fSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 13497a48dd6fSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 13507a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 13517a48dd6fSHong Zhang b->j = uj; 13527a48dd6fSHong Zhang b->i = ui; 13537a48dd6fSHong Zhang b->diag = 0; 13547a48dd6fSHong Zhang b->ilen = 0; 13557a48dd6fSHong Zhang b->imax = 0; 13567a48dd6fSHong Zhang b->row = perm; 13577a48dd6fSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 13587a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 13597a48dd6fSHong Zhang b->icol = perm; 13607a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 13617a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 136252e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 13637a48dd6fSHong Zhang b->maxnz = b->nz = ui[am]; 13647a48dd6fSHong Zhang 1365ed59401aSHong Zhang B->factor = FACTOR_CHOLESKY; 1366ed59401aSHong Zhang B->info.factor_mallocs = reallocs; 1367ed59401aSHong Zhang B->info.fill_ratio_given = fill; 13687a48dd6fSHong Zhang if (ai[am] != 0) { 1369ed59401aSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 13707a48dd6fSHong Zhang } else { 1371ed59401aSHong Zhang B->info.fill_ratio_needed = 0.0; 13727a48dd6fSHong Zhang } 137339e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1374622e2028SHong Zhang if (perm_identity){ 1375ed59401aSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1376ed59401aSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1377622e2028SHong Zhang } 1378a6175056SHong Zhang PetscFunctionReturn(0); 1379a6175056SHong Zhang } 1380d5d45c9bSBarry Smith 1381f76d2b81SHong Zhang #undef __FUNCT__ 1382f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1383dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1384f76d2b81SHong Zhang { 1385f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 138610c27e3dSHong Zhang Mat_SeqSBAIJ *b; 13872ed133dbSHong Zhang Mat B; 1388dfbe8321SBarry Smith PetscErrorCode ierr; 1389f76d2b81SHong Zhang PetscTruth perm_identity; 139010c27e3dSHong Zhang PetscReal fill = info->fill; 13911393bc97SHong Zhang PetscInt *rip,*riip,i,am=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow; 139210c27e3dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 13932ed133dbSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 139410c27e3dSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 139510c27e3dSHong Zhang PetscBT lnkbt; 13962ed133dbSHong Zhang IS iperm; 1397f76d2b81SHong Zhang 1398f76d2b81SHong Zhang PetscFunctionBegin; 13992ed133dbSHong Zhang /* check whether perm is the identity mapping */ 1400f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 14012ed133dbSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 14022ed133dbSHong Zhang 1403f76d2b81SHong Zhang if (!perm_identity){ 14042ed133dbSHong Zhang /* check if perm is symmetric! */ 14052ed133dbSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 14062ed133dbSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 14071393bc97SHong Zhang for (i=0; i<am; i++) { 14082ed133dbSHong Zhang if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation"); 14092ed133dbSHong Zhang } 14102ed133dbSHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 14112ed133dbSHong Zhang ierr = ISDestroy(iperm);CHKERRQ(ierr); 1412f76d2b81SHong Zhang } 141310c27e3dSHong Zhang 141410c27e3dSHong Zhang /* initialization */ 14151393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 141610c27e3dSHong Zhang ui[0] = 0; 141710c27e3dSHong Zhang 141810c27e3dSHong Zhang /* jl: linked list for storing indices of the pivot rows 14191393bc97SHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 14201393bc97SHong Zhang ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 14211393bc97SHong Zhang il = jl + am; 14221393bc97SHong Zhang cols = il + am; 14231393bc97SHong Zhang ui_ptr = (PetscInt**)(cols + am); 14241393bc97SHong Zhang for (i=0; i<am; i++){ 14251393bc97SHong Zhang jl[i] = am; il[i] = 0; 1426f76d2b81SHong Zhang } 142710c27e3dSHong Zhang 142810c27e3dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 14291393bc97SHong Zhang nlnk = am + 1; 14301393bc97SHong Zhang ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 143110c27e3dSHong Zhang 14321393bc97SHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 14331393bc97SHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 143410c27e3dSHong Zhang current_space = free_space; 143510c27e3dSHong Zhang 14361393bc97SHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 143710c27e3dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 143810c27e3dSHong Zhang nzk = 0; 14392ed133dbSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 14402ed133dbSHong Zhang ncols_upper = 0; 14412ed133dbSHong Zhang for (j=0; j<ncols; j++){ 1442622e2028SHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 14432ed133dbSHong Zhang if (i >= k){ /* only take upper triangular entry */ 14442ed133dbSHong Zhang cols[ncols_upper] = i; 14452ed133dbSHong Zhang ncols_upper++; 14462ed133dbSHong Zhang } 14472ed133dbSHong Zhang } 14481393bc97SHong Zhang ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 144910c27e3dSHong Zhang nzk += nlnk; 145010c27e3dSHong Zhang 145110c27e3dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 145210c27e3dSHong Zhang prow = jl[k]; /* 1st pivot row */ 145310c27e3dSHong Zhang 145410c27e3dSHong Zhang while (prow < k){ 145510c27e3dSHong Zhang nextprow = jl[prow]; 145610c27e3dSHong Zhang /* merge prow into k-th row */ 14571393bc97SHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 145810c27e3dSHong Zhang jmax = ui[prow+1]; 145910c27e3dSHong Zhang ncols = jmax-jmin; 14601393bc97SHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 14615a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 146210c27e3dSHong Zhang nzk += nlnk; 146310c27e3dSHong Zhang 146410c27e3dSHong Zhang /* update il and jl for prow */ 146510c27e3dSHong Zhang if (jmin < jmax){ 146610c27e3dSHong Zhang il[prow] = jmin; 14672ed133dbSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 146810c27e3dSHong Zhang } 146910c27e3dSHong Zhang prow = nextprow; 147010c27e3dSHong Zhang } 147110c27e3dSHong Zhang 147210c27e3dSHong Zhang /* if free space is not available, make more free space */ 147310c27e3dSHong Zhang if (current_space->local_remaining<nzk) { 14741393bc97SHong Zhang i = am - k + 1; /* num of unfactored rows */ 147510c27e3dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 147610c27e3dSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 147710c27e3dSHong Zhang reallocs++; 147810c27e3dSHong Zhang } 147910c27e3dSHong Zhang 148010c27e3dSHong Zhang /* copy data into free space, then initialize lnk */ 14811393bc97SHong Zhang ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 148210c27e3dSHong Zhang 148310c27e3dSHong Zhang /* add the k-th row into il and jl */ 148410c27e3dSHong Zhang if (nzk-1 > 0){ 14851393bc97SHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 148610c27e3dSHong Zhang jl[k] = jl[i]; jl[i] = k; 148710c27e3dSHong Zhang il[k] = ui[k] + 1; 148810c27e3dSHong Zhang } 14892ed133dbSHong Zhang ui_ptr[k] = current_space->array; 149010c27e3dSHong Zhang current_space->array += nzk; 149110c27e3dSHong Zhang current_space->local_used += nzk; 149210c27e3dSHong Zhang current_space->local_remaining -= nzk; 149310c27e3dSHong Zhang 149410c27e3dSHong Zhang ui[k+1] = ui[k] + nzk; 149510c27e3dSHong Zhang } 149610c27e3dSHong Zhang 149763ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 14981393bc97SHong Zhang if (ai[am] != 0) { 14991393bc97SHong Zhang PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]); 150063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr); 150163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr); 150263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr); 150310c27e3dSHong Zhang } else { 150463ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr); 150510c27e3dSHong Zhang } 150663ba0a88SBarry Smith #endif 150710c27e3dSHong Zhang 150810c27e3dSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 150910c27e3dSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 151010c27e3dSHong Zhang 151110c27e3dSHong Zhang /* destroy list of free space and other temporary array(s) */ 15121393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 151310c27e3dSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 151410c27e3dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 151510c27e3dSHong Zhang 151610c27e3dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 15171393bc97SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 15182ed133dbSHong Zhang B = *fact; 15192ed133dbSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 15202ed133dbSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 152110c27e3dSHong Zhang 15222ed133dbSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 152310c27e3dSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 152410c27e3dSHong Zhang b->singlemalloc = PETSC_FALSE; 152510c27e3dSHong Zhang /* the next line frees the default space generated by the Create() */ 152610c27e3dSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 152710c27e3dSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 15281393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 152910c27e3dSHong Zhang b->j = uj; 153010c27e3dSHong Zhang b->i = ui; 153110c27e3dSHong Zhang b->diag = 0; 153210c27e3dSHong Zhang b->ilen = 0; 153310c27e3dSHong Zhang b->imax = 0; 153410c27e3dSHong Zhang b->row = perm; 153510c27e3dSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 153610c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 153710c27e3dSHong Zhang b->icol = perm; 153810c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 15391393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 15401393bc97SHong Zhang ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 15411393bc97SHong Zhang b->maxnz = b->nz = ui[am]; 154210c27e3dSHong Zhang 15432ed133dbSHong Zhang B->factor = FACTOR_CHOLESKY; 15442ed133dbSHong Zhang B->info.factor_mallocs = reallocs; 15452ed133dbSHong Zhang B->info.fill_ratio_given = fill; 15461393bc97SHong Zhang if (ai[am] != 0) { 15471393bc97SHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 154810c27e3dSHong Zhang } else { 15492ed133dbSHong Zhang B->info.fill_ratio_needed = 0.0; 155010c27e3dSHong Zhang } 155139e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 15522ed133dbSHong Zhang if (perm_identity){ 155310c27e3dSHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 155410c27e3dSHong Zhang (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 15552ed133dbSHong Zhang } 1556f76d2b81SHong Zhang PetscFunctionReturn(0); 1557f76d2b81SHong Zhang } 1558