Actual source code: sbaij2.c

  1: /*$Id: sbaij2.c,v 1.32 2001/08/07 03:03:01 balay Exp $*/

 3:  #include src/mat/impls/baij/seq/baij.h
 4:  #include src/vec/vecimpl.h
 5:  #include src/inline/spops.h
 6:  #include src/inline/ilu.h
 7:  #include petscbt.h
 8:  #include src/mat/impls/sbaij/seq/sbaij.h

 10: int MatIncreaseOverlap_SeqSBAIJ(Mat A,int is_max,IS *is,int ov)
 11: {
 13:   SETERRQ(1,"Function not yet written for SBAIJ format");
 14:   /* return(0); */
 15: }

 17: int MatGetSubMatrix_SeqSBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
 18: {
 19:   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data,*c;
 20:   int          *smap,i,k,kstart,kend,ierr,oldcols = a->mbs,*lens;
 21:   int          row,mat_i,*mat_j,tcol,*mat_ilen;
 22:   int          *irow,nrows,*ssmap,bs=a->bs,bs2=a->bs2;
 23:   int          *aj = a->j,*ai = a->i;
 24:   MatScalar    *mat_a;
 25:   Mat          C;
 26:   PetscTruth   flag;

 29: 
 30:   if (isrow != iscol) SETERRQ(1,"MatGetSubmatrices_SeqSBAIJ: For symm. format, iscol must equal isro");
 31:   ISSorted(iscol,(PetscTruth*)&i);
 32:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");

 34:   ISGetIndices(isrow,&irow);
 35:   ISGetSize(isrow,&nrows);
 36: 
 37:   ierr  = PetscMalloc((1+oldcols)*sizeof(int),&smap);
 38:   ssmap = smap;
 39:   ierr  = PetscMalloc((1+nrows)*sizeof(int),&lens);
 40:   ierr  = PetscMemzero(smap,oldcols*sizeof(int));
 41:   for (i=0; i<nrows; i++) smap[irow[i]] = i+1; /* nrows = ncols */
 42:   /* determine lens of each row */
 43:   for (i=0; i<nrows; i++) {
 44:     kstart  = ai[irow[i]];
 45:     kend    = kstart + a->ilen[irow[i]];
 46:     lens[i] = 0;
 47:       for (k=kstart; k<kend; k++) {
 48:         if (ssmap[aj[k]]) {
 49:           lens[i]++;
 50:         }
 51:       }
 52:     }
 53:   /* Create and fill new matrix */
 54:   if (scall == MAT_REUSE_MATRIX) {
 55:     c = (Mat_SeqSBAIJ *)((*B)->data);

 57:     if (c->mbs!=nrows || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
 58:     PetscMemcmp(c->ilen,lens,c->mbs *sizeof(int),&flag);
 59:     if (flag == PETSC_FALSE) {
 60:       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
 61:     }
 62:     PetscMemzero(c->ilen,c->mbs*sizeof(int));
 63:     C = *B;
 64:   } else {
 65:     MatCreateSeqSBAIJ(A->comm,bs,nrows*bs,nrows*bs,0,lens,&C);
 66:   }
 67:   c = (Mat_SeqSBAIJ *)(C->data);
 68:   for (i=0; i<nrows; i++) {
 69:     row    = irow[i];
 70:     kstart = ai[row];
 71:     kend   = kstart + a->ilen[row];
 72:     mat_i  = c->i[i];
 73:     mat_j  = c->j + mat_i;
 74:     mat_a  = c->a + mat_i*bs2;
 75:     mat_ilen = c->ilen + i;
 76:     for (k=kstart; k<kend; k++) {
 77:       if ((tcol=ssmap[a->j[k]])) {
 78:         *mat_j++ = tcol - 1;
 79:         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));
 80:         mat_a   += bs2;
 81:         (*mat_ilen)++;
 82:       }
 83:     }
 84:   }
 85: 
 86:   /* Free work space */
 87:   PetscFree(smap);
 88:   PetscFree(lens);
 89:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 90:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 91: 
 92:   ISRestoreIndices(isrow,&irow);
 93:   *B = C;
 94:   return(0);
 95: }

 97: int MatGetSubMatrix_SeqSBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
 98: {
 99:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
100:   IS          is1;
101:   int         *vary,*iary,*irow,nrows,i,ierr,bs=a->bs,count;

104:   if (isrow != iscol) SETERRQ(1,"MatGetSubmatrices_SeqSBAIJ: For symm. format, iscol must equal isro");
105: 
106:   ISGetIndices(isrow,&irow);
107:   ISGetSize(isrow,&nrows);
108: 
109:   /* Verify if the indices corespond to each element in a block 
110:    and form the IS with compressed IS */
111:   PetscMalloc(2*(a->mbs+1)*sizeof(int),&vary);
112:   iary = vary + a->mbs;
113:   PetscMemzero(vary,(a->mbs)*sizeof(int));
114:   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
115: 
116:   count = 0;
117:   for (i=0; i<a->mbs; i++) {
118:     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(1,"Index set does not match blocks");
119:     if (vary[i]==bs) iary[count++] = i;
120:   }
121:   ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);
122: 
123:   ISRestoreIndices(isrow,&irow);
124:   PetscFree(vary);

126:   MatGetSubMatrix_SeqSBAIJ_Private(A,is1,is1,cs,scall,B);
127:   ISDestroy(is1);
128:   return(0);
129: }

131: int MatGetSubMatrices_SeqSBAIJ(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
132: {
133:   int ierr,i;

136:   if (scall == MAT_INITIAL_MATRIX) {
137:     PetscMalloc((n+1)*sizeof(Mat),B);
138:   }

140:   for (i=0; i<n; i++) {
141:     MatGetSubMatrix_SeqSBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
142:   }
143:   return(0);
144: }

146: /* -------------------------------------------------------*/
147: /* Should check that shapes of vectors and matrices match */
148: /* -------------------------------------------------------*/
149:  #include petscblaslapack.h

151: int MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz)
152: {
153:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
154:   PetscScalar     *x,*z,*xb,x1,zero=0.0;
155:   MatScalar       *v;
156:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

159:   VecSet(&zero,zz);
160:   VecGetArray(xx,&x);
161:   VecGetArray(zz,&z);

163:   v  = a->a;
164:   xb = x;
165: 
166:   for (i=0; i<mbs; i++) {
167:     n  = ai[1] - ai[0];  /* length of i_th row of A */
168:     x1 = xb[0];
169:     ib = aj + *ai;
170:     jmin = 0;
171:     if (*ib == i) {      /* (diag of A)*x */
172:       z[i] += *v++ * x[*ib++];
173:       jmin++;
174:     }
175:     for (j=jmin; j<n; j++) {
176:       cval    = *ib;
177:       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
178:       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
179:     }
180:     xb++; ai++;
181:   }

183:   VecRestoreArray(xx,&x);
184:   VecRestoreArray(zz,&z);
185:   PetscLogFlops(2*(a->s_nz*2 - A->m) - A->m);  /* s_nz = (nz+m)/2 */
186:   return(0);
187: }

