Actual source code: baijfact9.c

  1: #define PETSCMAT_DLL

  3: /*
  4:     Factorization code for BAIJ format. 
  5: */
 6:  #include ../src/mat/impls/baij/seq/baij.h
 7:  #include ../src/mat/blockinvert.h

  9: /* ------------------------------------------------------------*/
 10: /*
 11:       Version for when blocks are 5 by 5
 12: */
 15: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_inplace(Mat C,Mat A,const MatFactorInfo *info)
 16: {
 17:   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
 18:   IS              isrow = b->row,isicol = b->icol;
 19:   PetscErrorCode  ierr;
 20:   const PetscInt  *r,*ic,*bi = b->i,*bj = b->j,*ajtmpold,*ajtmp;
 21:   PetscInt        i,j,n = a->mbs,nz,row,idx,ipvt[5];
 22:   const PetscInt  *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
 23:   MatScalar       *w,*pv,*rtmp,*x,*pc;
 24:   const MatScalar *v,*aa = a->a;
 25:   MatScalar       p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
 26:   MatScalar       p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
 27:   MatScalar       x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
 28:   MatScalar       p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
 29:   MatScalar       m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
 30:   MatScalar       *ba = b->a,work[25];
 31:   PetscReal       shift = info->shiftamount;

 34:   ISGetIndices(isrow,&r);
 35:   ISGetIndices(isicol,&ic);
 36:   PetscMalloc(25*(n+1)*sizeof(MatScalar),&rtmp);

 38: #define PETSC_USE_MEMZERO 1
 39: #define PETSC_USE_MEMCPY 1

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

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

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

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

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

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

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

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

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

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

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

195:   PetscFree(rtmp);
196:   ISRestoreIndices(isicol,&ic);
197:   ISRestoreIndices(isrow,&r);
198:   C->ops->solve          = MatSolve_SeqBAIJ_5_inplace;
199:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_inplace;
200:   C->assembled = PETSC_TRUE;
201:   PetscLogFlops(1.333333333333*5*5*5*b->mbs); /* from inverting diagonal blocks */
202:   return(0);
203: }

205: /* MatLUFactorNumeric_SeqBAIJ_5 - 
206:      copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented 
207:        Kernel_A_gets_A_times_B()
208:        Kernel_A_gets_A_minus_B_times_C()
209:        Kernel_A_gets_inverse_A()
210: */

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

229:   ISGetIndices(isrow,&r);
230:   ISGetIndices(isicol,&ic);

232:   /* generate work space needed by the factorization */
233:   PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);
234:   PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));
235:   ics  = ic;

237:   for (i=0; i<n; i++){
238:     /* zero rtmp */
239:     /* L part */
240:     nz    = bi[i+1] - bi[i];
241:     bjtmp = bj + bi[i];
242:     for  (j=0; j<nz; j++){
243:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
244:     }

246:     /* U part */
247:     nz = bdiag[i] - bdiag[i+1];
248:     bjtmp = bj + bdiag[i+1]+1;
249:     for  (j=0; j<nz; j++){
250:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
251:     }
252: 
253:     /* load in initial (unfactored row) */
254:     nz    = ai[r[i]+1] - ai[r[i]];
255:     ajtmp = aj + ai[r[i]];
256:     v     = aa + bs2*ai[r[i]];
257:     for (j=0; j<nz; j++) {
258:       PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));
259:     }

261:     /* elimination */
262:     bjtmp = bj + bi[i];
263:     nzL   = bi[i+1] - bi[i];
264:     for(k=0;k < nzL;k++) {
265:       row = bjtmp[k];
266:       pc = rtmp + bs2*row;
267:       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
268:       if (flg) {
269:         pv = b->a + bs2*bdiag[row];
270:         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
271:         Kernel_A_gets_A_times_B_5(pc,pv,mwork);
272: 
273:         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
274:         pv = b->a + bs2*(bdiag[row+1]+1);
275:         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
276:         for (j=0; j<nz; j++) {
277:           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
278:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
279:           v    = rtmp + bs2*pj[j];
280:           Kernel_A_gets_A_minus_B_times_C_5(v,pc,pv);
281:           pv  += bs2;
282:         }
283:         PetscLogFlops(250*nz+225); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
284:       }
285:     }

