Actual source code: baijfact9.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 5 by 5
 11: */
 14: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_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;
 18:   PetscErrorCode  ierr;
 19:   const PetscInt  *r,*ic,*bi = b->i,*bj = b->j,*ajtmpold,*ajtmp;
 20:   PetscInt        i,j,n = a->mbs,nz,row,idx,ipvt[5];
 21:   const PetscInt  *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
 22:   MatScalar       *w,*pv,*rtmp,*x,*pc;
 23:   const MatScalar *v,*aa = a->a;
 24:   MatScalar       p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
 25:   MatScalar       p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
 26:   MatScalar       x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
 27:   MatScalar       p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
 28:   MatScalar       m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
 29:   MatScalar       *ba   = b->a,work[25];
 30:   PetscReal       shift = info->shiftamount;
 31:   PetscBool       allowzeropivot,zeropivotdetected;

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

 39: #define PETSC_USE_MEMZERO 1
 40: #define PETSC_USE_MEMCPY 1

 42:   for (i=0; i<n; i++) {
 43:     nz    = bi[i+1] - bi[i];
 44:     ajtmp = bj + bi[i];
 45:     for  (j=0; j<nz; j++) {
 46: #if defined(PETSC_USE_MEMZERO)
 47:       PetscMemzero(rtmp+25*ajtmp[j],25*sizeof(PetscScalar));
 48: #else
 49:       x     = rtmp+25*ajtmp[j];
 50:       x[0]  = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
 51:       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
 52:       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
 53: #endif
 54:     }
 55:     /* load in initial (unfactored row) */
 56:     idx      = r[i];
 57:     nz       = ai[idx+1] - ai[idx];
 58:     ajtmpold = aj + ai[idx];
 59:     v        = aa + 25*ai[idx];
 60:     for (j=0; j<nz; j++) {
 61: #if defined(PETSC_USE_MEMCPY)
 62:       PetscMemcpy(rtmp+25*ic[ajtmpold[j]],v,25*sizeof(PetscScalar));
 63: #else
 64:       x     = rtmp+25*ic[ajtmpold[j]];
 65:       x[0]  = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
 66:       x[4]  = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
 67:       x[9]  = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
 68:       x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17];
 69:       x[18] = v[18]; x[19] = v[19]; x[20] = v[20]; x[21] = v[21];
 70:       x[22] = v[22]; x[23] = v[23]; x[24] = v[24];
 71: #endif
 72:       v += 25;
 73:     }
 74:     row = *ajtmp++;
 75:     while (row < i) {
 76:       pc  = rtmp + 25*row;
 77:       p1  = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
 78:       p5  = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
 79:       p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
 80:       p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17]; p19 = pc[18];
 81:       p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
 82:       p25 = pc[24];
 83:       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
 84:           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
 85:           p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
 86:           || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 ||
 87:           p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 ||
 88:           p24 != 0.0 || p25 != 0.0) {
 89:         pv    = ba + 25*diag_offset[row];
 90:         pj    = bj + diag_offset[row] + 1;
 91:         x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
 92:         x5    = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
 93:         x10   = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
 94:         x15   = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17];
 95:         x19   = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21];
 96:         x23   = pv[22]; x24 = pv[23]; x25 = pv[24];
 97:         pc[0] = m1 = p1*x1 + p6*x2  + p11*x3 + p16*x4 + p21*x5;
 98:         pc[1] = m2 = p2*x1 + p7*x2  + p12*x3 + p17*x4 + p22*x5;
 99:         pc[2] = m3 = p3*x1 + p8*x2  + p13*x3 + p18*x4 + p23*x5;
100:         pc[3] = m4 = p4*x1 + p9*x2  + p14*x3 + p19*x4 + p24*x5;
101:         pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5;

