Actual source code: sbaijfact9.c
2: #include <../src/mat/impls/sbaij/seq/sbaij.h>
3: #include <petsc/private/kernels/blockinvert.h>
5: /* Version for when blocks are 6 by 6 */
6: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat C,Mat A,const MatFactorInfo *info)
7: {
8: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
9: IS perm = b->row;
10: const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
11: PetscInt i,j,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
12: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
13: MatScalar *u,*d,*w,*wp,u0,u1,u2,u3,u4,u5,u6,u7,u8,u9,u10,u11,u12;
14: MatScalar u13,u14,u15,u16,u17,u18,u19,u20,u21,u22,u23,u24,u25,u26,u27;
15: MatScalar u28,u29,u30,u31,u32,u33,u34,u35;
16: PetscReal shift = info->shiftamount;
17: PetscBool allowzeropivot,zeropivotdetected;
19: /* initialization */
20: allowzeropivot = PetscNot(A->erroriffailure);
21: PetscCalloc1(36*mbs,&w);
22: PetscMalloc2(mbs,&il,mbs,&jl);
23: il[0] = 0;
24: for (i=0; i<mbs; i++) jl[i] = mbs;
26: PetscMalloc2(36,&dk,36,&uik);
27: ISGetIndices(perm,&perm_ptr);
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: PetscMalloc1(36*ai[mbs],&aa);
35: PetscArraycpy(aa,a->a,36*ai[mbs]);
36: PetscMalloc1(ai[mbs],&a2anew);
37: PetscArraycpy(a2anew,a->a2anew,ai[mbs]);
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<36; k1++) {
45: dk[k1] = aa[k*36+k1];
46: aa[k*36+k1] = aa[j*36+k1];
47: aa[j*36+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*36; /* ptr to the beginning of j-th block of aa */
54: for (k=0; k<36; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
55: for (k=0; k<6; k++) { /* j-th block of aa <- dk^T */
56: for (k1=0; k1<6; k1++) *ap++ = dk[k + 6*k1];
57: }
58: }
59: }
60: }
61: PetscFree(a2anew);
62: }
64: /* for each row k */
65: for (k = 0; k<mbs; k++) {
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*36;
71: for (j = jmin; j < jmax; j++) {
72: vj = perm_ptr[aj[j]]; /* block col. index */
73: wp = w + vj*36;
74: for (i=0; i<36; i++) *wp++ = *ap++;
75: }
76: }
78: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
79: PetscArraycpy(dk,w+k*36,36);
80: i = jl[k]; /* first row to be added to k_th row */
82: while (i < mbs) {
83: nexti = jl[i]; /* next row to be added to k_th row */
85: /* compute multiplier */
86: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
88: /* uik = -inv(Di)*U_bar(i,k) */
89: d = ba + i*36;
90: u = ba + ili*36;
92: u0 = u[0]; u1 = u[1]; u2 = u[2]; u3 = u[3]; u4 = u[4]; u5 = u[5]; u6 = u[6];
93: u7 = u[7]; u8 = u[8]; u9 = u[9]; u10 = u[10]; u11 = u[11]; u12 = u[12]; u13 = u[13];
94: u14 = u[14]; u15 = u[15]; u16 = u[16]; u17 = u[17]; u18 = u[18]; u19 = u[19]; u20 = u[20];
95: u21 = u[21]; u22 = u[22]; u23 = u[23]; u24 = u[24]; u25 = u[25]; u26 = u[26]; u27 = u[27];
96: u28 = u[28]; u29 = u[29]; u30 = u[30]; u31 = u[31]; u32 = u[32]; u33 = u[33]; u34 = u[34];
97: u35 = u[35];
99: uik[0] = -(d[0]*u0 + d[6]*u1 + d[12]*u2 + d[18]*u3 + d[24]*u4 + d[30]*u5);
100: uik[1] = -(d[1]*u0 + d[7]*u1 + d[13]*u2 + d[19]*u3 + d[25]*u4 + d[31]*u5);
101: uik[2] = -(d[2]*u0 + d[8]*u1 + d[14]*u2 + d[20]*u3 + d[26]*u4 + d[32]*u5);
102: uik[3] = -(d[3]*u0 + d[9]*u1 + d[15]*u2 + d[21]*u3 + d[27]*u4 + d[33]*u5);
103: uik[4] = -(d[4]*u0+ d[10]*u1 + d[16]*u2 + d[22]*u3 + d[28]*u4 + d[34]*u5);
104: uik[5] = -(d[5]*u0+ d[11]*u1 + d[17]*u2 + d[23]*u3 + d[29]*u4 + d[35]*u5);
106: uik[6] = -(d[0]*u6 + d[6]*u7 + d[12]*u8 + d[18]*u9 + d[24]*u10 + d[30]*u11);
107: uik[7] = -(d[1]*u6 + d[7]*u7 + d[13]*u8 + d[19]*u9 + d[25]*u10 + d[31]*u11);
108: uik[8] = -(d[2]*u6 + d[8]*u7 + d[14]*u8 + d[20]*u9 + d[26]*u10 + d[32]*u11);
109: uik[9] = -(d[3]*u6 + d[9]*u7 + d[15]*u8 + d[21]*u9 + d[27]*u10 + d[33]*u11);
110: uik[10]= -(d[4]*u6+ d[10]*u7 + d[16]*u8 + d[22]*u9 + d[28]*u10 + d[34]*u11);
111: uik[11]= -(d[5]*u6+ d[11]*u7 + d[17]*u8 + d[23]*u9 + d[29]*u10 + d[35]*u11);
113: uik[12] = -(d[0]*u12 + d[6]*u13 + d[12]*u14 + d[18]*u15 + d[24]*u16 + d[30]*u17);
114: uik[13] = -(d[1]*u12 + d[7]*u13 + d[13]*u14 + d[19]*u15 + d[25]*u16 + d[31]*u17);
115: uik[14] = -(d[2]*u12 + d[8]*u13 + d[14]*u14 + d[20]*u15 + d[26]*u16 + d[32]*u17);
116: uik[15] = -(d[3]*u12 + d[9]*u13 + d[15]*u14 + d[21]*u15 + d[27]*u16 + d[33]*u17);
117: uik[16] = -(d[4]*u12+ d[10]*u13 + d[16]*u14 + d[22]*u15 + d[28]*u16 + d[34]*u17);
118: uik[17] = -(d[5]*u12+ d[11]*u13 + d[17]*u14 + d[23]*u15 + d[29]*u16 + d[35]*u17);
120: uik[18] = -(d[0]*u18 + d[6]*u19 + d[12]*u20 + d[18]*u21 + d[24]*u22 + d[30]*u23);
121: uik[19] = -(d[1]*u18 + d[7]*u19 + d[13]*u20 + d[19]*u21 + d[25]*u22 + d[31]*u23);
122: uik[20] = -(d[2]*u18 + d[8]*u19 + d[14]*u20 + d[20]*u21 + d[26]*u22 + d[32]*u23);
123: uik[21] = -(d[3]*u18 + d[9]*u19 + d[15]*u20 + d[21]*u21 + d[27]*u22 + d[33]*u23);
124: uik[22] = -(d[4]*u18+ d[10]*u19 + d[16]*u20 + d[22]*u21 + d[28]*u22 + d[34]*u23);
125: uik[23] = -(d[5]*u18+ d[11]*u19 + d[17]*u20 + d[23]*u21 + d[29]*u22 + d[35]*u23);
127: uik[24] = -(d[0]*u24 + d[6]*u25 + d[12]*u26 + d[18]*u27 + d[24]*u28 + d[30]*u29);
128: uik[25] = -(d[1]*u24 + d[7]*u25 + d[13]*u26 + d[19]*u27 + d[25]*u28 + d[31]*u29);
129: uik[26] = -(d[2]*u24 + d[8]*u25 + d[14]*u26 + d[20]*u27 + d[26]*u28 + d[32]*u29);
130: uik[27] = -(d[3]*u24 + d[9]*u25 + d[15]*u26 + d[21]*u27 + d[27]*u28 + d[33]*u29);
131: uik[28] = -(d[4]*u24+ d[10]*u25 + d[16]*u26 + d[22]*u27 + d[28]*u28 + d[34]*u29);
132: uik[29] = -(d[5]*u24+ d[11]*u25 + d[17]*u26 + d[23]*u27 + d[29]*u28 + d[35]*u29);
134: uik[30] = -(d[0]*u30 + d[6]*u31 + d[12]*u32 + d[18]*u33 + d[24]*u34 + d[30]*u35);
135: uik[31] = -(d[1]*u30 + d[7]*u31 + d[13]*u32 + d[19]*u33 + d[25]*u34 + d[31]*u35);
136: uik[32] = -(d[2]*u30 + d[8]*u31 + d[14]*u32 + d[20]*u33 + d[26]*u34 + d[32]*u35);
137: uik[33] = -(d[3]*u30 + d[9]*u31 + d[15]*u32 + d[21]*u33 + d[27]*u34 + d[33]*u35);
138: uik[34] = -(d[4]*u30+ d[10]*u31 + d[16]*u32 + d[22]*u33 + d[28]*u34 + d[34]*u35);
139: uik[35] = -(d[5]*u30+ d[11]*u31 + d[17]*u32 + d[23]*u33 + d[29]*u34 + d[35]*u35);
141: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
142: dk[0] += uik[0]*u0 + uik[1]*u1 + uik[2]*u2 + uik[3]*u3 + uik[4]*u4 + uik[5]*u5;
143: dk[1] += uik[6]*u0 + uik[7]*u1 + uik[8]*u2 + uik[9]*u3+ uik[10]*u4+ uik[11]*u5;
144: dk[2] += uik[12]*u0+ uik[13]*u1+ uik[14]*u2+ uik[15]*u3+ uik[16]*u4+ uik[17]*u5;
145: dk[3] += uik[18]*u0+ uik[19]*u1+ uik[20]*u2+ uik[21]*u3+ uik[22]*u4+ uik[23]*u5;
146: dk[4] += uik[24]*u0+ uik[25]*u1+ uik[26]*u2+ uik[27]*u3+ uik[28]*u4+ uik[29]*u5;
147: dk[5] += uik[30]*u0+ uik[31]*u1+ uik[32]*u2+ uik[33]*u3+ uik[34]*u4+ uik[35]*u5;
149: dk[6] += uik[0]*u6 + uik[1]*u7 + uik[2]*u8 + uik[3]*u9 + uik[4]*u10 + uik[5]*u11;
150: dk[7] += uik[6]*u6 + uik[7]*u7 + uik[8]*u8 + uik[9]*u9+ uik[10]*u10+ uik[11]*u11;
151: dk[8] += uik[12]*u6+ uik[13]*u7+ uik[14]*u8+ uik[15]*u9+ uik[16]*u10+ uik[17]*u11;
152: dk[9] += uik[18]*u6+ uik[19]*u7+ uik[20]*u8+ uik[21]*u9+ uik[22]*u10+ uik[23]*u11;
153: dk[10]+= uik[24]*u6+ uik[25]*u7+ uik[26]*u8+ uik[27]*u9+ uik[28]*u10+ uik[29]*u11;
154: dk[11]+= uik[30]*u6+ uik[31]*u7+ uik[32]*u8+ uik[33]*u9+ uik[34]*u10+ uik[35]*u11;
156: dk[12]+= uik[0]*u12 + uik[1]*u13 + uik[2]*u14 + uik[3]*u15 + uik[4]*u16 + uik[5]*u17;
157: dk[13]+= uik[6]*u12 + uik[7]*u13 + uik[8]*u14 + uik[9]*u15+ uik[10]*u16+ uik[11]*u17;
158: dk[14]+= uik[12]*u12+ uik[13]*u13+ uik[14]*u14+ uik[15]*u15+ uik[16]*u16+ uik[17]*u17;
159: dk[15]+= uik[18]*u12+ uik[19]*u13+ uik[20]*u14+ uik[21]*u15+ uik[22]*u16+ uik[23]*u17;
160: dk[16]+= uik[24]*u12+ uik[25]*u13+ uik[26]*u14+ uik[27]*u15+ uik[28]*u16+ uik[29]*u17;
161: dk[17]+= uik[30]*u12+ uik[31]*u13+ uik[32]*u14+ uik[33]*u15+ uik[34]*u16+ uik[35]*u17;
163: dk[18]+= uik[0]*u18 + uik[1]*u19 + uik[2]*u20 + uik[3]*u21 + uik[4]*u22 + uik[5]*u23;
164: dk[19]+= uik[6]*u18 + uik[7]*u19 + uik[8]*u20 + uik[9]*u21+ uik[10]*u22+ uik[11]*u23;
165: dk[20]+= uik[12]*u18+ uik[13]*u19+ uik[14]*u20+ uik[15]*u21+ uik[16]*u22+ uik[17]*u23;
166: dk[21]+= uik[18]*u18+ uik[19]*u19+ uik[20]*u20+ uik[21]*u21+ uik[22]*u22+ uik[23]*u23;
167: dk[22]+= uik[24]*u18+ uik[25]*u19+ uik[26]*u20+ uik[27]*u21+ uik[28]*u22+ uik[29]*u23;
168: dk[23]+= uik[30]*u18+ uik[31]*u19+ uik[32]*u20+ uik[33]*u21+ uik[34]*u22+ uik[35]*u23;
170: dk[24]+= uik[0]*u24 + uik[1]*u25 + uik[2]*u26 + uik[3]*u27 + uik[4]*u28 + uik[5]*u29;
171: dk[25]+= uik[6]*u24 + uik[7]*u25 + uik[8]*u26 + uik[9]*u27+ uik[10]*u28+ uik[11]*u29;
172: dk[26]+= uik[12]*u24+ uik[13]*u25+ uik[14]*u26+ uik[15]*u27+ uik[16]*u28+ uik[17]*u29;
173: dk[27]+= uik[18]*u24+ uik[19]*u25+ uik[20]*u26+ uik[21]*u27+ uik[22]*u28+ uik[23]*u29;
174: dk[28]+= uik[24]*u24+ uik[25]*u25+ uik[26]*u26+ uik[27]*u27+ uik[28]*u28+ uik[29]*u29;
175: dk[29]+= uik[30]*u24+ uik[31]*u25+ uik[32]*u26+ uik[33]*u27+ uik[34]*u28+ uik[35]*u29;
177: dk[30]+= uik[0]*u30 + uik[1]*u31 + uik[2]*u32 + uik[3]*u33 + uik[4]*u34 + uik[5]*u35;
178: dk[31]+= uik[6]*u30 + uik[7]*u31 + uik[8]*u32 + uik[9]*u33+ uik[10]*u34+ uik[11]*u35;
179: dk[32]+= uik[12]*u30+ uik[13]*u31+ uik[14]*u32+ uik[15]*u33+ uik[16]*u34+ uik[17]*u35;
180: dk[33]+= uik[18]*u30+ uik[19]*u31+ uik[20]*u32+ uik[21]*u33+ uik[22]*u34+ uik[23]*u35;
181: dk[34]+= uik[24]*u30+ uik[25]*u31+ uik[26]*u32+ uik[27]*u33+ uik[28]*u34+ uik[29]*u35;
182: dk[35]+= uik[30]*u30+ uik[31]*u31+ uik[32]*u32+ uik[33]*u33+ uik[34]*u34+ uik[35]*u35;
184: PetscLogFlops(216.0*4.0);
186: /* update -U(i,k) */
187: PetscArraycpy(ba+ili*36,uik,36);
189: /* add multiple of row i to k-th row ... */
190: jmin = ili + 1; jmax = bi[i+1];
191: if (jmin < jmax) {
192: for (j=jmin; j<jmax; j++) {
193: /* w += -U(i,k)^T * U_bar(i,j) */
194: wp = w + bj[j]*36;
195: u = ba + j*36;
197: u0 = u[0]; u1 = u[1]; u2 = u[2]; u3 = u[3]; u4 = u[4]; u5 = u[5]; u6 = u[6];
198: u7 = u[7]; u8 = u[8]; u9 = u[9]; u10 = u[10]; u11 = u[11]; u12 = u[12]; u13 = u[13];
199: u14 = u[14]; u15 = u[15]; u16 = u[16]; u17 = u[17]; u18 = u[18]; u19 = u[19]; u20 = u[20];
200: u21 = u[21]; u22 = u[22]; u23 = u[23]; u24 = u[24]; u25 = u[25]; u26 = u[26]; u27 = u[27];
201: u28 = u[28]; u29 = u[29]; u30 = u[30]; u31 = u[31]; u32 = u[32]; u33 = u[33]; u34 = u[34];
202: u35 = u[35];
204: wp[0] += uik[0]*u0 + uik[1]*u1 + uik[2]*u2 + uik[3]*u3 + uik[4]*u4 + uik[5]*u5;
205: wp[1] += uik[6]*u0 + uik[7]*u1 + uik[8]*u2 + uik[9]*u3+ uik[10]*u4+ uik[11]*u5;
206: wp[2] += uik[12]*u0+ uik[13]*u1+ uik[14]*u2+ uik[15]*u3+ uik[16]*u4+ uik[17]*u5;
207: wp[3] += uik[18]*u0+ uik[19]*u1+ uik[20]*u2+ uik[21]*u3+ uik[22]*u4+ uik[23]*u5;
208: wp[4] += uik[24]*u0+ uik[25]*u1+ uik[26]*u2+ uik[27]*u3+ uik[28]*u4+ uik[29]*u5;
209: wp[5] += uik[30]*u0+ uik[31]*u1+ uik[32]*u2+ uik[33]*u3+ uik[34]*u4+ uik[35]*u5;
211: wp[6] += uik[0]*u6 + uik[1]*u7 + uik[2]*u8 + uik[3]*u9 + uik[4]*u10 + uik[5]*u11;
212: wp[7] += uik[6]*u6 + uik[7]*u7 + uik[8]*u8 + uik[9]*u9+ uik[10]*u10+ uik[11]*u11;
213: wp[8] += uik[12]*u6+ uik[13]*u7+ uik[14]*u8+ uik[15]*u9+ uik[16]*u10+ uik[17]*u11;
214: wp[9] += uik[18]*u6+ uik[19]*u7+ uik[20]*u8+ uik[21]*u9+ uik[22]*u10+ uik[23]*u11;
215: wp[10]+= uik[24]*u6+ uik[25]*u7+ uik[26]*u8+ uik[27]*u9+ uik[28]*u10+ uik[29]*u11;
216: wp[11]+= uik[30]*u6+ uik[31]*u7+ uik[32]*u8+ uik[33]*u9+ uik[34]*u10+ uik[35]*u11;
218: wp[12]+= uik[0]*u12 + uik[1]*u13 + uik[2]*u14 + uik[3]*u15 + uik[4]*u16 + uik[5]*u17;
219: wp[13]+= uik[6]*u12 + uik[7]*u13 + uik[8]*u14 + uik[9]*u15+ uik[10]*u16+ uik[11]*u17;
220: wp[14]+= uik[12]*u12+ uik[13]*u13+ uik[14]*u14+ uik[15]*u15+ uik[16]*u16+ uik[17]*u17;
221: wp[15]+= uik[18]*u12+ uik[19]*u13+ uik[20]*u14+ uik[21]*u15+ uik[22]*u16+ uik[23]*u17;
222: wp[16]+= uik[24]*u12+ uik[25]*u13+ uik[26]*u14+ uik[27]*u15+ uik[28]*u16+ uik[29]*u17;
223: wp[17]+= uik[30]*u12+ uik[31]*u13+ uik[32]*u14+ uik[33]*u15+ uik[34]*u16+ uik[35]*u17;
225: wp[18]+= uik[0]*u18 + uik[1]*u19 + uik[2]*u20 + uik[3]*u21 + uik[4]*u22 + uik[5]*u23;
226: wp[19]+= uik[6]*u18 + uik[7]*u19 + uik[8]*u20 + uik[9]*u21+ uik[10]*u22+ uik[11]*u23;
227: wp[20]+= uik[12]*u18+ uik[13]*u19+ uik[14]*u20+ uik[15]*u21+ uik[16]*u22+ uik[17]*u23;
228: wp[21]+= uik[18]*u18+ uik[19]*u19+ uik[20]*u20+ uik[21]*u21+ uik[22]*u22+ uik[23]*u23;
229: wp[22]+= uik[24]*u18+ uik[25]*u19+ uik[26]*u20+ uik[27]*u21+ uik[28]*u22+ uik[29]*u23;
230: wp[23]+= uik[30]*u18+ uik[31]*u19+ uik[32]*u20+ uik[33]*u21+ uik[34]*u22+ uik[35]*u23;
232: wp[24]+= uik[0]*u24 + uik[1]*u25 + uik[2]*u26 + uik[3]*u27 + uik[4]*u28 + uik[5]*u29;
233: wp[25]+= uik[6]*u24 + uik[7]*u25 + uik[8]*u26 + uik[9]*u27+ uik[10]*u28+ uik[11]*u29;
234: wp[26]+= uik[12]*u24+ uik[13]*u25+ uik[14]*u26+ uik[15]*u27+ uik[16]*u28+ uik[17]*u29;
235: wp[27]+= uik[18]*u24+ uik[19]*u25+ uik[20]*u26+ uik[21]*u27+ uik[22]*u28+ uik[23]*u29;
236: wp[28]+= uik[24]*u24+ uik[25]*u25+ uik[26]*u26+ uik[27]*u27+ uik[28]*u28+ uik[29]*u29;
237: wp[29]+= uik[30]*u24+ uik[31]*u25+ uik[32]*u26+ uik[33]*u27+ uik[34]*u28+ uik[35]*u29;
239: wp[30]+= uik[0]*u30 + uik[1]*u31 + uik[2]*u32 + uik[3]*u33 + uik[4]*u34 + uik[5]*u35;
240: wp[31]+= uik[6]*u30 + uik[7]*u31 + uik[8]*u32 + uik[9]*u33+ uik[10]*u34+ uik[11]*u35;
241: wp[32]+= uik[12]*u30+ uik[13]*u31+ uik[14]*u32+ uik[15]*u33+ uik[16]*u34+ uik[17]*u35;
242: wp[33]+= uik[18]*u30+ uik[19]*u31+ uik[20]*u32+ uik[21]*u33+ uik[22]*u34+ uik[23]*u35;
243: wp[34]+= uik[24]*u30+ uik[25]*u31+ uik[26]*u32+ uik[27]*u33+ uik[28]*u34+ uik[29]*u35;
244: wp[35]+= uik[30]*u30+ uik[31]*u31+ uik[32]*u32+ uik[33]*u33+ uik[34]*u34+ uik[35]*u35;
245: }
246: PetscLogFlops(2.0*216.0*(jmax-jmin));
248: /* ... add i to row list for next nonzero entry */
249: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
250: j = bj[jmin];
251: jl[i] = jl[j]; jl[j] = i; /* update jl */
252: }
253: i = nexti;
254: }
256: /* save nonzero entries in k-th row of U ... */
258: /* invert diagonal block */
259: d = ba+k*36;
260: PetscArraycpy(d,dk,36);
261: PetscKernel_A_gets_inverse_A_6(d,shift,allowzeropivot,&zeropivotdetected);
262: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
264: jmin = bi[k]; jmax = bi[k+1];
265: if (jmin < jmax) {
266: for (j=jmin; j<jmax; j++) {
267: vj = bj[j]; /* block col. index of U */
268: u = ba + j*36;
269: wp = w + vj*36;
270: for (k1=0; k1<36; k1++) {
271: *u++ = *wp;
272: *wp++ = 0.0;
273: }
274: }
276: /* ... add k to row list for first nonzero entry in k-th row */
277: il[k] = jmin;
278: i = bj[jmin];
279: jl[k] = jl[i]; jl[i] = k;
280: }
281: }
283: PetscFree(w);
284: PetscFree2(il,jl);
285: PetscFree2(dk,uik);
286: if (a->permute) {
287: PetscFree(aa);
288: }
290: ISRestoreIndices(perm,&perm_ptr);
292: C->ops->solve = MatSolve_SeqSBAIJ_6_inplace;
293: C->ops->solvetranspose = MatSolve_SeqSBAIJ_6_inplace;
294: C->assembled = PETSC_TRUE;
295: C->preallocated = PETSC_TRUE;
297: PetscLogFlops(1.3333*216*b->mbs); /* from inverting diagonal blocks */
298: return 0;
299: }