1be1d678aSKris Buschelman #define PETSCMAT_DLL 2289bc588SBarry Smith 370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 41c4feecaSSatish Balay #include "src/inline/dot.h" 5ed480e8bSBarry Smith #include "src/inline/spops.h" 61393bc97SHong Zhang #include "petscbt.h" 71393bc97SHong Zhang #include "src/mat/utils/freespace.h" 83b2fbd54SBarry Smith 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ" 11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol) 123b2fbd54SBarry Smith { 133a40ed3dSBarry Smith PetscFunctionBegin; 143a40ed3dSBarry Smith 1529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not written"); 16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG) 17d64ed03dSBarry Smith PetscFunctionReturn(0); 18d64ed03dSBarry Smith #endif 193b2fbd54SBarry Smith } 203b2fbd54SBarry Smith 2186bacbd2SBarry Smith 22dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 2386bacbd2SBarry Smith 24e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE) 2538baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 269cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*); 2738baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*); 28e34fafa9SBarry Smith #endif 2986bacbd2SBarry Smith 304a2ae208SSatish Balay #undef __FUNCT__ 314a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 3286bacbd2SBarry Smith /* ------------------------------------------------------------ 3386bacbd2SBarry Smith 3486bacbd2SBarry Smith This interface was contribed by Tony Caola 3586bacbd2SBarry Smith 3686bacbd2SBarry Smith This routine is an interface to the pivoting drop-tolerance 375bd2ddc7SBarry Smith ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of 385bd2ddc7SBarry Smith SPARSEKIT2. 395bd2ddc7SBarry Smith 405bd2ddc7SBarry Smith The SPARSEKIT2 routines used here are covered by the GNU 415bd2ddc7SBarry Smith copyright; see the file gnu in this directory. 4286bacbd2SBarry Smith 4386bacbd2SBarry Smith Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their 4486bacbd2SBarry Smith help in getting this routine ironed out. 4586bacbd2SBarry Smith 465bd2ddc7SBarry Smith The major drawback to this routine is that if info->fill is 475bd2ddc7SBarry Smith not large enough it fails rather than allocating more space; 485bd2ddc7SBarry Smith this can be fixed by hacking/improving the f2c version of 495bd2ddc7SBarry Smith Yousef Saad's code. 5086bacbd2SBarry Smith 5186bacbd2SBarry Smith ------------------------------------------------------------ 5286bacbd2SBarry Smith */ 537529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 5486bacbd2SBarry Smith { 5560ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE) 5660ec86bdSBarry Smith PetscFunctionBegin; 57e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\ 5860ec86bdSBarry Smith You can obtain the drop tolerance routines by installing PETSc from\n\ 5960ec86bdSBarry Smith www.mcs.anl.gov/petsc\n"); 6060ec86bdSBarry Smith #else 6186bacbd2SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 6207d2056aSBarry Smith IS iscolf,isicol,isirow; 6386bacbd2SBarry Smith PetscTruth reorder; 649cc1f4e3SBarry Smith PetscErrorCode ierr,sierr; 659cc1f4e3SBarry Smith PetscInt *c,*r,*ic,i,n = A->m; 6638baddfdSBarry Smith PetscInt *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j; 6738baddfdSBarry Smith PetscInt *ordcol,*iwk,*iperm,*jw; 6838baddfdSBarry Smith PetscInt jmax,lfill,job,*o_i,*o_j; 6987828ca2SBarry Smith PetscScalar *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a; 7054a8161fSBarry Smith PetscReal af; 7186bacbd2SBarry Smith 7286bacbd2SBarry Smith PetscFunctionBegin; 7386bacbd2SBarry Smith 7486bacbd2SBarry Smith if (info->dt == PETSC_DEFAULT) info->dt = .005; 7538baddfdSBarry Smith if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax); 7686bacbd2SBarry Smith if (info->dtcol == PETSC_DEFAULT) info->dtcol = .01; 7733259db3SBarry Smith if (info->fill == PETSC_DEFAULT) info->fill = ((double)(n*(info->dtcount+1)))/a->nz; 7838baddfdSBarry Smith lfill = (PetscInt)(info->dtcount/2.0); 7938baddfdSBarry Smith jmax = (PetscInt)(info->fill*a->nz); 8086bacbd2SBarry Smith 8186bacbd2SBarry Smith 8286bacbd2SBarry Smith /* ------------------------------------------------------------ 8386bacbd2SBarry Smith If reorder=.TRUE., then the original matrix has to be 8486bacbd2SBarry Smith reordered to reflect the user selected ordering scheme, and 8586bacbd2SBarry Smith then de-reordered so it is in it's original format. 8686bacbd2SBarry Smith Because Saad's dperm() is NOT in place, we have to copy 8786bacbd2SBarry Smith the original matrix and allocate more storage. . . 8886bacbd2SBarry Smith ------------------------------------------------------------ 8986bacbd2SBarry Smith */ 9086bacbd2SBarry Smith 9186bacbd2SBarry Smith /* set reorder to true if either isrow or iscol is not identity */ 9286bacbd2SBarry Smith ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr); 9386bacbd2SBarry Smith if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);} 9486bacbd2SBarry Smith reorder = PetscNot(reorder); 9586bacbd2SBarry Smith 9686bacbd2SBarry Smith 9786bacbd2SBarry Smith /* storage for ilu factor */ 9838baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr); 9938baddfdSBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr); 10087828ca2SBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr); 10138baddfdSBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr); 10286bacbd2SBarry Smith 10386bacbd2SBarry Smith /* ------------------------------------------------------------ 10486bacbd2SBarry Smith Make sure that everything is Fortran formatted (1-Based) 10586bacbd2SBarry Smith ------------------------------------------------------------ 10686bacbd2SBarry Smith */ 107b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 108b19713cbSBarry Smith old_j[i]++; 10986bacbd2SBarry Smith } 110b19713cbSBarry Smith for(i=0;i<n+1;i++) { 111b19713cbSBarry Smith old_i[i]++; 112b19713cbSBarry Smith }; 113010a20caSHong Zhang 11486bacbd2SBarry Smith 11586bacbd2SBarry Smith if (reorder) { 116c0c2c81eSBarry Smith ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr); 117c0c2c81eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 118c0c2c81eSBarry Smith for(i=0;i<n;i++) { 119c0c2c81eSBarry Smith r[i] = r[i]+1; 120c0c2c81eSBarry Smith c[i] = c[i]+1; 121c0c2c81eSBarry Smith } 12238baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr); 12338baddfdSBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr); 12487828ca2SBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr); 1255bd2ddc7SBarry Smith job = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job); 126c0c2c81eSBarry Smith for (i=0;i<n;i++) { 127c0c2c81eSBarry Smith r[i] = r[i]-1; 128c0c2c81eSBarry Smith c[i] = c[i]-1; 129c0c2c81eSBarry Smith } 130c0c2c81eSBarry Smith ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 131c0c2c81eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 132b19713cbSBarry Smith o_a = old_a2; 133b19713cbSBarry Smith o_j = old_j2; 134b19713cbSBarry Smith o_i = old_i2; 135b19713cbSBarry Smith } else { 136b19713cbSBarry Smith o_a = old_a; 137b19713cbSBarry Smith o_j = old_j; 138b19713cbSBarry Smith o_i = old_i; 13986bacbd2SBarry Smith } 14086bacbd2SBarry Smith 14186bacbd2SBarry Smith /* ------------------------------------------------------------ 14286bacbd2SBarry Smith Call Saad's ilutp() routine to generate the factorization 14386bacbd2SBarry Smith ------------------------------------------------------------ 14486bacbd2SBarry Smith */ 14586bacbd2SBarry Smith 14638baddfdSBarry Smith ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr); 14738baddfdSBarry Smith ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr); 14887828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr); 14986bacbd2SBarry Smith 15054a8161fSBarry 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); 1516849ba73SBarry Smith if (sierr) { 1526849ba73SBarry Smith switch (sierr) { 15377431f27SBarry Smith case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax); 15477431f27SBarry Smith case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax); 155e005ede5SBarry Smith case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered"); 156e005ede5SBarry Smith case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong"); 15777431f27SBarry Smith case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax); 15877431f27SBarry Smith default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr); 15986bacbd2SBarry Smith } 16086bacbd2SBarry Smith } 16186bacbd2SBarry Smith 16286bacbd2SBarry Smith ierr = PetscFree(w);CHKERRQ(ierr); 16386bacbd2SBarry Smith ierr = PetscFree(jw);CHKERRQ(ierr); 16486bacbd2SBarry Smith 16586bacbd2SBarry Smith /* ------------------------------------------------------------ 16686bacbd2SBarry Smith Saad's routine gives the result in Modified Sparse Row (msr) 16786bacbd2SBarry Smith Convert to Compressed Sparse Row format (csr) 16886bacbd2SBarry Smith ------------------------------------------------------------ 16986bacbd2SBarry Smith */ 17086bacbd2SBarry Smith 17187828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr); 17238baddfdSBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr); 17386bacbd2SBarry Smith 1745bd2ddc7SBarry Smith SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk); 17586bacbd2SBarry Smith 17686bacbd2SBarry Smith ierr = PetscFree(iwk);CHKERRQ(ierr); 17786bacbd2SBarry Smith ierr = PetscFree(wk);CHKERRQ(ierr); 17886bacbd2SBarry Smith 17986bacbd2SBarry Smith if (reorder) { 18086bacbd2SBarry Smith ierr = PetscFree(old_a2);CHKERRQ(ierr); 18186bacbd2SBarry Smith ierr = PetscFree(old_j2);CHKERRQ(ierr); 18286bacbd2SBarry Smith ierr = PetscFree(old_i2);CHKERRQ(ierr); 183313ae024SBarry Smith } else { 184b19713cbSBarry Smith /* fix permutation of old_j that the factorization introduced */ 185f170005cSBarry Smith for (i=old_i[0]; i<old_i[n]; i++) { 186b19713cbSBarry Smith old_j[i-1] = iperm[old_j[i-1]-1]; 187b19713cbSBarry Smith } 188b19713cbSBarry Smith } 189b19713cbSBarry Smith 190b801d0f9SBarry Smith /* get rid of the shift to indices starting at 1 */ 19186bacbd2SBarry Smith for (i=0; i<n+1; i++) { 192b19713cbSBarry Smith old_i[i]--; 193b19713cbSBarry Smith } 194b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 195b19713cbSBarry Smith old_j[i]--; 196b19713cbSBarry Smith } 19786bacbd2SBarry Smith 198b19713cbSBarry Smith /* Make the factored matrix 0-based */ 19986bacbd2SBarry Smith for (i=0; i<n+1; i++) { 200b19713cbSBarry Smith new_i[i]--; 201b19713cbSBarry Smith } 202b19713cbSBarry Smith for (i=new_i[0];i<new_i[n];i++) { 203b19713cbSBarry Smith new_j[i]--; 204b19713cbSBarry Smith } 20586bacbd2SBarry Smith 20686bacbd2SBarry Smith /*-- due to the pivoting, we need to reorder iscol to correctly --*/ 20786bacbd2SBarry Smith /*-- permute the right-hand-side and solution vectors --*/ 2084c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2094c49b128SBarry Smith ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr); 210c0c2c81eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 21186bacbd2SBarry Smith for(i=0; i<n; i++) { 21286bacbd2SBarry Smith ordcol[i] = ic[iperm[i]-1]; 21386bacbd2SBarry Smith }; 21486bacbd2SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 21507d2056aSBarry Smith ierr = ISDestroy(isicol);CHKERRQ(ierr); 21686bacbd2SBarry Smith 217c0c2c81eSBarry Smith ierr = PetscFree(iperm);CHKERRQ(ierr); 218c0c2c81eSBarry Smith 219329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr); 22007d2056aSBarry Smith ierr = PetscFree(ordcol);CHKERRQ(ierr); 22186bacbd2SBarry Smith 22286bacbd2SBarry Smith /*----- put together the new matrix -----*/ 22386bacbd2SBarry Smith 224f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,fact);CHKERRQ(ierr); 225f69a0ea3SMatthew Knepley ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr); 226f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 227ab93d7beSBarry Smith ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 22886bacbd2SBarry Smith (*fact)->factor = FACTOR_LU; 22986bacbd2SBarry Smith (*fact)->assembled = PETSC_TRUE; 23086bacbd2SBarry Smith 23186bacbd2SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 232a96a251dSBarry Smith b->freedata = PETSC_TRUE; 23386bacbd2SBarry Smith b->sorted = PETSC_FALSE; 23407d2056aSBarry Smith b->singlemalloc = PETSC_FALSE; 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; 274a1a86e44SBarry Smith PetscFreeSpaceList 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; 301a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((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]]; 310*afcc9904SHong Zhang ierr = PetscLLAddPerm(nnz,ajtmp,ic,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) { 316d1170560SHong Zhang nzbd = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */ 317d1170560SHong Zhang ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 318d1170560SHong Zhang ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 3191393bc97SHong Zhang nzi += nlnk; 3201393bc97SHong Zhang row = lnk[row]; 321289bc588SBarry Smith } 3221393bc97SHong Zhang bi[i+1] = bi[i] + nzi; 3231393bc97SHong Zhang im[i] = nzi; 3241393bc97SHong Zhang 3251393bc97SHong Zhang /* mark bdiag */ 3261393bc97SHong Zhang nzbd = 0; 3271393bc97SHong Zhang nnz = nzi; 3281393bc97SHong Zhang k = lnk[n]; 3291393bc97SHong Zhang while (nnz-- && k < i){ 3301393bc97SHong Zhang nzbd++; 3311393bc97SHong Zhang k = lnk[k]; 3321393bc97SHong Zhang } 3331393bc97SHong Zhang bdiag[i] = bi[i] + nzbd; 3341393bc97SHong Zhang 3351393bc97SHong Zhang /* if free space is not available, make more free space */ 3361393bc97SHong Zhang if (current_space->local_remaining<nzi) { 3371393bc97SHong Zhang nnz = (n - i)*nzi; /* estimated and max additional space needed */ 338a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 3391393bc97SHong Zhang reallocs++; 3401393bc97SHong Zhang } 3411393bc97SHong Zhang 3421393bc97SHong Zhang /* copy data into free space, then initialize lnk */ 3431393bc97SHong Zhang ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 3441393bc97SHong Zhang bi_ptr[i] = current_space->array; 3451393bc97SHong Zhang current_space->array += nzi; 3461393bc97SHong Zhang current_space->local_used += nzi; 3471393bc97SHong Zhang current_space->local_remaining -= nzi; 348289bc588SBarry Smith } 34963ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 350464e8d44SSatish Balay if (ai[n] != 0) { 3511393bc97SHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 35263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr); 3530c451bc4SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %G or use \n",af));CHKERRQ(ierr); 35463ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af));CHKERRQ(ierr); 35563ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr); 35605bf559cSSatish Balay } else { 35763ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n"));CHKERRQ(ierr); 35805bf559cSSatish Balay } 35963ba0a88SBarry Smith #endif 3602fb3ed76SBarry Smith 361898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 362898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3631987afe7SBarry Smith 3641393bc97SHong Zhang /* destroy list of free space and other temporary array(s) */ 3651393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 366a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 3671393bc97SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 3681393bc97SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 369289bc588SBarry Smith 370289bc588SBarry Smith /* put together the new matrix */ 371f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,B);CHKERRQ(ierr); 372f69a0ea3SMatthew Knepley ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 373f204ca49SKris Buschelman ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr); 374ab93d7beSBarry Smith ierr = MatSeqAIJSetPreallocation_SeqAIJ(*B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 37552e6d16bSBarry Smith ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr); 376416022c9SBarry Smith b = (Mat_SeqAIJ*)(*B)->data; 377a96a251dSBarry Smith b->freedata = PETSC_TRUE; 3787c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 3791393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 3801393bc97SHong Zhang b->j = bj; 3811393bc97SHong Zhang b->i = bi; 3821393bc97SHong Zhang b->diag = bdiag; 383416022c9SBarry Smith b->ilen = 0; 384416022c9SBarry Smith b->imax = 0; 385416022c9SBarry Smith b->row = isrow; 386416022c9SBarry Smith b->col = iscol; 387c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 388c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 38982bf6240SBarry Smith b->icol = isicol; 39087828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 3911393bc97SHong Zhang 3921393bc97SHong Zhang /* In b structure: Free imax, ilen, old a, old j. Allocate solve_work, new a, new j */ 3931393bc97SHong Zhang ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 3941393bc97SHong Zhang b->maxnz = b->nz = bi[n] ; 3957fd98bd6SLois Curfman McInnes 396329f5518SBarry Smith (*B)->factor = FACTOR_LU; 397418422e8SSatish Balay (*B)->info.factor_mallocs = reallocs; 398ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 399703deb49SSatish Balay 400e93ccc4dSBarry Smith if (ai[n] != 0) { 4011393bc97SHong Zhang (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 402eea2913aSSatish Balay } else { 403eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 404eea2913aSSatish Balay } 4054846f1f5SKris Buschelman ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr); 40671c2f376SKris Buschelman (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 4073a40ed3dSBarry Smith PetscFunctionReturn(0); 408289bc588SBarry Smith } 4091393bc97SHong Zhang 4100a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 4114a2ae208SSatish Balay #undef __FUNCT__ 4124a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 413af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 414289bc588SBarry Smith { 415416022c9SBarry Smith Mat C=*B; 416416022c9SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 41782bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 4186849ba73SBarry Smith PetscErrorCode ierr; 419b3bf805bSHong Zhang PetscInt *r,*ic,i,j,n=A->m,*bi=b->i,*bj=b->j; 420b3bf805bSHong Zhang PetscInt *ajtmp,*bjtmp,nz,row,*ics; 421d2276718SHong Zhang PetscInt *diag_offset = b->diag,diag,*pj; 42287828ca2SBarry Smith PetscScalar *rtmp,*v,*pc,multiplier,*pv,*rtmps; 4236a7c8fc2SHong Zhang PetscScalar d; 4246a7c8fc2SHong Zhang PetscReal rs; 425b3bf805bSHong Zhang LUShift_Ctx sctx; 426d2276718SHong Zhang PetscInt newshift; 427289bc588SBarry Smith 4283a40ed3dSBarry Smith PetscFunctionBegin; 42978b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 43078b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 43187828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 43287828ca2SBarry Smith ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 433010a20caSHong Zhang rtmps = rtmp; ics = ic; 434289bc588SBarry Smith 4356cc28720Svictorle if (!a->diag) { 4360cf777b8SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 4370cf777b8SBarry Smith } 438118f5789SBarry Smith /* if both shift schemes are chosen by user, only use info->shiftpd */ 439118f5789SBarry Smith if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0; 4401a890ab1SHong Zhang if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 44138baddfdSBarry Smith PetscInt *aai = a->i,*ddiag = a->diag; 442118f5789SBarry Smith sctx.shift_top = 0; 4436cc28720Svictorle for (i=0; i<n; i++) { 4449f95998fSHong Zhang /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 4459f95998fSHong Zhang d = (a->a)[ddiag[i]]; 4461a890ab1SHong Zhang rs = -PetscAbsScalar(d) - PetscRealPart(d); 447010a20caSHong Zhang v = a->a+aai[i]; 448e24cfe64SHong Zhang nz = aai[i+1] - aai[i]; 449e24cfe64SHong Zhang for (j=0; j<nz; j++) 4501a890ab1SHong Zhang rs += PetscAbsScalar(v[j]); 4511a890ab1SHong Zhang if (rs>sctx.shift_top) sctx.shift_top = rs; 4526cc28720Svictorle } 453118f5789SBarry Smith if (sctx.shift_top == 0.0) sctx.shift_top += 1.e-12; 454118f5789SBarry Smith sctx.shift_top *= 1.1; 455d2276718SHong Zhang sctx.nshift_max = 5; 456d2276718SHong Zhang sctx.shift_lo = 0.; 457d2276718SHong Zhang sctx.shift_hi = 1.; 458a0a356efSHong Zhang } 459a0a356efSHong Zhang 460a0a356efSHong Zhang sctx.shift_amount = 0; 461a0a356efSHong Zhang sctx.nshift = 0; 462085256b3SBarry Smith do { 463d2276718SHong Zhang sctx.lushift = PETSC_FALSE; 464289bc588SBarry Smith for (i=0; i<n; i++){ 465b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 466b3bf805bSHong Zhang bjtmp = bj + bi[i]; 467b3bf805bSHong Zhang for (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0; 468289bc588SBarry Smith 469289bc588SBarry Smith /* load in initial (unfactored row) */ 470416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 471b3bf805bSHong Zhang ajtmp = a->j + a->i[r[i]]; 472010a20caSHong Zhang v = a->a + a->i[r[i]]; 473085256b3SBarry Smith for (j=0; j<nz; j++) { 474b3bf805bSHong Zhang rtmp[ics[ajtmp[j]]] = v[j]; 475085256b3SBarry Smith } 476d2276718SHong Zhang rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 477289bc588SBarry Smith 478b3bf805bSHong Zhang row = *bjtmp++; 479f2582111SSatish Balay while (row < i) { 4808c37ef55SBarry Smith pc = rtmp + row; 4818c37ef55SBarry Smith if (*pc != 0.0) { 482010a20caSHong Zhang pv = b->a + diag_offset[row]; 483010a20caSHong Zhang pj = b->j + diag_offset[row] + 1; 48435aab85fSBarry Smith multiplier = *pc / *pv++; 4858c37ef55SBarry Smith *pc = multiplier; 486b3bf805bSHong Zhang nz = bi[row+1] - diag_offset[row] - 1; 487f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 488efee365bSSatish Balay ierr = PetscLogFlops(2*nz);CHKERRQ(ierr); 489289bc588SBarry Smith } 490b3bf805bSHong Zhang row = *bjtmp++; 491289bc588SBarry Smith } 492416022c9SBarry Smith /* finished row so stick it into b->a */ 493b3bf805bSHong Zhang pv = b->a + bi[i] ; 494b3bf805bSHong Zhang pj = b->j + bi[i] ; 495b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 496b3bf805bSHong Zhang diag = diag_offset[i] - bi[i]; 497d3d32019SBarry Smith rs = 0.0; 498d3d32019SBarry Smith for (j=0; j<nz; j++) { 499d3d32019SBarry Smith pv[j] = rtmps[pj[j]]; 500d3d32019SBarry Smith if (j != diag) rs += PetscAbsScalar(pv[j]); 501d3d32019SBarry Smith } 502d2276718SHong Zhang 503b3bf805bSHong Zhang /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 504d2276718SHong Zhang sctx.rs = rs; 505d2276718SHong Zhang sctx.pv = pv[diag]; 5063e8c821fSHong Zhang ierr = MatLUCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 507d2276718SHong Zhang if (newshift == 1){ 508b3bf805bSHong Zhang break; /* sctx.shift_amount is updated */ 509d2276718SHong Zhang } else if (newshift == -1){ 510d2276718SHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs); 511d708749eSBarry Smith } 51235aab85fSBarry Smith } 513d2276718SHong Zhang 514118f5789SBarry Smith if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 5156cc28720Svictorle /* 5166c037d1bSvictorle * if no shift in this attempt & shifting & started shifting & can refine, 5176cc28720Svictorle * then try lower shift 5186cc28720Svictorle */ 519d2276718SHong Zhang sctx.shift_hi = info->shift_fraction; 520d2276718SHong Zhang info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 521d2276718SHong Zhang sctx.shift_amount = info->shift_fraction * sctx.shift_top; 522d2276718SHong Zhang sctx.lushift = PETSC_TRUE; 523d2276718SHong Zhang sctx.nshift++; 5246cc28720Svictorle } 525d2276718SHong Zhang } while (sctx.lushift); 5260f11f9aeSSatish Balay 52735aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 52835aab85fSBarry Smith for (i=0; i<n; i++) { 529010a20caSHong Zhang b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]]; 53035aab85fSBarry Smith } 53135aab85fSBarry Smith 532606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 53378b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 53478b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 535416022c9SBarry Smith C->factor = FACTOR_LU; 5367dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 537c456f294SBarry Smith C->assembled = PETSC_TRUE; 538efee365bSSatish Balay ierr = PetscLogFlops(C->n);CHKERRQ(ierr); 539d2276718SHong Zhang if (sctx.nshift){ 540118f5789SBarry Smith if (info->shiftnz) { 54163ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatLUFactorNumeric_SeqAIJ: number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 542118f5789SBarry Smith } else if (info->shiftpd) { 54363ba0a88SBarry 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); 5446cc28720Svictorle } 5450697b6caSHong Zhang } 5463a40ed3dSBarry Smith PetscFunctionReturn(0); 547289bc588SBarry Smith } 5480889b5dcSvictorle 5490889b5dcSvictorle #undef __FUNCT__ 5500889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ" 551dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A) 5520889b5dcSvictorle { 5530889b5dcSvictorle PetscFunctionBegin; 5540889b5dcSvictorle A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ; 5550889b5dcSvictorle A->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 5560889b5dcSvictorle PetscFunctionReturn(0); 5570889b5dcSvictorle } 5580889b5dcSvictorle 5590889b5dcSvictorle 5600a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5614a2ae208SSatish Balay #undef __FUNCT__ 5624a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 563dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 564da3a660dSBarry Smith { 565dfbe8321SBarry Smith PetscErrorCode ierr; 566416022c9SBarry Smith Mat C; 567416022c9SBarry Smith 5683a40ed3dSBarry Smith PetscFunctionBegin; 5699e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 570af281ebdSHong Zhang ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 571273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 57252e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr); 5733a40ed3dSBarry Smith PetscFunctionReturn(0); 574da3a660dSBarry Smith } 5750a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5764a2ae208SSatish Balay #undef __FUNCT__ 5774a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 578dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 5798c37ef55SBarry Smith { 580416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 581416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 5826849ba73SBarry Smith PetscErrorCode ierr; 58338baddfdSBarry Smith PetscInt *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 58438baddfdSBarry Smith PetscInt nz,*rout,*cout; 58587828ca2SBarry Smith PetscScalar *x,*b,*tmp,*tmps,*aa = a->a,sum,*v; 5868c37ef55SBarry Smith 5873a40ed3dSBarry Smith PetscFunctionBegin; 5883a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 58991bf9adeSBarry Smith 5901ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 5911ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 592416022c9SBarry Smith tmp = a->solve_work; 5938c37ef55SBarry Smith 5943b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 5953b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 5968c37ef55SBarry Smith 5978c37ef55SBarry Smith /* forward solve the lower triangular */ 5988c37ef55SBarry Smith tmp[0] = b[*r++]; 599010a20caSHong Zhang tmps = tmp; 6008c37ef55SBarry Smith for (i=1; i<n; i++) { 601010a20caSHong Zhang v = aa + ai[i] ; 602010a20caSHong Zhang vi = aj + ai[i] ; 60353e7425aSBarry Smith nz = a->diag[i] - ai[i]; 60453e7425aSBarry Smith sum = b[*r++]; 605ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 6068c37ef55SBarry Smith tmp[i] = sum; 6078c37ef55SBarry Smith } 6088c37ef55SBarry Smith 6098c37ef55SBarry Smith /* backward solve the upper triangular */ 6108c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 611010a20caSHong Zhang v = aa + a->diag[i] + 1; 612010a20caSHong Zhang vi = aj + a->diag[i] + 1; 613416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 6148c37ef55SBarry Smith sum = tmp[i]; 615ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 616010a20caSHong Zhang x[*c--] = tmp[i] = sum*aa[a->diag[i]]; 6178c37ef55SBarry Smith } 6188c37ef55SBarry Smith 6193b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6203b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 6211ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 6221ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 623efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr); 6243a40ed3dSBarry Smith PetscFunctionReturn(0); 6258c37ef55SBarry Smith } 626026e39d0SSatish Balay 627930ae53cSBarry Smith /* ----------------------------------------------------------- */ 6284a2ae208SSatish Balay #undef __FUNCT__ 6294a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 630dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 631930ae53cSBarry Smith { 632930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 6336849ba73SBarry Smith PetscErrorCode ierr; 63438baddfdSBarry Smith PetscInt n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag; 635362ced78SSatish Balay PetscScalar *x,*b,*aa = a->a; 636aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 63738baddfdSBarry Smith PetscInt adiag_i,i,*vi,nz,ai_i; 638362ced78SSatish Balay PetscScalar *v,sum; 639d85afc42SSatish Balay #endif 640930ae53cSBarry Smith 6413a40ed3dSBarry Smith PetscFunctionBegin; 6423a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 643930ae53cSBarry Smith 6441ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 6451ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 646930ae53cSBarry Smith 647aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 6481c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 6491c4feecaSSatish Balay #else 650930ae53cSBarry Smith /* forward solve the lower triangular */ 651930ae53cSBarry Smith x[0] = b[0]; 652930ae53cSBarry Smith for (i=1; i<n; i++) { 653930ae53cSBarry Smith ai_i = ai[i]; 654930ae53cSBarry Smith v = aa + ai_i; 655930ae53cSBarry Smith vi = aj + ai_i; 656930ae53cSBarry Smith nz = adiag[i] - ai_i; 657930ae53cSBarry Smith sum = b[i]; 658930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 659930ae53cSBarry Smith x[i] = sum; 660930ae53cSBarry Smith } 661930ae53cSBarry Smith 662930ae53cSBarry Smith /* backward solve the upper triangular */ 663930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 664930ae53cSBarry Smith adiag_i = adiag[i]; 665d708749eSBarry Smith v = aa + adiag_i + 1; 666d708749eSBarry Smith vi = aj + adiag_i + 1; 667930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 668930ae53cSBarry Smith sum = x[i]; 669930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 670930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 671930ae53cSBarry Smith } 6721c4feecaSSatish Balay #endif 673efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr); 6741ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 6751ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6763a40ed3dSBarry Smith PetscFunctionReturn(0); 677930ae53cSBarry Smith } 678930ae53cSBarry Smith 6794a2ae208SSatish Balay #undef __FUNCT__ 6804a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 681dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 682da3a660dSBarry Smith { 683416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 684416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 6856849ba73SBarry Smith PetscErrorCode ierr; 68638baddfdSBarry Smith PetscInt *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 68738baddfdSBarry Smith PetscInt nz,*rout,*cout; 68887828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,sum,*v; 689da3a660dSBarry Smith 6903a40ed3dSBarry Smith PetscFunctionBegin; 69178b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 692da3a660dSBarry Smith 6931ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 6941ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 695416022c9SBarry Smith tmp = a->solve_work; 696da3a660dSBarry Smith 6973b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6983b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 699da3a660dSBarry Smith 700da3a660dSBarry Smith /* forward solve the lower triangular */ 701da3a660dSBarry Smith tmp[0] = b[*r++]; 702da3a660dSBarry Smith for (i=1; i<n; i++) { 703010a20caSHong Zhang v = aa + ai[i] ; 704010a20caSHong Zhang vi = aj + ai[i] ; 705416022c9SBarry Smith nz = a->diag[i] - ai[i]; 706da3a660dSBarry Smith sum = b[*r++]; 707010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 708da3a660dSBarry Smith tmp[i] = sum; 709da3a660dSBarry Smith } 710da3a660dSBarry Smith 711da3a660dSBarry Smith /* backward solve the upper triangular */ 712da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 713010a20caSHong Zhang v = aa + a->diag[i] + 1; 714010a20caSHong Zhang vi = aj + a->diag[i] + 1; 715416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 716da3a660dSBarry Smith sum = tmp[i]; 717010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 718010a20caSHong Zhang tmp[i] = sum*aa[a->diag[i]]; 719da3a660dSBarry Smith x[*c--] += tmp[i]; 720da3a660dSBarry Smith } 721da3a660dSBarry Smith 7223b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7233b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7241ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 7251ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 726efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 727898183ecSLois Curfman McInnes 7283a40ed3dSBarry Smith PetscFunctionReturn(0); 729da3a660dSBarry Smith } 730da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 7314a2ae208SSatish Balay #undef __FUNCT__ 7324a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 733dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 734da3a660dSBarry Smith { 735416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7362235e667SBarry Smith IS iscol = a->col,isrow = a->row; 7376849ba73SBarry Smith PetscErrorCode ierr; 73838baddfdSBarry Smith PetscInt *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 73938baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 74087828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v,s1; 741da3a660dSBarry Smith 7423a40ed3dSBarry Smith PetscFunctionBegin; 7431ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 7441ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 745416022c9SBarry Smith tmp = a->solve_work; 746da3a660dSBarry Smith 7472235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7482235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 749da3a660dSBarry Smith 750da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 7512235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 752da3a660dSBarry Smith 753da3a660dSBarry Smith /* forward solve the U^T */ 754da3a660dSBarry Smith for (i=0; i<n; i++) { 755010a20caSHong Zhang v = aa + diag[i] ; 756010a20caSHong Zhang vi = aj + diag[i] + 1; 757f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 758c38d4ed2SBarry Smith s1 = tmp[i]; 759ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 760da3a660dSBarry Smith while (nz--) { 761010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*s1; 762da3a660dSBarry Smith } 763c38d4ed2SBarry Smith tmp[i] = s1; 764da3a660dSBarry Smith } 765da3a660dSBarry Smith 766da3a660dSBarry Smith /* backward solve the L^T */ 767da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 768010a20caSHong Zhang v = aa + diag[i] - 1 ; 769010a20caSHong Zhang vi = aj + diag[i] - 1 ; 770f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 771f1af5d2fSBarry Smith s1 = tmp[i]; 772da3a660dSBarry Smith while (nz--) { 773010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*s1; 774da3a660dSBarry Smith } 775da3a660dSBarry Smith } 776da3a660dSBarry Smith 777da3a660dSBarry Smith /* copy tmp into x according to permutation */ 778da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 779da3a660dSBarry Smith 7802235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7812235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7821ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 7831ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 784da3a660dSBarry Smith 785efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz-A->n);CHKERRQ(ierr); 7863a40ed3dSBarry Smith PetscFunctionReturn(0); 787da3a660dSBarry Smith } 788da3a660dSBarry Smith 7894a2ae208SSatish Balay #undef __FUNCT__ 7904a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 791dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 792da3a660dSBarry Smith { 793416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7942235e667SBarry Smith IS iscol = a->col,isrow = a->row; 7956849ba73SBarry Smith PetscErrorCode ierr; 79638baddfdSBarry Smith PetscInt *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 79738baddfdSBarry Smith PetscInt nz,*rout,*cout,*diag = a->diag; 79887828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v; 7996abc6512SBarry Smith 8003a40ed3dSBarry Smith PetscFunctionBegin; 80123bc6035SMatthew Knepley if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);} 8026abc6512SBarry Smith 8031ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 8041ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 805416022c9SBarry Smith tmp = a->solve_work; 8066abc6512SBarry Smith 8072235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8082235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 8096abc6512SBarry Smith 8106abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 8112235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 8126abc6512SBarry Smith 8136abc6512SBarry Smith /* forward solve the U^T */ 8146abc6512SBarry Smith for (i=0; i<n; i++) { 815010a20caSHong Zhang v = aa + diag[i] ; 816010a20caSHong Zhang vi = aj + diag[i] + 1; 817f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 8186abc6512SBarry Smith tmp[i] *= *v++; 8196abc6512SBarry Smith while (nz--) { 820010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*tmp[i]; 8216abc6512SBarry Smith } 8226abc6512SBarry Smith } 8236abc6512SBarry Smith 8246abc6512SBarry Smith /* backward solve the L^T */ 8256abc6512SBarry Smith for (i=n-1; i>=0; i--){ 826010a20caSHong Zhang v = aa + diag[i] - 1 ; 827010a20caSHong Zhang vi = aj + diag[i] - 1 ; 828f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 8296abc6512SBarry Smith while (nz--) { 830010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*tmp[i]; 8316abc6512SBarry Smith } 8326abc6512SBarry Smith } 8336abc6512SBarry Smith 8346abc6512SBarry Smith /* copy tmp into x according to permutation */ 8356abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 8366abc6512SBarry Smith 8372235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8382235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 8391ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 8401ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8416abc6512SBarry Smith 842efee365bSSatish Balay ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 8433a40ed3dSBarry Smith PetscFunctionReturn(0); 844da3a660dSBarry Smith } 8459e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 846dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat); 84708480c60SBarry Smith 8484a2ae208SSatish Balay #undef __FUNCT__ 8494a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 850dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 8519e25ed09SBarry Smith { 852416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 8539e25ed09SBarry Smith IS isicol; 8546849ba73SBarry Smith PetscErrorCode ierr; 8555a8e39fbSHong Zhang PetscInt *r,*ic,n=A->m,*ai=a->i,*aj=a->j; 8565a8e39fbSHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 8575a8e39fbSHong Zhang PetscInt *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0; 8585a8e39fbSHong Zhang PetscInt i,levels,diagonal_fill; 85977c4ece6SBarry Smith PetscTruth col_identity,row_identity; 860329f5518SBarry Smith PetscReal f; 8615a8e39fbSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 8625a8e39fbSHong Zhang PetscBT lnkbt; 8635a8e39fbSHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 864a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 865a1a86e44SBarry Smith PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 8669e25ed09SBarry Smith 8673a40ed3dSBarry Smith PetscFunctionBegin; 8686cf997b0SBarry Smith f = info->fill; 86938baddfdSBarry Smith levels = (PetscInt)info->levels; 87038baddfdSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 8714c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 87282bf6240SBarry Smith 87308480c60SBarry Smith /* special case that simply copies fill pattern */ 874be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 875be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 87686bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 8772e8a6d31SBarry Smith ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 87808480c60SBarry Smith (*fact)->factor = FACTOR_LU; 87908480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 88008480c60SBarry Smith if (!b->diag) { 8817c922b88SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 88208480c60SBarry Smith } 8837c922b88SBarry Smith ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 88408480c60SBarry Smith b->row = isrow; 88508480c60SBarry Smith b->col = iscol; 88682bf6240SBarry Smith b->icol = isicol; 88787828ca2SBarry Smith ierr = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 888f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 889c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 890c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 8913a40ed3dSBarry Smith PetscFunctionReturn(0); 89208480c60SBarry Smith } 89308480c60SBarry Smith 894898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 895898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 8969e25ed09SBarry Smith 8979e25ed09SBarry Smith /* get new row pointers */ 8985a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 8995a8e39fbSHong Zhang bi[0] = 0; 9005a8e39fbSHong Zhang /* bdiag is location of diagonal in factor */ 9015a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 9025a8e39fbSHong Zhang bdiag[0] = 0; 9036cf997b0SBarry Smith 9045a8e39fbSHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr); 9055a8e39fbSHong Zhang bjlvl_ptr = (PetscInt**)(bj_ptr + n); 9065a8e39fbSHong Zhang 9075a8e39fbSHong Zhang /* create a linked list for storing column indices of the active row */ 9085a8e39fbSHong Zhang nlnk = n + 1; 9095a8e39fbSHong Zhang ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 9105a8e39fbSHong Zhang 9115a8e39fbSHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 912a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 9135a8e39fbSHong Zhang current_space = free_space; 914a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 9155a8e39fbSHong Zhang current_space_lvl = free_space_lvl; 9165a8e39fbSHong Zhang 9175a8e39fbSHong Zhang for (i=0; i<n; i++) { 9185a8e39fbSHong Zhang nzi = 0; 9196cf997b0SBarry Smith /* copy current row into linked list */ 9205a8e39fbSHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 9215a8e39fbSHong Zhang if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 9225a8e39fbSHong Zhang cols = aj + ai[r[i]]; 9235a8e39fbSHong Zhang lnk[i] = -1; /* marker to indicate if diagonal exists */ 9245a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 9255a8e39fbSHong Zhang nzi += nlnk; 9266cf997b0SBarry Smith 9276cf997b0SBarry Smith /* make sure diagonal entry is included */ 9285a8e39fbSHong Zhang if (diagonal_fill && lnk[i] == -1) { 9296cf997b0SBarry Smith fm = n; 9305a8e39fbSHong Zhang while (lnk[fm] < i) fm = lnk[fm]; 9315a8e39fbSHong Zhang lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 9325a8e39fbSHong Zhang lnk[fm] = i; 9335a8e39fbSHong Zhang lnk_lvl[i] = 0; 9345a8e39fbSHong Zhang nzi++; dcount++; 9356cf997b0SBarry Smith } 9366cf997b0SBarry Smith 9375a8e39fbSHong Zhang /* add pivot rows into the active row */ 9385a8e39fbSHong Zhang nzbd = 0; 9395a8e39fbSHong Zhang prow = lnk[n]; 9405a8e39fbSHong Zhang while (prow < i) { 9415a8e39fbSHong Zhang nnz = bdiag[prow]; 9425a8e39fbSHong Zhang cols = bj_ptr[prow] + nnz + 1; 9435a8e39fbSHong Zhang cols_lvl = bjlvl_ptr[prow] + nnz + 1; 9445a8e39fbSHong Zhang nnz = bi[prow+1] - bi[prow] - nnz - 1; 9455a8e39fbSHong Zhang ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 9465a8e39fbSHong Zhang nzi += nlnk; 9475a8e39fbSHong Zhang prow = lnk[prow]; 9485a8e39fbSHong Zhang nzbd++; 94956beaf04SBarry Smith } 9505a8e39fbSHong Zhang bdiag[i] = nzbd; 9515a8e39fbSHong Zhang bi[i+1] = bi[i] + nzi; 952ecf371e4SBarry Smith 9535a8e39fbSHong Zhang /* if free space is not available, make more free space */ 9545a8e39fbSHong Zhang if (current_space->local_remaining<nzi) { 9555a8e39fbSHong Zhang nnz = nzi*(n - i); /* estimated and max additional space needed */ 956a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 957a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 9585a8e39fbSHong Zhang reallocs++; 9595a8e39fbSHong Zhang } 960ecf371e4SBarry Smith 9615a8e39fbSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 9625a8e39fbSHong Zhang ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 9635a8e39fbSHong Zhang bj_ptr[i] = current_space->array; 9645a8e39fbSHong Zhang bjlvl_ptr[i] = current_space_lvl->array; 9655a8e39fbSHong Zhang 9665a8e39fbSHong Zhang /* make sure the active row i has diagonal entry */ 9675a8e39fbSHong Zhang if (*(bj_ptr[i]+bdiag[i]) != i) { 96877431f27SBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 9695a8e39fbSHong Zhang try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i); 9706cf997b0SBarry Smith } 9715a8e39fbSHong Zhang 9725a8e39fbSHong Zhang current_space->array += nzi; 9735a8e39fbSHong Zhang current_space->local_used += nzi; 9745a8e39fbSHong Zhang current_space->local_remaining -= nzi; 9755a8e39fbSHong Zhang current_space_lvl->array += nzi; 9765a8e39fbSHong Zhang current_space_lvl->local_used += nzi; 9775a8e39fbSHong Zhang current_space_lvl->local_remaining -= nzi; 9789e25ed09SBarry Smith } 9795a8e39fbSHong Zhang 980898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 981898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 9825a8e39fbSHong Zhang 9835a8e39fbSHong Zhang /* destroy list of free space and other temporary arrays */ 9845a8e39fbSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 985a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 9865a8e39fbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 987a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 9885a8e39fbSHong Zhang ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 9899e25ed09SBarry Smith 99063ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 991f64065bbSBarry Smith { 9925a8e39fbSHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 99363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr); 99463ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af));CHKERRQ(ierr); 99563ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af));CHKERRQ(ierr); 99663ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr); 997335d9088SBarry Smith if (diagonal_fill) { 99863ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount));CHKERRQ(ierr); 999335d9088SBarry Smith } 1000f64065bbSBarry Smith } 100163ba0a88SBarry Smith #endif 1002bbb6d6a8SBarry Smith 10039e25ed09SBarry Smith /* put together the new matrix */ 1004f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,fact);CHKERRQ(ierr); 1005f69a0ea3SMatthew Knepley ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr); 1006f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 1007ab93d7beSBarry Smith ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 100852e6d16bSBarry Smith ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr); 1009416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 1010a96a251dSBarry Smith b->freedata = PETSC_TRUE; 10117c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 10125a8e39fbSHong Zhang len = (bi[n] )*sizeof(PetscScalar); 1013b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 10145a8e39fbSHong Zhang b->j = bj; 10155a8e39fbSHong Zhang b->i = bi; 10165a8e39fbSHong Zhang for (i=0; i<n; i++) bdiag[i] += bi[i]; 10175a8e39fbSHong Zhang b->diag = bdiag; 1018416022c9SBarry Smith b->ilen = 0; 1019416022c9SBarry Smith b->imax = 0; 1020416022c9SBarry Smith b->row = isrow; 1021416022c9SBarry Smith b->col = iscol; 1022c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1023c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 102482bf6240SBarry Smith b->icol = isicol; 102587828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1026416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 10275a8e39fbSHong Zhang Allocate bdiag, solve_work, new a, new j */ 10285a8e39fbSHong Zhang ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 10295a8e39fbSHong Zhang b->maxnz = b->nz = bi[n] ; 10309e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 1031418422e8SSatish Balay (*fact)->info.factor_mallocs = reallocs; 1032ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 10335a8e39fbSHong Zhang (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 103471c2f376SKris Buschelman 103554e71613SHong Zhang ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr); 103671c2f376SKris Buschelman (*fact)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 10374846f1f5SKris Buschelman 10383a40ed3dSBarry Smith PetscFunctionReturn(0); 10399e25ed09SBarry Smith } 1040d5d45c9bSBarry Smith 10413c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1042a6175056SHong Zhang #undef __FUNCT__ 1043a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1044af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B) 1045a6175056SHong Zhang { 10462ed133dbSHong Zhang Mat C = *B; 10470968510aSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 10482ed133dbSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 10492ed133dbSHong Zhang IS ip=b->row; 10502ed133dbSHong Zhang PetscErrorCode ierr; 10512ed133dbSHong Zhang PetscInt *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol; 10522ed133dbSHong Zhang PetscInt *ai=a->i,*aj=a->j; 10532ed133dbSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1054622e2028SHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1055540703acSHong Zhang PetscReal zeropivot,rs,shiftnz; 1056fbf22428SSatish Balay PetscReal shiftpd; 1057540703acSHong Zhang ChShift_Ctx sctx; 1058540703acSHong Zhang PetscInt newshift; 1059d5d45c9bSBarry Smith 1060a6175056SHong Zhang PetscFunctionBegin; 1061540703acSHong Zhang shiftnz = info->shiftnz; 1062540703acSHong Zhang shiftpd = info->shiftpd; 1063ee45ca4aSHong Zhang zeropivot = info->zeropivot; 1064ee45ca4aSHong Zhang 10652ed133dbSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 10662ed133dbSHong Zhang 10672ed133dbSHong Zhang /* initialization */ 10682ed133dbSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 10692ed133dbSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 10702ed133dbSHong Zhang jl = il + mbs; 10712ed133dbSHong Zhang rtmp = (MatScalar*)(jl + mbs); 10722ed133dbSHong Zhang 1073540703acSHong Zhang sctx.shift_amount = 0; 1074540703acSHong Zhang sctx.nshift = 0; 10752ed133dbSHong Zhang do { 1076540703acSHong Zhang sctx.chshift = PETSC_FALSE; 10772ed133dbSHong Zhang for (i=0; i<mbs; i++) { 10782ed133dbSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1079a6175056SHong Zhang } 10802ed133dbSHong Zhang 10812ed133dbSHong Zhang for (k = 0; k<mbs; k++){ 1082c05c3958SHong Zhang bval = ba + bi[k]; 10832ed133dbSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 10842ed133dbSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 10852ed133dbSHong Zhang for (j = jmin; j < jmax; j++){ 10862ed133dbSHong Zhang col = rip[aj[j]]; 10872ed133dbSHong Zhang if (col >= k){ /* only take upper triangular entry */ 10882ed133dbSHong Zhang rtmp[col] = aa[j]; 1089c05c3958SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 10902ed133dbSHong Zhang } 10912ed133dbSHong Zhang } 109239e3d630SHong Zhang /* shift the diagonal of the matrix */ 1093540703acSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 10942ed133dbSHong Zhang 10952ed133dbSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 10962ed133dbSHong Zhang dk = rtmp[k]; 10972ed133dbSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 10982ed133dbSHong Zhang 10992ed133dbSHong Zhang while (i < k){ 11002ed133dbSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 11012ed133dbSHong Zhang 11022ed133dbSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 11032ed133dbSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 11042ed133dbSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 11052ed133dbSHong Zhang dk += uikdi*ba[ili]; 11062ed133dbSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 11072ed133dbSHong Zhang 11082ed133dbSHong Zhang /* add multiple of row i to k-th row */ 11092ed133dbSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 11102ed133dbSHong Zhang if (jmin < jmax){ 11112ed133dbSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 11122ed133dbSHong Zhang /* update il and jl for row i */ 11132ed133dbSHong Zhang il[i] = jmin; 11142ed133dbSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 11152ed133dbSHong Zhang } 11162ed133dbSHong Zhang i = nexti; 11172ed133dbSHong Zhang } 11182ed133dbSHong Zhang 11190697b6caSHong Zhang /* shift the diagonals when zero pivot is detected */ 11200697b6caSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 11210697b6caSHong Zhang rs = 0.0; 11223c9e8598SHong Zhang jmin = bi[k]+1; 11233c9e8598SHong Zhang nz = bi[k+1] - jmin; 11243c9e8598SHong Zhang if (nz){ 11253c9e8598SHong Zhang bcol = bj + jmin; 11263c9e8598SHong Zhang while (nz--){ 112789f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 112889f845c8SHong Zhang bcol++; 11293c9e8598SHong Zhang } 11303c9e8598SHong Zhang } 1131540703acSHong Zhang 1132540703acSHong Zhang sctx.rs = rs; 1133540703acSHong Zhang sctx.pv = dk; 11343e8c821fSHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 1135540703acSHong Zhang if (newshift == 1){ 1136540703acSHong Zhang break; /* sctx.shift_amount is updated */ 1137540703acSHong Zhang } else if (newshift == -1){ 11380697b6caSHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 11393c9e8598SHong Zhang } 11403c9e8598SHong Zhang 11413c9e8598SHong Zhang /* copy data into U(k,:) */ 114239e3d630SHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 114339e3d630SHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 114439e3d630SHong Zhang if (jmin < jmax) { 114539e3d630SHong Zhang for (j=jmin; j<jmax; j++){ 114639e3d630SHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 11473c9e8598SHong Zhang } 114839e3d630SHong Zhang /* add the k-th row into il and jl */ 11493c9e8598SHong Zhang il[k] = jmin; 11503c9e8598SHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 11513c9e8598SHong Zhang } 11523c9e8598SHong Zhang } 1153540703acSHong Zhang } while (sctx.chshift); 11543c9e8598SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 11553c9e8598SHong Zhang 115639e3d630SHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 11573c9e8598SHong Zhang C->factor = FACTOR_CHOLESKY; 11583c9e8598SHong Zhang C->assembled = PETSC_TRUE; 11593c9e8598SHong Zhang C->preallocated = PETSC_TRUE; 1160efee365bSSatish Balay ierr = PetscLogFlops(C->m);CHKERRQ(ierr); 1161540703acSHong Zhang if (sctx.nshift){ 1162540703acSHong Zhang if (shiftnz) { 116363ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 1164540703acSHong Zhang } else if (shiftpd) { 116563ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 11660697b6caSHong Zhang } 11673c9e8598SHong Zhang } 11683c9e8598SHong Zhang PetscFunctionReturn(0); 11693c9e8598SHong Zhang } 1170a6175056SHong Zhang 1171a6175056SHong Zhang #undef __FUNCT__ 1172a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 1173dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1174a6175056SHong Zhang { 11750968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1176ed59401aSHong Zhang Mat_SeqSBAIJ *b; 1177ed59401aSHong Zhang Mat B; 1178dfbe8321SBarry Smith PetscErrorCode ierr; 1179653a6975SHong Zhang PetscTruth perm_identity; 1180622e2028SHong Zhang PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui; 1181ed59401aSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 1182391eac55SHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 11835a8e39fbSHong Zhang PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 1184ed59401aSHong Zhang PetscReal fill=info->fill,levels=info->levels; 1185a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1186a1a86e44SBarry Smith PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 11877a48dd6fSHong Zhang PetscBT lnkbt; 1188a6175056SHong Zhang 1189a6175056SHong Zhang PetscFunctionBegin; 1190653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 11917a48dd6fSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 11927a48dd6fSHong Zhang 119339e3d630SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 119439e3d630SHong Zhang ui[0] = 0; 119539e3d630SHong Zhang 1196622e2028SHong Zhang /* special case that simply copies fill pattern */ 1197622e2028SHong Zhang if (!levels && perm_identity) { 1198622e2028SHong Zhang ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 1199ed59401aSHong Zhang for (i=0; i<am; i++) { 120039e3d630SHong Zhang ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 1201ed59401aSHong Zhang } 120239e3d630SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 120339e3d630SHong Zhang cols = uj; 1204ed59401aSHong Zhang for (i=0; i<am; i++) { 1205ed59401aSHong Zhang aj = a->j + a->diag[i]; 120639e3d630SHong Zhang ncols = ui[i+1] - ui[i]; 120739e3d630SHong Zhang for (j=0; j<ncols; j++) *cols++ = *aj++; 1208ed59401aSHong Zhang } 120939e3d630SHong Zhang } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 12107a48dd6fSHong Zhang /* initialization */ 12115a8e39fbSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 12127a48dd6fSHong Zhang 12137a48dd6fSHong Zhang /* jl: linked list for storing indices of the pivot rows 12147a48dd6fSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 12151393bc97SHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 12167a48dd6fSHong Zhang il = jl + am; 12177a48dd6fSHong Zhang uj_ptr = (PetscInt**)(il + am); 12187a48dd6fSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 12197a48dd6fSHong Zhang for (i=0; i<am; i++){ 12207a48dd6fSHong Zhang jl[i] = am; il[i] = 0; 12217a48dd6fSHong Zhang } 12227a48dd6fSHong Zhang 12237a48dd6fSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 12247a48dd6fSHong Zhang nlnk = am + 1; 12252ed133dbSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12267a48dd6fSHong Zhang 12277a48dd6fSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 1228a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 12297a48dd6fSHong Zhang current_space = free_space; 1230a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 12317a48dd6fSHong Zhang current_space_lvl = free_space_lvl; 12327a48dd6fSHong Zhang 12337a48dd6fSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 12347a48dd6fSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 12357a48dd6fSHong Zhang nzk = 0; 1236622e2028SHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 1237622e2028SHong Zhang ncols_upper = 0; 1238622e2028SHong Zhang for (j=0; j<ncols; j++){ 12395a8e39fbSHong Zhang i = *(aj + ai[rip[k]] + j); 12405a8e39fbSHong Zhang if (rip[i] >= k){ /* only take upper triangular entry */ 12415a8e39fbSHong Zhang ajtmp[ncols_upper] = i; 1242622e2028SHong Zhang ncols_upper++; 1243622e2028SHong Zhang } 1244622e2028SHong Zhang } 12455a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 12467a48dd6fSHong Zhang nzk += nlnk; 12477a48dd6fSHong Zhang 12487a48dd6fSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 12497a48dd6fSHong Zhang prow = jl[k]; /* 1st pivot row */ 12507a48dd6fSHong Zhang 12517a48dd6fSHong Zhang while (prow < k){ 12527a48dd6fSHong Zhang nextprow = jl[prow]; 12537a48dd6fSHong Zhang 12547a48dd6fSHong Zhang /* merge prow into k-th row */ 12557a48dd6fSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 12567a48dd6fSHong Zhang jmax = ui[prow+1]; 12577a48dd6fSHong Zhang ncols = jmax-jmin; 1258ed59401aSHong Zhang i = jmin - ui[prow]; 1259ed59401aSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 12605a8e39fbSHong Zhang uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 12615a8e39fbSHong Zhang j = *(uj - 1); 12625a8e39fbSHong Zhang ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 12637a48dd6fSHong Zhang nzk += nlnk; 12647a48dd6fSHong Zhang 12657a48dd6fSHong Zhang /* update il and jl for prow */ 12667a48dd6fSHong Zhang if (jmin < jmax){ 12677a48dd6fSHong Zhang il[prow] = jmin; 12687a48dd6fSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 12697a48dd6fSHong Zhang } 12707a48dd6fSHong Zhang prow = nextprow; 12717a48dd6fSHong Zhang } 12727a48dd6fSHong Zhang 12737a48dd6fSHong Zhang /* if free space is not available, make more free space */ 12747a48dd6fSHong Zhang if (current_space->local_remaining<nzk) { 12757a48dd6fSHong Zhang i = am - k + 1; /* num of unfactored rows */ 12767a48dd6fSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1277a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1278a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 12797a48dd6fSHong Zhang reallocs++; 12807a48dd6fSHong Zhang } 12817a48dd6fSHong Zhang 12827a48dd6fSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 12832ed133dbSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 12847a48dd6fSHong Zhang 12857a48dd6fSHong Zhang /* add the k-th row into il and jl */ 128639e3d630SHong Zhang if (nzk > 1){ 12877a48dd6fSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 12887a48dd6fSHong Zhang jl[k] = jl[i]; jl[i] = k; 12897a48dd6fSHong Zhang il[k] = ui[k] + 1; 12907a48dd6fSHong Zhang } 12917a48dd6fSHong Zhang uj_ptr[k] = current_space->array; 12927a48dd6fSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 12937a48dd6fSHong Zhang 12947a48dd6fSHong Zhang current_space->array += nzk; 12957a48dd6fSHong Zhang current_space->local_used += nzk; 12967a48dd6fSHong Zhang current_space->local_remaining -= nzk; 12977a48dd6fSHong Zhang 12987a48dd6fSHong Zhang current_space_lvl->array += nzk; 12997a48dd6fSHong Zhang current_space_lvl->local_used += nzk; 13007a48dd6fSHong Zhang current_space_lvl->local_remaining -= nzk; 13017a48dd6fSHong Zhang 13027a48dd6fSHong Zhang ui[k+1] = ui[k] + nzk; 13037a48dd6fSHong Zhang } 13047a48dd6fSHong Zhang 130563ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 13067a48dd6fSHong Zhang if (ai[am] != 0) { 130739e3d630SHong Zhang PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 130863ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr); 130963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr); 131063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr); 13117a48dd6fSHong Zhang } else { 131263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr); 13137a48dd6fSHong Zhang } 131463ba0a88SBarry Smith #endif 13157a48dd6fSHong Zhang 13167a48dd6fSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 13177a48dd6fSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 13185a8e39fbSHong Zhang ierr = PetscFree(ajtmp);CHKERRQ(ierr); 13197a48dd6fSHong Zhang 13207a48dd6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 13217a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1322a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 13232ed133dbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1324a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 13257a48dd6fSHong Zhang 132639e3d630SHong Zhang } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 132739e3d630SHong Zhang 13287a48dd6fSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 1329f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr); 1330f69a0ea3SMatthew Knepley ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr); 1331ed59401aSHong Zhang B = *fact; 1332ed59401aSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 1333ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 13347a48dd6fSHong Zhang 1335ed59401aSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 13367a48dd6fSHong Zhang b->singlemalloc = PETSC_FALSE; 13377a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 13387a48dd6fSHong Zhang b->j = uj; 13397a48dd6fSHong Zhang b->i = ui; 13407a48dd6fSHong Zhang b->diag = 0; 13417a48dd6fSHong Zhang b->ilen = 0; 13427a48dd6fSHong Zhang b->imax = 0; 13437a48dd6fSHong Zhang b->row = perm; 13447a48dd6fSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 13457a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 13467a48dd6fSHong Zhang b->icol = perm; 13477a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 13487a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 134952e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 13507a48dd6fSHong Zhang b->maxnz = b->nz = ui[am]; 13517a48dd6fSHong Zhang 1352ed59401aSHong Zhang B->factor = FACTOR_CHOLESKY; 1353ed59401aSHong Zhang B->info.factor_mallocs = reallocs; 1354ed59401aSHong Zhang B->info.fill_ratio_given = fill; 13557a48dd6fSHong Zhang if (ai[am] != 0) { 1356ed59401aSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 13577a48dd6fSHong Zhang } else { 1358ed59401aSHong Zhang B->info.fill_ratio_needed = 0.0; 13597a48dd6fSHong Zhang } 136039e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1361622e2028SHong Zhang if (perm_identity){ 1362ed59401aSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1363ed59401aSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1364622e2028SHong Zhang } 1365a6175056SHong Zhang PetscFunctionReturn(0); 1366a6175056SHong Zhang } 1367d5d45c9bSBarry Smith 1368f76d2b81SHong Zhang #undef __FUNCT__ 1369f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1370dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1371f76d2b81SHong Zhang { 1372f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 137310c27e3dSHong Zhang Mat_SeqSBAIJ *b; 13742ed133dbSHong Zhang Mat B; 1375dfbe8321SBarry Smith PetscErrorCode ierr; 1376f76d2b81SHong Zhang PetscTruth perm_identity; 137710c27e3dSHong Zhang PetscReal fill = info->fill; 13781393bc97SHong Zhang PetscInt *rip,*riip,i,am=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow; 137910c27e3dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 13802ed133dbSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 1381a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 138210c27e3dSHong Zhang PetscBT lnkbt; 13832ed133dbSHong Zhang IS iperm; 1384f76d2b81SHong Zhang 1385f76d2b81SHong Zhang PetscFunctionBegin; 13862ed133dbSHong Zhang /* check whether perm is the identity mapping */ 1387f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 13882ed133dbSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 13892ed133dbSHong Zhang 1390f76d2b81SHong Zhang if (!perm_identity){ 13912ed133dbSHong Zhang /* check if perm is symmetric! */ 13922ed133dbSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 13932ed133dbSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 13941393bc97SHong Zhang for (i=0; i<am; i++) { 13952ed133dbSHong Zhang if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation"); 13962ed133dbSHong Zhang } 13972ed133dbSHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 13982ed133dbSHong Zhang ierr = ISDestroy(iperm);CHKERRQ(ierr); 1399f76d2b81SHong Zhang } 140010c27e3dSHong Zhang 140110c27e3dSHong Zhang /* initialization */ 14021393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 140310c27e3dSHong Zhang ui[0] = 0; 140410c27e3dSHong Zhang 140510c27e3dSHong Zhang /* jl: linked list for storing indices of the pivot rows 14061393bc97SHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 14071393bc97SHong Zhang ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 14081393bc97SHong Zhang il = jl + am; 14091393bc97SHong Zhang cols = il + am; 14101393bc97SHong Zhang ui_ptr = (PetscInt**)(cols + am); 14111393bc97SHong Zhang for (i=0; i<am; i++){ 14121393bc97SHong Zhang jl[i] = am; il[i] = 0; 1413f76d2b81SHong Zhang } 141410c27e3dSHong Zhang 141510c27e3dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 14161393bc97SHong Zhang nlnk = am + 1; 14171393bc97SHong Zhang ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 141810c27e3dSHong Zhang 14191393bc97SHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 1420a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 142110c27e3dSHong Zhang current_space = free_space; 142210c27e3dSHong Zhang 14231393bc97SHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 142410c27e3dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 142510c27e3dSHong Zhang nzk = 0; 14262ed133dbSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 14272ed133dbSHong Zhang ncols_upper = 0; 14282ed133dbSHong Zhang for (j=0; j<ncols; j++){ 1429622e2028SHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 14302ed133dbSHong Zhang if (i >= k){ /* only take upper triangular entry */ 14312ed133dbSHong Zhang cols[ncols_upper] = i; 14322ed133dbSHong Zhang ncols_upper++; 14332ed133dbSHong Zhang } 14342ed133dbSHong Zhang } 14351393bc97SHong Zhang ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 143610c27e3dSHong Zhang nzk += nlnk; 143710c27e3dSHong Zhang 143810c27e3dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 143910c27e3dSHong Zhang prow = jl[k]; /* 1st pivot row */ 144010c27e3dSHong Zhang 144110c27e3dSHong Zhang while (prow < k){ 144210c27e3dSHong Zhang nextprow = jl[prow]; 144310c27e3dSHong Zhang /* merge prow into k-th row */ 14441393bc97SHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 144510c27e3dSHong Zhang jmax = ui[prow+1]; 144610c27e3dSHong Zhang ncols = jmax-jmin; 14471393bc97SHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 14485a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 144910c27e3dSHong Zhang nzk += nlnk; 145010c27e3dSHong Zhang 145110c27e3dSHong Zhang /* update il and jl for prow */ 145210c27e3dSHong Zhang if (jmin < jmax){ 145310c27e3dSHong Zhang il[prow] = jmin; 14542ed133dbSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 145510c27e3dSHong Zhang } 145610c27e3dSHong Zhang prow = nextprow; 145710c27e3dSHong Zhang } 145810c27e3dSHong Zhang 145910c27e3dSHong Zhang /* if free space is not available, make more free space */ 146010c27e3dSHong Zhang if (current_space->local_remaining<nzk) { 14611393bc97SHong Zhang i = am - k + 1; /* num of unfactored rows */ 146210c27e3dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1463a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 146410c27e3dSHong Zhang reallocs++; 146510c27e3dSHong Zhang } 146610c27e3dSHong Zhang 146710c27e3dSHong Zhang /* copy data into free space, then initialize lnk */ 14681393bc97SHong Zhang ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 146910c27e3dSHong Zhang 147010c27e3dSHong Zhang /* add the k-th row into il and jl */ 147110c27e3dSHong Zhang if (nzk-1 > 0){ 14721393bc97SHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 147310c27e3dSHong Zhang jl[k] = jl[i]; jl[i] = k; 147410c27e3dSHong Zhang il[k] = ui[k] + 1; 147510c27e3dSHong Zhang } 14762ed133dbSHong Zhang ui_ptr[k] = current_space->array; 147710c27e3dSHong Zhang current_space->array += nzk; 147810c27e3dSHong Zhang current_space->local_used += nzk; 147910c27e3dSHong Zhang current_space->local_remaining -= nzk; 148010c27e3dSHong Zhang 148110c27e3dSHong Zhang ui[k+1] = ui[k] + nzk; 148210c27e3dSHong Zhang } 148310c27e3dSHong Zhang 148463ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 14851393bc97SHong Zhang if (ai[am] != 0) { 14861393bc97SHong Zhang PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]); 148763ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr); 148863ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr); 148963ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr); 149010c27e3dSHong Zhang } else { 149163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr); 149210c27e3dSHong Zhang } 149363ba0a88SBarry Smith #endif 149410c27e3dSHong Zhang 149510c27e3dSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 149610c27e3dSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 149710c27e3dSHong Zhang 149810c27e3dSHong Zhang /* destroy list of free space and other temporary array(s) */ 14991393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1500a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 150110c27e3dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 150210c27e3dSHong Zhang 150310c27e3dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 1504f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr); 1505f69a0ea3SMatthew Knepley ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr); 15062ed133dbSHong Zhang B = *fact; 15072ed133dbSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 1508ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 150910c27e3dSHong Zhang 15102ed133dbSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 151110c27e3dSHong Zhang b->singlemalloc = PETSC_FALSE; 15121393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 151310c27e3dSHong Zhang b->j = uj; 151410c27e3dSHong Zhang b->i = ui; 151510c27e3dSHong Zhang b->diag = 0; 151610c27e3dSHong Zhang b->ilen = 0; 151710c27e3dSHong Zhang b->imax = 0; 151810c27e3dSHong Zhang b->row = perm; 151910c27e3dSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 152010c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 152110c27e3dSHong Zhang b->icol = perm; 152210c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 15231393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 15241393bc97SHong Zhang ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 15251393bc97SHong Zhang b->maxnz = b->nz = ui[am]; 152610c27e3dSHong Zhang 15272ed133dbSHong Zhang B->factor = FACTOR_CHOLESKY; 15282ed133dbSHong Zhang B->info.factor_mallocs = reallocs; 15292ed133dbSHong Zhang B->info.fill_ratio_given = fill; 15301393bc97SHong Zhang if (ai[am] != 0) { 15311393bc97SHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 153210c27e3dSHong Zhang } else { 15332ed133dbSHong Zhang B->info.fill_ratio_needed = 0.0; 153410c27e3dSHong Zhang } 153539e3d630SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 15362ed133dbSHong Zhang if (perm_identity){ 153710c27e3dSHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 153810c27e3dSHong Zhang (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 15392ed133dbSHong Zhang } 1540f76d2b81SHong Zhang PetscFunctionReturn(0); 1541f76d2b81SHong Zhang } 1542