Actual source code: baijfact7.c


  2: /*
  3:     Factorization code for BAIJ format.
  4: */
  5: #include <../src/mat/impls/baij/seq/baij.h>
  6: #include <petsc/private/kernels/blockinvert.h>

  8: /* ------------------------------------------------------------*/
  9: /*
 10:       Version for when blocks are 6 by 6
 11: */
 12: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_inplace(Mat C,Mat A,const MatFactorInfo *info)
 13: {
 14:   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
 15:   IS             isrow = b->row,isicol = b->icol;
 16:   const PetscInt *ajtmpold,*ajtmp,*diag_offset = b->diag,*r,*ic,*bi = b->i,*bj = b->j,*ai=a->i,*aj=a->j,*pj;
 17:   PetscInt       nz,row,i,j,n = a->mbs,idx;
 18:   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
 19:   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
 20:   MatScalar      p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
 21:   MatScalar      x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
 22:   MatScalar      p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
 23:   MatScalar      m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
 24:   MatScalar      p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
 25:   MatScalar      x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
 26:   MatScalar      m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
 27:   MatScalar      *ba   = b->a,*aa = a->a;
 28:   PetscReal      shift = info->shiftamount;
 29:   PetscBool      allowzeropivot,zeropivotdetected;

 31:   allowzeropivot = PetscNot(A->erroriffailure);
 32:   ISGetIndices(isrow,&r);
 33:   ISGetIndices(isicol,&ic);
 34:   PetscMalloc1(36*(n+1),&rtmp);

 36:   for (i=0; i<n; i++) {
 37:     nz    = bi[i+1] - bi[i];
 38:     ajtmp = bj + bi[i];
 39:     for  (j=0; j<nz; j++) {
 40:       x     = rtmp+36*ajtmp[j];
 41:       x[0]  = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
 42:       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
 43:       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
 44:       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
 45:       x[34] = x[35] = 0.0;
 46:     }
 47:     /* load in initial (unfactored row) */
 48:     idx      = r[i];
 49:     nz       = ai[idx+1] - ai[idx];
 50:     ajtmpold = aj + ai[idx];
 51:     v        = aa + 36*ai[idx];
 52:     for (j=0; j<nz; j++) {
 53:       x     = rtmp+36*ic[ajtmpold[j]];
 54:       x[0]  =  v[0];  x[1] =  v[1];  x[2] =  v[2];  x[3] =  v[3];
 55:       x[4]  =  v[4];  x[5] =  v[5];  x[6] =  v[6];  x[7] =  v[7];
 56:       x[8]  =  v[8];  x[9] =  v[9];  x[10] = v[10]; x[11] = v[11];
 57:       x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
 58:       x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
 59:       x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
 60:       x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
 61:       x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
 62:       x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
 63:       v    += 36;
 64:     }
 65:     row = *ajtmp++;
 66:     while (row < i) {
 67:       pc  =  rtmp + 36*row;
 68:       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
 69:       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];
 70:       p9  = pc[8];  p10 = pc[9];  p11 = pc[10]; p12 = pc[11];
 71:       p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
 72:       p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
 73:       p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
 74:       p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
 75:       p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
 76:       p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
 77:       if (p1  != 0.0 || p2  != 0.0 || p3  != 0.0 || p4  != 0.0 ||
 78:           p5  != 0.0 || p6  != 0.0 || p7  != 0.0 || p8  != 0.0 ||
 79:           p9  != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
 80:           p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
 81:           p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
 82:           p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
 83:           p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
 84:           p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
 85:           p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
 86:         pv    = ba + 36*diag_offset[row];
 87:         pj    = bj + diag_offset[row] + 1;
 88:         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
 89:         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
 90:         x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
 91:         x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
 92:         x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
 93:         x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
 94:         x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
 95:         x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
 96:         x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
 97:         pc[0] = m1  = p1*x1  + p7*x2   + p13*x3  + p19*x4  + p25*x5  + p31*x6;
 98:         pc[1] = m2  = p2*x1  + p8*x2   + p14*x3  + p20*x4  + p26*x5  + p32*x6;
 99:         pc[2] = m3  = p3*x1  + p9*x2   + p15*x3  + p21*x4  + p27*x5  + p33*x6;
100:         pc[3] = m4  = p4*x1  + p10*x2  + p16*x3  + p22*x4  + p28*x5  + p34*x6;
101:         pc[4] = m5  = p5*x1  + p11*x2  + p17*x3  + p23*x4  + p29*x5  + p35*x6;
102:         pc[5] = m6  = p6*x1  + p12*x2  + p18*x3  + p24*x4  + p30*x5  + p36*x6;

104:         pc[6]  = m7  = p1*x7  + p7*x8   + p13*x9  + p19*x10 + p25*x11 + p31*x12;
105:         pc[7]  = m8  = p2*x7  + p8*x8   + p14*x9  + p20*x10 + p26*x11 + p32*x12;
106:         pc[8]  = m9  = p3*x7  + p9*x8   + p15*x9  + p21*x10 + p27*x11 + p33*x12;
107:         pc[9]  = m10 = p4*x7  + p10*x8  + p16*x9  + p22*x10 + p28*x11 + p34*x12;
108:         pc[10] = m11 = p5*x7  + p11*x8  + p17*x9  + p23*x10 + p29*x11 + p35*x12;
109:         pc[11] = m12 = p6*x7  + p12*x8  + p18*x9  + p24*x10 + p30*x11 + p36*x12;

111:         pc[12] = m13 = p1*x13 + p7*x14  + p13*x15 + p19*x16 + p25*x17 + p31*x18;
112:         pc[13] = m14 = p2*x13 + p8*x14  + p14*x15 + p20*x16 + p26*x17 + p32*x18;
113:         pc[14] = m15 = p3*x13 + p9*x14  + p15*x15 + p21*x16 + p27*x17 + p33*x18;
114:         pc[15] = m16 = p4*x13 + p10*x14 + p16*x15 + p22*x16 + p28*x17 + p34*x18;
115:         pc[16] = m17 = p5*x13 + p11*x14 + p17*x15 + p23*x16 + p29*x17 + p35*x18;
116:         pc[17] = m18 = p6*x13 + p12*x14 + p18*x15 + p24*x16 + p30*x17 + p36*x18;

118:         pc[18] = m19 = p1*x19 + p7*x20  + p13*x21 + p19*x22 + p25*x23 + p31*x24;
119:         pc[19] = m20 = p2*x19 + p8*x20  + p14*x21 + p20*x22 + p26*x23 + p32*x24;
120:         pc[20] = m21 = p3*x19 + p9*x20  + p15*x21 + p21*x22 + p27*x23 + p33*x24;
121:         pc[21] = m22 = p4*x19 + p10*x20 + p16*x21 + p22*x22 + p28*x23 + p34*x24;
122:         pc[22] = m23 = p5*x19 + p11*x20 + p17*x21 + p23*x22 + p29*x23 + p35*x24;
123:         pc[23] = m24 = p6*x19 + p12*x20 + p18*x21 + p24*x22 + p30*x23 + p36*x24;

125:         pc[24] = m25 = p1*x25 + p7*x26  + p13*x27 + p19*x28 + p25*x29 + p31*x30;
126:         pc[25] = m26 = p2*x25 + p8*x26  + p14*x27 + p20*x28 + p26*x29 + p32*x30;
127:         pc[26] = m27 = p3*x25 + p9*x26  + p15*x27 + p21*x28 + p27*x29 + p33*x30;
128:         pc[27] = m28 = p4*x25 + p10*x26 + p16*x27 + p22*x28 + p28*x29 + p34*x30;
129:         pc[28] = m29 = p5*x25 + p11*x26 + p17*x27 + p23*x28 + p29*x29 + p35*x30;
130:         pc[29] = m30 = p6*x25 + p12*x26 + p18*x27 + p24*x28 + p30*x29 + p36*x30;

132:         pc[30] = m31 = p1*x31 + p7*x32  + p13*x33 + p19*x34 + p25*x35 + p31*x36;
133:         pc[31] = m32 = p2*x31 + p8*x32  + p14*x33 + p20*x34 + p26*x35 + p32*x36;
134:         pc[32] = m33 = p3*x31 + p9*x32  + p15*x33 + p21*x34 + p27*x35 + p33*x36;
135:         pc[33] = m34 = p4*x31 + p10*x32 + p16*x33 + p22*x34 + p28*x35 + p34*x36;
136:         pc[34] = m35 = p5*x31 + p11*x32 + p17*x33 + p23*x34 + p29*x35 + p35*x36;
137:         pc[35] = m36 = p6*x31 + p12*x32 + p18*x33 + p24*x34 + p30*x35 + p36*x36;

139:         nz  = bi[row+1] - diag_offset[row] - 1;
140:         pv += 36;
141:         for (j=0; j<nz; j++) {
142:           x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
143:           x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
144:           x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
145:           x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
146:           x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
147:           x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
148:           x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
149:           x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
150:           x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
151:           x     = rtmp + 36*pj[j];
152:           x[0] -= m1*x1  + m7*x2   + m13*x3  + m19*x4  + m25*x5  + m31*x6;
153:           x[1] -= m2*x1  + m8*x2   + m14*x3  + m20*x4  + m26*x5  + m32*x6;
154:           x[2] -= m3*x1  + m9*x2   + m15*x3  + m21*x4  + m27*x5  + m33*x6;
155:           x[3] -= m4*x1  + m10*x2  + m16*x3  + m22*x4  + m28*x5  + m34*x6;
156:           x[4] -= m5*x1  + m11*x2  + m17*x3  + m23*x4  + m29*x5  + m35*x6;
157:           x[5] -= m6*x1  + m12*x2  + m18*x3  + m24*x4  + m30*x5  + m36*x6;

159:           x[6]  -= m1*x7  + m7*x8   + m13*x9  + m19*x10 + m25*x11 + m31*x12;
160:           x[7]  -= m2*x7  + m8*x8   + m14*x9  + m20*x10 + m26*x11 + m32*x12;
161:           x[8]  -= m3*x7  + m9*x8   + m15*x9  + m21*x10 + m27*x11 + m33*x12;
162:           x[9]  -= m4*x7  + m10*x8  + m16*x9  + m22*x10 + m28*x11 + m34*x12;
163:           x[10] -= m5*x7  + m11*x8  + m17*x9  + m23*x10 + m29*x11 + m35*x12;
164:           x[11] -= m6*x7  + m12*x8  + m18*x9  + m24*x10 + m30*x11 + m36*x12;

166:           x[12] -= m1*x13 + m7*x14  + m13*x15 + m19*x16 + m25*x17 + m31*x18;
167:           x[13] -= m2*x13 + m8*x14  + m14*x15 + m20*x16 + m26*x17 + m32*x18;
168:           x[14] -= m3*x13 + m9*x14  + m15*x15 + m21*x16 + m27*x17 + m33*x18;
169:           x[15] -= m4*x13 + m10*x14 + m16*x15 + m22*x16 + m28*x17 + m34*x18;
170:           x[16] -= m5*x13 + m11*x14 + m17*x15 + m23*x16 + m29*x17 + m35*x18;
171:           x[17] -= m6*x13 + m12*x14 + m18*x15 + m24*x16 + m30*x17 + m36*x18;

173:           x[18] -= m1*x19 + m7*x20  + m13*x21 + m19*x22 + m25*x23 + m31*x24;
174:           x[19] -= m2*x19 + m8*x20  + m14*x21 + m20*x22 + m26*x23 + m32*x24;
175:           x[20] -= m3*x19 + m9*x20  + m15*x21 + m21*x22 + m27*x23 + m33*x24;
176:           x[21] -= m4*x19 + m10*x20 + m16*x21 + m22*x22 + m28*x23 + m34*x24;
177:           x[22] -= m5*x19 + m11*x20 + m17*x21 + m23*x22 + m29*x23 + m35*x24;
178:           x[23] -= m6*x19 + m12*x20 + m18*x21 + m24*x22 + m30*x23 + m36*x24;

180:           x[24] -= m1*x25 + m7*x26  + m13*x27 + m19*x28 + m25*x29 + m31*x30;
181:           x[25] -= m2*x25 + m8*x26  + m14*x27 + m20*x28 + m26*x29 + m32*x30;
182:           x[26] -= m3*x25 + m9*x26  + m15*x27 + m21*x28 + m27*x29 + m33*x30;
183:           x[27] -= m4*x25 + m10*x26 + m16*x27 + m22*x28 + m28*x29 + m34*x30;
184:           x[28] -= m5*x25 + m11*x26 + m17*x27 + m23*x28 + m29*x29 + m35*x30;
185:           x[29] -= m6*x25 + m12*x26 + m18*x27 + m24*x28 + m30*x29 + m36*x30;

187:           x[30] -= m1*x31 + m7*x32  + m13*x33 + m19*x34 + m25*x35 + m31*x36;
188:           x[31] -= m2*x31 + m8*x32  + m14*x33 + m20*x34 + m26*x35 + m32*x36;
189:           x[32] -= m3*x31 + m9*x32  + m15*x33 + m21*x34 + m27*x35 + m33*x36;
190:           x[33] -= m4*x31 + m10*x32 + m16*x33 + m22*x34 + m28*x35 + m34*x36;
191:           x[34] -= m5*x31 + m11*x32 + m17*x33 + m23*x34 + m29*x35 + m35*x36;
192:           x[35] -= m6*x31 + m12*x32 + m18*x33 + m24*x34 + m30*x35 + m36*x36;

194:           pv += 36;
195:         }
196:         PetscLogFlops(432.0*nz+396.0);
197:       }
198:       row = *ajtmp++;
199:     }
200:     /* finished row so stick it into b->a */
201:     pv = ba + 36*bi[i];
202:     pj = bj + bi[i];
203:     nz = bi[i+1] - bi[i];
204:     for (j=0; j<nz; j++) {
205:       x      = rtmp+36*pj[j];
206:       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
207:       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7];
208:       pv[8]  = x[8];  pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11];
209:       pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
210:       pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
211:       pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
212:       pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
213:       pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
214:       pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
215:       pv    += 36;
216:     }
217:     /* invert diagonal block */
218:     w    = ba + 36*diag_offset[i];
219:     PetscKernel_A_gets_inverse_A_6(w,shift,allowzeropivot,&zeropivotdetected);
220:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
221:   }