103:         pc[5] = m6 = p1*x6 + p6*x7  + p11*x8 + p16*x9 + p21*x10;
104:         pc[6] = m7 = p2*x6 + p7*x7  + p12*x8 + p17*x9 + p22*x10;
105:         pc[7] = m8 = p3*x6 + p8*x7  + p13*x8 + p18*x9 + p23*x10;
106:         pc[8] = m9 = p4*x6 + p9*x7  + p14*x8 + p19*x9 + p24*x10;
107:         pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10;

109:         pc[10] = m11 = p1*x11 + p6*x12  + p11*x13 + p16*x14 + p21*x15;
110:         pc[11] = m12 = p2*x11 + p7*x12  + p12*x13 + p17*x14 + p22*x15;
111:         pc[12] = m13 = p3*x11 + p8*x12  + p13*x13 + p18*x14 + p23*x15;
112:         pc[13] = m14 = p4*x11 + p9*x12  + p14*x13 + p19*x14 + p24*x15;
113:         pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15;

115:         pc[15] = m16 = p1*x16 + p6*x17  + p11*x18 + p16*x19 + p21*x20;
116:         pc[16] = m17 = p2*x16 + p7*x17  + p12*x18 + p17*x19 + p22*x20;
117:         pc[17] = m18 = p3*x16 + p8*x17  + p13*x18 + p18*x19 + p23*x20;
118:         pc[18] = m19 = p4*x16 + p9*x17  + p14*x18 + p19*x19 + p24*x20;
119:         pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20;

121:         pc[20] = m21 = p1*x21 + p6*x22  + p11*x23 + p16*x24 + p21*x25;
122:         pc[21] = m22 = p2*x21 + p7*x22  + p12*x23 + p17*x24 + p22*x25;
123:         pc[22] = m23 = p3*x21 + p8*x22  + p13*x23 + p18*x24 + p23*x25;
124:         pc[23] = m24 = p4*x21 + p9*x22  + p14*x23 + p19*x24 + p24*x25;
125:         pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25;

127:         nz  = bi[row+1] - diag_offset[row] - 1;
128:         pv += 25;
129:         for (j=0; j<nz; j++) {
130:           x1    = pv[0];  x2 = pv[1];   x3  = pv[2];  x4  = pv[3];
131:           x5    = pv[4];  x6 = pv[5];   x7  = pv[6];  x8  = pv[7]; x9 = pv[8];
132:           x10   = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
133:           x14   = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16];
134:           x18   = pv[17]; x19 = pv[18]; x20 = pv[19]; x21 = pv[20];
135:           x22   = pv[21]; x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
136:           x     = rtmp + 25*pj[j];
137:           x[0] -= m1*x1 + m6*x2  + m11*x3 + m16*x4 + m21*x5;
138:           x[1] -= m2*x1 + m7*x2  + m12*x3 + m17*x4 + m22*x5;
139:           x[2] -= m3*x1 + m8*x2  + m13*x3 + m18*x4 + m23*x5;
140:           x[3] -= m4*x1 + m9*x2  + m14*x3 + m19*x4 + m24*x5;
141:           x[4] -= m5*x1 + m10*x2 + m15*x3 + m20*x4 + m25*x5;

143:           x[5] -= m1*x6 + m6*x7  + m11*x8 + m16*x9 + m21*x10;
144:           x[6] -= m2*x6 + m7*x7  + m12*x8 + m17*x9 + m22*x10;
145:           x[7] -= m3*x6 + m8*x7  + m13*x8 + m18*x9 + m23*x10;
146:           x[8] -= m4*x6 + m9*x7  + m14*x8 + m19*x9 + m24*x10;
147:           x[9] -= m5*x6 + m10*x7 + m15*x8 + m20*x9 + m25*x10;

149:           x[10] -= m1*x11 + m6*x12  + m11*x13 + m16*x14 + m21*x15;
150:           x[11] -= m2*x11 + m7*x12  + m12*x13 + m17*x14 + m22*x15;
151:           x[12] -= m3*x11 + m8*x12  + m13*x13 + m18*x14 + m23*x15;
152:           x[13] -= m4*x11 + m9*x12  + m14*x13 + m19*x14 + m24*x15;
153:           x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15;

