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