189: int MatMult_SeqSBAIJ_2(Mat A,Vec xx,Vec zz)
190: {
191:   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
192:   PetscScalar      *x,*z,*xb,x1,x2,zero=0.0;
193:   MatScalar        *v;
194:   int              mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;


198:   VecSet(&zero,zz);
199:   VecGetArray(xx,&x);
200:   VecGetArray(zz,&z);
201: 
202:   v     = a->a;
203:   xb = x;

205:   for (i=0; i<mbs; i++) {
206:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
207:     x1 = xb[0]; x2 = xb[1];
208:     ib = aj + *ai;
209:     jmin = 0;
210:     if (*ib == i){     /* (diag of A)*x */
211:       z[2*i]   += v[0]*x1 + v[2]*x2;
212:       z[2*i+1] += v[2]*x1 + v[3]*x2;
213:       v += 4; jmin++;
214:     }
215:     for (j=jmin; j<n; j++) {
216:       /* (strict lower triangular part of A)*x  */
217:       cval       = ib[j]*2;
218:       z[cval]     += v[0]*x1 + v[1]*x2;
219:       z[cval+1]   += v[2]*x1 + v[3]*x2;
220:       /* (strict upper triangular part of A)*x  */
221:       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
222:       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
223:       v  += 4;
224:     }
225:     xb +=2; ai++;
226:   }

228:   VecRestoreArray(xx,&x);
229:   VecRestoreArray(zz,&z);
230:   PetscLogFlops(8*(a->s_nz*2 - A->m) - A->m);
231:   return(0);
232: }

234: int MatMult_SeqSBAIJ_3(Mat A,Vec xx,Vec zz)
235: {
236:   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
237:   PetscScalar   *x,*z,*xb,x1,x2,x3,zero=0.0;
238:   MatScalar     *v;
239:   int           mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;


243:   VecSet(&zero,zz);
244:   VecGetArray(xx,&x);
245:   VecGetArray(zz,&z);
246: 
247:   v     = a->a;
248:   xb = x;

250:   for (i=0; i<mbs; i++) {
251:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
252:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
253:     ib = aj + *ai;
254:     jmin = 0;
255:     if (*ib == i){     /* (diag of A)*x */
256:       z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
257:       z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
258:       z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
259:       v += 9; jmin++;
260:     }
261:     for (j=jmin; j<n; j++) {
262:       /* (strict lower triangular part of A)*x  */
263:       cval       = ib[j]*3;
264:       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
265:       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
266:       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
267:       /* (strict upper triangular part of A)*x  */
268:       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
269:       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
270:       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
271:       v  += 9;
272:     }
273:     xb +=3; ai++;
274:   }

276:   VecRestoreArray(xx,&x);
277:   VecRestoreArray(zz,&z);
278:   PetscLogFlops(18*(a->s_nz*2 - A->m) - A->m);
279:   return(0);
280: }

282: int MatMult_SeqSBAIJ_4(Mat A,Vec xx,Vec zz)
283: {
284:   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
285:   PetscScalar      *x,*z,*xb,x1,x2,x3,x4,zero=0.0;
286:   MatScalar        *v;
287:   int              mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

290:   VecSet(&zero,zz);
291:   VecGetArray(xx,&x);
292:   VecGetArray(zz,&z);
293: 
294:   v     = a->a;
295:   xb = x;

297:   for (i=0; i<mbs; i++) {
298:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
299:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
300:     ib = aj + *ai;
301:     jmin = 0;
302:     if (*ib == i){     /* (diag of A)*x */
303:       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
304:       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
305:       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
306:       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
307:       v += 16; jmin++;
308:     }
309:     for (j=jmin; j<n; j++) {
310:       /* (strict lower triangular part of A)*x  */
311:       cval       = ib[j]*4;
312:       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
313:       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
314:       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
315:       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
316:       /* (strict upper triangular part of A)*x  */
317:       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
318:       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
319:       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
320:       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
321:       v  += 16;
322:     }
323:     xb +=4; ai++;
324:   }

326:   VecRestoreArray(xx,&x);
327:   VecRestoreArray(zz,&z);
328:   PetscLogFlops(32*(a->s_nz*2 - A->m) - A->m);
329:   return(0);
330: }

332: int MatMult_SeqSBAIJ_5(Mat A,Vec xx,Vec zz)
333: {
334:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
335:   PetscScalar     *x,*z,*xb,x1,x2,x3,x4,x5,zero=0.0;
336:   MatScalar       *v;
337:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

340:   VecSet(&zero,zz);
341:   VecGetArray(xx,&x);
342:   VecGetArray(zz,&z);
343: 
344:   v     = a->a;
345:   xb = x;

347:   for (i=0; i<mbs; i++) {
348:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
349:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
350:     ib = aj + *ai;
351:     jmin = 0;
352:     if (*ib == i){      /* (diag of A)*x */
353:       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
354:       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
355:       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
356:       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
357:       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
358:       v += 25; jmin++;
359:     }
360:     for (j=jmin; j<n; j++) {
361:       /* (strict lower triangular part of A)*x  */
362:       cval       = ib[j]*5;
363:       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
364:       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
365:       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
366:       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
367:       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
368:       /* (strict upper triangular part of A)*x  */
369:       z[5*i]   +=v[0]*x[cval]+v[5]*x[cval+1]+v[10]*x[cval+2]+v[15]*x[cval+3]+v[20]*x[cval+4];
370:       z[5*i+1] +=v[1]*x[cval]+v[6]*x[cval+1]+v[11]*x[cval+2]+v[16]*x[cval+3]+v[21]*x[cval+4];
371:       z[5*i+2] +=v[2]*x[cval]+v[7]*x[cval+1]+v[12]*x[cval+2]+v[17]*x[cval+3]+v[22]*x[cval+4];
372:       z[5*i+3] +=v[3]*x[cval]+v[8]*x[cval+1]+v[13]*x[cval+2]+v[18]*x[cval+3]+v[23]*x[cval+4];
373:       z[5*i+4] +=v[4]*x[cval]+v[9]*x[cval+1]+v[14]*x[cval+2]+v[19]*x[cval+3]+v[24]*x[cval+4];
374:       v  += 25;
375:     }
376:     xb +=5; ai++;
377:   }

379:   VecRestoreArray(xx,&x);
380:   VecRestoreArray(zz,&z);
381:   PetscLogFlops(50*(a->s_nz*2 - A->m) - A->m);
382:   return(0);
383: }