223:   PetscFree(rtmp);
224:   ISRestoreIndices(isicol,&ic);
225:   ISRestoreIndices(isrow,&r);

227:   C->ops->solve          = MatSolve_SeqBAIJ_6_inplace;
228:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_inplace;
229:   C->assembled           = PETSC_TRUE;

231:   PetscLogFlops(1.333333333333*6*6*6*b->mbs); /* from inverting diagonal blocks */
232:   return 0;
233: }

235: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat B,Mat A,const MatFactorInfo *info)
236: {
237:   Mat            C     = B;
238:   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
239:   IS             isrow = b->row,isicol = b->icol;
240:   const PetscInt *r,*ic;
241:   PetscInt       i,j,k,nz,nzL,row;
242:   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
243:   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
244:   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
245:   PetscInt       flg;
246:   PetscReal      shift = info->shiftamount;
247:   PetscBool      allowzeropivot,zeropivotdetected;

249:   allowzeropivot = PetscNot(A->erroriffailure);
250:   ISGetIndices(isrow,&r);
251:   ISGetIndices(isicol,&ic);

253:   /* generate work space needed by the factorization */
254:   PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
255:   PetscArrayzero(rtmp,bs2*n);

257:   for (i=0; i<n; i++) {
258:     /* zero rtmp */
259:     /* L part */
260:     nz    = bi[i+1] - bi[i];
261:     bjtmp = bj + bi[i];
262:     for  (j=0; j<nz; j++) {
263:       PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
264:     }

266:     /* U part */
267:     nz    = bdiag[i] - bdiag[i+1];
268:     bjtmp = bj + bdiag[i+1]+1;
269:     for  (j=0; j<nz; j++) {
270:       PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
271:     }

273:     /* load in initial (unfactored row) */
274:     nz    = ai[r[i]+1] - ai[r[i]];
275:     ajtmp = aj + ai[r[i]];
276:     v     = aa + bs2*ai[r[i]];
277:     for (j=0; j<nz; j++) {
278:       PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2);
279:     }

281:     /* elimination */
282:     bjtmp = bj + bi[i];
283:     nzL   = bi[i+1] - bi[i];
284:     for (k=0; k < nzL; k++) {
285:       row = bjtmp[k];
286:       pc  = rtmp + bs2*row;
287:       for (flg=0,j=0; j<bs2; j++) {
288:         if (pc[j]!=0.0) {
289:           flg = 1;
290:           break;
291:         }
292:       }
293:       if (flg) {
294:         pv = b->a + bs2*bdiag[row];
295:         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
296:         PetscKernel_A_gets_A_times_B_6(pc,pv,mwork);

298:         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
299:         pv = b->a + bs2*(bdiag[row+1]+1);
300:         nz = bdiag[row] - bdiag[row+1] -  1; /* num of entries inU(row,:), excluding diag */
301:         for (j=0; j<nz; j++) {
302:           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
303:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
304:           v    = rtmp + bs2*pj[j];
305:           PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv);
306:           pv  += bs2;
307:         }
308:         PetscLogFlops(432.0*nz+396); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
309:       }
310:     }

