Actual source code: sbaijfact10.c
1: /*$Id: sbaijfact10.c,v 1.5 2001/03/23 23:22:21 balay Exp $*/
2: #include sbaij.h
3: #include src/inline/ilu.h
5: /*
6: Version for when blocks are 6 by 6 Using natural ordering
7: */
8: int MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat A,Mat *B)
9: {
10: Mat C = *B;
11: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
12: int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
13: int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
14: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
15: MatScalar *u,*d,*w,*wp,u0,u1,u2,u3,u4,u5,u6,u7,u8,u9,u10,u11,u12;
16: MatScalar u13,u14,u15,u16,u17,u18,u19,u20,u21,u22,u23,u24,u25,u26,u27;
17: MatScalar u28,u29,u30,u31,u32,u33,u34,u35;
18:
20:
21: /* initialization */
22: PetscMalloc(36*mbs*sizeof(MatScalar),&w);
23: PetscMemzero(w,36*mbs*sizeof(MatScalar));
24: PetscMalloc(2*mbs*sizeof(int),&il);
25: jl = il + mbs;
26: for (i=0; i<mbs; i++) {
27: jl[i] = mbs; il[0] = 0;
28: }
29: PetscMalloc(72*sizeof(MatScalar),&dk);
30: uik = dk + 36;
32: ai = a->i; aj = a->j; aa = a->a;
34: /* for each row k */
35: for (k = 0; k<mbs; k++){
37: /*initialize k-th row with elements nonzero in row k of A */
38: jmin = ai[k]; jmax = ai[k+1];
39: if (jmin < jmax) {
40: ap = aa + jmin*36;
41: for (j = jmin; j < jmax; j++){
42: vj = aj[j]; /* block col. index */
43: wp = w + vj*36;
44: for (i=0; i<36; i++) *wp++ = *ap++;
45: }
46: }
48: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
49: PetscMemcpy(dk,w+k*36,36*sizeof(MatScalar));
50: i = jl[k]; /* first row to be added to k_th row */
52: while (i < mbs){
53: nexti = jl[i]; /* next row to be added to k_th row */
55: /* compute multiplier */
56: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
58: /* uik = -inv(Di)*U_bar(i,k) */
59: d = ba + i*36;
60: u = ba + ili*36;
62: u0 = u[0]; u1 = u[1]; u2 = u[2]; u3 = u[3]; u4 = u[4]; u5 = u[5]; u6 = u[6];
63: u7 = u[7]; u8 = u[8]; u9 = u[9]; u10 = u[10]; u11 = u[11]; u12 = u[12]; u13 = u[13];
64: u14 = u[14]; u15 = u[15]; u16 = u[16]; u17 = u[17]; u18 = u[18]; u19 = u[19]; u20 = u[20];
65: u21 = u[21]; u22 = u[22]; u23 = u[23]; u24 = u[24]; u25 = u[25]; u26 = u[26]; u27 = u[27];
66: u28 = u[28]; u29 = u[29]; u30 = u[30]; u31 = u[31]; u32 = u[32]; u33 = u[33]; u34 = u[34];
67: u35 = u[35];
69: uik[0] = -(d[0]*u0 + d[6]*u1 + d[12]*u2 + d[18]*u3 + d[24]*u4 + d[30]*u5);
70: uik[1] = -(d[1]*u0 + d[7]*u1 + d[13]*u2 + d[19]*u3 + d[25]*u4 + d[31]*u5);
71: uik[2] = -(d[2]*u0 + d[8]*u1 + d[14]*u2 + d[20]*u3 + d[26]*u4 + d[32]*u5);
72: uik[3] = -(d[3]*u0 + d[9]*u1 + d[15]*u2 + d[21]*u3 + d[27]*u4 + d[33]*u5);
73: uik[4] = -(d[4]*u0+ d[10]*u1 + d[16]*u2 + d[22]*u3 + d[28]*u4 + d[34]*u5);
74: uik[5] = -(d[5]*u0+ d[11]*u1 + d[17]*u2 + d[23]*u3 + d[29]*u4 + d[35]*u5);
76: uik[6] = -(d[0]*u6 + d[6]*u7 + d[12]*u8 + d[18]*u9 + d[24]*u10 + d[30]*u11);
77: uik[7] = -(d[1]*u6 + d[7]*u7 + d[13]*u8 + d[19]*u9 + d[25]*u10 + d[31]*u11);
78: uik[8] = -(d[2]*u6 + d[8]*u7 + d[14]*u8 + d[20]*u9 + d[26]*u10 + d[32]*u11);
79: uik[9] = -(d[3]*u6 + d[9]*u7 + d[15]*u8 + d[21]*u9 + d[27]*u10 + d[33]*u11);
80: uik[10]= -(d[4]*u6+ d[10]*u7 + d[16]*u8 + d[22]*u9 + d[28]*u10 + d[34]*u11);
81: uik[11]= -(d[5]*u6+ d[11]*u7 + d[17]*u8 + d[23]*u9 + d[29]*u10 + d[35]*u11);
83: uik[12] = -(d[0]*u12 + d[6]*u13 + d[12]*u14 + d[18]*u15 + d[24]*u16 + d[30]*u17);
84: uik[13] = -(d[1]*u12 + d[7]*u13 + d[13]*u14 + d[19]*u15 + d[25]*u16 + d[31]*u17);
85: uik[14] = -(d[2]*u12 + d[8]*u13 + d[14]*u14 + d[20]*u15 + d[26]*u16 + d[32]*u17);
86: uik[15] = -(d[3]*u12 + d[9]*u13 + d[15]*u14 + d[21]*u15 + d[27]*u16 + d[33]*u17);
87: uik[16] = -(d[4]*u12+ d[10]*u13 + d[16]*u14 + d[22]*u15 + d[28]*u16 + d[34]*u17);
88: uik[17] = -(d[5]*u12+ d[11]*u13 + d[17]*u14 + d[23]*u15 + d[29]*u16 + d[35]*u17);
90: uik[18] = -(d[0]*u18 + d[6]*u19 + d[12]*u20 + d[18]*u21 + d[24]*u22 + d[30]*u23);
91: uik[19] = -(d[1]*u18 + d[7]*u19 + d[13]*u20 + d[19]*u21 + d[25]*u22 + d[31]*u23);
92: uik[20] = -(d[2]*u18 + d[8]*u19 + d[14]*u20 + d[20]*u21 + d[26]*u22 + d[32]*u23);
93: uik[21] = -(d[3]*u18 + d[9]*u19 + d[15]*u20 + d[21]*u21 + d[27]*u22 + d[33]*u23);
94: uik[22] = -(d[4]*u18+ d[10]*u19 + d[16]*u20 + d[22]*u21 + d[28]*u22 + d[34]*u23);
95: uik[23] = -(d[5]*u18+ d[11]*u19 + d[17]*u20 + d[23]*u21 + d[29]*u22 + d[35]*u23);
97: uik[24] = -(d[0]*u24 + d[6]*u25 + d[12]*u26 + d[18]*u27 + d[24]*u28 + d[30]*u29);
98: uik[25] = -(d[1]*u24 + d[7]*u25 + d[13]*u26 + d[19]*u27 + d[25]*u28 + d[31]*u29);
99: uik[26] = -(d[2]*u24 + d[8]*u25 + d[14]*u26 + d[20]*u27 + d[26]*u28 + d[32]*u29);
100: uik[27] = -(d[3]*u24 + d[9]*u25 + d[15]*u26 + d[21]*u27 + d[27]*u28 + d[33]*u29);
101: uik[28] = -(d[4]*u24+ d[10]*u25 + d[16]*u26 + d[22]*u27 + d[28]*u28 + d[34]*u29);
102: uik[29] = -(d[5]*u24+ d[11]*u25 + d[17]*u26 + d[23]*u27 + d[29]*u28 + d[35]*u29);
104: uik[30] = -(d[0]*u30 + d[6]*u31 + d[12]*u32 + d[18]*u33 + d[24]*u34 + d[30]*u35);
105: uik[31] = -(d[1]*u30 + d[7]*u31 + d[13]*u32 + d[19]*u33 + d[25]*u34 + d[31]*u35);
106: uik[32] = -(d[2]*u30 + d[8]*u31 + d[14]*u32 + d[20]*u33 + d[26]*u34 + d[32]*u35);
107: uik[33] = -(d[3]*u30 + d[9]*u31 + d[15]*u32 + d[21]*u33 + d[27]*u34 + d[33]*u35);
108: uik[34] = -(d[4]*u30+ d[10]*u31 + d[16]*u32 + d[22]*u33 + d[28]*u34 + d[34]*u35);
109: uik[35] = -(d[5]*u30+ d[11]*u31 + d[17]*u32 + d[23]*u33 + d[29]*u34 + d[35]*u35);
111: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
112: dk[0] += uik[0]*u0 + uik[1]*u1 + uik[2]*u2 + uik[3]*u3 + uik[4]*u4 + uik[5]*u5;
113: dk[1] += uik[6]*u0 + uik[7]*u1 + uik[8]*u2 + uik[9]*u3+ uik[10]*u4+ uik[11]*u5;
114: dk[2] += uik[12]*u0+ uik[13]*u1+ uik[14]*u2+ uik[15]*u3+ uik[16]*u4+ uik[17]*u5;
115: dk[3] += uik[18]*u0+ uik[19]*u1+ uik[20]*u2+ uik[21]*u3+ uik[22]*u4+ uik[23]*u5;
116: dk[4] += uik[24]*u0+ uik[25]*u1+ uik[26]*u2+ uik[27]*u3+ uik[28]*u4+ uik[29]*u5;
117: dk[5] += uik[30]*u0+ uik[31]*u1+ uik[32]*u2+ uik[33]*u3+ uik[34]*u4+ uik[35]*u5;
119: dk[6] += uik[0]*u6 + uik[1]*u7 + uik[2]*u8 + uik[3]*u9 + uik[4]*u10 + uik[5]*u11;
120: dk[7] += uik[6]*u6 + uik[7]*u7 + uik[8]*u8 + uik[9]*u9+ uik[10]*u10+ uik[11]*u11;
121: dk[8] += uik[12]*u6+ uik[13]*u7+ uik[14]*u8+ uik[15]*u9+ uik[16]*u10+ uik[17]*u11;
122: dk[9] += uik[18]*u6+ uik[19]*u7+ uik[20]*u8+ uik[21]*u9+ uik[22]*u10+ uik[23]*u11;
123: dk[10]+= uik[24]*u6+ uik[25]*u7+ uik[26]*u8+ uik[27]*u9+ uik[28]*u10+ uik[29]*u11;
124: dk[11]+= uik[30]*u6+ uik[31]*u7+ uik[32]*u8+ uik[33]*u9+ uik[34]*u10+ uik[35]*u11;
126: dk[12]+= uik[0]*u12 + uik[1]*u13 + uik[2]*u14 + uik[3]*u15 + uik[4]*u16 + uik[5]*u17;
127: dk[13]+= uik[6]*u12 + uik[7]*u13 + uik[8]*u14 + uik[9]*u15+ uik[10]*u16+ uik[11]*u17;
128: dk[14]+= uik[12]*u12+ uik[13]*u13+ uik[14]*u14+ uik[15]*u15+ uik[16]*u16+ uik[17]*u17;
129: dk[15]+= uik[18]*u12+ uik[19]*u13+ uik[20]*u14+ uik[21]*u15+ uik[22]*u16+ uik[23]*u17;
130: dk[16]+= uik[24]*u12+ uik[25]*u13+ uik[26]*u14+ uik[27]*u15+ uik[28]*u16+ uik[29]*u17;
131: dk[17]+= uik[30]*u12+ uik[31]*u13+ uik[32]*u14+ uik[33]*u15+ uik[34]*u16+ uik[35]*u17;
133: dk[18]+= uik[0]*u18 + uik[1]*u19 + uik[2]*u20 + uik[3]*u21 + uik[4]*u22 + uik[5]*u23;
134: dk[19]+= uik[6]*u18 + uik[7]*u19 + uik[8]*u20 + uik[9]*u21+ uik[10]*u22+ uik[11]*u23;
135: dk[20]+= uik[12]*u18+ uik[13]*u19+ uik[14]*u20+ uik[15]*u21+ uik[16]*u22+ uik[17]*u23;
136: dk[21]+= uik[18]*u18+ uik[19]*u19+ uik[20]*u20+ uik[21]*u21+ uik[22]*u22+ uik[23]*u23;
137: dk[22]+= uik[24]*u18+ uik[25]*u19+ uik[26]*u20+ uik[27]*u21+ uik[28]*u22+ uik[29]*u23;
138: dk[23]+= uik[30]*u18+ uik[31]*u19+ uik[32]*u20+ uik[33]*u21+ uik[34]*u22+ uik[35]*u23;
140: dk[24]+= uik[0]*u24 + uik[1]*u25 + uik[2]*u26 + uik[3]*u27 + uik[4]*u28 + uik[5]*u29;
141: dk[25]+= uik[6]*u24 + uik[7]*u25 + uik[8]*u26 + uik[9]*u27+ uik[10]*u28+ uik[11]*u29;
142: dk[26]+= uik[12]*u24+ uik[13]*u25+ uik[14]*u26+ uik[15]*u27+ uik[16]*u28+ uik[17]*u29;
143: dk[27]+= uik[18]*u24+ uik[19]*u25+ uik[20]*u26+ uik[21]*u27+ uik[22]*u28+ uik[23]*u29;
144: dk[28]+= uik[24]*u24+ uik[25]*u25+ uik[26]*u26+ uik[27]*u27+ uik[28]*u28+ uik[29]*u29;
145: dk[29]+= uik[30]*u24+ uik[31]*u25+ uik[32]*u26+ uik[33]*u27+ uik[34]*u28+ uik[35]*u29;
147: dk[30]+= uik[0]*u30 + uik[1]*u31 + uik[2]*u32 + uik[3]*u33 + uik[4]*u34 + uik[5]*u35;
148: dk[31]+= uik[6]*u30 + uik[7]*u31 + uik[8]*u32 + uik[9]*u33+ uik[10]*u34+ uik[11]*u35;
149: dk[32]+= uik[12]*u30+ uik[13]*u31+ uik[14]*u32+ uik[15]*u33+ uik[16]*u34+ uik[17]*u35;
150: dk[33]+= uik[18]*u30+ uik[19]*u31+ uik[20]*u32+ uik[21]*u33+ uik[22]*u34+ uik[23]*u35;
151: dk[34]+= uik[24]*u30+ uik[25]*u31+ uik[26]*u32+ uik[27]*u33+ uik[28]*u34+ uik[29]*u35;
152: dk[35]+= uik[30]*u30+ uik[31]*u31+ uik[32]*u32+ uik[33]*u33+ uik[34]*u34+ uik[35]*u35;
153:
154: /* update -U(i,k) */
155: PetscMemcpy(ba+ili*36,uik,36*sizeof(MatScalar));
157: /* add multiple of row i to k-th row ... */
158: jmin = ili + 1; jmax = bi[i+1];
159: if (jmin < jmax){
160: for (j=jmin; j<jmax; j++) {
161: /* w += -U(i,k)^T * U_bar(i,j) */
162: wp = w + bj[j]*36;
163: u = ba + j*36;
165: u0 = u[0]; u1 = u[1]; u2 = u[2]; u3 = u[3]; u4 = u[4]; u5 = u[5]; u6 = u[6];
166: u7 = u[7]; u8 = u[8]; u9 = u[9]; u10 = u[10]; u11 = u[11]; u12 = u[12]; u13 = u[13];
167: u14 = u[14]; u15 = u[15]; u16 = u[16]; u17 = u[17]; u18 = u[18]; u19 = u[19]; u20 = u[20];
168: u21 = u[21]; u22 = u[22]; u23 = u[23]; u24 = u[24]; u25 = u[25]; u26 = u[26]; u27 = u[27];
169: u28 = u[28]; u29 = u[29]; u30 = u[30]; u31 = u[31]; u32 = u[32]; u33 = u[33]; u34 = u[34];
170: u35 = u[35];
172: wp[0] += uik[0]*u0 + uik[1]*u1 + uik[2]*u2 + uik[3]*u3 + uik[4]*u4 + uik[5]*u5;
173: wp[1] += uik[6]*u0 + uik[7]*u1 + uik[8]*u2 + uik[9]*u3+ uik[10]*u4+ uik[11]*u5;
174: wp[2] += uik[12]*u0+ uik[13]*u1+ uik[14]*u2+ uik[15]*u3+ uik[16]*u4+ uik[17]*u5;
175: wp[3] += uik[18]*u0+ uik[19]*u1+ uik[20]*u2+ uik[21]*u3+ uik[22]*u4+ uik[23]*u5;
176: wp[4] += uik[24]*u0+ uik[25]*u1+ uik[26]*u2+ uik[27]*u3+ uik[28]*u4+ uik[29]*u5;
177: wp[5] += uik[30]*u0+ uik[31]*u1+ uik[32]*u2+ uik[33]*u3+ uik[34]*u4+ uik[35]*u5;
179: wp[6] += uik[0]*u6 + uik[1]*u7 + uik[2]*u8 + uik[3]*u9 + uik[4]*u10 + uik[5]*u11;
180: wp[7] += uik[6]*u6 + uik[7]*u7 + uik[8]*u8 + uik[9]*u9+ uik[10]*u10+ uik[11]*u11;
181: wp[8] += uik[12]*u6+ uik[13]*u7+ uik[14]*u8+ uik[15]*u9+ uik[16]*u10+ uik[17]*u11;
182: wp[9] += uik[18]*u6+ uik[19]*u7+ uik[20]*u8+ uik[21]*u9+ uik[22]*u10+ uik[23]*u11;
183: wp[10]+= uik[24]*u6+ uik[25]*u7+ uik[26]*u8+ uik[27]*u9+ uik[28]*u10+ uik[29]*u11;
184: wp[11]+= uik[30]*u6+ uik[31]*u7+ uik[32]*u8+ uik[33]*u9+ uik[34]*u10+ uik[35]*u11;
186: wp[12]+= uik[0]*u12 + uik[1]*u13 + uik[2]*u14 + uik[3]*u15 + uik[4]*u16 + uik[5]*u17;
187: wp[13]+= uik[6]*u12 + uik[7]*u13 + uik[8]*u14 + uik[9]*u15+ uik[10]*u16+ uik[11]*u17;
188: wp[14]+= uik[12]*u12+ uik[13]*u13+ uik[14]*u14+ uik[15]*u15+ uik[16]*u16+ uik[17]*u17;
189: wp[15]+= uik[18]*u12+ uik[19]*u13+ uik[20]*u14+ uik[21]*u15+ uik[22]*u16+ uik[23]*u17;
190: wp[16]+= uik[24]*u12+ uik[25]*u13+ uik[26]*u14+ uik[27]*u15+ uik[28]*u16+ uik[29]*u17;
191: wp[17]+= uik[30]*u12+ uik[31]*u13+ uik[32]*u14+ uik[33]*u15+ uik[34]*u16+ uik[35]*u17;
193: wp[18]+= uik[0]*u18 + uik[1]*u19 + uik[2]*u20 + uik[3]*u21 + uik[4]*u22 + uik[5]*u23;
194: wp[19]+= uik[6]*u18 + uik[7]*u19 + uik[8]*u20 + uik[9]*u21+ uik[10]*u22+ uik[11]*u23;
195: wp[20]+= uik[12]*u18+ uik[13]*u19+ uik[14]*u20+ uik[15]*u21+ uik[16]*u22+ uik[17]*u23;
196: wp[21]+= uik[18]*u18+ uik[19]*u19+ uik[20]*u20+ uik[21]*u21+ uik[22]*u22+ uik[23]*u23;
197: wp[22]+= uik[24]*u18+ uik[25]*u19+ uik[26]*u20+ uik[27]*u21+ uik[28]*u22+ uik[29]*u23;
198: wp[23]+= uik[30]*u18+ uik[31]*u19+ uik[32]*u20+ uik[33]*u21+ uik[34]*u22+ uik[35]*u23;
200: wp[24]+= uik[0]*u24 + uik[1]*u25 + uik[2]*u26 + uik[3]*u27 + uik[4]*u28 + uik[5]*u29;
201: wp[25]+= uik[6]*u24 + uik[7]*u25 + uik[8]*u26 + uik[9]*u27+ uik[10]*u28+ uik[11]*u29;
202: wp[26]+= uik[12]*u24+ uik[13]*u25+ uik[14]*u26+ uik[15]*u27+ uik[16]*u28+ uik[17]*u29;
203: wp[27]+= uik[18]*u24+ uik[19]*u25+ uik[20]*u26+ uik[21]*u27+ uik[22]*u28+ uik[23]*u29;
204: wp[28]+= uik[24]*u24+ uik[25]*u25+ uik[26]*u26+ uik[27]*u27+ uik[28]*u28+ uik[29]*u29;
205: wp[29]+= uik[30]*u24+ uik[31]*u25+ uik[32]*u26+ uik[33]*u27+ uik[34]*u28+ uik[35]*u29;
207: wp[30]+= uik[0]*u30 + uik[1]*u31 + uik[2]*u32 + uik[3]*u33 + uik[4]*u34 + uik[5]*u35;
208: wp[31]+= uik[6]*u30 + uik[7]*u31 + uik[8]*u32 + uik[9]*u33+ uik[10]*u34+ uik[11]*u35;
209: wp[32]+= uik[12]*u30+ uik[13]*u31+ uik[14]*u32+ uik[15]*u33+ uik[16]*u34+ uik[17]*u35;
210: wp[33]+= uik[18]*u30+ uik[19]*u31+ uik[20]*u32+ uik[21]*u33+ uik[22]*u34+ uik[23]*u35;
211: wp[34]+= uik[24]*u30+ uik[25]*u31+ uik[26]*u32+ uik[27]*u33+ uik[28]*u34+ uik[29]*u35;
212: wp[35]+= uik[30]*u30+ uik[31]*u31+ uik[32]*u32+ uik[33]*u33+ uik[34]*u34+ uik[35]*u35;
213: }
214:
215: /* ... add i to row list for next nonzero entry */
216: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
217: j = bj[jmin];
218: jl[i] = jl[j]; jl[j] = i; /* update jl */
219: }
220: i = nexti;
221: }
223: /* save nonzero entries in k-th row of U ... */
225: /* invert diagonal block */
226: d = ba+k*36;
227: PetscMemcpy(d,dk,36*sizeof(MatScalar));
228: Kernel_A_gets_inverse_A_6(d);
229:
230: jmin = bi[k]; jmax = bi[k+1];
231: if (jmin < jmax) {
232: for (j=jmin; j<jmax; j++){
233: vj = bj[j]; /* block col. index of U */
234: u = ba + j*36;
235: wp = w + vj*36;
236: for (k1=0; k1<36; k1++){
237: *u++ = *wp;
238: *wp++ = 0.0;
239: }
240: }
241:
242: /* ... add k to row list for first nonzero entry in k-th row */
243: il[k] = jmin;
244: i = bj[jmin];
245: jl[k] = jl[i]; jl[i] = k;
246: }
247: }
249: PetscFree(w);
250: PetscFree(il);
251: PetscFree(dk);
253: C->factor = FACTOR_CHOLESKY;
254: C->assembled = PETSC_TRUE;
255: C->preallocated = PETSC_TRUE;
256: PetscLogFlops(1.3333*216*b->mbs); /* from inverting diagonal blocks */
257: return(0);
258: }