386: int MatMult_SeqSBAIJ_6(Mat A,Vec xx,Vec zz)
387: {
388:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
389:   PetscScalar     *x,*z,*xb,x1,x2,x3,x4,x5,x6,zero=0.0;
390:   MatScalar       *v;
391:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

394:   VecSet(&zero,zz);
395:   VecGetArray(xx,&x);
396:   VecGetArray(zz,&z);
397: 
398:   v     = a->a;
399:   xb = x;

401:   for (i=0; i<mbs; i++) {
402:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
403:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
404:     ib = aj + *ai;
405:     jmin = 0;
406:     if (*ib == i){      /* (diag of A)*x */
407:       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
408:       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
409:       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
410:       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
411:       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
412:       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
413:       v += 36; jmin++;
414:     }
415:     for (j=jmin; j<n; j++) {
416:       /* (strict lower triangular part of A)*x  */
417:       cval       = ib[j]*6;
418:       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
419:       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
420:       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
421:       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
422:       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
423:       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
424:       /* (strict upper triangular part of A)*x  */
425:       z[6*i]   +=v[0]*x[cval]+v[6]*x[cval+1]+v[12]*x[cval+2]+v[18]*x[cval+3]+v[24]*x[cval+4]+v[30]*x[cval+5];
426:       z[6*i+1] +=v[1]*x[cval]+v[7]*x[cval+1]+v[13]*x[cval+2]+v[19]*x[cval+3]+v[25]*x[cval+4]+v[31]*x[cval+5];
427:       z[6*i+2] +=v[2]*x[cval]+v[8]*x[cval+1]+v[14]*x[cval+2]+v[20]*x[cval+3]+v[26]*x[cval+4]+v[32]*x[cval+5];
428:       z[6*i+3] +=v[3]*x[cval]+v[9]*x[cval+1]+v[15]*x[cval+2]+v[21]*x[cval+3]+v[27]*x[cval+4]+v[33]*x[cval+5];
429:       z[6*i+4] +=v[4]*x[cval]+v[10]*x[cval+1]+v[16]*x[cval+2]+v[22]*x[cval+3]+v[28]*x[cval+4]+v[34]*x[cval+5];
430:       z[6*i+5] +=v[5]*x[cval]+v[11]*x[cval+1]+v[17]*x[cval+2]+v[23]*x[cval+3]+v[29]*x[cval+4]+v[35]*x[cval+5];
431:       v  += 36;
432:     }
433:     xb +=6; ai++;
434:   }

436:   VecRestoreArray(xx,&x);
437:   VecRestoreArray(zz,&z);
438:   PetscLogFlops(72*(a->s_nz*2 - A->m) - A->m);
439:   return(0);
440: }
441: int MatMult_SeqSBAIJ_7(Mat A,Vec xx,Vec zz)
442: {
443:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
444:   PetscScalar     *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7,zero=0.0;
445:   MatScalar       *v;
446:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

449:   VecSet(&zero,zz);
450:   VecGetArray(xx,&x);
451:   VecGetArray(zz,&z);
452: 
453:   v     = a->a;
454:   xb = x;

456:   for (i=0; i<mbs; i++) {
457:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
458:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
459:     ib = aj + *ai;
460:     jmin = 0;
461:     if (*ib == i){      /* (diag of A)*x */
462:       z[7*i]   += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4+ v[28]*x5 + v[35]*x6+ v[42]*x7;
463:       z[7*i+1] += v[7]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4+ v[29]*x5 + v[36]*x6+ v[43]*x7;
464:       z[7*i+2] += v[14]*x1+ v[15]*x2 +v[16]*x3 + v[23]*x4+ v[30]*x5 + v[37]*x6+ v[44]*x7;
465:       z[7*i+3] += v[21]*x1+ v[22]*x2 +v[23]*x3 + v[24]*x4+ v[31]*x5 + v[38]*x6+ v[45]*x7;
466:       z[7*i+4] += v[28]*x1+ v[29]*x2 +v[30]*x3 + v[31]*x4+ v[32]*x5 + v[39]*x6+ v[46]*x7;
467:       z[7*i+5] += v[35]*x1+ v[36]*x2 +v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[47]*x7;
468:       z[7*i+6] += v[42]*x1+ v[43]*x2 +v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
469:       v += 49; jmin++;
470:     }
471:     for (j=jmin; j<n; j++) {
472:       /* (strict lower triangular part of A)*x  */
473:       cval       = ib[j]*7;
474:       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
475:       z[cval+1] += v[7]*x1  + v[8]*x2 + v[9]*x3 + v[10]*x4+ v[11]*x5 + v[12]*x6+ v[13]*x7;
476:       z[cval+2] += v[14]*x1  + v[15]*x2 + v[16]*x3 + v[17]*x4+ v[18]*x5 + v[19]*x6+ v[20]*x7;
477:       z[cval+3] += v[21]*x1  + v[22]*x2 + v[23]*x3 + v[24]*x4+ v[25]*x5 + v[26]*x6+ v[27]*x7;
478:       z[cval+4] += v[28]*x1  + v[29]*x2 + v[30]*x3 + v[31]*x4+ v[32]*x5 + v[33]*x6+ v[34]*x7;
479:       z[cval+5] += v[35]*x1  + v[36]*x2 + v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[41]*x7;
480:       z[cval+6] += v[42]*x1  + v[43]*x2 + v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
481:       /* (strict upper triangular part of A)*x  */
482:       z[7*i]  +=v[0]*x[cval]+v[7]*x[cval+1]+v[14]*x[cval+2]+v[21]*x[cval+3]+v[28]*x[cval+4]+v[35]*x[cval+5]+v[42]*x[cval+6];
483:       z[7*i+1]+=v[1]*x[cval]+v[8]*x[cval+1]+v[15]*x[cval+2]+v[22]*x[cval+3]+v[29]*x[cval+4]+v[36]*x[cval+5]+v[43]*x[cval+6];
484:       z[7*i+2]+=v[2]*x[cval]+v[9]*x[cval+1]+v[16]*x[cval+2]+v[23]*x[cval+3]+v[30]*x[cval+4]+v[37]*x[cval+5]+v[44]*x[cval+6];
485:       z[7*i+3]+=v[3]*x[cval]+v[10]*x[cval+1]+v[17]*x[cval+2]+v[24]*x[cval+3]+v[31]*x[cval+4]+v[38]*x[cval+5]+v[45]*x[cval+6];
486:       z[7*i+4]+=v[4]*x[cval]+v[11]*x[cval+1]+v[18]*x[cval+2]+v[25]*x[cval+3]+v[32]*x[cval+4]+v[39]*x[cval+5]+v[46]*x[cval+6];
487:       z[7*i+5]+=v[5]*x[cval]+v[12]*x[cval+1]+v[19]*x[cval+2]+v[26]*x[cval+3]+v[33]*x[cval+4]+v[40]*x[cval+5]+v[47]*x[cval+6];
488:       z[7*i+6]+=v[6]*x[cval]+v[13]*x[cval+1]+v[20]*x[cval+2]+v[27]*x[cval+3]+v[34]*x[cval+4]+v[41]*x[cval+5]+v[48]*x[cval+6];
489:       v  += 49;
490:     }
491:     xb +=7; ai++;
492:   }
493:   VecRestoreArray(xx,&x);
494:   VecRestoreArray(zz,&z);
495:   PetscLogFlops(98*(a->s_nz*2 - A->m) - A->m);
496:   return(0);
497: }