312:     /* finished row so stick it into b->a */
313:     /* L part */
314:     pv = b->a + bs2*bi[i];
315:     pj = b->j + bi[i];
316:     nz = bi[i+1] - bi[i];
317:     for (j=0; j<nz; j++) {
318:       PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
319:     }

321:     /* Mark diagonal and invert diagonal for simpler triangular solves */
322:     pv   = b->a + bs2*bdiag[i];
323:     pj   = b->j + bdiag[i];
324:     PetscArraycpy(pv,rtmp+bs2*pj[0],bs2);
325:     PetscKernel_A_gets_inverse_A_6(pv,shift,allowzeropivot,&zeropivotdetected);
326:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

328:     /* U part */
329:     pv = b->a + bs2*(bdiag[i+1]+1);
330:     pj = b->j + bdiag[i+1]+1;
331:     nz = bdiag[i] - bdiag[i+1] - 1;
332:     for (j=0; j<nz; j++) {
333:       PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
334:     }
335:   }

337:   PetscFree2(rtmp,mwork);
338:   ISRestoreIndices(isicol,&ic);
339:   ISRestoreIndices(isrow,&r);

341:   C->ops->solve          = MatSolve_SeqBAIJ_6;
342:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6;
343:   C->assembled           = PETSC_TRUE;

