Actual source code: baijfact7.c

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

 35:   allowzeropivot = PetscNot(A->erroriffailure);
 36:   ISGetIndices(isrow,&r);
 37:   ISGetIndices(isicol,&ic);
 38:   PetscMalloc1(36*(n+1),&rtmp);

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

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

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

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

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

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

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

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

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

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

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

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

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

227:   PetscFree(rtmp);
228:   ISRestoreIndices(isicol,&ic);
229:   ISRestoreIndices(isrow,&r);

231:   C->ops->solve          = MatSolve_SeqBAIJ_6_inplace;
232:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_inplace;
233:   C->assembled           = PETSC_TRUE;

235:   PetscLogFlops(1.333333333333*6*6*6*b->mbs); /* from inverting diagonal blocks */
236:   return(0);
237: }

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

257:   allowzeropivot = PetscNot(A->erroriffailure);
258:   ISGetIndices(isrow,&r);
259:   ISGetIndices(isicol,&ic);

261:   /* generate work space needed by the factorization */
262:   PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
263:   PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));

265:   for (i=0; i<n; i++) {
266:     /* zero rtmp */
267:     /* L part */
268:     nz    = bi[i+1] - bi[i];
269:     bjtmp = bj + bi[i];
270:     for  (j=0; j<nz; j++) {
271:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
272:     }

274:     /* U part */
275:     nz    = bdiag[i] - bdiag[i+1];
276:     bjtmp = bj + bdiag[i+1]+1;
277:     for  (j=0; j<nz; j++) {
278:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
279:     }

281:     /* load in initial (unfactored row) */
282:     nz    = ai[r[i]+1] - ai[r[i]];
283:     ajtmp = aj + ai[r[i]];
284:     v     = aa + bs2*ai[r[i]];
285:     for (j=0; j<nz; j++) {
286:       PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));
287:     }

289:     /* elimination */
290:     bjtmp = bj + bi[i];
291:     nzL   = bi[i+1] - bi[i];
292:     for (k=0; k < nzL; k++) {
293:       row = bjtmp[k];
294:       pc  = rtmp + bs2*row;
295:       for (flg=0,j=0; j<bs2; j++) {
296:         if (pc[j]!=0.0) {
297:           flg = 1;
298:           break;
299:         }
300:       }
301:       if (flg) {
302:         pv = b->a + bs2*bdiag[row];
303:         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
304:         PetscKernel_A_gets_A_times_B_6(pc,pv,mwork);

306:         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
307:         pv = b->a + bs2*(bdiag[row+1]+1);
308:         nz = bdiag[row] - bdiag[row+1] -  1; /* num of entries inU(row,:), excluding diag */
309:         for (j=0; j<nz; j++) {
310:           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
311:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
312:           v    = rtmp + bs2*pj[j];
313:           PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv);
314:           pv  += bs2;
315:         }
316:         PetscLogFlops(432*nz+396); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
317:       }
318:     }

320:     /* finished row so stick it into b->a */
321:     /* L part */
322:     pv = b->a + bs2*bi[i];
323:     pj = b->j + bi[i];
324:     nz = bi[i+1] - bi[i];
325:     for (j=0; j<nz; j++) {
326:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
327:     }

329:     /* Mark diagonal and invert diagonal for simplier triangular solves */
330:     pv   = b->a + bs2*bdiag[i];
331:     pj   = b->j + bdiag[i];
332:     PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
333:     PetscKernel_A_gets_inverse_A_6(pv,shift,allowzeropivot,&zeropivotdetected);
334:     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

336:     /* U part */
337:     pv = b->a + bs2*(bdiag[i+1]+1);
338:     pj = b->j + bdiag[i+1]+1;
339:     nz = bdiag[i] - bdiag[i+1] - 1;
340:     for (j=0; j<nz; j++) {
341:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
342:     }
343:   }

345:   PetscFree2(rtmp,mwork);
346:   ISRestoreIndices(isicol,&ic);
347:   ISRestoreIndices(isrow,&r);