155:           x[15] -= m1*x16 + m6*x17  + m11*x18 + m16*x19 + m21*x20;
156:           x[16] -= m2*x16 + m7*x17  + m12*x18 + m17*x19 + m22*x20;
157:           x[17] -= m3*x16 + m8*x17  + m13*x18 + m18*x19 + m23*x20;
158:           x[18] -= m4*x16 + m9*x17  + m14*x18 + m19*x19 + m24*x20;
159:           x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20;

161:           x[20] -= m1*x21 + m6*x22  + m11*x23 + m16*x24 + m21*x25;
162:           x[21] -= m2*x21 + m7*x22  + m12*x23 + m17*x24 + m22*x25;
163:           x[22] -= m3*x21 + m8*x22  + m13*x23 + m18*x24 + m23*x25;
164:           x[23] -= m4*x21 + m9*x22  + m14*x23 + m19*x24 + m24*x25;
165:           x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25;

167:           pv += 25;
168:         }
169:         PetscLogFlops(250.0*nz+225.0);
170:       }
171:       row = *ajtmp++;
172:     }
173:     /* finished row so stick it into b->a */
174:     pv = ba + 25*bi[i];
175:     pj = bj + bi[i];
176:     nz = bi[i+1] - bi[i];
177:     for (j=0; j<nz; j++) {
178: #if defined(PETSC_USE_MEMCPY)
179:       PetscMemcpy(pv,rtmp+25*pj[j],25*sizeof(PetscScalar));
180: #else
181:       x      = rtmp+25*pj[j];
182:       pv[0]  = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
183:       pv[4]  = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
184:       pv[9]  = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
185:       pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16];
186:       pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20];
187:       pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23]; pv[24] = x[24];
188: #endif
189:       pv += 25;
190:     }
191:     /* invert diagonal block */
192:     w    = ba + 25*diag_offset[i];
193:     PetscKernel_A_gets_inverse_A_5(w,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
194:     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
195:   }

197:   PetscFree(rtmp);
198:   ISRestoreIndices(isicol,&ic);
199:   ISRestoreIndices(isrow,&r);

201:   C->ops->solve          = MatSolve_SeqBAIJ_5_inplace;
202:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_inplace;
203:   C->assembled           = PETSC_TRUE;

205:   PetscLogFlops(1.333333333333*5*5*5*b->mbs); /* from inverting diagonal blocks */
206:   return(0);
207: }

209: /* MatLUFactorNumeric_SeqBAIJ_5 -
210:      copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented
211:        PetscKernel_A_gets_A_times_B()
212:        PetscKernel_A_gets_A_minus_B_times_C()
213:        PetscKernel_A_gets_inverse_A()
214: */

218: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat B,Mat A,const MatFactorInfo *info)
219: {
220:   Mat            C     =B;
221:   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
222:   IS             isrow = b->row,isicol = b->icol;
224:   const PetscInt *r,*ic;
225:   PetscInt       i,j,k,nz,nzL,row;
226:   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
227:   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
228:   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a,work[25];
229:   PetscInt       flg,ipvt[5];
230:   PetscReal      shift = info->shiftamount;
231:   PetscBool      allowzeropivot,zeropivotdetected;

234:   allowzeropivot = PetscNot(A->erroriffailure);
235:   ISGetIndices(isrow,&r);
236:   ISGetIndices(isicol,&ic);

238:   /* generate work space needed by the factorization */
239:   PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
240:   PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));

242:   for (i=0; i<n; i++) {
243:     /* zero rtmp */
244:     /* L part */
245:     nz    = bi[i+1] - bi[i];
246:     bjtmp = bj + bi[i];
247:     for  (j=0; j<nz; j++) {
248:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
249:     }

251:     /* U part */
252:     nz    = bdiag[i] - bdiag[i+1];
253:     bjtmp = bj + bdiag[i+1]+1;
254:     for  (j=0; j<nz; j++) {
255:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
256:     }

258:     /* load in initial (unfactored row) */
259:     nz    = ai[r[i]+1] - ai[r[i]];
260:     ajtmp = aj + ai[r[i]];
261:     v     = aa + bs2*ai[r[i]];
262:     for (j=0; j<nz; j++) {
263:       PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));
264:     }