499: /*
500:     This will not work with MatScalar == float because it calls the BLAS
501: */
502: int MatMult_SeqSBAIJ_N(Mat A,Vec xx,Vec zz)
503: {
504:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
505:   PetscScalar     *x,*x_ptr,*z,*z_ptr,*xb,*zb,*work,*workt,zero=0.0;
506:   MatScalar       *v;
507:   int             ierr,mbs=a->mbs,i,*idx,*aj,*ii,bs=a->bs,j,n,bs2=a->bs2;
508:   int             ncols,k;

511:   VecSet(&zero,zz);
512:   VecGetArray(xx,&x); x_ptr=x;
513:   VecGetArray(zz,&z); z_ptr=z;

515:   aj   = a->j;
516:   v    = a->a;
517:   ii   = a->i;

519:   if (!a->mult_work) {
520:     PetscMalloc((A->m+1)*sizeof(PetscScalar),&a->mult_work);
521:   }
522:   work = a->mult_work;
523: 
524:   for (i=0; i<mbs; i++) {
525:     n     = ii[1] - ii[0]; ncols = n*bs;
526:     workt = work; idx=aj+ii[0];

528:     /* upper triangular part */
529:     for (j=0; j<n; j++) {
530:       xb = x_ptr + bs*(*idx++);
531:       for (k=0; k<bs; k++) workt[k] = xb[k];
532:       workt += bs;
533:     }
534:     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
535:     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
536: 
537:     /* strict lower triangular part */
538:     idx = aj+ii[0];
539:     if (*idx == i){
540:       ncols -= bs; v += bs2; idx++; n--;
541:     }
542: 
543:     if (ncols > 0){
544:       workt = work;
545:       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));
546:       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
547:       for (j=0; j<n; j++) {
548:         zb = z_ptr + bs*(*idx++);
549:         for (k=0; k<bs; k++) zb[k] += workt[k] ;
550:         workt += bs;
551:       }
552:     }
553:     x += bs; v += n*bs2; z += bs; ii++;
554:   }
555: 
556:   VecRestoreArray(xx,&x);
557:   VecRestoreArray(zz,&z);
558:   PetscLogFlops(2*(a->s_nz*2 - A->m)*bs2 - A->m);
559:   return(0);
560: }

562: int MatMultAdd_SeqSBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
563: {
564:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
565:   PetscScalar     *x,*y,*z,*xb,x1;
566:   MatScalar       *v;
567:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

570:   VecGetArray(xx,&x);
571:   if (yy != xx) {
572:     VecGetArray(yy,&y);
573:   } else {
574:     y = x;
575:   }
576:   if (zz != yy) {
577:     VecCopy(yy,zz);
578:     VecGetArray(zz,&z);
579:   } else {
580:     z = y;
581:   }

583:   v  = a->a;
584:   xb = x;

586:   for (i=0; i<mbs; i++) {
587:     n  = ai[1] - ai[0];  /* length of i_th row of A */
588:     x1 = xb[0];
589:     ib = aj + *ai;
590:     jmin = 0;
591:     if (*ib == i) {            /* (diag of A)*x */
592:       z[i] += *v++ * x[*ib++]; jmin++;
593:     }
594:     for (j=jmin; j<n; j++) {
595:       cval    = *ib;
596:       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
597:       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
598:     }
599:     xb++; ai++;
600:   }

602:   VecRestoreArray(xx,&x);
603:   if (yy != xx) VecRestoreArray(yy,&y);
604:   if (zz != yy) VecRestoreArray(zz,&z);
605: 
606:   PetscLogFlops(2*(a->s_nz*2 - A->m));
607:   return(0);
608: }

610: int MatMultAdd_SeqSBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
611: {
612:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
613:   PetscScalar     *x,*y,*z,*xb,x1,x2;
614:   MatScalar       *v;
615:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

618:   VecGetArray(xx,&x);
619:   if (yy != xx) {
620:     VecGetArray(yy,&y);
621:   } else {
622:     y = x;
623:   }
624:   if (zz != yy) {
625:     VecCopy(yy,zz);
626:     VecGetArray(zz,&z);
627:   } else {
628:     z = y;
629:   }

631:   v     = a->a;
632:   xb = x;

634:   for (i=0; i<mbs; i++) {
635:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
636:     x1 = xb[0]; x2 = xb[1];
637:     ib = aj + *ai;
638:     jmin = 0;
639:     if (*ib == i){      /* (diag of A)*x */
640:       z[2*i]   += v[0]*x1 + v[2]*x2;
641:       z[2*i+1] += v[2]*x1 + v[3]*x2;
642:       v += 4; jmin++;
643:     }
644:     for (j=jmin; j<n; j++) {
645:       /* (strict lower triangular part of A)*x  */
646:       cval       = ib[j]*2;
647:       z[cval]     += v[0]*x1 + v[1]*x2;
648:       z[cval+1]   += v[2]*x1 + v[3]*x2;
649:       /* (strict upper triangular part of A)*x  */
650:       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
651:       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
652:       v  += 4;
653:     }
654:     xb +=2; ai++;
655:   }

657:   VecRestoreArray(xx,&x);
658:   if (yy != xx) VecRestoreArray(yy,&y);
659:   if (zz != yy) VecRestoreArray(zz,&z);

661:   PetscLogFlops(4*(a->s_nz*2 - A->m));
662:   return(0);
663: }

665: int MatMultAdd_SeqSBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
666: {
667:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
668:   PetscScalar     *x,*y,*z,*xb,x1,x2,x3;
669:   MatScalar       *v;
670:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

673:   VecGetArray(xx,&x);
674:   if (yy != xx) {
675:     VecGetArray(yy,&y);
676:   } else {
677:     y = x;
678:   }
679:   if (zz != yy) {
680:     VecCopy(yy,zz);
681:     VecGetArray(zz,&z);
682:   } else {
683:     z = y;
684:   }

686:   v     = a->a;
687:   xb = x;

689:   for (i=0; i<mbs; i++) {
690:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
691:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
692:     ib = aj + *ai;
693:     jmin = 0;
694:     if (*ib == i){     /* (diag of A)*x */
695:      z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
696:      z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
697:      z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
698:      v += 9; jmin++;
699:     }
700:     for (j=jmin; j<n; j++) {
701:       /* (strict lower triangular part of A)*x  */
702:       cval       = ib[j]*3;
703:       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
704:       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
705:       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
706:       /* (strict upper triangular part of A)*x  */
707:       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
708:       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
709:       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
710:       v  += 9;
711:     }
712:     xb +=3; ai++;
713:   }

715:   VecRestoreArray(xx,&x);
716:   if (yy != xx) VecRestoreArray(yy,&y);
717:   if (zz != yy) VecRestoreArray(zz,&z);

719:   PetscLogFlops(18*(a->s_nz*2 - A->m));
720:   return(0);
721: }

