xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact7.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 
2 #include <../src/mat/impls/sbaij/seq/sbaij.h>
3 #include <petsc/private/kernels/blockinvert.h>
4 
5 /* Version for when blocks are 5 by 5  */
6 #undef __FUNCT__
7 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_5"
8 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat C,Mat A,const MatFactorInfo *info)
9 {
10   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
11   IS             perm = b->row;
12   PetscErrorCode ierr;
13   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
14   PetscInt       i,j,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili,ipvt[5];
15   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
16   MatScalar      *u,*d,*rtmp,*rtmp_ptr,work[25];
17   PetscReal      shift = info->shiftamount;
18 
19   PetscFunctionBegin;
20   /* initialization */
21   ierr = PetscCalloc1(25*mbs,&rtmp);CHKERRQ(ierr);
22   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
23   for (i=0; i<mbs; i++) {
24     jl[i] = mbs; il[0] = 0;
25   }
26   ierr = PetscMalloc2(25,&dk,25,&uik);CHKERRQ(ierr);
27   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
28 
29   /* check permutation */
30   if (!a->permute) {
31     ai = a->i; aj = a->j; aa = a->a;
32   } else {
33     ai   = a->inew; aj = a->jnew;
34     ierr = PetscMalloc1(25*ai[mbs],&aa);CHKERRQ(ierr);
35     ierr = PetscMemcpy(aa,a->a,25*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
36     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
37     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
38 
39     for (i=0; i<mbs; i++) {
40       jmin = ai[i]; jmax = ai[i+1];
41       for (j=jmin; j<jmax; j++) {
42         while (a2anew[j] != j) {
43           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
44           for (k1=0; k1<25; k1++) {
45             dk[k1]      = aa[k*25+k1];
46             aa[k*25+k1] = aa[j*25+k1];
47             aa[j*25+k1] = dk[k1];
48           }
49         }
50         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
51         if (i > aj[j]) {
52           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
53           ap = aa + j*25;                     /* ptr to the beginning of j-th block of aa */
54           for (k=0; k<25; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
55           for (k=0; k<5; k++) {               /* j-th block of aa <- dk^T */
56             for (k1=0; k1<5; k1++) *ap++ = dk[k + 5*k1];
57           }
58         }
59       }
60     }
61     ierr = PetscFree(a2anew);CHKERRQ(ierr);
62   }
63 
64   /* for each row k */
65   for (k = 0; k<mbs; k++) {
66 
67     /*initialize k-th row with elements nonzero in row perm(k) of A */
68     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
69     if (jmin < jmax) {
70       ap = aa + jmin*25;
71       for (j = jmin; j < jmax; j++) {
72         vj       = perm_ptr[aj[j]];   /* block col. index */
73         rtmp_ptr = rtmp + vj*25;
74         for (i=0; i<25; i++) *rtmp_ptr++ = *ap++;
75       }
76     }
77 
78     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
79     ierr = PetscMemcpy(dk,rtmp+k*25,25*sizeof(MatScalar));CHKERRQ(ierr);
80     i    = jl[k]; /* first row to be added to k_th row  */
81 
82     while (i < mbs) {
83       nexti = jl[i]; /* next row to be added to k_th row */
84 
85       /* compute multiplier */
86       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
87 
88       /* uik = -inv(Di)*U_bar(i,k) */
89       d = ba + i*25;
90       u = ba + ili*25;
91 
92       uik[0] = -(d[0]*u[0] + d[5]*u[1] + d[10]*u[2] + d[15]*u[3] + d[20]*u[4]);
93       uik[1] = -(d[1]*u[0] + d[6]*u[1] + d[11]*u[2] + d[16]*u[3] + d[21]*u[4]);
94       uik[2] = -(d[2]*u[0] + d[7]*u[1] + d[12]*u[2] + d[17]*u[3] + d[22]*u[4]);
95       uik[3] = -(d[3]*u[0] + d[8]*u[1] + d[13]*u[2] + d[18]*u[3] + d[23]*u[4]);
96       uik[4] = -(d[4]*u[0] + d[9]*u[1] + d[14]*u[2] + d[19]*u[3] + d[24]*u[4]);
97 
98       uik[5] = -(d[0]*u[5] + d[5]*u[6] + d[10]*u[7] + d[15]*u[8] + d[20]*u[9]);
99       uik[6] = -(d[1]*u[5] + d[6]*u[6] + d[11]*u[7] + d[16]*u[8] + d[21]*u[9]);
100       uik[7] = -(d[2]*u[5] + d[7]*u[6] + d[12]*u[7] + d[17]*u[8] + d[22]*u[9]);
101       uik[8] = -(d[3]*u[5] + d[8]*u[6] + d[13]*u[7] + d[18]*u[8] + d[23]*u[9]);
102       uik[9] = -(d[4]*u[5] + d[9]*u[6] + d[14]*u[7] + d[19]*u[8] + d[24]*u[9]);
103 
104       uik[10]= -(d[0]*u[10] + d[5]*u[11] + d[10]*u[12] + d[15]*u[13] + d[20]*u[14]);
105       uik[11]= -(d[1]*u[10] + d[6]*u[11] + d[11]*u[12] + d[16]*u[13] + d[21]*u[14]);
106       uik[12]= -(d[2]*u[10] + d[7]*u[11] + d[12]*u[12] + d[17]*u[13] + d[22]*u[14]);
107       uik[13]= -(d[3]*u[10] + d[8]*u[11] + d[13]*u[12] + d[18]*u[13] + d[23]*u[14]);
108       uik[14]= -(d[4]*u[10] + d[9]*u[11] + d[14]*u[12] + d[19]*u[13] + d[24]*u[14]);
109 
110       uik[15]= -(d[0]*u[15] + d[5]*u[16] + d[10]*u[17] + d[15]*u[18] + d[20]*u[19]);
111       uik[16]= -(d[1]*u[15] + d[6]*u[16] + d[11]*u[17] + d[16]*u[18] + d[21]*u[19]);
112       uik[17]= -(d[2]*u[15] + d[7]*u[16] + d[12]*u[17] + d[17]*u[18] + d[22]*u[19]);
113       uik[18]= -(d[3]*u[15] + d[8]*u[16] + d[13]*u[17] + d[18]*u[18] + d[23]*u[19]);
114       uik[19]= -(d[4]*u[15] + d[9]*u[16] + d[14]*u[17] + d[19]*u[18] + d[24]*u[19]);
115 
116       uik[20]= -(d[0]*u[20] + d[5]*u[21] + d[10]*u[22] + d[15]*u[23] + d[20]*u[24]);
117       uik[21]= -(d[1]*u[20] + d[6]*u[21] + d[11]*u[22] + d[16]*u[23] + d[21]*u[24]);
118       uik[22]= -(d[2]*u[20] + d[7]*u[21] + d[12]*u[22] + d[17]*u[23] + d[22]*u[24]);
119       uik[23]= -(d[3]*u[20] + d[8]*u[21] + d[13]*u[22] + d[18]*u[23] + d[23]*u[24]);
120       uik[24]= -(d[4]*u[20] + d[9]*u[21] + d[14]*u[22] + d[19]*u[23] + d[24]*u[24]);
121 
122 
123       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
124       dk[0] +=  uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3] + uik[4]*u[4];
125       dk[1] +=  uik[5]*u[0] + uik[6]*u[1] + uik[7]*u[2] + uik[8]*u[3] + uik[9]*u[4];
126       dk[2] += uik[10]*u[0]+ uik[11]*u[1]+ uik[12]*u[2]+ uik[13]*u[3]+ uik[14]*u[4];
127       dk[3] += uik[15]*u[0]+ uik[16]*u[1]+ uik[17]*u[2]+ uik[18]*u[3]+ uik[19]*u[4];
128       dk[4] += uik[20]*u[0]+ uik[21]*u[1]+ uik[22]*u[2]+ uik[23]*u[3]+ uik[24]*u[4];
129 
130       dk[5] +=  uik[0]*u[5] + uik[1]*u[6] + uik[2]*u[7] + uik[3]*u[8] + uik[4]*u[9];
131       dk[6] +=  uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8] + uik[9]*u[9];
132       dk[7] += uik[10]*u[5]+ uik[11]*u[6]+ uik[12]*u[7]+ uik[13]*u[8]+ uik[14]*u[9];
133       dk[8] += uik[15]*u[5]+ uik[16]*u[6]+ uik[17]*u[7]+ uik[18]*u[8]+ uik[19]*u[9];
134       dk[9] += uik[20]*u[5]+ uik[21]*u[6]+ uik[22]*u[7]+ uik[23]*u[8]+ uik[24]*u[9];
135 
136       dk[10] +=  uik[0]*u[10] + uik[1]*u[11] + uik[2]*u[12] + uik[3]*u[13] + uik[4]*u[14];
137       dk[11] +=  uik[5]*u[10] + uik[6]*u[11] + uik[7]*u[12] + uik[8]*u[13] + uik[9]*u[14];
138       dk[12] += uik[10]*u[10]+ uik[11]*u[11]+ uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14];
139       dk[13] += uik[15]*u[10]+ uik[16]*u[11]+ uik[17]*u[12]+ uik[18]*u[13]+ uik[19]*u[14];
140       dk[14] += uik[20]*u[10]+ uik[21]*u[11]+ uik[22]*u[12]+ uik[23]*u[13]+ uik[24]*u[14];
141 
142       dk[15] +=  uik[0]*u[15] + uik[1]*u[16] + uik[2]*u[17] + uik[3]*u[18] + uik[4]*u[19];
143       dk[16] +=  uik[5]*u[15] + uik[6]*u[16] + uik[7]*u[17] + uik[8]*u[18] + uik[9]*u[19];
144       dk[17] += uik[10]*u[15]+ uik[11]*u[16]+ uik[12]*u[17]+ uik[13]*u[18]+ uik[14]*u[19];
145       dk[18] += uik[15]*u[15]+ uik[16]*u[16]+ uik[17]*u[17]+ uik[18]*u[18]+ uik[19]*u[19];
146       dk[19] += uik[20]*u[15]+ uik[21]*u[16]+ uik[22]*u[17]+ uik[23]*u[18]+ uik[24]*u[19];
147 
148       dk[20] +=  uik[0]*u[20] + uik[1]*u[21] + uik[2]*u[22] + uik[3]*u[23] + uik[4]*u[24];
149       dk[21] +=  uik[5]*u[20] + uik[6]*u[21] + uik[7]*u[22] + uik[8]*u[23] + uik[9]*u[24];
150       dk[22] += uik[10]*u[20]+ uik[11]*u[21]+ uik[12]*u[22]+ uik[13]*u[23]+ uik[14]*u[24];
151       dk[23] += uik[15]*u[20]+ uik[16]*u[21]+ uik[17]*u[22]+ uik[18]*u[23]+ uik[19]*u[24];
152       dk[24] += uik[20]*u[20]+ uik[21]*u[21]+ uik[22]*u[22]+ uik[23]*u[23]+ uik[24]*u[24];
153 
154       ierr = PetscLogFlops(125.0*4.0);CHKERRQ(ierr);
155 
156       /* update -U(i,k) */
157       ierr = PetscMemcpy(ba+ili*25,uik,25*sizeof(MatScalar));CHKERRQ(ierr);
158 
159       /* add multiple of row i to k-th row ... */
160       jmin = ili + 1; jmax = bi[i+1];
161       if (jmin < jmax) {
162         for (j=jmin; j<jmax; j++) {
163           /* rtmp += -U(i,k)^T * U_bar(i,j) */
164           rtmp_ptr     = rtmp + bj[j]*25;
165           u            = ba + j*25;
166           rtmp_ptr[0] +=  uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3] + uik[4]*u[4];
167           rtmp_ptr[1] +=  uik[5]*u[0] + uik[6]*u[1] + uik[7]*u[2] + uik[8]*u[3] + uik[9]*u[4];
168           rtmp_ptr[2] += uik[10]*u[0]+ uik[11]*u[1]+ uik[12]*u[2]+ uik[13]*u[3]+ uik[14]*u[4];
169           rtmp_ptr[3] += uik[15]*u[0]+ uik[16]*u[1]+ uik[17]*u[2]+ uik[18]*u[3]+ uik[19]*u[4];
170           rtmp_ptr[4] += uik[20]*u[0]+ uik[21]*u[1]+ uik[22]*u[2]+ uik[23]*u[3]+ uik[24]*u[4];
171 
172           rtmp_ptr[5] +=  uik[0]*u[5] + uik[1]*u[6] + uik[2]*u[7] + uik[3]*u[8] + uik[4]*u[9];
173           rtmp_ptr[6] +=  uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8] + uik[9]*u[9];
174           rtmp_ptr[7] += uik[10]*u[5]+ uik[11]*u[6]+ uik[12]*u[7]+ uik[13]*u[8]+ uik[14]*u[9];
175           rtmp_ptr[8] += uik[15]*u[5]+ uik[16]*u[6]+ uik[17]*u[7]+ uik[18]*u[8]+ uik[19]*u[9];
176           rtmp_ptr[9] += uik[20]*u[5]+ uik[21]*u[6]+ uik[22]*u[7]+ uik[23]*u[8]+ uik[24]*u[9];
177 
178           rtmp_ptr[10] +=  uik[0]*u[10] + uik[1]*u[11] + uik[2]*u[12] + uik[3]*u[13] + uik[4]*u[14];
179           rtmp_ptr[11] +=  uik[5]*u[10] + uik[6]*u[11] + uik[7]*u[12] + uik[8]*u[13] + uik[9]*u[14];
180           rtmp_ptr[12] += uik[10]*u[10]+ uik[11]*u[11]+ uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14];
181           rtmp_ptr[13] += uik[15]*u[10]+ uik[16]*u[11]+ uik[17]*u[12]+ uik[18]*u[13]+ uik[19]*u[14];
182           rtmp_ptr[14] += uik[20]*u[10]+ uik[21]*u[11]+ uik[22]*u[12]+ uik[23]*u[13]+ uik[24]*u[14];
183 
184           rtmp_ptr[15] +=  uik[0]*u[15] + uik[1]*u[16] + uik[2]*u[17] + uik[3]*u[18] + uik[4]*u[19];
185           rtmp_ptr[16] +=  uik[5]*u[15] + uik[6]*u[16] + uik[7]*u[17] + uik[8]*u[18] + uik[9]*u[19];
186           rtmp_ptr[17] += uik[10]*u[15]+ uik[11]*u[16]+ uik[12]*u[17]+ uik[13]*u[18]+ uik[14]*u[19];
187           rtmp_ptr[18] += uik[15]*u[15]+ uik[16]*u[16]+ uik[17]*u[17]+ uik[18]*u[18]+ uik[19]*u[19];
188           rtmp_ptr[19] += uik[20]*u[15]+ uik[21]*u[16]+ uik[22]*u[17]+ uik[23]*u[18]+ uik[24]*u[19];
189 
190           rtmp_ptr[20] +=  uik[0]*u[20] + uik[1]*u[21] + uik[2]*u[22] + uik[3]*u[23] + uik[4]*u[24];
191           rtmp_ptr[21] +=  uik[5]*u[20] + uik[6]*u[21] + uik[7]*u[22] + uik[8]*u[23] + uik[9]*u[24];
192           rtmp_ptr[22] += uik[10]*u[20]+ uik[11]*u[21]+ uik[12]*u[22]+ uik[13]*u[23]+ uik[14]*u[24];
193           rtmp_ptr[23] += uik[15]*u[20]+ uik[16]*u[21]+ uik[17]*u[22]+ uik[18]*u[23]+ uik[19]*u[24];
194           rtmp_ptr[24] += uik[20]*u[20]+ uik[21]*u[21]+ uik[22]*u[22]+ uik[23]*u[23]+ uik[24]*u[24];
195         }
196         ierr = PetscLogFlops(2.0*125.0*(jmax-jmin));CHKERRQ(ierr);
197 
198         /* ... add i to row list for next nonzero entry */
199         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
200         j     = bj[jmin];
201         jl[i] = jl[j]; jl[j] = i; /* update jl */
202       }
203       i = nexti;
204     }
205 
206     /* save nonzero entries in k-th row of U ... */
207 
208     /* invert diagonal block */
209     d    = ba+k*25;
210     ierr = PetscMemcpy(d,dk,25*sizeof(MatScalar));CHKERRQ(ierr);
211     ierr = PetscKernel_A_gets_inverse_A_5(d,ipvt,work,shift);CHKERRQ(ierr);
212 
213     jmin = bi[k]; jmax = bi[k+1];
214     if (jmin < jmax) {
215       for (j=jmin; j<jmax; j++) {
216         vj       = bj[j];      /* block col. index of U */
217         u        = ba + j*25;
218         rtmp_ptr = rtmp + vj*25;
219         for (k1=0; k1<25; k1++) {
220           *u++        = *rtmp_ptr;
221           *rtmp_ptr++ = 0.0;
222         }
223       }
224 
225       /* ... add k to row list for first nonzero entry in k-th row */
226       il[k] = jmin;
227       i     = bj[jmin];
228       jl[k] = jl[i]; jl[i] = k;
229     }
230   }
231 
232   ierr = PetscFree(rtmp);CHKERRQ(ierr);
233   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
234   ierr = PetscFree2(dk,uik);CHKERRQ(ierr);
235   if (a->permute) {
236     ierr = PetscFree(aa);CHKERRQ(ierr);
237   }
238 
239   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
240 
241   C->ops->solve          = MatSolve_SeqSBAIJ_5_inplace;
242   C->ops->solvetranspose = MatSolve_SeqSBAIJ_5_inplace;
243   C->assembled           = PETSC_TRUE;
244   C->preallocated        = PETSC_TRUE;
245 
246   ierr = PetscLogFlops(1.3333*125*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
247   PetscFunctionReturn(0);
248 }
249