266:     /* elimination */
267:     bjtmp = bj + bi[i];
268:     nzL   = bi[i+1] - bi[i];
269:     for (k=0; k < nzL; k++) {
270:       row = bjtmp[k];
271:       pc  = rtmp + bs2*row;
272:       for (flg=0,j=0; j<bs2; j++) {
273:         if (pc[j]!=0.0) {
274:           flg = 1;
275:           break;
276:         }
277:       }
278:       if (flg) {
279:         pv = b->a + bs2*bdiag[row];
280:         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
281:         PetscKernel_A_gets_A_times_B_5(pc,pv,mwork);

283:         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
284:         pv = b->a + bs2*(bdiag[row+1]+1);
285:         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
286:         for (j=0; j<nz; j++) {
287:           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
288:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
289:           v    = rtmp + bs2*pj[j];
290:           PetscKernel_A_gets_A_minus_B_times_C_5(v,pc,pv);
291:           pv  += bs2;
292:         }
293:         PetscLogFlops(250*nz+225); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
294:       }
295:     }

297:     /* finished row so stick it into b->a */
298:     /* L part */
299:     pv = b->a + bs2*bi[i];
300:     pj = b->j + bi[i];
301:     nz = bi[i+1] - bi[i];
302:     for (j=0; j<nz; j++) {
303:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
304:     }

306:     /* Mark diagonal and invert diagonal for simplier triangular solves */
307:     pv   = b->a + bs2*bdiag[i];
308:     pj   = b->j + bdiag[i];
309:     PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
310:     PetscKernel_A_gets_inverse_A_5(pv,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
311:     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

313:     /* U part */
314:     pv = b->a + bs2*(bdiag[i+1]+1);
315:     pj = b->j + bdiag[i+1]+1;
316:     nz = bdiag[i] - bdiag[i+1] - 1;
317:     for (j=0; j<nz; j++) {
318:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
319:     }
320:   }

322:   PetscFree2(rtmp,mwork);
323:   ISRestoreIndices(isicol,&ic);
324:   ISRestoreIndices(isrow,&r);

326:   C->ops->solve          = MatSolve_SeqBAIJ_5;
327:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5;
328:   C->assembled           = PETSC_TRUE;

330:   PetscLogFlops(1.333333333333*5*5*5*n); /* from inverting diagonal blocks */
331:   return(0);
332: }

