xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact4.c (revision a4efd8ea4bb0600a5aef92b2c23cbfa314827193)
1 #define PETSCMAT_DLL
2 
3 #include "src/mat/impls/sbaij/seq/sbaij.h"
4 #include "src/inline/ilu.h"
5 
6 /*
7       Version for when blocks are 3 by 3 Using natural ordering
8 */
9 #undef __FUNCT__
10 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering"
11 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
12 {
13   Mat            C = *B;
14   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
15   PetscErrorCode ierr;
16   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
17   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
18   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
19   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
20   PetscReal      shift = info->shiftinblocks;
21 
22   PetscFunctionBegin;
23 
24   /* initialization */
25   ierr = PetscMalloc(9*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
26   ierr = PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));CHKERRQ(ierr);
27   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
28   jl   = il + mbs;
29   for (i=0; i<mbs; i++) {
30     jl[i] = mbs; il[0] = 0;
31   }
32   ierr = PetscMalloc(18*sizeof(MatScalar),&dk);CHKERRQ(ierr);
33   uik  = dk + 9;
34 
35   ai   = a->i; aj = a->j; aa = a->a;
36 
37   /* for each row k */
38   for (k = 0; k<mbs; k++){
39 
40     /*initialize k-th row with elements nonzero in row k of A */
41     jmin = ai[k]; jmax = ai[k+1];
42     if (jmin < jmax) {
43       ap = aa + jmin*9;
44       for (j = jmin; j < jmax; j++){
45         vj = aj[j];         /* block col. index */
46         rtmp_ptr = rtmp + vj*9;
47         for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
48       }
49     }
50 
51     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
52     ierr = PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));CHKERRQ(ierr);
53     i = jl[k]; /* first row to be added to k_th row  */
54 
55     while (i < mbs){
56       nexti = jl[i]; /* next row to be added to k_th row */
57 
58       /* compute multiplier */
59       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
60 
61       /* uik = -inv(Di)*U_bar(i,k) */
62       diag = ba + i*9;
63       u    = ba + ili*9;
64 
65       uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
66       uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
67       uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);
68 
69       uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
70       uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
71       uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);
72 
73       uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
74       uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
75       uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);
76 
77       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
78       dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
79       dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
80       dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
81 
82       dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
83       dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
84       dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
85 
86       dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
87       dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
88       dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
89 
90       ierr = PetscLogFlops(27*4);CHKERRQ(ierr);
91 
92       /* update -U(i,k) */
93       ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr);
94 
95       /* add multiple of row i to k-th row ... */
96       jmin = ili + 1; jmax = bi[i+1];
97       if (jmin < jmax){
98         for (j=jmin; j<jmax; j++) {
99           /* rtmp += -U(i,k)^T * U_bar(i,j) */
100           rtmp_ptr = rtmp + bj[j]*9;
101           u = ba + j*9;
102           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
103           rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
104           rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
105 
106           rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
107           rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
108           rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
109 
110           rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
111           rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
112           rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
113         }
114         ierr = PetscLogFlops(2*27*(jmax-jmin));CHKERRQ(ierr);
115 
116         /* ... add i to row list for next nonzero entry */
117         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
118         j     = bj[jmin];
119         jl[i] = jl[j]; jl[j] = i; /* update jl */
120       }
121       i = nexti;
122     }
123 
124     /* save nonzero entries in k-th row of U ... */
125 
126     /* invert diagonal block */
127     diag = ba+k*9;
128     ierr = PetscMemcpy(diag,dk,9*sizeof(MatScalar));CHKERRQ(ierr);
129     ierr = Kernel_A_gets_inverse_A_3(diag,shift);CHKERRQ(ierr);
130 
131     jmin = bi[k]; jmax = bi[k+1];
132     if (jmin < jmax) {
133       for (j=jmin; j<jmax; j++){
134          vj = bj[j];           /* block col. index of U */
135          u   = ba + j*9;
136          rtmp_ptr = rtmp + vj*9;
137          for (k1=0; k1<9; k1++){
138            *u++        = *rtmp_ptr;
139            *rtmp_ptr++ = 0.0;
140          }
141       }
142 
143       /* ... add k to row list for first nonzero entry in k-th row */
144       il[k] = jmin;
145       i     = bj[jmin];
146       jl[k] = jl[i]; jl[i] = k;
147     }
148   }
149 
150   ierr = PetscFree(rtmp);CHKERRQ(ierr);
151   ierr = PetscFree(il);CHKERRQ(ierr);
152   ierr = PetscFree(dk);CHKERRQ(ierr);
153 
154   C->ops->solve          = MatSolve_SeqSBAIJ_3_NaturalOrdering;
155   C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering;
156   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering;
157   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering;
158 
159   C->assembled = PETSC_TRUE;
160   C->preallocated = PETSC_TRUE;
161   ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
162   PetscFunctionReturn(0);
163 }
164