345:   PetscLogFlops(1.333333333333*6*6*6*n); /* from inverting diagonal blocks */
346:   return 0;
347: }

349: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
350: {
351:   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
352:   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
353:   PetscInt       *ajtmpold,*ajtmp,nz,row;
354:   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
355:   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
356:   MatScalar      x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15;
357:   MatScalar      x16,x17,x18,x19,x20,x21,x22,x23,x24,x25;
358:   MatScalar      p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15;
359:   MatScalar      p16,p17,p18,p19,p20,p21,p22,p23,p24,p25;
360:   MatScalar      m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15;
361:   MatScalar      m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
362:   MatScalar      p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
363:   MatScalar      x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
364:   MatScalar      m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
365:   MatScalar      *ba   = b->a,*aa = a->a;
366:   PetscReal      shift = info->shiftamount;
367:   PetscBool      allowzeropivot,zeropivotdetected;

369:   allowzeropivot = PetscNot(A->erroriffailure);
370:   PetscMalloc1(36*(n+1),&rtmp);
371:   for (i=0; i<n; i++) {
372:     nz    = bi[i+1] - bi[i];
373:     ajtmp = bj + bi[i];
374:     for  (j=0; j<nz; j++) {
375:       x     = rtmp+36*ajtmp[j];
376:       x[0]  = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
377:       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
378:       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
379:       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
380:       x[34] = x[35] = 0.0;
381:     }
382:     /* load in initial (unfactored row) */
383:     nz       = ai[i+1] - ai[i];
384:     ajtmpold = aj + ai[i];
385:     v        = aa + 36*ai[i];
386:     for (j=0; j<nz; j++) {
387:       x     = rtmp+36*ajtmpold[j];
388:       x[0]  =  v[0];  x[1] =  v[1];  x[2] =  v[2];  x[3] =  v[3];
389:       x[4]  =  v[4];  x[5] =  v[5];  x[6] =  v[6];  x[7] =  v[7];
390:       x[8]  =  v[8];  x[9] =  v[9];  x[10] = v[10]; x[11] = v[11];
391:       x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
392:       x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
393:       x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
394:       x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
395:       x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
396:       x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
397:       v    += 36;
398:     }
399:     row = *ajtmp++;
400:     while (row < i) {
401:       pc  = rtmp + 36*row;
402:       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
403:       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];
404:       p9  = pc[8];  p10 = pc[9];  p11 = pc[10]; p12 = pc[11];
405:       p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
406:       p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
407:       p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
408:       p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
409:       p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
410:       p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
411:       if (p1  != 0.0 || p2  != 0.0 || p3  != 0.0 || p4  != 0.0 ||
412:           p5  != 0.0 || p6  != 0.0 || p7  != 0.0 || p8  != 0.0 ||
413:           p9  != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
414:           p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
415:           p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
416:           p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
417:           p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
418:           p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
419:           p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
420:         pv    = ba + 36*diag_offset[row];
421:         pj    = bj + diag_offset[row] + 1;
422:         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
423:         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
424:         x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
425:         x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
426:         x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
427:         x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
428:         x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
429:         x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
430:         x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
431:         pc[0] = m1  = p1*x1  + p7*x2   + p13*x3  + p19*x4  + p25*x5  + p31*x6;
432:         pc[1] = m2  = p2*x1  + p8*x2   + p14*x3  + p20*x4  + p26*x5  + p32*x6;
433:         pc[2] = m3  = p3*x1  + p9*x2   + p15*x3  + p21*x4  + p27*x5  + p33*x6;
434:         pc[3] = m4  = p4*x1  + p10*x2  + p16*x3  + p22*x4  + p28*x5  + p34*x6;
435:         pc[4] = m5  = p5*x1  + p11*x2  + p17*x3  + p23*x4  + p29*x5  + p35*x6;
436:         pc[5] = m6  = p6*x1  + p12*x2  + p18*x3  + p24*x4  + p30*x5  + p36*x6;

438:         pc[6]  = m7  = p1*x7  + p7*x8   + p13*x9  + p19*x10 + p25*x11 + p31*x12;
439:         pc[7]  = m8  = p2*x7  + p8*x8   + p14*x9  + p20*x10 + p26*x11 + p32*x12;
440:         pc[8]  = m9  = p3*x7  + p9*x8   + p15*x9  + p21*x10 + p27*x11 + p33*x12;
441:         pc[9]  = m10 = p4*x7  + p10*x8  + p16*x9  + p22*x10 + p28*x11 + p34*x12;
442:         pc[10] = m11 = p5*x7  + p11*x8  + p17*x9  + p23*x10 + p29*x11 + p35*x12;
443:         pc[11] = m12 = p6*x7  + p12*x8  + p18*x9  + p24*x10 + p30*x11 + p36*x12;

445:         pc[12] = m13 = p1*x13 + p7*x14  + p13*x15 + p19*x16 + p25*x17 + p31*x18;
446:         pc[13] = m14 = p2*x13 + p8*x14  + p14*x15 + p20*x16 + p26*x17 + p32*x18;
447:         pc[14] = m15 = p3*x13 + p9*x14  + p15*x15 + p21*x16 + p27*x17 + p33*x18;
448:         pc[15] = m16 = p4*x13 + p10*x14 + p16*x15 + p22*x16 + p28*x17 + p34*x18;
449:         pc[16] = m17 = p5*x13 + p11*x14 + p17*x15 + p23*x16 + p29*x17 + p35*x18;
450:         pc[17] = m18 = p6*x13 + p12*x14 + p18*x15 + p24*x16 + p30*x17 + p36*x18;

452:         pc[18] = m19 = p1*x19 + p7*x20  + p13*x21 + p19*x22 + p25*x23 + p31*x24;
453:         pc[19] = m20 = p2*x19 + p8*x20  + p14*x21 + p20*x22 + p26*x23 + p32*x24;
454:         pc[20] = m21 = p3*x19 + p9*x20  + p15*x21 + p21*x22 + p27*x23 + p33*x24;
455:         pc[21] = m22 = p4*x19 + p10*x20 + p16*x21 + p22*x22 + p28*x23 + p34*x24;
456:         pc[22] = m23 = p5*x19 + p11*x20 + p17*x21 + p23*x22 + p29*x23 + p35*x24;
457:         pc[23] = m24 = p6*x19 + p12*x20 + p18*x21 + p24*x22 + p30*x23 + p36*x24;

459:         pc[24] = m25 = p1*x25 + p7*x26  + p13*x27 + p19*x28 + p25*x29 + p31*x30;
460:         pc[25] = m26 = p2*x25 + p8*x26  + p14*x27 + p20*x28 + p26*x29 + p32*x30;
461:         pc[26] = m27 = p3*x25 + p9*x26  + p15*x27 + p21*x28 + p27*x29 + p33*x30;
462:         pc[27] = m28 = p4*x25 + p10*x26 + p16*x27 + p22*x28 + p28*x29 + p34*x30;
463:         pc[28] = m29 = p5*x25 + p11*x26 + p17*x27 + p23*x28 + p29*x29 + p35*x30;
464:         pc[29] = m30 = p6*x25 + p12*x26 + p18*x27 + p24*x28 + p30*x29 + p36*x30;

466:         pc[30] = m31 = p1*x31 + p7*x32  + p13*x33 + p19*x34 + p25*x35 + p31*x36;
467:         pc[31] = m32 = p2*x31 + p8*x32  + p14*x33 + p20*x34 + p26*x35 + p32*x36;
468:         pc[32] = m33 = p3*x31 + p9*x32  + p15*x33 + p21*x34 + p27*x35 + p33*x36;
469:         pc[33] = m34 = p4*x31 + p10*x32 + p16*x33 + p22*x34 + p28*x35 + p34*x36;
470:         pc[34] = m35 = p5*x31 + p11*x32 + p17*x33 + p23*x34 + p29*x35 + p35*x36;
471:         pc[35] = m36 = p6*x31 + p12*x32 + p18*x33 + p24*x34 + p30*x35 + p36*x36;

473:         nz  = bi[row+1] - diag_offset[row] - 1;
474:         pv += 36;
475:         for (j=0; j<nz; j++) {
476:           x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
477:           x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
478:           x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
479:           x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
480:           x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
481:           x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
482:           x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
483:           x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
484:           x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
485:           x     = rtmp + 36*pj[j];
486:           x[0] -= m1*x1  + m7*x2   + m13*x3  + m19*x4  + m25*x5  + m31*x6;
487:           x[1] -= m2*x1  + m8*x2   + m14*x3  + m20*x4  + m26*x5  + m32*x6;
488:           x[2] -= m3*x1  + m9*x2   + m15*x3  + m21*x4  + m27*x5  + m33*x6;
489:           x[3] -= m4*x1  + m10*x2  + m16*x3  + m22*x4  + m28*x5  + m34*x6;
490:           x[4] -= m5*x1  + m11*x2  + m17*x3  + m23*x4  + m29*x5  + m35*x6;
491:           x[5] -= m6*x1  + m12*x2  + m18*x3  + m24*x4  + m30*x5  + m36*x6;

493:           x[6]  -= m1*x7  + m7*x8   + m13*x9  + m19*x10 + m25*x11 + m31*x12;
494:           x[7]  -= m2*x7  + m8*x8   + m14*x9  + m20*x10 + m26*x11 + m32*x12;
495:           x[8]  -= m3*x7  + m9*x8   + m15*x9  + m21*x10 + m27*x11 + m33*x12;
496:           x[9]  -= m4*x7  + m10*x8  + m16*x9  + m22*x10 + m28*x11 + m34*x12;
497:           x[10] -= m5*x7  + m11*x8  + m17*x9  + m23*x10 + m29*x11 + m35*x12;
498:           x[11] -= m6*x7  + m12*x8  + m18*x9  + m24*x10 + m30*x11 + m36*x12;

500:           x[12] -= m1*x13 + m7*x14  + m13*x15 + m19*x16 + m25*x17 + m31*x18;
501:           x[13] -= m2*x13 + m8*x14  + m14*x15 + m20*x16 + m26*x17 + m32*x18;
502:           x[14] -= m3*x13 + m9*x14  + m15*x15 + m21*x16 + m27*x17 + m33*x18;
503:           x[15] -= m4*x13 + m10*x14 + m16*x15 + m22*x16 + m28*x17 + m34*x18;
504:           x[16] -= m5*x13 + m11*x14 + m17*x15 + m23*x16 + m29*x17 + m35*x18;
505:           x[17] -= m6*x13 + m12*x14 + m18*x15 + m24*x16 + m30*x17 + m36*x18;

507:           x[18] -= m1*x19 + m7*x20  + m13*x21 + m19*x22 + m25*x23 + m31*x24;
508:           x[19] -= m2*x19 + m8*x20  + m14*x21 + m20*x22 + m26*x23 + m32*x24;
509:           x[20] -= m3*x19 + m9*x20  + m15*x21 + m21*x22 + m27*x23 + m33*x24;
510:           x[21] -= m4*x19 + m10*x20 + m16*x21 + m22*x22 + m28*x23 + m34*x24;
511:           x[22] -= m5*x19 + m11*x20 + m17*x21 + m23*x22 + m29*x23 + m35*x24;
512:           x[23] -= m6*x19 + m12*x20 + m18*x21 + m24*x22 + m30*x23 + m36*x24;

514:           x[24] -= m1*x25 + m7*x26  + m13*x27 + m19*x28 + m25*x29 + m31*x30;
515:           x[25] -= m2*x25 + m8*x26  + m14*x27 + m20*x28 + m26*x29 + m32*x30;
516:           x[26] -= m3*x25 + m9*x26  + m15*x27 + m21*x28 + m27*x29 + m33*x30;
517:           x[27] -= m4*x25 + m10*x26 + m16*x27 + m22*x28 + m28*x29 + m34*x30;
518:           x[28] -= m5*x25 + m11*x26 + m17*x27 + m23*x28 + m29*x29 + m35*x30;
519:           x[29] -= m6*x25 + m12*x26 + m18*x27 + m24*x28 + m30*x29 + m36*x30;

521:           x[30] -= m1*x31 + m7*x32  + m13*x33 + m19*x34 + m25*x35 + m31*x36;
522:           x[31] -= m2*x31 + m8*x32  + m14*x33 + m20*x34 + m26*x35 + m32*x36;
523:           x[32] -= m3*x31 + m9*x32  + m15*x33 + m21*x34 + m27*x35 + m33*x36;
524:           x[33] -= m4*x31 + m10*x32 + m16*x33 + m22*x34 + m28*x35 + m34*x36;
525:           x[34] -= m5*x31 + m11*x32 + m17*x33 + m23*x34 + m29*x35 + m35*x36;
526:           x[35] -= m6*x31 + m12*x32 + m18*x33 + m24*x34 + m30*x35 + m36*x36;

528:           pv += 36;
529:         }
530:         PetscLogFlops(432.0*nz+396.0);
531:       }
532:       row = *ajtmp++;
533:     }
534:     /* finished row so stick it into b->a */
535:     pv = ba + 36*bi[i];
536:     pj = bj + bi[i];
537:     nz = bi[i+1] - bi[i];
538:     for (j=0; j<nz; j++) {
539:       x      = rtmp+36*pj[j];
540:       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
541:       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7];
542:       pv[8]  = x[8];  pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11];
543:       pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
544:       pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
545:       pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
546:       pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
547:       pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
548:       pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
549:       pv    += 36;
550:     }
551:     /* invert diagonal block */
552:     w    = ba + 36*diag_offset[i];
553:     PetscKernel_A_gets_inverse_A_6(w,shift,allowzeropivot,&zeropivotdetected);
554:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
555:   }