334: /*
335:       Version for when blocks are 5 by 5 Using natural ordering
336: */
339: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
340: {
341:   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
343:   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j,ipvt[5];
344:   PetscInt       *ajtmpold,*ajtmp,nz,row;
345:   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
346:   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
347:   MatScalar      x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15;
348:   MatScalar      x16,x17,x18,x19,x20,x21,x22,x23,x24,x25;
349:   MatScalar      p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15;
350:   MatScalar      p16,p17,p18,p19,p20,p21,p22,p23,p24,p25;
351:   MatScalar      m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15;
352:   MatScalar      m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
353:   MatScalar      *ba   = b->a,*aa = a->a,work[25];
354:   PetscReal      shift = info->shiftamount;
355:   PetscBool      allowzeropivot,zeropivotdetected;

358:   allowzeropivot = PetscNot(A->erroriffailure);
359:   PetscMalloc1(25*(n+1),&rtmp);
360:   for (i=0; i<n; i++) {
361:     nz    = bi[i+1] - bi[i];
362:     ajtmp = bj + bi[i];
363:     for  (j=0; j<nz; j++) {
364:       x     = rtmp+25*ajtmp[j];
365:       x[0]  = x[1]  = x[2]  = x[3]  = x[4]  = x[5]  = x[6] = x[7] = x[8] = x[9] = 0.0;
366:       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0;
367:       x[16] = x[17] = x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
368:     }
369:     /* load in initial (unfactored row) */
370:     nz       = ai[i+1] - ai[i];
371:     ajtmpold = aj + ai[i];
372:     v        = aa + 25*ai[i];
373:     for (j=0; j<nz; j++) {
374:       x     = rtmp+25*ajtmpold[j];
375:       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
376:       x[4]  = v[4];  x[5]  = v[5];  x[6]  = v[6];  x[7]  = v[7];  x[8]  = v[8];
377:       x[9]  = v[9];  x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
378:       x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17]; x[18] = v[18];
379:       x[19] = v[19]; x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
380:       x[24] = v[24];
381:       v    += 25;
382:     }
383:     row = *ajtmp++;
384:     while (row < i) {
385:       pc  = rtmp + 25*row;
386:       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
387:       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];  p9  = pc[8];
388:       p10 = pc[9];  p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
389:       p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17];
390:       p19 = pc[18]; p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22];
391:       p24 = pc[23]; p25 = pc[24];
392:       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
393:           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
394:           p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
395:           || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0
396:           || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0) {
397:         pv    = ba + 25*diag_offset[row];
398:         pj    = bj + diag_offset[row] + 1;
399:         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
400:         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];  x9  = pv[8];
401:         x10   = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
402:         x15   = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17]; x19 = pv[18];
403:         x20   = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
404:         x25   = pv[24];
405:         pc[0] = m1 = p1*x1 + p6*x2  + p11*x3 + p16*x4 + p21*x5;
406:         pc[1] = m2 = p2*x1 + p7*x2  + p12*x3 + p17*x4 + p22*x5;
407:         pc[2] = m3 = p3*x1 + p8*x2  + p13*x3 + p18*x4 + p23*x5;
408:         pc[3] = m4 = p4*x1 + p9*x2  + p14*x3 + p19*x4 + p24*x5;
409:         pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5;

411:         pc[5] = m6  = p1*x6 + p6*x7  + p11*x8 + p16*x9 + p21*x10;
412:         pc[6] = m7  = p2*x6 + p7*x7  + p12*x8 + p17*x9 + p22*x10;
413:         pc[7] = m8  = p3*x6 + p8*x7  + p13*x8 + p18*x9 + p23*x10;
414:         pc[8] = m9  = p4*x6 + p9*x7  + p14*x8 + p19*x9 + p24*x10;
415:         pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10;

417:         pc[10] = m11 = p1*x11 + p6*x12  + p11*x13 + p16*x14 + p21*x15;
418:         pc[11] = m12 = p2*x11 + p7*x12  + p12*x13 + p17*x14 + p22*x15;
419:         pc[12] = m13 = p3*x11 + p8*x12  + p13*x13 + p18*x14 + p23*x15;
420:         pc[13] = m14 = p4*x11 + p9*x12  + p14*x13 + p19*x14 + p24*x15;
421:         pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15;

423:         pc[15] = m16 = p1*x16 + p6*x17  + p11*x18 + p16*x19 + p21*x20;
424:         pc[16] = m17 = p2*x16 + p7*x17  + p12*x18 + p17*x19 + p22*x20;
425:         pc[17] = m18 = p3*x16 + p8*x17  + p13*x18 + p18*x19 + p23*x20;
426:         pc[18] = m19 = p4*x16 + p9*x17  + p14*x18 + p19*x19 + p24*x20;
427:         pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20;

429:         pc[20] = m21 = p1*x21 + p6*x22  + p11*x23 + p16*x24 + p21*x25;
430:         pc[21] = m22 = p2*x21 + p7*x22  + p12*x23 + p17*x24 + p22*x25;
431:         pc[22] = m23 = p3*x21 + p8*x22  + p13*x23 + p18*x24 + p23*x25;
432:         pc[23] = m24 = p4*x21 + p9*x22  + p14*x23 + p19*x24 + p24*x25;
433:         pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25;