723: int MatMultAdd_SeqSBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
724: {
725:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
726:   PetscScalar     *x,*y,*z,*xb,x1,x2,x3,x4;
727:   MatScalar       *v;
728:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

731:   VecGetArray(xx,&x);
732:   if (yy != xx) {
733:     VecGetArray(yy,&y);
734:   } else {
735:     y = x;
736:   }
737:   if (zz != yy) {
738:     VecCopy(yy,zz);
739:     VecGetArray(zz,&z);
740:   } else {
741:     z = y;
742:   }

744:   v     = a->a;
745:   xb = x;

747:   for (i=0; i<mbs; i++) {
748:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
749:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
750:     ib = aj + *ai;
751:     jmin = 0;
752:     if (*ib == i){      /* (diag of A)*x */
753:       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
754:       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
755:       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
756:       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
757:       v += 16; jmin++;
758:     }
759:     for (j=jmin; j<n; j++) {
760:       /* (strict lower triangular part of A)*x  */
761:       cval       = ib[j]*4;
762:       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
763:       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
764:       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
765:       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
766:       /* (strict upper triangular part of A)*x  */
767:       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
768:       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
769:       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
770:       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
771:       v  += 16;
772:     }
773:     xb +=4; ai++;
774:   }

776:   VecRestoreArray(xx,&x);
777:   if (yy != xx) VecRestoreArray(yy,&y);
778:   if (zz != yy) VecRestoreArray(zz,&z);

780:   PetscLogFlops(32*(a->s_nz*2 - A->m));
781:   return(0);
782: }

784: int MatMultAdd_SeqSBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
785: {
786:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
787:   PetscScalar     *x,*y,*z,*xb,x1,x2,x3,x4,x5;
788:   MatScalar       *v;
789:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

792:   VecGetArray(xx,&x);
793:   if (yy != xx) {
794:     VecGetArray(yy,&y);
795:   } else {
796:     y = x;
797:   }
798:   if (zz != yy) {
799:     VecCopy(yy,zz);
800:     VecGetArray(zz,&z);
801:   } else {
802:     z = y;
803:   }

805:   v     = a->a;
806:   xb = x;

808:   for (i=0; i<mbs; i++) {
809:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
810:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
811:     ib = aj + *ai;
812:     jmin = 0;
813:     if (*ib == i){      /* (diag of A)*x */
814:       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
815:       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
816:       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
817:       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
818:       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
819:       v += 25; jmin++;
820:     }
821:     for (j=jmin; j<n; j++) {
822:       /* (strict lower triangular part of A)*x  */
823:       cval       = ib[j]*5;
824:       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
825:       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
826:       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
827:       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
828:       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
829:       /* (strict upper triangular part of A)*x  */
830:       z[5*i]   +=v[0]*x[cval]+v[5]*x[cval+1]+v[10]*x[cval+2]+v[15]*x[cval+3]+v[20]*x[cval+4];
831:       z[5*i+1] +=v[1]*x[cval]+v[6]*x[cval+1]+v[11]*x[cval+2]+v[16]*x[cval+3]+v[21]*x[cval+4];
832:       z[5*i+2] +=v[2]*x[cval]+v[7]*x[cval+1]+v[12]*x[cval+2]+v[17]*x[cval+3]+v[22]*x[cval+4];
833:       z[5*i+3] +=v[3]*x[cval]+v[8]*x[cval+1]+v[13]*x[cval+2]+v[18]*x[cval+3]+v[23]*x[cval+4];
834:       z[5*i+4] +=v[4]*x[cval]+v[9]*x[cval+1]+v[14]*x[cval+2]+v[19]*x[cval+3]+v[24]*x[cval+4];
835:       v  += 25;
836:     }
837:     xb +=5; ai++;
838:   }

840:   VecRestoreArray(xx,&x);
841:   if (yy != xx) VecRestoreArray(yy,&y);
842:   if (zz != yy) VecRestoreArray(zz,&z);

844:   PetscLogFlops(50*(a->s_nz*2 - A->m));
845:   return(0);
846: }
847: int MatMultAdd_SeqSBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
848: {
849:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
850:   PetscScalar     *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6;
851:   MatScalar       *v;
852:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

855:   VecGetArray(xx,&x);
856:   if (yy != xx) {
857:     VecGetArray(yy,&y);
858:   } else {
859:     y = x;
860:   }
861:   if (zz != yy) {
862:     VecCopy(yy,zz);
863:     VecGetArray(zz,&z);
864:   } else {
865:     z = y;
866:   }

868:   v     = a->a;
869:   xb = x;

871:   for (i=0; i<mbs; i++) {
872:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
873:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
874:     ib = aj + *ai;
875:     jmin = 0;
876:     if (*ib == i){     /* (diag of A)*x */
877:       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
878:       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
879:       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
880:       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
881:       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
882:       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
883:       v += 36; jmin++;
884:     }
885:     for (j=jmin; j<n; j++) {
886:       /* (strict lower triangular part of A)*x  */
887:       cval       = ib[j]*6;
888:       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
889:       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
890:       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
891:       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
892:       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
893:       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
894:       /* (strict upper triangular part of A)*x  */
895:       z[6*i]   +=v[0]*x[cval]+v[6]*x[cval+1]+v[12]*x[cval+2]+v[18]*x[cval+3]+v[24]*x[cval+4]+v[30]*x[cval+5];
896:       z[6*i+1] +=v[1]*x[cval]+v[7]*x[cval+1]+v[13]*x[cval+2]+v[19]*x[cval+3]+v[25]*x[cval+4]+v[31]*x[cval+5];
897:       z[6*i+2] +=v[2]*x[cval]+v[8]*x[cval+1]+v[14]*x[cval+2]+v[20]*x[cval+3]+v[26]*x[cval+4]+v[32]*x[cval+5];
898:       z[6*i+3] +=v[3]*x[cval]+v[9]*x[cval+1]+v[15]*x[cval+2]+v[21]*x[cval+3]+v[27]*x[cval+4]+v[33]*x[cval+5];
899:       z[6*i+4] +=v[4]*x[cval]+v[10]*x[cval+1]+v[16]*x[cval+2]+v[22]*x[cval+3]+v[28]*x[cval+4]+v[34]*x[cval+5];
900:       z[6*i+5] +=v[5]*x[cval]+v[11]*x[cval+1]+v[17]*x[cval+2]+v[23]*x[cval+3]+v[29]*x[cval+4]+v[35]*x[cval+5];
901:       v  += 36;
902:     }
903:     xb +=6; ai++;
904:   }

906:   VecRestoreArray(xx,&x);
907:   if (yy != xx) VecRestoreArray(yy,&y);
908:   if (zz != yy) VecRestoreArray(zz,&z);

910:   PetscLogFlops(72*(a->s_nz*2 - A->m));
911:   return(0);
912: }

