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