435:         nz  = bi[row+1] - diag_offset[row] - 1;
436:         pv += 25;
437:         for (j=0; j<nz; j++) {
438:           x1    = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
439:           x5    = pv[4];  x6  = pv[5];   x7 = pv[6];  x8  = pv[7]; x9 = pv[8];
440:           x10   = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
441:           x14   = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17];
442:           x19   = pv[18];  x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22];
443:           x24   = pv[23];  x25 = pv[24];
444:           x     = rtmp + 25*pj[j];
445:           x[0] -= m1*x1 + m6*x2   + m11*x3  + m16*x4 + m21*x5;
446:           x[1] -= m2*x1 + m7*x2   + m12*x3  + m17*x4 + m22*x5;
447:           x[2] -= m3*x1 + m8*x2   + m13*x3  + m18*x4 + m23*x5;
448:           x[3] -= m4*x1 + m9*x2   + m14*x3  + m19*x4 + m24*x5;
449:           x[4] -= m5*x1 + m10*x2  + m15*x3  + m20*x4 + m25*x5;

451:           x[5] -= m1*x6 + m6*x7   + m11*x8  + m16*x9 + m21*x10;
452:           x[6] -= m2*x6 + m7*x7   + m12*x8  + m17*x9 + m22*x10;
453:           x[7] -= m3*x6 + m8*x7   + m13*x8  + m18*x9 + m23*x10;
454:           x[8] -= m4*x6 + m9*x7   + m14*x8  + m19*x9 + m24*x10;
455:           x[9] -= m5*x6 + m10*x7  + m15*x8  + m20*x9 + m25*x10;

457:           x[10] -= m1*x11 + m6*x12  + m11*x13 + m16*x14 + m21*x15;
458:           x[11] -= m2*x11 + m7*x12  + m12*x13 + m17*x14 + m22*x15;
459:           x[12] -= m3*x11 + m8*x12  + m13*x13 + m18*x14 + m23*x15;
460:           x[13] -= m4*x11 + m9*x12  + m14*x13 + m19*x14 + m24*x15;
461:           x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15;

463:           x[15] -= m1*x16 + m6*x17  + m11*x18 + m16*x19 + m21*x20;
464:           x[16] -= m2*x16 + m7*x17  + m12*x18 + m17*x19 + m22*x20;
465:           x[17] -= m3*x16 + m8*x17  + m13*x18 + m18*x19 + m23*x20;
466:           x[18] -= m4*x16 + m9*x17  + m14*x18 + m19*x19 + m24*x20;
467:           x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20;

469:           x[20] -= m1*x21 + m6*x22  + m11*x23 + m16*x24 + m21*x25;
470:           x[21] -= m2*x21 + m7*x22  + m12*x23 + m17*x24 + m22*x25;
471:           x[22] -= m3*x21 + m8*x22  + m13*x23 + m18*x24 + m23*x25;
472:           x[23] -= m4*x21 + m9*x22  + m14*x23 + m19*x24 + m24*x25;
473:           x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25;
474:           pv    += 25;
475:         }
476:         PetscLogFlops(250.0*nz+225.0);
477:       }
478:       row = *ajtmp++;
479:     }
480:     /* finished row so stick it into b->a */
481:     pv = ba + 25*bi[i];
482:     pj = bj + bi[i];
483:     nz = bi[i+1] - bi[i];
484:     for (j=0; j<nz; j++) {
485:       x      = rtmp+25*pj[j];
486:       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
487:       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7]; pv[8] = x[8];
488:       pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
489:       pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16]; pv[17] = x[17];
490:       pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22];
491:       pv[23] = x[23]; pv[24] = x[24];
492:       pv    += 25;
493:     }
494:     /* invert diagonal block */
495:     w    = ba + 25*diag_offset[i];
496:     PetscKernel_A_gets_inverse_A_5(w,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
497:     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
498:   }