287:     /* finished row so stick it into b->a */
288:     /* L part */
289:     pv   = b->a + bs2*bi[i] ;
290:     pj   = b->j + bi[i] ;
291:     nz   = bi[i+1] - bi[i];
292:     for (j=0; j<nz; j++) {
293:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
294:     }

296:     /* Mark diagonal and invert diagonal for simplier triangular solves */
297:     pv   = b->a + bs2*bdiag[i];
298:     pj   = b->j + bdiag[i];
299:     PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
300:     /* Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work); */
301:     Kernel_A_gets_inverse_A_5(pv,ipvt,work,shift);
302: 
303:     /* U part */
304:     pv = b->a + bs2*(bdiag[i+1]+1);
305:     pj = b->j + bdiag[i+1]+1;
306:     nz = bdiag[i] - bdiag[i+1] - 1;
307:     for (j=0; j<nz; j++){
308:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
309:     }
310:   }

312:   PetscFree2(rtmp,mwork);
313:   ISRestoreIndices(isicol,&ic);
314:   ISRestoreIndices(isrow,&r);
315:   C->ops->solve          = MatSolve_SeqBAIJ_5;
316:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5;
317:   C->assembled = PETSC_TRUE;
318:   PetscLogFlops(1.333333333333*5*5*5*n); /* from inverting diagonal blocks */
319:   return(0);
320: }

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

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

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

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

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

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

421:         nz = bi[row+1] - diag_offset[row] - 1;
422:         pv += 25;
423:         for (j=0; j<nz; j++) {
424:           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
425:           x5   = pv[4];  x6  = pv[5];   x7 = pv[6];  x8  = pv[7]; x9 = pv[8];
426:           x10  = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
427:           x14  = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17];
428:           x19 = pv[18];  x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22];
429:           x24 = pv[23];  x25 = pv[24];
430:           x    = rtmp + 25*pj[j];
431:           x[0] -= m1*x1 + m6*x2   + m11*x3  + m16*x4 + m21*x5;
432:           x[1] -= m2*x1 + m7*x2   + m12*x3  + m17*x4 + m22*x5;
433:           x[2] -= m3*x1 + m8*x2   + m13*x3  + m18*x4 + m23*x5;
434:           x[3] -= m4*x1 + m9*x2   + m14*x3  + m19*x4 + m24*x5;
435:           x[4] -= m5*x1 + m10*x2  + m15*x3  + m20*x4 + m25*x5;

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

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

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

455:           x[20] -= m1*x21 + m6*x22  + m11*x23 + m16*x24 + m21*x25;
456:           x[21] -= m2*x21 + m7*x22  + m12*x23 + m17*x24 + m22*x25;
457:           x[22] -= m3*x21 + m8*x22  + m13*x23 + m18*x24 + m23*x25;
458:           x[23] -= m4*x21 + m9*x22  + m14*x23 + m19*x24 + m24*x25;
459:           x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25;
460:           pv   += 25;
461:         }
462:         PetscLogFlops(250.0*nz+225.0);
463:       }
464:       row = *ajtmp++;
465:     }
466:     /* finished row so stick it into b->a */
467:     pv = ba + 25*bi[i];
468:     pj = bj + bi[i];
469:     nz = bi[i+1] - bi[i];
470:     for (j=0; j<nz; j++) {
471:       x      = rtmp+25*pj[j];
472:       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
473:       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7]; pv[8] = x[8];
474:       pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
475:       pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16]; pv[17] = x[17];
476:       pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22];
477:       pv[23] = x[23]; pv[24] = x[24];
478:       pv   += 25;
479:     }
480:     /* invert diagonal block */
481:     w = ba + 25*diag_offset[i];
482:     Kernel_A_gets_inverse_A_5(w,ipvt,work,shift);
483:   }

