173f4d377SMatthew Knepley /*$Id: sbaijfact.c,v 1.61 2001/08/06 21:15:47 bsmith Exp $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h" 43a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h" 549b5e25fSSatish Balay #include "src/vec/vecimpl.h" 649b5e25fSSatish Balay #include "src/inline/ilu.h" 749a6740bSHong Zhang #include "include/petscis.h" 849b5e25fSSatish Balay 98dc9108eSHong Zhang #if !defined(PETSC_USE_COMPLEX) 105f9f512dSHong Zhang /* 115f9f512dSHong Zhang input: 12*c037c3f7SHong Zhang F -- numeric factor 135f9f512dSHong Zhang output: 14*c037c3f7SHong Zhang nneg, nzero, npos: matrix inertia 155f9f512dSHong Zhang */ 165f9f512dSHong Zhang 175f9f512dSHong Zhang #undef __FUNCT__ 185f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ" 19*c037c3f7SHong Zhang int MatGetInertia_SeqSBAIJ(Mat F,int *nneig,int *nzero,int *npos) 205f9f512dSHong Zhang { 215f9f512dSHong Zhang PetscScalar *dd; 225f9f512dSHong Zhang Mat_SeqSBAIJ *fact_ptr; 23*c037c3f7SHong Zhang int m = F->m,i; 245f9f512dSHong Zhang 255f9f512dSHong Zhang PetscFunctionBegin; 26*c037c3f7SHong Zhang if (!F->assembled) { 27*c037c3f7SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Numeric factor mat is not assembled"); 28*c037c3f7SHong Zhang } 29*c037c3f7SHong Zhang fact_ptr = (Mat_SeqSBAIJ*)F->data; 305f9f512dSHong Zhang dd = fact_ptr->a; 315f9f512dSHong Zhang *nneig = 0; 325f9f512dSHong Zhang for (i=0; i<m; i++){ 33b36a9721SBarry Smith if (PetscRealPart(dd[i]) < 0.0) (*nneig)++; 345f9f512dSHong Zhang } 355f9f512dSHong Zhang 365f9f512dSHong Zhang PetscFunctionReturn(0); 375f9f512dSHong Zhang } 388dc9108eSHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 395f9f512dSHong Zhang 405f9f512dSHong Zhang /* Using Modified Sparse Row (MSR) storage. 415f9f512dSHong Zhang See page 85, "Iterative Methods ..." by Saad. */ 425f9f512dSHong Zhang /* 435f9f512dSHong Zhang Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. 445f9f512dSHong Zhang */ 4587fe1012SHong Zhang /* Use Modified Sparse Row storage for u and ju, see Sasd pp.85 */ 464a2ae208SSatish Balay #undef __FUNCT__ 474a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" 4877f73120SHong Zhang int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,PetscReal f,Mat *B) 4949b5e25fSSatish Balay { 5049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; 51cb718733SHong Zhang int *rip,ierr,i,mbs = a->mbs,*ai,*aj; 52066653e3SSatish Balay int *jutmp,bs = a->bs,bs2=a->bs2; 5387fe1012SHong Zhang int m,realloc = 0,prow; 54d60e5362SHong Zhang int *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd; 55671cb588SHong Zhang PetscTruth perm_identity; 5649b5e25fSSatish Balay 5749b5e25fSSatish Balay PetscFunctionBegin; 582d007305SHong Zhang 59cb718733SHong Zhang /* check whether perm is the identity mapping */ 60671cb588SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 61671cb588SHong Zhang if (!perm_identity) a->permute = PETSC_TRUE; 622d007305SHong Zhang 63671cb588SHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 64671cb588SHong Zhang 65671cb588SHong Zhang if (perm_identity){ /* without permutation */ 662d007305SHong Zhang ai = a->i; aj = a->j; 672d007305SHong Zhang } else { /* non-trivial permutation */ 68323b833fSBarry Smith ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 692d007305SHong Zhang ai = a->inew; aj = a->jnew; 702d007305SHong Zhang } 712d007305SHong Zhang 7249b5e25fSSatish Balay /* initialization */ 73b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr); 7449b5e25fSSatish Balay umax = (int)(f*ai[mbs] + 1); umax += mbs + 1; 7582502324SSatish Balay ierr = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr); 7649b5e25fSSatish Balay iu[0] = mbs+1; 7787fe1012SHong Zhang juidx = mbs + 1; /* index for ju */ 7887fe1012SHong Zhang ierr = PetscMalloc(2*mbs*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for pivot row */ 7987fe1012SHong Zhang q = jl + mbs; /* linked list for col index */ 8049b5e25fSSatish Balay for (i=0; i<mbs; i++){ 8187fe1012SHong Zhang jl[i] = mbs; 8287fe1012SHong Zhang q[i] = 0; 8349b5e25fSSatish Balay } 8449b5e25fSSatish Balay 8549b5e25fSSatish Balay /* for each row k */ 8649b5e25fSSatish Balay for (k=0; k<mbs; k++){ 8749b5e25fSSatish Balay nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */ 8849b5e25fSSatish Balay q[k] = mbs; 8949b5e25fSSatish Balay /* initialize nonzero structure of k-th row to row rip[k] of A */ 9049b5e25fSSatish Balay jmin = ai[rip[k]]; 9149b5e25fSSatish Balay jmax = ai[rip[k]+1]; 9249b5e25fSSatish Balay for (j=jmin; j<jmax; j++){ 93cb718733SHong Zhang vj = rip[aj[j]]; /* col. value */ 9449b5e25fSSatish Balay if(vj > k){ 9549b5e25fSSatish Balay qm = k; 9649b5e25fSSatish Balay do { 9749b5e25fSSatish Balay m = qm; qm = q[m]; 9849b5e25fSSatish Balay } while(qm < vj); 9949b5e25fSSatish Balay if (qm == vj) { 1003078e140SBarry Smith SETERRQ(1," error: duplicate entry in A\n"); 10149b5e25fSSatish Balay } 10249b5e25fSSatish Balay nzk++; 10349b5e25fSSatish Balay q[m] = vj; 10449b5e25fSSatish Balay q[vj] = qm; 10549b5e25fSSatish Balay } /* if(vj > k) */ 10649b5e25fSSatish Balay } /* for (j=jmin; j<jmax; j++) */ 10749b5e25fSSatish Balay 10849b5e25fSSatish Balay /* modify nonzero structure of k-th row by computing fill-in 10949b5e25fSSatish Balay for each row i to be merged in */ 11087fe1012SHong Zhang prow = k; 11187fe1012SHong Zhang prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */ 11287fe1012SHong Zhang 11387fe1012SHong Zhang while (prow < k){ 11487fe1012SHong Zhang /* merge row prow into k-th row */ 11587fe1012SHong Zhang jmin = iu[prow] + 1; jmax = iu[prow+1]; 11649b5e25fSSatish Balay qm = k; 11787fe1012SHong Zhang for (j=jmin; j<jmax; j++){ 11849b5e25fSSatish Balay vj = ju[j]; 11949b5e25fSSatish Balay do { 12049b5e25fSSatish Balay m = qm; qm = q[m]; 12149b5e25fSSatish Balay } while (qm < vj); 12249b5e25fSSatish Balay if (qm != vj){ 12349b5e25fSSatish Balay nzk++; q[m] = vj; q[vj] = qm; qm = vj; 12449b5e25fSSatish Balay } 12549b5e25fSSatish Balay } 12687fe1012SHong Zhang prow = jl[prow]; /* next pivot row */ 12749b5e25fSSatish Balay } 12849b5e25fSSatish Balay 12949b5e25fSSatish Balay /* add k to row list for first nonzero element in k-th row */ 13049b5e25fSSatish Balay if (nzk > 0){ 13149b5e25fSSatish Balay i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */ 13249b5e25fSSatish Balay jl[k] = jl[i]; jl[i] = k; 13349b5e25fSSatish Balay } 13487fe1012SHong Zhang iu[k+1] = iu[k] + nzk; 13549b5e25fSSatish Balay 13649b5e25fSSatish Balay /* allocate more space to ju if needed */ 1373078e140SBarry Smith if (iu[k+1] > umax) { 13849b5e25fSSatish Balay /* estimate how much additional space we will need */ 13949b5e25fSSatish Balay /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 14049b5e25fSSatish Balay /* just double the memory each time */ 14149b5e25fSSatish Balay maxadd = umax; 14249b5e25fSSatish Balay if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; 14349b5e25fSSatish Balay umax += maxadd; 14449b5e25fSSatish Balay 1459f9cb213SHong Zhang /* allocate a longer ju */ 14682502324SSatish Balay ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr); 14749b5e25fSSatish Balay ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr); 14849b5e25fSSatish Balay ierr = PetscFree(ju);CHKERRQ(ierr); 1499f9cb213SHong Zhang ju = jutmp; 15049b5e25fSSatish Balay realloc++; /* count how many times we realloc */ 15149b5e25fSSatish Balay } 15249b5e25fSSatish Balay 15349b5e25fSSatish Balay /* save nonzero structure of k-th row in ju */ 15449b5e25fSSatish Balay i=k; 15587fe1012SHong Zhang while (nzk --) { 15649b5e25fSSatish Balay i = q[i]; 15787fe1012SHong Zhang ju[juidx++] = i; 15849b5e25fSSatish Balay } 15977f73120SHong Zhang } 16049b5e25fSSatish Balay 16149b5e25fSSatish Balay if (ai[mbs] != 0) { 16249b5e25fSSatish Balay PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 163b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 1643078e140SBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 1653078e140SBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af); 166b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n"); 16749b5e25fSSatish Balay } else { 168b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"); 16949b5e25fSSatish Balay } 17049b5e25fSSatish Balay 17177f73120SHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 17287fe1012SHong Zhang /* ierr = PetscFree(q);CHKERRQ(ierr); */ 17349b5e25fSSatish Balay ierr = PetscFree(jl);CHKERRQ(ierr); 17449b5e25fSSatish Balay 17549b5e25fSSatish Balay /* put together the new matrix */ 17649b5e25fSSatish Balay ierr = MatCreateSeqSBAIJ(A->comm,bs,bs*mbs,bs*mbs,0,PETSC_NULL,B);CHKERRQ(ierr); 177b0a32e0cSBarry Smith /* PetscLogObjectParent(*B,iperm); */ 17849b5e25fSSatish Balay b = (Mat_SeqSBAIJ*)(*B)->data; 17949b5e25fSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 18049b5e25fSSatish Balay b->singlemalloc = PETSC_FALSE; 18149b5e25fSSatish Balay /* the next line frees the default space generated by the Create() */ 18249b5e25fSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 18349b5e25fSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 18482502324SSatish Balay ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 18549b5e25fSSatish Balay b->j = ju; 18649b5e25fSSatish Balay b->i = iu; 18749b5e25fSSatish Balay b->diag = 0; 18849b5e25fSSatish Balay b->ilen = 0; 18949b5e25fSSatish Balay b->imax = 0; 19077f73120SHong Zhang b->row = perm; 191bcd9e38bSBarry Smith b->pivotinblocks = PETSC_FALSE; /* need to get from MatCholeskyInfo */ 19277f73120SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 193cb718733SHong Zhang b->icol = perm; 194cb718733SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 19587828ca2SBarry Smith ierr = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 19649b5e25fSSatish Balay /* In b structure: Free imax, ilen, old a, old j. 19749b5e25fSSatish Balay Allocate idnew, solve_work, new a, new j */ 198b0a32e0cSBarry Smith PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar))); 19949b5e25fSSatish Balay b->s_maxnz = b->s_nz = iu[mbs]; 20049b5e25fSSatish Balay 2015c0bcdfcSHong Zhang (*B)->factor = FACTOR_CHOLESKY; 20249b5e25fSSatish Balay (*B)->info.factor_mallocs = realloc; 20349b5e25fSSatish Balay (*B)->info.fill_ratio_given = f; 20449b5e25fSSatish Balay if (ai[mbs] != 0) { 20549b5e25fSSatish Balay (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 20649b5e25fSSatish Balay } else { 20749b5e25fSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 20849b5e25fSSatish Balay } 2090aa0ce06SHong Zhang 210671cb588SHong Zhang if (perm_identity){ 211671cb588SHong Zhang switch (bs) { 212671cb588SHong Zhang case 1: 213671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 214671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 2154d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n"); 216671cb588SHong Zhang break; 217671cb588SHong Zhang case 2: 218671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 219671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 2204d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 221671cb588SHong Zhang break; 222671cb588SHong Zhang case 3: 223671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 224671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 2254d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n"); 226671cb588SHong Zhang break; 227671cb588SHong Zhang case 4: 228671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 229671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 2304d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 231671cb588SHong Zhang break; 232671cb588SHong Zhang case 5: 233671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 234671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 2354d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 236671cb588SHong Zhang break; 237671cb588SHong Zhang case 6: 238671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 239671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 2404d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 241671cb588SHong Zhang break; 242671cb588SHong Zhang case 7: 243671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 244671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 2454d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 246671cb588SHong Zhang break; 247671cb588SHong Zhang default: 248671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 249671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 2504d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n"); 251671cb588SHong Zhang break; 252671cb588SHong Zhang } 253671cb588SHong Zhang } 254671cb588SHong Zhang 25549b5e25fSSatish Balay PetscFunctionReturn(0); 25649b5e25fSSatish Balay } 25749b5e25fSSatish Balay 25835d8aa7fSBarry Smith 2594a2ae208SSatish Balay #undef __FUNCT__ 2604a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" 2616f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B) 26249b5e25fSSatish Balay { 26349b5e25fSSatish Balay Mat C = *B; 2644c16a6a6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 2654c16a6a6SHong Zhang IS perm = b->row; 2664c16a6a6SHong Zhang int *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 2674c16a6a6SHong Zhang int *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 2684c16a6a6SHong Zhang int bs=a->bs,bs2 = a->bs2; 2694c16a6a6SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 2704c16a6a6SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 27128de702eSHong Zhang MatScalar *work; 2724c16a6a6SHong Zhang int *pivots; 2734c16a6a6SHong Zhang 2744c16a6a6SHong Zhang PetscFunctionBegin; 2759706f043SHong Zhang 2764c16a6a6SHong Zhang /* initialization */ 27782502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 2784c16a6a6SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 27982502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 28028de702eSHong Zhang jl = il + mbs; 2814c16a6a6SHong Zhang for (i=0; i<mbs; i++) { 2824c16a6a6SHong Zhang jl[i] = mbs; il[0] = 0; 2834c16a6a6SHong Zhang } 284b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 28528de702eSHong Zhang uik = dk + bs2; 28628de702eSHong Zhang work = uik + bs2; 28782502324SSatish Balay ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr); 2884c16a6a6SHong Zhang 2894c16a6a6SHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 2904c16a6a6SHong Zhang 2914c16a6a6SHong Zhang /* check permutation */ 2924c16a6a6SHong Zhang if (!a->permute){ 2934c16a6a6SHong Zhang ai = a->i; aj = a->j; aa = a->a; 2944c16a6a6SHong Zhang } else { 2954c16a6a6SHong Zhang ai = a->inew; aj = a->jnew; 29682502324SSatish Balay ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 2974c16a6a6SHong Zhang ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 29882502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr); 2994c16a6a6SHong Zhang ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 3004c16a6a6SHong Zhang 3014c16a6a6SHong Zhang for (i=0; i<mbs; i++){ 3024c16a6a6SHong Zhang jmin = ai[i]; jmax = ai[i+1]; 3034c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 3044c16a6a6SHong Zhang while (a2anew[j] != j){ 3054c16a6a6SHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 3064c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 3074c16a6a6SHong Zhang dk[k1] = aa[k*bs2+k1]; 3084c16a6a6SHong Zhang aa[k*bs2+k1] = aa[j*bs2+k1]; 3094c16a6a6SHong Zhang aa[j*bs2+k1] = dk[k1]; 3104c16a6a6SHong Zhang } 3114c16a6a6SHong Zhang } 3124c16a6a6SHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 3134c16a6a6SHong Zhang if (i > aj[j]){ 3144c16a6a6SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 3154c16a6a6SHong Zhang ap = aa + j*bs2; /* ptr to the beginning of j-th block of aa */ 3164c16a6a6SHong Zhang for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 3174c16a6a6SHong Zhang for (k=0; k<bs; k++){ /* j-th block of aa <- dk^T */ 3184c16a6a6SHong Zhang for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1]; 3194c16a6a6SHong Zhang } 3204c16a6a6SHong Zhang } 3214c16a6a6SHong Zhang } 3224c16a6a6SHong Zhang } 323323b833fSBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 3244c16a6a6SHong Zhang } 3254c16a6a6SHong Zhang 3264c16a6a6SHong Zhang /* for each row k */ 3274c16a6a6SHong Zhang for (k = 0; k<mbs; k++){ 3284c16a6a6SHong Zhang 3294c16a6a6SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 3304c16a6a6SHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 331057f5ba7SHong Zhang 3324c16a6a6SHong Zhang ap = aa + jmin*bs2; 3334c16a6a6SHong Zhang for (j = jmin; j < jmax; j++){ 3344c16a6a6SHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 3354c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 3364c16a6a6SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 3374c16a6a6SHong Zhang } 3384c16a6a6SHong Zhang 3394c16a6a6SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 3404c16a6a6SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 3414c16a6a6SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 3424c16a6a6SHong Zhang 343057f5ba7SHong Zhang while (i < k){ 3444c16a6a6SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 3454c16a6a6SHong Zhang 3464c16a6a6SHong Zhang /* compute multiplier */ 3474c16a6a6SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 3484c16a6a6SHong Zhang 3494c16a6a6SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 3504c16a6a6SHong Zhang diag = ba + i*bs2; 3514c16a6a6SHong Zhang u = ba + ili*bs2; 3524c16a6a6SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 3534c16a6a6SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 3544c16a6a6SHong Zhang 3554c16a6a6SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 3564c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 3574c16a6a6SHong Zhang 3584c16a6a6SHong Zhang /* update -U(i,k) */ 3594c16a6a6SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 3604c16a6a6SHong Zhang 3614c16a6a6SHong Zhang /* add multiple of row i to k-th row ... */ 3624c16a6a6SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 3634c16a6a6SHong Zhang if (jmin < jmax){ 3644c16a6a6SHong Zhang for (j=jmin; j<jmax; j++) { 3654c16a6a6SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 3664c16a6a6SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 3674c16a6a6SHong Zhang u = ba + j*bs2; 3684c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 3694c16a6a6SHong Zhang } 3704c16a6a6SHong Zhang 3714c16a6a6SHong Zhang /* ... add i to row list for next nonzero entry */ 3724c16a6a6SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 3734c16a6a6SHong Zhang j = bj[jmin]; 3744c16a6a6SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 3754c16a6a6SHong Zhang } 3764c16a6a6SHong Zhang i = nexti; 3774c16a6a6SHong Zhang } 3784c16a6a6SHong Zhang 3794c16a6a6SHong Zhang /* save nonzero entries in k-th row of U ... */ 3804c16a6a6SHong Zhang 3814c16a6a6SHong Zhang /* invert diagonal block */ 3824c16a6a6SHong Zhang diag = ba+k*bs2; 3834c16a6a6SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 3844c16a6a6SHong Zhang Kernel_A_gets_inverse_A(bs,diag,pivots,work); 3854c16a6a6SHong Zhang 3864c16a6a6SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 3874c16a6a6SHong Zhang if (jmin < jmax) { 3884c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 3894c16a6a6SHong Zhang vj = bj[j]; /* block col. index of U */ 3904c16a6a6SHong Zhang u = ba + j*bs2; 3914c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 3924c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 3934c16a6a6SHong Zhang *u++ = *rtmp_ptr; 3944c16a6a6SHong Zhang *rtmp_ptr++ = 0.0; 3954c16a6a6SHong Zhang } 3964c16a6a6SHong Zhang } 3974c16a6a6SHong Zhang 3984c16a6a6SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 3994c16a6a6SHong Zhang il[k] = jmin; 4004c16a6a6SHong Zhang i = bj[jmin]; 4014c16a6a6SHong Zhang jl[k] = jl[i]; jl[i] = k; 4024c16a6a6SHong Zhang } 4034c16a6a6SHong Zhang } 4044c16a6a6SHong Zhang 4054c16a6a6SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 4064c16a6a6SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 4074c16a6a6SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 4084c16a6a6SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 4094c16a6a6SHong Zhang if (a->permute){ 4104c16a6a6SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 4114c16a6a6SHong Zhang } 4124c16a6a6SHong Zhang 4134c16a6a6SHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 4144c16a6a6SHong Zhang C->factor = FACTOR_CHOLESKY; 4154c16a6a6SHong Zhang C->assembled = PETSC_TRUE; 4164c16a6a6SHong Zhang C->preallocated = PETSC_TRUE; 417b0a32e0cSBarry Smith PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ 4184c16a6a6SHong Zhang PetscFunctionReturn(0); 4194c16a6a6SHong Zhang } 420d4edadd4SHong Zhang 4214a2ae208SSatish Balay #undef __FUNCT__ 4224a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" 423671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B) 424671cb588SHong Zhang { 425671cb588SHong Zhang Mat C = *B; 426671cb588SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 427671cb588SHong Zhang int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 4281ad1de41SHong Zhang int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 429671cb588SHong Zhang int bs=a->bs,bs2 = a->bs2; 430671cb588SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 431671cb588SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 43228de702eSHong Zhang MatScalar *work; 433671cb588SHong Zhang int *pivots; 434671cb588SHong Zhang 435671cb588SHong Zhang PetscFunctionBegin; 436671cb588SHong Zhang 437671cb588SHong Zhang /* initialization */ 438671cb588SHong Zhang 43982502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 440671cb588SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 44182502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 44228de702eSHong Zhang jl = il + mbs; 443671cb588SHong Zhang for (i=0; i<mbs; i++) { 444671cb588SHong Zhang jl[i] = mbs; il[0] = 0; 445671cb588SHong Zhang } 446b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 44728de702eSHong Zhang uik = dk + bs2; 44828de702eSHong Zhang work = uik + bs2; 44982502324SSatish Balay ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr); 450671cb588SHong Zhang 451671cb588SHong Zhang ai = a->i; aj = a->j; aa = a->a; 452671cb588SHong Zhang 453671cb588SHong Zhang /* for each row k */ 454671cb588SHong Zhang for (k = 0; k<mbs; k++){ 455671cb588SHong Zhang 456671cb588SHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 457671cb588SHong Zhang jmin = ai[k]; jmax = ai[k+1]; 458671cb588SHong Zhang ap = aa + jmin*bs2; 459671cb588SHong Zhang for (j = jmin; j < jmax; j++){ 460671cb588SHong Zhang vj = aj[j]; /* block col. index */ 461671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 462671cb588SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 463671cb588SHong Zhang } 464671cb588SHong Zhang 465671cb588SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 466671cb588SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 467671cb588SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 468671cb588SHong Zhang 469057f5ba7SHong Zhang while (i < k){ 470671cb588SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 471671cb588SHong Zhang 472671cb588SHong Zhang /* compute multiplier */ 473671cb588SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 474671cb588SHong Zhang 475671cb588SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 476671cb588SHong Zhang diag = ba + i*bs2; 477671cb588SHong Zhang u = ba + ili*bs2; 478671cb588SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 479671cb588SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 480671cb588SHong Zhang 481671cb588SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 482671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 483671cb588SHong Zhang 484671cb588SHong Zhang /* update -U(i,k) */ 485671cb588SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 486671cb588SHong Zhang 487671cb588SHong Zhang /* add multiple of row i to k-th row ... */ 488671cb588SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 489671cb588SHong Zhang if (jmin < jmax){ 490671cb588SHong Zhang for (j=jmin; j<jmax; j++) { 491671cb588SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 492671cb588SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 493671cb588SHong Zhang u = ba + j*bs2; 494671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 495671cb588SHong Zhang } 496671cb588SHong Zhang 497671cb588SHong Zhang /* ... add i to row list for next nonzero entry */ 498671cb588SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 499671cb588SHong Zhang j = bj[jmin]; 500671cb588SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 501671cb588SHong Zhang } 502671cb588SHong Zhang i = nexti; 503671cb588SHong Zhang } 504671cb588SHong Zhang 505671cb588SHong Zhang /* save nonzero entries in k-th row of U ... */ 506671cb588SHong Zhang 507671cb588SHong Zhang /* invert diagonal block */ 508671cb588SHong Zhang diag = ba+k*bs2; 509671cb588SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 510671cb588SHong Zhang Kernel_A_gets_inverse_A(bs,diag,pivots,work); 511671cb588SHong Zhang 512671cb588SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 513671cb588SHong Zhang if (jmin < jmax) { 514671cb588SHong Zhang for (j=jmin; j<jmax; j++){ 515671cb588SHong Zhang vj = bj[j]; /* block col. index of U */ 516671cb588SHong Zhang u = ba + j*bs2; 517671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 518671cb588SHong Zhang for (k1=0; k1<bs2; k1++){ 519671cb588SHong Zhang *u++ = *rtmp_ptr; 520671cb588SHong Zhang *rtmp_ptr++ = 0.0; 521671cb588SHong Zhang } 522671cb588SHong Zhang } 523671cb588SHong Zhang 524671cb588SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 525671cb588SHong Zhang il[k] = jmin; 526671cb588SHong Zhang i = bj[jmin]; 527671cb588SHong Zhang jl[k] = jl[i]; jl[i] = k; 528671cb588SHong Zhang } 529671cb588SHong Zhang } 530671cb588SHong Zhang 531671cb588SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 532671cb588SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 533671cb588SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 534671cb588SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 535671cb588SHong Zhang 536671cb588SHong Zhang C->factor = FACTOR_CHOLESKY; 537671cb588SHong Zhang C->assembled = PETSC_TRUE; 538671cb588SHong Zhang C->preallocated = PETSC_TRUE; 539b0a32e0cSBarry Smith PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ 540671cb588SHong Zhang PetscFunctionReturn(0); 541671cb588SHong Zhang } 542671cb588SHong Zhang 54349b5e25fSSatish Balay /* 544fcf159c0SHong Zhang Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 545cc0c071aSHong Zhang Version for blocks 2 by 2. 54649b5e25fSSatish Balay */ 5474a2ae208SSatish Balay #undef __FUNCT__ 5484a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2" 5496f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B) 55049b5e25fSSatish Balay { 55149b5e25fSSatish Balay Mat C = *B; 55291602c66SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 553cc0c071aSHong Zhang IS perm = b->row; 554cc0c071aSHong Zhang int *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 555cc0c071aSHong Zhang int *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 556a1723e09SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 557cc0c071aSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 55849b5e25fSSatish Balay 55949b5e25fSSatish Balay PetscFunctionBegin; 560671cb588SHong Zhang 56191602c66SHong Zhang /* initialization */ 56291602c66SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 56391602c66SHong Zhang window U(0:k, k:mbs-1). 56491602c66SHong Zhang jl: list of rows to be added to uneliminated rows 56591602c66SHong Zhang i>= k: jl(i) is the first row to be added to row i 56691602c66SHong Zhang i< k: jl(i) is the row following row i in some list of rows 56791602c66SHong Zhang jl(i) = mbs indicates the end of a list 56891602c66SHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 56991602c66SHong Zhang row i of U */ 57082502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 571cc0c071aSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 57282502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 57328de702eSHong Zhang jl = il + mbs; 57491602c66SHong Zhang for (i=0; i<mbs; i++) { 5753845f261SHong Zhang jl[i] = mbs; il[0] = 0; 57691602c66SHong Zhang } 57782502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 57828de702eSHong Zhang uik = dk + 4; 579cc0c071aSHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 580cc0c071aSHong Zhang 581cc0c071aSHong Zhang /* check permutation */ 582cc0c071aSHong Zhang if (!a->permute){ 583cc0c071aSHong Zhang ai = a->i; aj = a->j; aa = a->a; 584cc0c071aSHong Zhang } else { 585cc0c071aSHong Zhang ai = a->inew; aj = a->jnew; 58682502324SSatish Balay ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 587cc0c071aSHong Zhang ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 58882502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr); 589cc0c071aSHong Zhang ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 590cc0c071aSHong Zhang 591cc0c071aSHong Zhang for (i=0; i<mbs; i++){ 592cc0c071aSHong Zhang jmin = ai[i]; jmax = ai[i+1]; 593cc0c071aSHong Zhang for (j=jmin; j<jmax; j++){ 594cc0c071aSHong Zhang while (a2anew[j] != j){ 595cc0c071aSHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 596cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 597cc0c071aSHong Zhang dk[k1] = aa[k*4+k1]; 598cc0c071aSHong Zhang aa[k*4+k1] = aa[j*4+k1]; 599cc0c071aSHong Zhang aa[j*4+k1] = dk[k1]; 600cc0c071aSHong Zhang } 601cc0c071aSHong Zhang } 602cc0c071aSHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 603cc0c071aSHong Zhang if (i > aj[j]){ 604a1723e09SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 605cc0c071aSHong Zhang ap = aa + j*4; /* ptr to the beginning of the block */ 606cc0c071aSHong Zhang dk[1] = ap[1]; /* swap ap[1] and ap[2] */ 607cc0c071aSHong Zhang ap[1] = ap[2]; 608cc0c071aSHong Zhang ap[2] = dk[1]; 609cc0c071aSHong Zhang } 610cc0c071aSHong Zhang } 611cc0c071aSHong Zhang } 612ac355199SBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 613cc0c071aSHong Zhang } 6143845f261SHong Zhang 61591602c66SHong Zhang /* for each row k */ 61691602c66SHong Zhang for (k = 0; k<mbs; k++){ 61791602c66SHong Zhang 61891602c66SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 619cc0c071aSHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 620cc0c071aSHong Zhang ap = aa + jmin*4; 62191602c66SHong Zhang for (j = jmin; j < jmax; j++){ 622cc0c071aSHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 623cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 624cc0c071aSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 62591602c66SHong Zhang } 62691602c66SHong Zhang 62791602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 628cc0c071aSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 62991602c66SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 63091602c66SHong Zhang 631057f5ba7SHong Zhang while (i < k){ 63291602c66SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 63391602c66SHong Zhang 6343845f261SHong Zhang /* compute multiplier */ 63591602c66SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 6363845f261SHong Zhang 6373845f261SHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 638cc0c071aSHong Zhang diag = ba + i*4; 639cc0c071aSHong Zhang u = ba + ili*4; 640cc0c071aSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 641cc0c071aSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 642cc0c071aSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 643cc0c071aSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 6443845f261SHong Zhang 6453845f261SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 646cc0c071aSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 647cc0c071aSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 648cc0c071aSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 649cc0c071aSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 6503845f261SHong Zhang 6513845f261SHong Zhang /* update -U(i,k): ba[ili] = uik */ 652cc0c071aSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 65391602c66SHong Zhang 65491602c66SHong Zhang /* add multiple of row i to k-th row ... */ 65591602c66SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 65691602c66SHong Zhang if (jmin < jmax){ 6573845f261SHong Zhang for (j=jmin; j<jmax; j++) { 6583845f261SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 659cc0c071aSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 660cc0c071aSHong Zhang u = ba + j*4; 661cc0c071aSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 662cc0c071aSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 663cc0c071aSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 664cc0c071aSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 6653845f261SHong Zhang } 6663845f261SHong Zhang 66791602c66SHong Zhang /* ... add i to row list for next nonzero entry */ 66891602c66SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 66991602c66SHong Zhang j = bj[jmin]; 67091602c66SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 67191602c66SHong Zhang } 672a1723e09SHong Zhang i = nexti; 67391602c66SHong Zhang } 674cc0c071aSHong Zhang 67591602c66SHong Zhang /* save nonzero entries in k-th row of U ... */ 6763845f261SHong Zhang 677cc0c071aSHong Zhang /* invert diagonal block */ 678cc0c071aSHong Zhang diag = ba+k*4; 679cc0c071aSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 6803845f261SHong Zhang ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 6813845f261SHong Zhang 68291602c66SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 68391602c66SHong Zhang if (jmin < jmax) { 68491602c66SHong Zhang for (j=jmin; j<jmax; j++){ 685cc0c071aSHong Zhang vj = bj[j]; /* block col. index of U */ 686cc0c071aSHong Zhang u = ba + j*4; 687cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 688cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 689cc0c071aSHong Zhang *u++ = *rtmp_ptr; 690cc0c071aSHong Zhang *rtmp_ptr++ = 0.0; 691cc0c071aSHong Zhang } 692cc0c071aSHong Zhang } 6933845f261SHong Zhang 69491602c66SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 69591602c66SHong Zhang il[k] = jmin; 69691602c66SHong Zhang i = bj[jmin]; 69791602c66SHong Zhang jl[k] = jl[i]; jl[i] = k; 69891602c66SHong Zhang } 69991602c66SHong Zhang } 7003845f261SHong Zhang 70149b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 70291602c66SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 7033845f261SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 70491602c66SHong Zhang if (a->permute) { 70591602c66SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 70691602c66SHong Zhang } 707cc0c071aSHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 70891602c66SHong Zhang C->factor = FACTOR_CHOLESKY; 70949b5e25fSSatish Balay C->assembled = PETSC_TRUE; 7105c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 711b0a32e0cSBarry Smith PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 71249b5e25fSSatish Balay PetscFunctionReturn(0); 71349b5e25fSSatish Balay } 71491602c66SHong Zhang 71549b5e25fSSatish Balay /* 71649b5e25fSSatish Balay Version for when blocks are 2 by 2 Using natural ordering 71749b5e25fSSatish Balay */ 7184a2ae208SSatish Balay #undef __FUNCT__ 7194a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" 7206f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B) 72149b5e25fSSatish Balay { 72249b5e25fSSatish Balay Mat C = *B; 723ab27746eSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 724ab27746eSHong Zhang int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 72571149678SHong Zhang int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 726ab27746eSHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 727ab27746eSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 72849b5e25fSSatish Balay 72949b5e25fSSatish Balay PetscFunctionBegin; 730671cb588SHong Zhang 731ab27746eSHong Zhang /* initialization */ 732ab27746eSHong Zhang /* il and jl record the first nonzero element in each row of the accessing 733ab27746eSHong Zhang window U(0:k, k:mbs-1). 734ab27746eSHong Zhang jl: list of rows to be added to uneliminated rows 735ab27746eSHong Zhang i>= k: jl(i) is the first row to be added to row i 736ab27746eSHong Zhang i< k: jl(i) is the row following row i in some list of rows 737ab27746eSHong Zhang jl(i) = mbs indicates the end of a list 738ab27746eSHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 739ab27746eSHong Zhang row i of U */ 74082502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 741ab27746eSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 74282502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 74328de702eSHong Zhang jl = il + mbs; 744ab27746eSHong Zhang for (i=0; i<mbs; i++) { 745ab27746eSHong Zhang jl[i] = mbs; il[0] = 0; 74649b5e25fSSatish Balay } 74782502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 74828de702eSHong Zhang uik = dk + 4; 749ab27746eSHong Zhang 750ab27746eSHong Zhang ai = a->i; aj = a->j; aa = a->a; 751ab27746eSHong Zhang 752ab27746eSHong Zhang /* for each row k */ 753ab27746eSHong Zhang for (k = 0; k<mbs; k++){ 754ab27746eSHong Zhang 755ab27746eSHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 756ab27746eSHong Zhang jmin = ai[k]; jmax = ai[k+1]; 757ab27746eSHong Zhang ap = aa + jmin*4; 758ab27746eSHong Zhang for (j = jmin; j < jmax; j++){ 759ab27746eSHong Zhang vj = aj[j]; /* block col. index */ 760ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 761ab27746eSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 76249b5e25fSSatish Balay } 763ab27746eSHong Zhang 764ab27746eSHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 765ab27746eSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 766ab27746eSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 767ab27746eSHong Zhang 768057f5ba7SHong Zhang while (i < k){ 769ab27746eSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 770ab27746eSHong Zhang 771ab27746eSHong Zhang /* compute multiplier */ 772ab27746eSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 773ab27746eSHong Zhang 774ab27746eSHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 775ab27746eSHong Zhang diag = ba + i*4; 776ab27746eSHong Zhang u = ba + ili*4; 777ab27746eSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 778ab27746eSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 779ab27746eSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 780ab27746eSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 781ab27746eSHong Zhang 782ab27746eSHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 783ab27746eSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 784ab27746eSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 785ab27746eSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 786ab27746eSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 787ab27746eSHong Zhang 788ab27746eSHong Zhang /* update -U(i,k): ba[ili] = uik */ 789ab27746eSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 790ab27746eSHong Zhang 791ab27746eSHong Zhang /* add multiple of row i to k-th row ... */ 792ab27746eSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 793ab27746eSHong Zhang if (jmin < jmax){ 794ab27746eSHong Zhang for (j=jmin; j<jmax; j++) { 795ab27746eSHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 796ab27746eSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 797ab27746eSHong Zhang u = ba + j*4; 798ab27746eSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 799ab27746eSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 800ab27746eSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 801ab27746eSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 80249b5e25fSSatish Balay } 803ab27746eSHong Zhang 804ab27746eSHong Zhang /* ... add i to row list for next nonzero entry */ 805ab27746eSHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 806ab27746eSHong Zhang j = bj[jmin]; 807ab27746eSHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 80849b5e25fSSatish Balay } 809ab27746eSHong Zhang i = nexti; 81049b5e25fSSatish Balay } 811ab27746eSHong Zhang 812ab27746eSHong Zhang /* save nonzero entries in k-th row of U ... */ 813ab27746eSHong Zhang 81449b5e25fSSatish Balay /* invert diagonal block */ 815ab27746eSHong Zhang diag = ba+k*4; 816ab27746eSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 817ab27746eSHong Zhang ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 818ab27746eSHong Zhang 819ab27746eSHong Zhang jmin = bi[k]; jmax = bi[k+1]; 820ab27746eSHong Zhang if (jmin < jmax) { 821ab27746eSHong Zhang for (j=jmin; j<jmax; j++){ 822ab27746eSHong Zhang vj = bj[j]; /* block col. index of U */ 823ab27746eSHong Zhang u = ba + j*4; 824ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 825ab27746eSHong Zhang for (k1=0; k1<4; k1++){ 826ab27746eSHong Zhang *u++ = *rtmp_ptr; 827ab27746eSHong Zhang *rtmp_ptr++ = 0.0; 828ab27746eSHong Zhang } 829ab27746eSHong Zhang } 830ab27746eSHong Zhang 831ab27746eSHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 832ab27746eSHong Zhang il[k] = jmin; 833ab27746eSHong Zhang i = bj[jmin]; 834ab27746eSHong Zhang jl[k] = jl[i]; jl[i] = k; 835ab27746eSHong Zhang } 83649b5e25fSSatish Balay } 83749b5e25fSSatish Balay 83849b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 839ab27746eSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 840ab27746eSHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 841ab27746eSHong Zhang 842ab27746eSHong Zhang C->factor = FACTOR_CHOLESKY; 84349b5e25fSSatish Balay C->assembled = PETSC_TRUE; 8445c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 845b0a32e0cSBarry Smith PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 84649b5e25fSSatish Balay PetscFunctionReturn(0); 84749b5e25fSSatish Balay } 84849b5e25fSSatish Balay 84949b5e25fSSatish Balay /* 8505c0bcdfcSHong Zhang Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 85191602c66SHong Zhang Version for blocks are 1 by 1. 85249b5e25fSSatish Balay */ 8534a2ae208SSatish Balay #undef __FUNCT__ 8544a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1" 8556f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B) 85649b5e25fSSatish Balay { 85749b5e25fSSatish Balay Mat C = *B; 85849b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 85949b5e25fSSatish Balay IS ip = b->row; 860351d0355SHong Zhang int *rip,ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j; 861cb718733SHong Zhang int *ai,*aj,*r; 862ab27746eSHong Zhang int k,jmin,jmax,*jl,*il,vj,nexti,ili; 863066653e3SSatish Balay MatScalar *rtmp; 8642d007305SHong Zhang MatScalar *ba = b->a,*aa,ak; 86549b5e25fSSatish Balay MatScalar dk,uikdi; 86649b5e25fSSatish Balay 86749b5e25fSSatish Balay PetscFunctionBegin; 86849b5e25fSSatish Balay ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 869cb718733SHong Zhang if (!a->permute){ 8702d007305SHong Zhang ai = a->i; aj = a->j; aa = a->a; 8712d007305SHong Zhang } else { 8722d007305SHong Zhang ai = a->inew; aj = a->jnew; 87382502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 874cb718733SHong Zhang ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 87582502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr); 8762d007305SHong Zhang ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 8772d007305SHong Zhang 8782d007305SHong Zhang jmin = ai[0]; jmax = ai[mbs]; 8792d007305SHong Zhang for (j=jmin; j<jmax; j++){ 880cb718733SHong Zhang while (r[j] != j){ 8812d007305SHong Zhang k = r[j]; r[j] = r[k]; r[k] = k; 8822d007305SHong Zhang ak = aa[k]; aa[k] = aa[j]; aa[j] = ak; 8832d007305SHong Zhang } 8842d007305SHong Zhang } 885ac355199SBarry Smith ierr = PetscFree(r);CHKERRQ(ierr); 8862d007305SHong Zhang } 8872d007305SHong Zhang 88891602c66SHong Zhang /* initialization */ 88949b5e25fSSatish Balay /* il and jl record the first nonzero element in each row of the accessing 89049b5e25fSSatish Balay window U(0:k, k:mbs-1). 89149b5e25fSSatish Balay jl: list of rows to be added to uneliminated rows 89249b5e25fSSatish Balay i>= k: jl(i) is the first row to be added to row i 89349b5e25fSSatish Balay i< k: jl(i) is the row following row i in some list of rows 89449b5e25fSSatish Balay jl(i) = mbs indicates the end of a list 89549b5e25fSSatish Balay il(i): points to the first nonzero element in columns k,...,mbs-1 of 89649b5e25fSSatish Balay row i of U */ 89782502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 89882502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 89928de702eSHong Zhang jl = il + mbs; 90049b5e25fSSatish Balay for (i=0; i<mbs; i++) { 90149b5e25fSSatish Balay rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 90249b5e25fSSatish Balay } 90349b5e25fSSatish Balay 90491602c66SHong Zhang /* for each row k */ 90549b5e25fSSatish Balay for (k = 0; k<mbs; k++){ 90649b5e25fSSatish Balay 90791602c66SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 90849b5e25fSSatish Balay jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 909057f5ba7SHong Zhang 91049b5e25fSSatish Balay for (j = jmin; j < jmax; j++){ 911351d0355SHong Zhang vj = rip[aj[j]]; 912ab27746eSHong Zhang rtmp[vj] = aa[j]; 91349b5e25fSSatish Balay } 91449b5e25fSSatish Balay 91591602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 91649b5e25fSSatish Balay dk = rtmp[k]; 91749b5e25fSSatish Balay i = jl[k]; /* first row to be added to k_th row */ 91849b5e25fSSatish Balay 919057f5ba7SHong Zhang while (i < k){ 92049b5e25fSSatish Balay nexti = jl[i]; /* next row to be added to k_th row */ 92149b5e25fSSatish Balay 92291602c66SHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 92349b5e25fSSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 92449b5e25fSSatish Balay uikdi = - ba[ili]*ba[i]; 92549b5e25fSSatish Balay dk += uikdi*ba[ili]; 926658e7b3eSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 92749b5e25fSSatish Balay 92891602c66SHong Zhang /* add multiple of row i to k-th row ... */ 92949b5e25fSSatish Balay jmin = ili + 1; jmax = bi[i+1]; 93049b5e25fSSatish Balay if (jmin < jmax){ 93149b5e25fSSatish Balay for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 93291602c66SHong Zhang /* ... add i to row list for next nonzero entry */ 93349b5e25fSSatish Balay il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 93449b5e25fSSatish Balay j = bj[jmin]; 93549b5e25fSSatish Balay jl[i] = jl[j]; jl[j] = i; /* update jl */ 93649b5e25fSSatish Balay } 937ab27746eSHong Zhang i = nexti; 93849b5e25fSSatish Balay } 93949b5e25fSSatish Balay 94091602c66SHong Zhang /* check for zero pivot and save diagoanl element */ 94149b5e25fSSatish Balay if (dk == 0.0){ 94229bbc08cSBarry Smith SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 943671cb588SHong Zhang /* 9441223c01bSHong Zhang } else if (PetscRealPart(dk) < 0.0){ 9451223c01bSHong Zhang SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Negative pivot: d[%d] = %g\n",k,dk); 946671cb588SHong Zhang */ 94749b5e25fSSatish Balay } 94849b5e25fSSatish Balay 94991602c66SHong Zhang /* save nonzero entries in k-th row of U ... */ 950f3dd7b2dSKris Buschelman ba[k] = 1.0/dk; 95149b5e25fSSatish Balay jmin = bi[k]; jmax = bi[k+1]; 95249b5e25fSSatish Balay if (jmin < jmax) { 95349b5e25fSSatish Balay for (j=jmin; j<jmax; j++){ 954cc0c071aSHong Zhang vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0; 95549b5e25fSSatish Balay } 95691602c66SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 95749b5e25fSSatish Balay il[k] = jmin; 95849b5e25fSSatish Balay i = bj[jmin]; 95949b5e25fSSatish Balay jl[k] = jl[i]; jl[i] = k; 96049b5e25fSSatish Balay } 96149b5e25fSSatish Balay } 96249b5e25fSSatish Balay 96349b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 96449b5e25fSSatish Balay ierr = PetscFree(il);CHKERRQ(ierr); 965cb718733SHong Zhang if (a->permute){ 966cb718733SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 967cb718733SHong Zhang } 96849b5e25fSSatish Balay 96949b5e25fSSatish Balay ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 9708be8c0c7SHong Zhang C->factor = FACTOR_CHOLESKY; 97149b5e25fSSatish Balay C->assembled = PETSC_TRUE; 9725c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 973b0a32e0cSBarry Smith PetscLogFlops(b->mbs); 97449b5e25fSSatish Balay PetscFunctionReturn(0); 97549b5e25fSSatish Balay } 97649b5e25fSSatish Balay 977671cb588SHong Zhang /* 978671cb588SHong Zhang Version for when blocks are 1 by 1 Using natural ordering 979671cb588SHong Zhang */ 9804a2ae208SSatish Balay #undef __FUNCT__ 9814a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" 982671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B) 983671cb588SHong Zhang { 984671cb588SHong Zhang Mat C = *B; 985671cb588SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 986671cb588SHong Zhang int ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j; 9871ad1de41SHong Zhang int *ai,*aj; 988671cb588SHong Zhang int k,jmin,jmax,*jl,*il,vj,nexti,ili; 9891ad1de41SHong Zhang MatScalar *rtmp,*ba = b->a,*aa,dk,uikdi; 990671cb588SHong Zhang 991671cb588SHong Zhang PetscFunctionBegin; 992671cb588SHong Zhang 993671cb588SHong Zhang /* initialization */ 994671cb588SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 995671cb588SHong Zhang window U(0:k, k:mbs-1). 996671cb588SHong Zhang jl: list of rows to be added to uneliminated rows 997671cb588SHong Zhang i>= k: jl(i) is the first row to be added to row i 998671cb588SHong Zhang i< k: jl(i) is the row following row i in some list of rows 999671cb588SHong Zhang jl(i) = mbs indicates the end of a list 1000671cb588SHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 1001671cb588SHong Zhang row i of U */ 100282502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 100382502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 100428de702eSHong Zhang jl = il + mbs; 1005671cb588SHong Zhang for (i=0; i<mbs; i++) { 1006671cb588SHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1007671cb588SHong Zhang } 1008671cb588SHong Zhang 1009671cb588SHong Zhang ai = a->i; aj = a->j; aa = a->a; 1010671cb588SHong Zhang 1011671cb588SHong Zhang /* for each row k */ 1012671cb588SHong Zhang for (k = 0; k<mbs; k++){ 1013671cb588SHong Zhang 1014671cb588SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 1015671cb588SHong Zhang jmin = ai[k]; jmax = ai[k+1]; 1016057f5ba7SHong Zhang 1017671cb588SHong Zhang for (j = jmin; j < jmax; j++){ 1018671cb588SHong Zhang vj = aj[j]; 1019671cb588SHong Zhang rtmp[vj] = aa[j]; 1020671cb588SHong Zhang } 1021671cb588SHong Zhang 1022671cb588SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 1023671cb588SHong Zhang dk = rtmp[k]; 1024671cb588SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 1025671cb588SHong Zhang 1026057f5ba7SHong Zhang while (i < k){ 1027671cb588SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 1028671cb588SHong Zhang 1029671cb588SHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 1030671cb588SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1031671cb588SHong Zhang uikdi = - ba[ili]*ba[i]; 1032671cb588SHong Zhang dk += uikdi*ba[ili]; 1033671cb588SHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 1034671cb588SHong Zhang 1035671cb588SHong Zhang /* add multiple of row i to k-th row ... */ 1036671cb588SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 1037671cb588SHong Zhang if (jmin < jmax){ 1038671cb588SHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1039671cb588SHong Zhang /* ... add i to row list for next nonzero entry */ 1040671cb588SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 1041671cb588SHong Zhang j = bj[jmin]; 1042671cb588SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 1043671cb588SHong Zhang } 1044671cb588SHong Zhang i = nexti; 1045671cb588SHong Zhang } 1046671cb588SHong Zhang 1047671cb588SHong Zhang /* check for zero pivot and save diagoanl element */ 1048671cb588SHong Zhang if (dk == 0.0){ 1049671cb588SHong Zhang SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 1050671cb588SHong Zhang /* 1051671cb588SHong Zhang } else if (PetscRealPart(dk) < 0){ 10521223c01bSHong Zhang SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Negative pivot: d[%d] = %g\n",k,dk); 1053671cb588SHong Zhang */ 1054671cb588SHong Zhang } 1055671cb588SHong Zhang 1056671cb588SHong Zhang /* save nonzero entries in k-th row of U ... */ 1057671cb588SHong Zhang ba[k] = 1.0/dk; 1058671cb588SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 1059671cb588SHong Zhang if (jmin < jmax) { 1060671cb588SHong Zhang for (j=jmin; j<jmax; j++){ 1061671cb588SHong Zhang vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0; 1062671cb588SHong Zhang } 1063671cb588SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 1064671cb588SHong Zhang il[k] = jmin; 1065671cb588SHong Zhang i = bj[jmin]; 1066671cb588SHong Zhang jl[k] = jl[i]; jl[i] = k; 1067671cb588SHong Zhang } 1068671cb588SHong Zhang } 1069671cb588SHong Zhang 1070671cb588SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 1071671cb588SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1072671cb588SHong Zhang 1073671cb588SHong Zhang C->factor = FACTOR_CHOLESKY; 1074671cb588SHong Zhang C->assembled = PETSC_TRUE; 1075671cb588SHong Zhang C->preallocated = PETSC_TRUE; 1076b0a32e0cSBarry Smith PetscLogFlops(b->mbs); 1077671cb588SHong Zhang PetscFunctionReturn(0); 1078671cb588SHong Zhang } 1079671cb588SHong Zhang 10804a2ae208SSatish Balay #undef __FUNCT__ 10814a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" 1082b8b23757SHong Zhang int MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,PetscReal f) 108349b5e25fSSatish Balay { 10844445ddedSHong Zhang int ierr; 108549b5e25fSSatish Balay Mat C; 108649b5e25fSSatish Balay 108749b5e25fSSatish Balay PetscFunctionBegin; 1088b8b23757SHong Zhang ierr = MatCholeskyFactorSymbolic(A,perm,f,&C);CHKERRQ(ierr); 1089a4ada70bSHong Zhang ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr); 10904445ddedSHong Zhang ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 109149b5e25fSSatish Balay PetscFunctionReturn(0); 109249b5e25fSSatish Balay } 109349b5e25fSSatish Balay 109449b5e25fSSatish Balay 1095