500:   PetscFree(rtmp);

502:   C->ops->solve          = MatSolve_SeqBAIJ_5_NaturalOrdering_inplace;
503:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace;
504:   C->assembled           = PETSC_TRUE;

506:   PetscLogFlops(1.333333333333*5*5*5*b->mbs); /* from inverting diagonal blocks */
507:   return(0);
508: }

512: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
513: {
514:   Mat            C =B;
515:   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
517:   PetscInt       i,j,k,nz,nzL,row;
518:   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
519:   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
520:   MatScalar      *rtmp,*pc,*mwork,*v,*vv,*pv,*aa=a->a,work[25];
521:   PetscInt       flg,ipvt[5];
522:   PetscReal      shift = info->shiftamount;
523:   PetscBool      allowzeropivot,zeropivotdetected;

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

528:   /* generate work space needed by the factorization */
529:   PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
530:   PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));

532:   for (i=0; i<n; i++) {
533:     /* zero rtmp */
534:     /* L part */
535:     nz    = bi[i+1] - bi[i];
536:     bjtmp = bj + bi[i];
537:     for  (j=0; j<nz; j++) {
538:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
539:     }

541:     /* U part */
542:     nz    = bdiag[i] - bdiag[i+1];
543:     bjtmp = bj + bdiag[i+1]+1;
544:     for  (j=0; j<nz; j++) {
545:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
546:     }

548:     /* load in initial (unfactored row) */
549:     nz    = ai[i+1] - ai[i];
550:     ajtmp = aj + ai[i];
551:     v     = aa + bs2*ai[i];
552:     for (j=0; j<nz; j++) {
553:       PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));
554:     }

556:     /* elimination */
557:     bjtmp = bj + bi[i];
558:     nzL   = bi[i+1] - bi[i];
559:     for (k=0; k < nzL; k++) {
560:       row = bjtmp[k];
561:       pc  = rtmp + bs2*row;
562:       for (flg=0,j=0; j<bs2; j++) {
563:         if (pc[j]!=0.0) {
564:           flg = 1;
565:           break;
566:         }
567:       }
568:       if (flg) {
569:         pv = b->a + bs2*bdiag[row];
570:         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
571:         PetscKernel_A_gets_A_times_B_5(pc,pv,mwork);

573:         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
574:         pv = b->a + bs2*(bdiag[row+1]+1);
575:         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
576:         for (j=0; j<nz; j++) {
577:           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
578:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
579:           vv   = rtmp + bs2*pj[j];
580:           PetscKernel_A_gets_A_minus_B_times_C_5(vv,pc,pv);
581:           pv  += bs2;
582:         }
583:         PetscLogFlops(250*nz+225); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
584:       }
585:     }

587:     /* finished row so stick it into b->a */
588:     /* L part */
589:     pv = b->a + bs2*bi[i];
590:     pj = b->j + bi[i];
591:     nz = bi[i+1] - bi[i];
592:     for (j=0; j<nz; j++) {
593:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
594:     }

596:     /* Mark diagonal and invert diagonal for simplier triangular solves */
597:     pv   = b->a + bs2*bdiag[i];
598:     pj   = b->j + bdiag[i];
599:     PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
600:     PetscKernel_A_gets_inverse_A_5(pv,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
601:     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

603:     /* U part */
604:     pv = b->a + bs2*(bdiag[i+1]+1);
605:     pj = b->j + bdiag[i+1]+1;
606:     nz = bdiag[i] - bdiag[i+1] - 1;
607:     for (j=0; j<nz; j++) {
608:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
609:     }
610:   }
611:   PetscFree2(rtmp,mwork);

613:   C->ops->solve          = MatSolve_SeqBAIJ_5_NaturalOrdering;
614:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering;
615:   C->assembled           = PETSC_TRUE;

617:   PetscLogFlops(1.333333333333*5*5*5*n); /* from inverting diagonal blocks */
618:   return(0);
619: }