914: int MatMultAdd_SeqSBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
915: {
916:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
917:   PetscScalar     *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6,x7;
918:   MatScalar       *v;
919:   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;

922:   VecGetArray(xx,&x);
923:   if (yy != xx) {
924:     VecGetArray(yy,&y);
925:   } else {
926:     y = x;
927:   }
928:   if (zz != yy) {
929:     VecCopy(yy,zz);
930:     VecGetArray(zz,&z);
931:   } else {
932:     z = y;
933:   }

935:   v     = a->a;
936:   xb = x;

938:   for (i=0; i<mbs; i++) {
939:     n  = ai[1] - ai[0]; /* length of i_th block row of A */
940:     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
941:     ib = aj + *ai;
942:     jmin = 0;
943:     if (*ib == i){     /* (diag of A)*x */
944:       z[7*i]   += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4+ v[28]*x5 + v[35]*x6+ v[42]*x7;
945:       z[7*i+1] += v[7]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4+ v[29]*x5 + v[36]*x6+ v[43]*x7;
946:       z[7*i+2] += v[14]*x1+ v[15]*x2 +v[16]*x3 + v[23]*x4+ v[30]*x5 + v[37]*x6+ v[44]*x7;
947:       z[7*i+3] += v[21]*x1+ v[22]*x2 +v[23]*x3 + v[24]*x4+ v[31]*x5 + v[38]*x6+ v[45]*x7;
948:       z[7*i+4] += v[28]*x1+ v[29]*x2 +v[30]*x3 + v[31]*x4+ v[32]*x5 + v[39]*x6+ v[46]*x7;
949:       z[7*i+5] += v[35]*x1+ v[36]*x2 +v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[47]*x7;
950:       z[7*i+6] += v[42]*x1+ v[43]*x2 +v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
951:       v += 49; jmin++;
952:     }
953:     for (j=jmin; j<n; j++) {
954:       /* (strict lower triangular part of A)*x  */
955:       cval       = ib[j]*7;
956:       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
957:       z[cval+1] += v[7]*x1  + v[8]*x2 + v[9]*x3 + v[10]*x4+ v[11]*x5 + v[12]*x6+ v[13]*x7;
958:       z[cval+2] += v[14]*x1  + v[15]*x2 + v[16]*x3 + v[17]*x4+ v[18]*x5 + v[19]*x6+ v[20]*x7;
959:       z[cval+3] += v[21]*x1  + v[22]*x2 + v[23]*x3 + v[24]*x4+ v[25]*x5 + v[26]*x6+ v[27]*x7;
960:       z[cval+4] += v[28]*x1  + v[29]*x2 + v[30]*x3 + v[31]*x4+ v[32]*x5 + v[33]*x6+ v[34]*x7;
961:       z[cval+5] += v[35]*x1  + v[36]*x2 + v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[41]*x7;
962:       z[cval+6] += v[42]*x1  + v[43]*x2 + v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
963:       /* (strict upper triangular part of A)*x  */
964:       z[7*i]  +=v[0]*x[cval]+v[7]*x[cval+1]+v[14]*x[cval+2]+v[21]*x[cval+3]+v[28]*x[cval+4]+v[35]*x[cval+5]+v[42]*x[cval+6];
965:       z[7*i+1]+=v[1]*x[cval]+v[8]*x[cval+1]+v[15]*x[cval+2]+v[22]*x[cval+3]+v[29]*x[cval+4]+v[36]*x[cval+5]+v[43]*x[cval+6];
966:       z[7*i+2]+=v[2]*x[cval]+v[9]*x[cval+1]+v[16]*x[cval+2]+v[23]*x[cval+3]+v[30]*x[cval+4]+v[37]*x[cval+5]+v[44]*x[cval+6];
967:       z[7*i+3]+=v[3]*x[cval]+v[10]*x[cval+1]+v[17]*x[cval+2]+v[24]*x[cval+3]+v[31]*x[cval+4]+v[38]*x[cval+5]+v[45]*x[cval+6];
968:       z[7*i+4]+=v[4]*x[cval]+v[11]*x[cval+1]+v[18]*x[cval+2]+v[25]*x[cval+3]+v[32]*x[cval+4]+v[39]*x[cval+5]+v[46]*x[cval+6];
969:       z[7*i+5]+=v[5]*x[cval]+v[12]*x[cval+1]+v[19]*x[cval+2]+v[26]*x[cval+3]+v[33]*x[cval+4]+v[40]*x[cval+5]+v[47]*x[cval+6];
970:       z[7*i+6]+=v[6]*x[cval]+v[13]*x[cval+1]+v[20]*x[cval+2]+v[27]*x[cval+3]+v[34]*x[cval+4]+v[41]*x[cval+5]+v[48]*x[cval+6];
971:       v  += 49;
972:     }
973:     xb +=7; ai++;
974:   }

976:   VecRestoreArray(xx,&x);
977:   if (yy != xx) VecRestoreArray(yy,&y);
978:   if (zz != yy) VecRestoreArray(zz,&z);

980:   PetscLogFlops(98*(a->s_nz*2 - A->m));
981:   return(0);
982: }

984: int MatMultAdd_SeqSBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
985: {
986:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
987:   PetscScalar     *x,*x_ptr,*y,*z,*z_ptr=0,*xb,*zb,*work,*workt;
988:   MatScalar       *v;
989:   int             ierr,mbs=a->mbs,i,*idx,*aj,*ii,bs=a->bs,j,n,bs2=a->bs2;
990:   int             ncols,k;

993:   VecGetArray(xx,&x); x_ptr=x;
994:   if (yy != xx) {
995:     VecGetArray(yy,&y);
996:   } else {
997:     y = x;
998:   }
999:   if (zz != yy) {
1000:     VecCopy(yy,zz);
1001:     VecGetArray(zz,&z); z_ptr=z;
1002:   } else {
1003:     z = y;
1004:   }

1006:   aj   = a->j;
1007:   v    = a->a;
1008:   ii   = a->i;

1010:   if (!a->mult_work) {
1011:     PetscMalloc((A->m+1)*sizeof(PetscScalar),&a->mult_work);
1012:   }
1013:   work = a->mult_work;
1014: 
1015: 
1016:   for (i=0; i<mbs; i++) {
1017:     n     = ii[1] - ii[0]; ncols = n*bs;
1018:     workt = work; idx=aj+ii[0];

1020:     /* upper triangular part */
1021:     for (j=0; j<n; j++) {
1022:       xb = x_ptr + bs*(*idx++);
1023:       for (k=0; k<bs; k++) workt[k] = xb[k];
1024:       workt += bs;
1025:     }
1026:     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
1027:     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);

1029:     /* strict lower triangular part */
1030:     idx = aj+ii[0];
1031:     if (*idx == i){
1032:       ncols -= bs; v += bs2; idx++; n--;
1033:     }
1034:     if (ncols > 0){
1035:       workt = work;
1036:       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));
1037:       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
1038:       for (j=0; j<n; j++) {
1039:         zb = z_ptr + bs*(*idx++);
1040:         /* idx++; */
1041:         for (k=0; k<bs; k++) zb[k] += workt[k] ;
1042:         workt += bs;
1043:       }
1044:     }

