1 /* 2 Factorization code for BAIJ format. 3 */ 4 #include "src/mat/impls/baij/seq/baij.h" 5 #include "src/inline/ilu.h" 6 7 /* ----------------------------------------------------------- */ 8 #undef __FUNCT__ 9 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_N" 10 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat A,Mat *B) 11 { 12 Mat C = *B; 13 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 14 IS isrow = b->row,isicol = b->icol; 15 PetscErrorCode ierr; 16 int *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 17 int *ajtmpold,*ajtmp,nz,row,bslog,*ai=a->i,*aj=a->j,k,flg; 18 int *diag_offset=b->diag,diag,bs=a->bs,bs2 = a->bs2,*v_pivots,*pj; 19 MatScalar *ba = b->a,*aa = a->a,*pv,*v,*rtmp,*multiplier,*v_work,*pc,*w; 20 21 PetscFunctionBegin; 22 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 23 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 24 ierr = PetscMalloc(bs2*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 25 ierr = PetscMemzero(rtmp,bs2*(n+1)*sizeof(MatScalar));CHKERRQ(ierr); 26 /* generate work space needed by dense LU factorization */ 27 ierr = PetscMalloc(bs*sizeof(int) + (bs+bs2)*sizeof(MatScalar),&v_work);CHKERRQ(ierr); 28 multiplier = v_work + bs; 29 v_pivots = (int*)(multiplier + bs2); 30 31 /* flops in while loop */ 32 bslog = 2*bs*bs2; 33 34 for (i=0; i<n; i++) { 35 nz = bi[i+1] - bi[i]; 36 ajtmp = bj + bi[i]; 37 for (j=0; j<nz; j++) { 38 ierr = PetscMemzero(rtmp+bs2*ajtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 39 } 40 /* load in initial (unfactored row) */ 41 nz = ai[r[i]+1] - ai[r[i]]; 42 ajtmpold = aj + ai[r[i]]; 43 v = aa + bs2*ai[r[i]]; 44 for (j=0; j<nz; j++) { 45 ierr = PetscMemcpy(rtmp+bs2*ic[ajtmpold[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 46 } 47 row = *ajtmp++; 48 while (row < i) { 49 pc = rtmp + bs2*row; 50 /* if (*pc) { */ 51 for (flg=0,k=0; k<bs2; k++) { if (pc[k]!=0.0) { flg = 1; break; }} 52 if (flg) { 53 pv = ba + bs2*diag_offset[row]; 54 pj = bj + diag_offset[row] + 1; 55 Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); 56 nz = bi[row+1] - diag_offset[row] - 1; 57 pv += bs2; 58 for (j=0; j<nz; j++) { 59 Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 60 } 61 PetscLogFlops(bslog*(nz+1)-bs); 62 } 63 row = *ajtmp++; 64 } 65 /* finished row so stick it into b->a */ 66 pv = ba + bs2*bi[i]; 67 pj = bj + bi[i]; 68 nz = bi[i+1] - bi[i]; 69 for (j=0; j<nz; j++) { 70 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 71 } 72 diag = diag_offset[i] - bi[i]; 73 /* invert diagonal block */ 74 w = pv + bs2*diag; 75 Kernel_A_gets_inverse_A(bs,w,v_pivots,v_work); 76 } 77 78 ierr = PetscFree(rtmp);CHKERRQ(ierr); 79 ierr = PetscFree(v_work);CHKERRQ(ierr); 80 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 81 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 82 C->factor = FACTOR_LU; 83 C->assembled = PETSC_TRUE; 84 PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ 85 PetscFunctionReturn(0); 86 } 87