557:   PetscFree(rtmp);

559:   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering_inplace;
560:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace;
561:   C->assembled           = PETSC_TRUE;

563:   PetscLogFlops(1.333333333333*6*6*6*b->mbs); /* from inverting diagonal blocks */
564:   return 0;
565: }

567: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
568: {
569:   Mat            C =B;
570:   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
571:   PetscInt       i,j,k,nz,nzL,row;
572:   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
573:   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
574:   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
575:   PetscInt       flg;
576:   PetscReal      shift = info->shiftamount;
577:   PetscBool      allowzeropivot,zeropivotdetected;

579:   allowzeropivot = PetscNot(A->erroriffailure);

581:   /* generate work space needed by the factorization */
582:   PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
583:   PetscArrayzero(rtmp,bs2*n);

585:   for (i=0; i<n; i++) {
586:     /* zero rtmp */
587:     /* L part */
588:     nz    = bi[i+1] - bi[i];
589:     bjtmp = bj + bi[i];
590:     for  (j=0; j<nz; j++) {
591:       PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
592:     }

594:     /* U part */
595:     nz    = bdiag[i] - bdiag[i+1];
596:     bjtmp = bj + bdiag[i+1]+1;
597:     for  (j=0; j<nz; j++) {
598:       PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
599:     }

601:     /* load in initial (unfactored row) */
602:     nz    = ai[i+1] - ai[i];
603:     ajtmp = aj + ai[i];
604:     v     = aa + bs2*ai[i];
605:     for (j=0; j<nz; j++) {
606:       PetscArraycpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2);
607:     }

609:     /* elimination */
610:     bjtmp = bj + bi[i];
611:     nzL   = bi[i+1] - bi[i];
612:     for (k=0; k < nzL; k++) {
613:       row = bjtmp[k];
614:       pc  = rtmp + bs2*row;
615:       for (flg=0,j=0; j<bs2; j++) {
616:         if (pc[j]!=0.0) {
617:           flg = 1;
618:           break;
619:         }
620:       }
621:       if (flg) {
622:         pv = b->a + bs2*bdiag[row];
623:         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
624:         PetscKernel_A_gets_A_times_B_6(pc,pv,mwork);

626:         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
627:         pv = b->a + bs2*(bdiag[row+1]+1);
628:         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
629:         for (j=0; j<nz; j++) {
630:           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
631:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
632:           v    = rtmp + bs2*pj[j];
633:           PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv);
634:           pv  += bs2;
635:         }
636:         PetscLogFlops(432.0*nz+396); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
637:       }
638:     }

640:     /* finished row so stick it into b->a */
641:     /* L part */
642:     pv = b->a + bs2*bi[i];
643:     pj = b->j + bi[i];
644:     nz = bi[i+1] - bi[i];
645:     for (j=0; j<nz; j++) {
646:       PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
647:     }

649:     /* Mark diagonal and invert diagonal for simpler triangular solves */
650:     pv   = b->a + bs2*bdiag[i];
651:     pj   = b->j + bdiag[i];
652:     PetscArraycpy(pv,rtmp+bs2*pj[0],bs2);
653:     PetscKernel_A_gets_inverse_A_6(pv,shift,allowzeropivot,&zeropivotdetected);
654:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

656:     /* U part */
657:     pv = b->a + bs2*(bdiag[i+1]+1);
658:     pj = b->j + bdiag[i+1]+1;
659:     nz = bdiag[i] - bdiag[i+1] - 1;
660:     for (j=0; j<nz; j++) {
661:       PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
662:     }
663:   }
664:   PetscFree2(rtmp,mwork);

666:   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering;
667:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
668:   C->assembled           = PETSC_TRUE;

670:   PetscLogFlops(1.333333333333*6*6*6*n); /* from inverting diagonal blocks */
671:   return 0;
672: }