1046:     x += bs; v += n*bs2; z += bs; ii++;
1047:   }

1049:   VecRestoreArray(xx,&x);
1050:   if (yy != xx) VecRestoreArray(yy,&y);
1051:   if (zz != yy) VecRestoreArray(zz,&z);

1053:   PetscLogFlops(2*(a->s_nz*2 - A->m));
1054:   return(0);
1055: }

1057: int MatMultTranspose_SeqSBAIJ(Mat A,Vec xx,Vec zz)
1058: {
1060:   SETERRQ(1,"Matrix is symmetric. Call MatMult().");
1061:   /* return(0); */
1062: }

1064: int MatMultTransposeAdd_SeqSBAIJ(Mat A,Vec xx,Vec yy,Vec zz)

1066: {
1068:   SETERRQ(1,"Matrix is symmetric. Call MatMultAdd().");
1069:   /* return(0); */
1070: }

1072: int MatScale_SeqSBAIJ(PetscScalar *alpha,Mat inA)
1073: {
1074:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
1075:   int         one = 1,totalnz = a->bs2*a->s_nz;

1078:   BLscal_(&totalnz,alpha,a->a,&one);
1079:   PetscLogFlops(totalnz);
1080:   return(0);
1081: }

1083: int MatNorm_SeqSBAIJ(Mat A,NormType type,PetscReal *norm)
1084: {
1085:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1086:   MatScalar   *v = a->a;
1087:   PetscReal   sum_diag = 0.0, sum_off = 0.0, *sum;
1088:   int         i,j,k,bs = a->bs,bs2=a->bs2,k1,mbs=a->mbs,*aj=a->j;
1089:   int         *jl,*il,jmin,jmax,ierr,nexti,ik,*col;
1090: 
1092:   if (type == NORM_FROBENIUS) {
1093:     for (k=0; k<mbs; k++){
1094:       jmin = a->i[k]; jmax = a->i[k+1];
1095:       col  = aj + jmin;
1096:       if (*col == k){         /* diagonal block */
1097:         for (i=0; i<bs2; i++){
1098: #if defined(PETSC_USE_COMPLEX)
1099:           sum_diag += PetscRealPart(PetscConj(*v)*(*v)); v++;
1100: #else
1101:           sum_diag += (*v)*(*v); v++;
1102: #endif
1103:         }
1104:         jmin++;
1105:       }
1106:       for (j=jmin; j<jmax; j++){  /* off-diagonal blocks */
1107:         for (i=0; i<bs2; i++){
1108: #if defined(PETSC_USE_COMPLEX)
1109:           sum_off += PetscRealPart(PetscConj(*v)*(*v)); v++;
1110: #else
1111:           sum_off += (*v)*(*v); v++;
1112: #endif  
1113:         }
1114:       }
1115:     }
1116:     *norm = sqrt(sum_diag + 2*sum_off);

1118:   }  else if (type == NORM_INFINITY) { /* maximum row sum */
1119:     PetscMalloc(mbs*sizeof(int),&il);
1120:     PetscMalloc(mbs*sizeof(int),&jl);
1121:     PetscMalloc(bs*sizeof(PetscReal),&sum);
1122:     for (i=0; i<mbs; i++) {
1123:       jl[i] = mbs; il[0] = 0;
1124:     }

1126:     *norm = 0.0;
1127:     for (k=0; k<mbs; k++) { /* k_th block row */
1128:       for (j=0; j<bs; j++) sum[j]=0.0;

1130:       /*-- col sum --*/
1131:       i = jl[k]; /* first |A(i,k)| to be added */
1132:       /* jl[k]=i: first nozero element in row i for submatrix A(1:k,k:n) (active window)
1133:                   at step k */
1134:       while (i<mbs){
1135:         nexti = jl[i];  /* next block row to be added */
1136:         ik    = il[i];  /* block index of A(i,k) in the array a */
1137:         for (j=0; j<bs; j++){
1138:           v = a->a + ik*bs2 + j*bs;
1139:           for (k1=0; k1<bs; k1++) {
1140:             sum[j] += PetscAbsScalar(*v); v++;
1141:           }
1142:         }
1143:         /* update il, jl */
1144:         jmin = ik + 1; /* block index of array a: points to the next nonzero of A in row i */
1145:         jmax = a->i[i+1];
1146:         if (jmin < jmax){
1147:           il[i] = jmin;
1148:           j   = a->j[jmin];
1149:           jl[i] = jl[j]; jl[j]=i;
1150:         }
1151:         i = nexti;
1152:       }
1153: 
1154:       /*-- row sum --*/
1155:       jmin = a->i[k]; jmax = a->i[k+1];
1156:       for (i=jmin; i<jmax; i++) {
1157:         for (j=0; j<bs; j++){
1158:           v = a->a + i*bs2 + j;
1159:           for (k1=0; k1<bs; k1++){
1160:             sum[j] += PetscAbsScalar(*v);
1161:             v   += bs;
1162:           }
1163:         }
1164:       }
1165:       /* add k_th block row to il, jl */
1166:       col = aj+jmin;
1167:       if (*col == k) jmin++;
1168:       if (jmin < jmax){
1169:         il[k] = jmin;
1170:         j   = a->j[jmin];
1171:         jl[k] = jl[j]; jl[j] = k;
1172:       }
1173:       for (j=0; j<bs; j++){
1174:         if (sum[j] > *norm) *norm = sum[j];
1175:       }
1176:     }
1177:     PetscFree(il);
1178:     PetscFree(jl);
1179:     PetscFree(sum);
1180:   } else {
1181:     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
1182:   }
1183:   return(0);
1184: }

1186: int MatEqual_SeqSBAIJ(Mat A,Mat B,PetscTruth* flg)
1187: {
1188:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data,*b = (Mat_SeqSBAIJ *)B->data;
1189:   int          ierr;
1190:   PetscTruth   flag;

1193:   PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&flag);
1194:   if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");

1196:   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1197:   if ((A->m != B->m) || (A->n != B->n) || (a->bs != b->bs)|| (a->s_nz != b->s_nz)) {
1198:     *flg = PETSC_FALSE;
1199:     return(0);
1200:   }
1201: 
1202:   /* if the a->i are the same */
1203:   PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(int),flg);
1204:   if (*flg == PETSC_FALSE) {
1205:     return(0);
1206:   }
1207: 
1208:   /* if a->j are the same */
1209:   PetscMemcmp(a->j,b->j,(a->s_nz)*sizeof(int),flg);
1210:   if (*flg == PETSC_FALSE) {
1211:     return(0);
1212:   }
1213:   /* if a->a are the same */
1214:   PetscMemcmp(a->a,b->a,(a->s_nz)*(a->bs)*(a->bs)*sizeof(PetscScalar),flg);
1215:   return(0);
1216: 
1217: }