485:   PetscFree(rtmp);
486:   C->ops->solve          = MatSolve_SeqBAIJ_5_NaturalOrdering_inplace;
487:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace;
488:   C->assembled = PETSC_TRUE;
489:   PetscLogFlops(1.333333333333*5*5*5*b->mbs); /* from inverting diagonal blocks */
490:   return(0);
491: }

495: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
496: {
497:   Mat            C=B;
498:   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
500:   PetscInt       i,j,k,nz,nzL,row;
501:   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
502:   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
503:   MatScalar      *rtmp,*pc,*mwork,*v,*vv,*pv,*aa=a->a,work[25];
504:   PetscInt       flg,ipvt[5];
505:   PetscReal      shift = info->shiftamount;

508:   /* generate work space needed by the factorization */
509:   PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);
510:   PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));

512:   for (i=0; i<n; i++){
513:     /* zero rtmp */
514:     /* L part */
515:     nz    = bi[i+1] - bi[i];
516:     bjtmp = bj + bi[i];
517:     for  (j=0; j<nz; j++){
518:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
519:     }

521:     /* U part */
522:     nz = bdiag[i] - bdiag[i+1];
523:     bjtmp = bj + bdiag[i+1]+1;
524:     for  (j=0; j<nz; j++){
525:       PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
526:     }
527: 
528:     /* load in initial (unfactored row) */
529:     nz    = ai[i+1] - ai[i];
530:     ajtmp = aj + ai[i];
531:     v     = aa + bs2*ai[i];
532:     for (j=0; j<nz; j++) {
533:       PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));
534:     }

536:     /* elimination */
537:     bjtmp = bj + bi[i];
538:     nzL   = bi[i+1] - bi[i];
539:     for(k=0;k < nzL;k++) {
540:       row = bjtmp[k];
541:       pc = rtmp + bs2*row;
542:       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
543:       if (flg) {
544:         pv = b->a + bs2*bdiag[row];
545:         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
546:         Kernel_A_gets_A_times_B_5(pc,pv,mwork);
547: 
548:         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
549:         pv = b->a + bs2*(bdiag[row+1]+1);
550:         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
551:         for (j=0; j<nz; j++) {
552:           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
553:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
554:           vv    = rtmp + bs2*pj[j];
555:           Kernel_A_gets_A_minus_B_times_C_5(vv,pc,pv);
556:           pv  += bs2;
557:         }
558:         PetscLogFlops(250*nz+225); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
559:       }
560:     }

562:     /* finished row so stick it into b->a */
563:     /* L part */
564:     pv   = b->a + bs2*bi[i] ;
565:     pj   = b->j + bi[i] ;
566:     nz   = bi[i+1] - bi[i];
567:     for (j=0; j<nz; j++) {
568:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
569:     }
570: 
571:     /* Mark diagonal and invert diagonal for simplier triangular solves */
572:     pv   = b->a + bs2*bdiag[i];
573:     pj   = b->j + bdiag[i];
574:     PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
575:     /* Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work); */
576:     Kernel_A_gets_inverse_A_5(pv,ipvt,work,shift);
577: 
578:     /* U part */
579:     pv = b->a + bs2*(bdiag[i+1]+1);
580:     pj = b->j + bdiag[i+1]+1;
581:     nz = bdiag[i] - bdiag[i+1] - 1;
582:     for (j=0; j<nz; j++){
583:       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
584:     }
585:   }
586:   PetscFree2(rtmp,mwork);
587:   C->ops->solve          = MatSolve_SeqBAIJ_5_NaturalOrdering;
588:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering;
589:   C->assembled = PETSC_TRUE;
590:   PetscLogFlops(1.333333333333*5*5*5*n); /* from inverting diagonal blocks */
591:   return(0);
592: }