349:   C->ops->solve          = MatSolve_SeqBAIJ_6;
350:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6;
351:   C->assembled           = PETSC_TRUE;

353:   PetscLogFlops(1.333333333333*6*6*6*n); /* from inverting diagonal blocks */
354:   return(0);
355: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

569:   PetscFree(rtmp);

571:   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering_inplace;
572:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace;
573:   C->assembled           = PETSC_TRUE;

575:   PetscLogFlops(1.333333333333*6*6*6*b->mbs); /* from inverting diagonal blocks */
576:   return(0);
577: }

581: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
582: {
583:   Mat            C =B;
584:   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
586:   PetscInt       i,j,k,nz,nzL,row;
587:   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
588:   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
589:   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
590:   PetscInt       flg;
591:   PetscReal      shift = info->shiftamount;
592:   PetscBool      allowzeropivot,zeropivotdetected;

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

597:   /* generate work space needed by the factorization */
598:   PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
599:   PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));

601:   for (i=0; i<n; i++) {
602:     /* zero rtmp */
603:     /* L part */
604:     nz    = bi[i+1] - bi[i];
605:     bjtmp = bj + bi[i];
606:     for  (j=0; j<nz; j++) {
607:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
608:     }

610:     /* U part */
611:     nz    = bdiag[i] - bdiag[i+1];
612:     bjtmp = bj + bdiag[i+1]+1;
613:     for  (j=0; j<nz; j++) {
614:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
615:     }

617:     /* load in initial (unfactored row) */
618:     nz    = ai[i+1] - ai[i];
619:     ajtmp = aj + ai[i];
620:     v     = aa + bs2*ai[i];
621:     for (j=0; j<nz; j++) {
622:       PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));
623:     }

625:     /* elimination */
626:     bjtmp = bj + bi[i];
627:     nzL   = bi[i+1] - bi[i];
628:     for (k=0; k < nzL; k++) {
629:       row = bjtmp[k];
630:       pc  = rtmp + bs2*row;
631:       for (flg=0,j=0; j<bs2; j++) {
632:         if (pc[j]!=0.0) {
633:           flg = 1;
634:           break;
635:         }
636:       }
637:       if (flg) {
638:         pv = b->a + bs2*bdiag[row];
639:         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
640:         PetscKernel_A_gets_A_times_B_6(pc,pv,mwork);

642:         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
643:         pv = b->a + bs2*(bdiag[row+1]+1);
644:         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
645:         for (j=0; j<nz; j++) {
646:           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
647:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
648:           v    = rtmp + bs2*pj[j];
649:           PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv);
650:           pv  += bs2;
651:         }
652:         PetscLogFlops(432*nz+396); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
653:       }
654:     }

656:     /* finished row so stick it into b->a */
657:     /* L part */
658:     pv = b->a + bs2*bi[i];
659:     pj = b->j + bi[i];
660:     nz = bi[i+1] - bi[i];
661:     for (j=0; j<nz; j++) {
662:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
663:     }

665:     /* Mark diagonal and invert diagonal for simplier triangular solves */
666:     pv   = b->a + bs2*bdiag[i];
667:     pj   = b->j + bdiag[i];
668:     PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
669:     PetscKernel_A_gets_inverse_A_6(pv,shift,allowzeropivot,&zeropivotdetected);
670:     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

672:     /* U part */
673:     pv = b->a + bs2*(bdiag[i+1]+1);
674:     pj = b->j + bdiag[i+1]+1;
675:     nz = bdiag[i] - bdiag[i+1] - 1;
676:     for (j=0; j<nz; j++) {
677:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
678:     }
679:   }
680:   PetscFree2(rtmp,mwork);

682:   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering;
683:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
684:   C->assembled           = PETSC_TRUE;

686:   PetscLogFlops(1.333333333333*6*6*6*n); /* from inverting diagonal blocks */
687:   return(0);
688: }