1219: int MatGetDiagonal_SeqSBAIJ(Mat A,Vec v)
1220: {
1221:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1222:   int          ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
1223:   PetscScalar  *x,zero = 0.0;
1224:   MatScalar    *aa,*aa_j;

1227:   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1228:   bs   = a->bs;
1229:   aa   = a->a;
1230:   ai   = a->i;
1231:   aj   = a->j;
1232:   ambs = a->mbs;
1233:   bs2  = a->bs2;

1235:   VecSet(&zero,v);
1236:   VecGetArray(v,&x);
1237:   VecGetLocalSize(v,&n);
1238:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1239:   for (i=0; i<ambs; i++) {
1240:     j=ai[i];
1241:     if (aj[j] == i) {             /* if this is a diagonal element */
1242:       row  = i*bs;
1243:       aa_j = aa + j*bs2;
1244:       for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
1245:     }
1246:   }
1247:   VecRestoreArray(v,&x);
1248:   return(0);
1249: }

1251: int MatDiagonalScale_SeqSBAIJ(Mat A,Vec ll,Vec rr)
1252: {
1253:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1254:   PetscScalar  *l,*r,x,*li,*ri;
1255:   MatScalar    *aa,*v;
1256:   int          ierr,i,j,k,lm,rn,M,m,*ai,*aj,mbs,tmp,bs,bs2;

1259:   ai  = a->i;
1260:   aj  = a->j;
1261:   aa  = a->a;
1262:   m   = A->m;
1263:   bs  = a->bs;
1264:   mbs = a->mbs;
1265:   bs2 = a->bs2;

1267:   if (ll != rr) {
1268:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be samen");
1269:   }
1270:   if (ll) {
1271:     VecGetArray(ll,&l);
1272:     VecGetLocalSize(ll,&lm);
1273:     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1274:     for (i=0; i<mbs; i++) { /* for each block row */
1275:       M  = ai[i+1] - ai[i];
1276:       li = l + i*bs;
1277:       v  = aa + bs2*ai[i];
1278:       for (j=0; j<M; j++) { /* for each block */
1279:         for (k=0; k<bs2; k++) {
1280:           (*v++) *= li[k%bs];
1281:         }
1282: #ifdef CONT
1283:         /* will be used to replace the above loop */
1284:         ri = l + bs*aj[ai[i]+j];
1285:         for (k=0; k<bs; k++) { /* column value */
1286:           x = ri[k];
1287:           for (tmp=0; tmp<bs; tmp++) (*v++) *= li[tmp]*x;
1288:         }
1289: #endif

1291:       }
1292:     }
1293:     VecRestoreArray(ll,&l);
1294:     PetscLogFlops(2*a->s_nz);
1295:   }
1296:   /* will be deleted */
1297:   if (rr) {
1298:     VecGetArray(rr,&r);
1299:     VecGetLocalSize(rr,&rn);
1300:     if (rn != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1301:     for (i=0; i<mbs; i++) { /* for each block row */
1302:       M  = ai[i+1] - ai[i];
1303:       v  = aa + bs2*ai[i];
1304:       for (j=0; j<M; j++) { /* for each block */
1305:         ri = r + bs*aj[ai[i]+j];
1306:         for (k=0; k<bs; k++) {
1307:           x = ri[k];
1308:           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
1309:         }
1310:       }
1311:     }
1312:     VecRestoreArray(rr,&r);
1313:     PetscLogFlops(a->s_nz);
1314:   }
1315:   return(0);
1316: }

1318: int MatGetInfo_SeqSBAIJ(Mat A,MatInfoType flag,MatInfo *info)
1319: {
1320:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

1323:   info->rows_global    = (double)A->m;
1324:   info->columns_global = (double)A->m;
1325:   info->rows_local     = (double)A->m;
1326:   info->columns_local  = (double)A->m;
1327:   info->block_size     = a->bs2;
1328:   info->nz_allocated   = a->s_maxnz; /*num. of nonzeros in upper triangular part */
1329:   info->nz_used        = a->bs2*a->s_nz; /*num. of nonzeros in upper triangular part */
1330:   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
1331:   info->assemblies   = A->num_ass;
1332:   info->mallocs      = a->reallocs;
1333:   info->memory       = A->mem;
1334:   if (A->factor) {
1335:     info->fill_ratio_given  = A->info.fill_ratio_given;
1336:     info->fill_ratio_needed = A->info.fill_ratio_needed;
1337:     info->factor_mallocs    = A->info.factor_mallocs;
1338:   } else {
1339:     info->fill_ratio_given  = 0;
1340:     info->fill_ratio_needed = 0;
1341:     info->factor_mallocs    = 0;
1342:   }
1343:   return(0);
1344: }


1347: int MatZeroEntries_SeqSBAIJ(Mat A)
1348: {
1349:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1350:   int         ierr;

1353:   PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));
1354:   return(0);
1355: }

1357: int MatGetRowMax_SeqSBAIJ(Mat A,Vec v)
1358: {
1359:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1360:   int          ierr,i,j,n,row,col,bs,*ai,*aj,mbs;
1361:   PetscReal    atmp;
1362:   MatScalar    *aa;
1363:   PetscScalar  zero = 0.0,*x;
1364:   int          ncols,brow,bcol,krow,kcol;

1367:   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1368:   bs   = a->bs;
1369:   aa   = a->a;
1370:   ai   = a->i;
1371:   aj   = a->j;
1372:   mbs = a->mbs;

1374:   VecSet(&zero,v);
1375:   VecGetArray(v,&x);
1376:   VecGetLocalSize(v,&n);
1377:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1378:   for (i=0; i<mbs; i++) {
1379:     ncols = ai[1] - ai[0]; ai++;
1380:     brow  = bs*i;
1381:     for (j=0; j<ncols; j++){
1382:       bcol = bs*(*aj);
1383:       for (kcol=0; kcol<bs; kcol++){
1384:         col = bcol + kcol;      /* col index */
1385:         for (krow=0; krow<bs; krow++){
1386:           atmp = PetscAbsScalar(*aa); aa++;
1387:           row = brow + krow;    /* row index */
1388:           /* printf("val[%d,%d]: %gn",row,col,atmp); */
1389:           if (PetscRealPart(x[row]) < atmp) x[row] = atmp;
1390:           if (*aj > i && PetscRealPart(x[col]) < atmp) x[col] = atmp;
1391:         }
1392:       }
1393:       aj++;
1394:     }
1395:   }
1396:   VecRestoreArray(v,&x);
1397:   return(0);
1398: }