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