xref: /petsc/src/mat/impls/maij/maij.c (revision 7529efd4b68aedd45db7642b3ca14a621048ee72)
1 /*
2     Defines the basic matrix operations for the MAIJ  matrix storage format.
3   This format is used for restriction and interpolation operations for
4   multicomponent problems. It interpolates each component the same way
5   independently.
6 
7      We provide:
8          MatMult()
9          MatMultTranspose()
10          MatMultTransposeAdd()
11          MatMultAdd()
12           and
13          MatCreateMAIJ(Mat,dof,Mat*)
14 
15      This single directory handles both the sequential and parallel codes
16 */
17 
18 #include "src/mat/impls/maij/maij.h"
19 #include "src/mat/utils/freespace.h"
20 #include "vecimpl.h"
21 
22 #undef __FUNCT__
23 #define __FUNCT__ "MatMAIJGetAIJ"
24 PetscErrorCode MatMAIJGetAIJ(Mat A,Mat *B)
25 {
26   PetscErrorCode ierr;
27   PetscTruth     ismpimaij,isseqmaij;
28 
29   PetscFunctionBegin;
30   ierr = PetscTypeCompare((PetscObject)A,MATMPIMAIJ,&ismpimaij);CHKERRQ(ierr);
31   ierr = PetscTypeCompare((PetscObject)A,MATSEQMAIJ,&isseqmaij);CHKERRQ(ierr);
32   if (ismpimaij) {
33     Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
34 
35     *B = b->A;
36   } else if (isseqmaij) {
37     Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
38 
39     *B = b->AIJ;
40   } else {
41     *B = A;
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 #undef __FUNCT__
47 #define __FUNCT__ "MatMAIJRedimension"
48 PetscErrorCode MatMAIJRedimension(Mat A,PetscInt dof,Mat *B)
49 {
50   PetscErrorCode ierr;
51   Mat            Aij;
52 
53   PetscFunctionBegin;
54   ierr = MatMAIJGetAIJ(A,&Aij);CHKERRQ(ierr);
55   ierr = MatCreateMAIJ(Aij,dof,B);CHKERRQ(ierr);
56   PetscFunctionReturn(0);
57 }
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "MatDestroy_SeqMAIJ"
61 PetscErrorCode MatDestroy_SeqMAIJ(Mat A)
62 {
63   PetscErrorCode ierr;
64   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
65 
66   PetscFunctionBegin;
67   if (b->AIJ) {
68     ierr = MatDestroy(b->AIJ);CHKERRQ(ierr);
69   }
70   ierr = PetscFree(b);CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "MatView_SeqMAIJ"
76 PetscErrorCode MatView_SeqMAIJ(Mat A,PetscViewer viewer)
77 {
78   PetscErrorCode ierr;
79   Mat            B;
80 
81   PetscFunctionBegin;
82   ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
83   ierr = MatView(B,viewer);CHKERRQ(ierr);
84   ierr = MatDestroy(B);CHKERRQ(ierr);
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "MatView_MPIMAIJ"
90 PetscErrorCode MatView_MPIMAIJ(Mat A,PetscViewer viewer)
91 {
92   PetscErrorCode ierr;
93   Mat            B;
94 
95   PetscFunctionBegin;
96   ierr = MatConvert(A,MATMPIAIJ,MAT_INITIAL_MATRIX,&B);
97   ierr = MatView(B,viewer);CHKERRQ(ierr);
98   ierr = MatDestroy(B);CHKERRQ(ierr);
99   PetscFunctionReturn(0);
100 }
101 
102 #undef __FUNCT__
103 #define __FUNCT__ "MatDestroy_MPIMAIJ"
104 PetscErrorCode MatDestroy_MPIMAIJ(Mat A)
105 {
106   PetscErrorCode ierr;
107   Mat_MPIMAIJ    *b = (Mat_MPIMAIJ*)A->data;
108 
109   PetscFunctionBegin;
110   if (b->AIJ) {
111     ierr = MatDestroy(b->AIJ);CHKERRQ(ierr);
112   }
113   if (b->OAIJ) {
114     ierr = MatDestroy(b->OAIJ);CHKERRQ(ierr);
115   }
116   if (b->A) {
117     ierr = MatDestroy(b->A);CHKERRQ(ierr);
118   }
119   if (b->ctx) {
120     ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr);
121   }
122   if (b->w) {
123     ierr = VecDestroy(b->w);CHKERRQ(ierr);
124   }
125   ierr = PetscFree(b);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 /*MC
130   MATMAIJ - MATMAIJ = "maij" - A matrix type to be used for restriction and interpolation operations for
131   multicomponent problems, interpolating or restricting each component the same way independently.
132   The matrix type is based on MATSEQAIJ for sequential matrices, and MATMPIAIJ for distributed matrices.
133 
134   Operations provided:
135 . MatMult
136 . MatMultTranspose
137 . MatMultAdd
138 . MatMultTransposeAdd
139 
140   Level: advanced
141 
142 .seealso: MatCreateSeqDense
143 M*/
144 
145 EXTERN_C_BEGIN
146 #undef __FUNCT__
147 #define __FUNCT__ "MatCreate_MAIJ"
148 PetscErrorCode MatCreate_MAIJ(Mat A)
149 {
150   PetscErrorCode ierr;
151   Mat_MPIMAIJ    *b;
152 
153   PetscFunctionBegin;
154   ierr     = PetscNew(Mat_MPIMAIJ,&b);CHKERRQ(ierr);
155   A->data  = (void*)b;
156   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
157   A->factor           = 0;
158   A->mapping          = 0;
159 
160   b->AIJ  = 0;
161   b->dof  = 0;
162   b->OAIJ = 0;
163   b->ctx  = 0;
164   b->w    = 0;
165   PetscFunctionReturn(0);
166 }
167 EXTERN_C_END
168 
169 /* --------------------------------------------------------------------------------------*/
170 #undef __FUNCT__
171 #define __FUNCT__ "MatMult_SeqMAIJ_2"
172 PetscErrorCode MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
173 {
174   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
175   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
176   PetscScalar    *x,*y,*v,sum1, sum2;
177   PetscErrorCode ierr;
178   PetscInt       m = b->AIJ->m,*idx,*ii;
179   PetscInt       n,i,jrow,j;
180 
181   PetscFunctionBegin;
182   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
183   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
184   idx  = a->j;
185   v    = a->a;
186   ii   = a->i;
187 
188   for (i=0; i<m; i++) {
189     jrow = ii[i];
190     n    = ii[i+1] - jrow;
191     sum1  = 0.0;
192     sum2  = 0.0;
193     for (j=0; j<n; j++) {
194       sum1 += v[jrow]*x[2*idx[jrow]];
195       sum2 += v[jrow]*x[2*idx[jrow]+1];
196       jrow++;
197      }
198     y[2*i]   = sum1;
199     y[2*i+1] = sum2;
200   }
201 
202   ierr = PetscLogFlops(4*a->nz - 2*m);CHKERRQ(ierr);
203   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
204   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
205   PetscFunctionReturn(0);
206 }
207 
208 #undef __FUNCT__
209 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_2"
210 PetscErrorCode MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
211 {
212   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
213   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
214   PetscScalar    *x,*y,*v,alpha1,alpha2,zero = 0.0;
215   PetscErrorCode ierr;
216   PetscInt       m = b->AIJ->m,n,i,*idx;
217 
218   PetscFunctionBegin;
219   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
220   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
221   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
222 
223   for (i=0; i<m; i++) {
224     idx    = a->j + a->i[i] ;
225     v      = a->a + a->i[i] ;
226     n      = a->i[i+1] - a->i[i];
227     alpha1 = x[2*i];
228     alpha2 = x[2*i+1];
229     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
230   }
231   ierr = PetscLogFlops(4*a->nz - 2*b->AIJ->n);CHKERRQ(ierr);
232   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
233   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "MatMultAdd_SeqMAIJ_2"
239 PetscErrorCode MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
240 {
241   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
242   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
243   PetscScalar    *x,*y,*v,sum1, sum2;
244   PetscErrorCode ierr;
245   PetscInt       m = b->AIJ->m,*idx,*ii;
246   PetscInt       n,i,jrow,j;
247 
248   PetscFunctionBegin;
249   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
250   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
251   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
252   idx  = a->j;
253   v    = a->a;
254   ii   = a->i;
255 
256   for (i=0; i<m; i++) {
257     jrow = ii[i];
258     n    = ii[i+1] - jrow;
259     sum1  = 0.0;
260     sum2  = 0.0;
261     for (j=0; j<n; j++) {
262       sum1 += v[jrow]*x[2*idx[jrow]];
263       sum2 += v[jrow]*x[2*idx[jrow]+1];
264       jrow++;
265      }
266     y[2*i]   += sum1;
267     y[2*i+1] += sum2;
268   }
269 
270   ierr = PetscLogFlops(4*a->nz - 2*m);CHKERRQ(ierr);
271   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
272   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
273   PetscFunctionReturn(0);
274 }
275 #undef __FUNCT__
276 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_2"
277 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
278 {
279   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
280   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
281   PetscScalar    *x,*y,*v,alpha1,alpha2;
282   PetscErrorCode ierr;
283   PetscInt       m = b->AIJ->m,n,i,*idx;
284 
285   PetscFunctionBegin;
286   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
287   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
288   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
289 
290   for (i=0; i<m; i++) {
291     idx   = a->j + a->i[i] ;
292     v     = a->a + a->i[i] ;
293     n     = a->i[i+1] - a->i[i];
294     alpha1 = x[2*i];
295     alpha2 = x[2*i+1];
296     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
297   }
298   ierr = PetscLogFlops(4*a->nz - 2*b->AIJ->n);CHKERRQ(ierr);
299   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
300   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
301   PetscFunctionReturn(0);
302 }
303 /* --------------------------------------------------------------------------------------*/
304 #undef __FUNCT__
305 #define __FUNCT__ "MatMult_SeqMAIJ_3"
306 PetscErrorCode MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
307 {
308   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
309   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
310   PetscScalar    *x,*y,*v,sum1, sum2, sum3;
311   PetscErrorCode ierr;
312   PetscInt       m = b->AIJ->m,*idx,*ii;
313   PetscInt       n,i,jrow,j;
314 
315   PetscFunctionBegin;
316   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
317   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
318   idx  = a->j;
319   v    = a->a;
320   ii   = a->i;
321 
322   for (i=0; i<m; i++) {
323     jrow = ii[i];
324     n    = ii[i+1] - jrow;
325     sum1  = 0.0;
326     sum2  = 0.0;
327     sum3  = 0.0;
328     for (j=0; j<n; j++) {
329       sum1 += v[jrow]*x[3*idx[jrow]];
330       sum2 += v[jrow]*x[3*idx[jrow]+1];
331       sum3 += v[jrow]*x[3*idx[jrow]+2];
332       jrow++;
333      }
334     y[3*i]   = sum1;
335     y[3*i+1] = sum2;
336     y[3*i+2] = sum3;
337   }
338 
339   ierr = PetscLogFlops(6*a->nz - 3*m);CHKERRQ(ierr);
340   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
341   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
342   PetscFunctionReturn(0);
343 }
344 
345 #undef __FUNCT__
346 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_3"
347 PetscErrorCode MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
348 {
349   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
350   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
351   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0;
352   PetscErrorCode ierr;
353   PetscInt       m = b->AIJ->m,n,i,*idx;
354 
355   PetscFunctionBegin;
356   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
357   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
358   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
359 
360   for (i=0; i<m; i++) {
361     idx    = a->j + a->i[i];
362     v      = a->a + a->i[i];
363     n      = a->i[i+1] - a->i[i];
364     alpha1 = x[3*i];
365     alpha2 = x[3*i+1];
366     alpha3 = x[3*i+2];
367     while (n-->0) {
368       y[3*(*idx)]   += alpha1*(*v);
369       y[3*(*idx)+1] += alpha2*(*v);
370       y[3*(*idx)+2] += alpha3*(*v);
371       idx++; v++;
372     }
373   }
374   ierr = PetscLogFlops(6*a->nz - 3*b->AIJ->n);CHKERRQ(ierr);
375   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
376   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
377   PetscFunctionReturn(0);
378 }
379 
380 #undef __FUNCT__
381 #define __FUNCT__ "MatMultAdd_SeqMAIJ_3"
382 PetscErrorCode MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
383 {
384   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
385   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
386   PetscScalar    *x,*y,*v,sum1, sum2, sum3;
387   PetscErrorCode ierr;
388   PetscInt       m = b->AIJ->m,*idx,*ii;
389   PetscInt       n,i,jrow,j;
390 
391   PetscFunctionBegin;
392   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
393   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
394   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
395   idx  = a->j;
396   v    = a->a;
397   ii   = a->i;
398 
399   for (i=0; i<m; i++) {
400     jrow = ii[i];
401     n    = ii[i+1] - jrow;
402     sum1  = 0.0;
403     sum2  = 0.0;
404     sum3  = 0.0;
405     for (j=0; j<n; j++) {
406       sum1 += v[jrow]*x[3*idx[jrow]];
407       sum2 += v[jrow]*x[3*idx[jrow]+1];
408       sum3 += v[jrow]*x[3*idx[jrow]+2];
409       jrow++;
410      }
411     y[3*i]   += sum1;
412     y[3*i+1] += sum2;
413     y[3*i+2] += sum3;
414   }
415 
416   ierr = PetscLogFlops(6*a->nz);CHKERRQ(ierr);
417   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
418   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
419   PetscFunctionReturn(0);
420 }
421 #undef __FUNCT__
422 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_3"
423 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
424 {
425   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
426   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
427   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3;
428   PetscErrorCode ierr;
429   PetscInt       m = b->AIJ->m,n,i,*idx;
430 
431   PetscFunctionBegin;
432   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
433   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
434   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
435   for (i=0; i<m; i++) {
436     idx    = a->j + a->i[i] ;
437     v      = a->a + a->i[i] ;
438     n      = a->i[i+1] - a->i[i];
439     alpha1 = x[3*i];
440     alpha2 = x[3*i+1];
441     alpha3 = x[3*i+2];
442     while (n-->0) {
443       y[3*(*idx)]   += alpha1*(*v);
444       y[3*(*idx)+1] += alpha2*(*v);
445       y[3*(*idx)+2] += alpha3*(*v);
446       idx++; v++;
447     }
448   }
449   ierr = PetscLogFlops(6*a->nz);CHKERRQ(ierr);
450   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
451   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
452   PetscFunctionReturn(0);
453 }
454 
455 /* ------------------------------------------------------------------------------*/
456 #undef __FUNCT__
457 #define __FUNCT__ "MatMult_SeqMAIJ_4"
458 PetscErrorCode MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
459 {
460   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
461   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
462   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4;
463   PetscErrorCode ierr;
464   PetscInt       m = b->AIJ->m,*idx,*ii;
465   PetscInt       n,i,jrow,j;
466 
467   PetscFunctionBegin;
468   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
469   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
470   idx  = a->j;
471   v    = a->a;
472   ii   = a->i;
473 
474   for (i=0; i<m; i++) {
475     jrow = ii[i];
476     n    = ii[i+1] - jrow;
477     sum1  = 0.0;
478     sum2  = 0.0;
479     sum3  = 0.0;
480     sum4  = 0.0;
481     for (j=0; j<n; j++) {
482       sum1 += v[jrow]*x[4*idx[jrow]];
483       sum2 += v[jrow]*x[4*idx[jrow]+1];
484       sum3 += v[jrow]*x[4*idx[jrow]+2];
485       sum4 += v[jrow]*x[4*idx[jrow]+3];
486       jrow++;
487      }
488     y[4*i]   = sum1;
489     y[4*i+1] = sum2;
490     y[4*i+2] = sum3;
491     y[4*i+3] = sum4;
492   }
493 
494   ierr = PetscLogFlops(8*a->nz - 4*m);CHKERRQ(ierr);
495   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
496   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
497   PetscFunctionReturn(0);
498 }
499 
500 #undef __FUNCT__
501 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_4"
502 PetscErrorCode MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
503 {
504   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
505   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
506   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0;
507   PetscErrorCode ierr;
508   PetscInt       m = b->AIJ->m,n,i,*idx;
509 
510   PetscFunctionBegin;
511   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
512   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
513   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
514   for (i=0; i<m; i++) {
515     idx    = a->j + a->i[i] ;
516     v      = a->a + a->i[i] ;
517     n      = a->i[i+1] - a->i[i];
518     alpha1 = x[4*i];
519     alpha2 = x[4*i+1];
520     alpha3 = x[4*i+2];
521     alpha4 = x[4*i+3];
522     while (n-->0) {
523       y[4*(*idx)]   += alpha1*(*v);
524       y[4*(*idx)+1] += alpha2*(*v);
525       y[4*(*idx)+2] += alpha3*(*v);
526       y[4*(*idx)+3] += alpha4*(*v);
527       idx++; v++;
528     }
529   }
530   ierr = PetscLogFlops(8*a->nz - 4*b->AIJ->n);CHKERRQ(ierr);
531   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
532   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
533   PetscFunctionReturn(0);
534 }
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "MatMultAdd_SeqMAIJ_4"
538 PetscErrorCode MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
539 {
540   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
541   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
542   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4;
543   PetscErrorCode ierr;
544   PetscInt       m = b->AIJ->m,*idx,*ii;
545   PetscInt       n,i,jrow,j;
546 
547   PetscFunctionBegin;
548   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
549   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
550   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
551   idx  = a->j;
552   v    = a->a;
553   ii   = a->i;
554 
555   for (i=0; i<m; i++) {
556     jrow = ii[i];
557     n    = ii[i+1] - jrow;
558     sum1  = 0.0;
559     sum2  = 0.0;
560     sum3  = 0.0;
561     sum4  = 0.0;
562     for (j=0; j<n; j++) {
563       sum1 += v[jrow]*x[4*idx[jrow]];
564       sum2 += v[jrow]*x[4*idx[jrow]+1];
565       sum3 += v[jrow]*x[4*idx[jrow]+2];
566       sum4 += v[jrow]*x[4*idx[jrow]+3];
567       jrow++;
568      }
569     y[4*i]   += sum1;
570     y[4*i+1] += sum2;
571     y[4*i+2] += sum3;
572     y[4*i+3] += sum4;
573   }
574 
575   ierr = PetscLogFlops(8*a->nz - 4*m);CHKERRQ(ierr);
576   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
577   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
578   PetscFunctionReturn(0);
579 }
580 #undef __FUNCT__
581 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_4"
582 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
583 {
584   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
585   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
586   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4;
587   PetscErrorCode ierr;
588   PetscInt       m = b->AIJ->m,n,i,*idx;
589 
590   PetscFunctionBegin;
591   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
592   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
593   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
594 
595   for (i=0; i<m; i++) {
596     idx    = a->j + a->i[i] ;
597     v      = a->a + a->i[i] ;
598     n      = a->i[i+1] - a->i[i];
599     alpha1 = x[4*i];
600     alpha2 = x[4*i+1];
601     alpha3 = x[4*i+2];
602     alpha4 = x[4*i+3];
603     while (n-->0) {
604       y[4*(*idx)]   += alpha1*(*v);
605       y[4*(*idx)+1] += alpha2*(*v);
606       y[4*(*idx)+2] += alpha3*(*v);
607       y[4*(*idx)+3] += alpha4*(*v);
608       idx++; v++;
609     }
610   }
611   ierr = PetscLogFlops(8*a->nz - 4*b->AIJ->n);CHKERRQ(ierr);
612   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
613   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
614   PetscFunctionReturn(0);
615 }
616 /* ------------------------------------------------------------------------------*/
617 
618 #undef __FUNCT__
619 #define __FUNCT__ "MatMult_SeqMAIJ_5"
620 PetscErrorCode MatMult_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
621 {
622   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
623   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
624   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
625   PetscErrorCode ierr;
626   PetscInt       m = b->AIJ->m,*idx,*ii;
627   PetscInt       n,i,jrow,j;
628 
629   PetscFunctionBegin;
630   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
631   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
632   idx  = a->j;
633   v    = a->a;
634   ii   = a->i;
635 
636   for (i=0; i<m; i++) {
637     jrow = ii[i];
638     n    = ii[i+1] - jrow;
639     sum1  = 0.0;
640     sum2  = 0.0;
641     sum3  = 0.0;
642     sum4  = 0.0;
643     sum5  = 0.0;
644     for (j=0; j<n; j++) {
645       sum1 += v[jrow]*x[5*idx[jrow]];
646       sum2 += v[jrow]*x[5*idx[jrow]+1];
647       sum3 += v[jrow]*x[5*idx[jrow]+2];
648       sum4 += v[jrow]*x[5*idx[jrow]+3];
649       sum5 += v[jrow]*x[5*idx[jrow]+4];
650       jrow++;
651      }
652     y[5*i]   = sum1;
653     y[5*i+1] = sum2;
654     y[5*i+2] = sum3;
655     y[5*i+3] = sum4;
656     y[5*i+4] = sum5;
657   }
658 
659   ierr = PetscLogFlops(10*a->nz - 5*m);CHKERRQ(ierr);
660   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
661   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
662   PetscFunctionReturn(0);
663 }
664 
665 #undef __FUNCT__
666 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_5"
667 PetscErrorCode MatMultTranspose_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
668 {
669   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
670   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
671   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,zero = 0.0;
672   PetscErrorCode ierr;
673   PetscInt       m = b->AIJ->m,n,i,*idx;
674 
675   PetscFunctionBegin;
676   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
677   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
678   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
679 
680   for (i=0; i<m; i++) {
681     idx    = a->j + a->i[i] ;
682     v      = a->a + a->i[i] ;
683     n      = a->i[i+1] - a->i[i];
684     alpha1 = x[5*i];
685     alpha2 = x[5*i+1];
686     alpha3 = x[5*i+2];
687     alpha4 = x[5*i+3];
688     alpha5 = x[5*i+4];
689     while (n-->0) {
690       y[5*(*idx)]   += alpha1*(*v);
691       y[5*(*idx)+1] += alpha2*(*v);
692       y[5*(*idx)+2] += alpha3*(*v);
693       y[5*(*idx)+3] += alpha4*(*v);
694       y[5*(*idx)+4] += alpha5*(*v);
695       idx++; v++;
696     }
697   }
698   ierr = PetscLogFlops(10*a->nz - 5*b->AIJ->n);CHKERRQ(ierr);
699   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
700   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
701   PetscFunctionReturn(0);
702 }
703 
704 #undef __FUNCT__
705 #define __FUNCT__ "MatMultAdd_SeqMAIJ_5"
706 PetscErrorCode MatMultAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
707 {
708   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
709   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
710   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
711   PetscErrorCode ierr;
712   PetscInt       m = b->AIJ->m,*idx,*ii;
713   PetscInt       n,i,jrow,j;
714 
715   PetscFunctionBegin;
716   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
717   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
718   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
719   idx  = a->j;
720   v    = a->a;
721   ii   = a->i;
722 
723   for (i=0; i<m; i++) {
724     jrow = ii[i];
725     n    = ii[i+1] - jrow;
726     sum1  = 0.0;
727     sum2  = 0.0;
728     sum3  = 0.0;
729     sum4  = 0.0;
730     sum5  = 0.0;
731     for (j=0; j<n; j++) {
732       sum1 += v[jrow]*x[5*idx[jrow]];
733       sum2 += v[jrow]*x[5*idx[jrow]+1];
734       sum3 += v[jrow]*x[5*idx[jrow]+2];
735       sum4 += v[jrow]*x[5*idx[jrow]+3];
736       sum5 += v[jrow]*x[5*idx[jrow]+4];
737       jrow++;
738      }
739     y[5*i]   += sum1;
740     y[5*i+1] += sum2;
741     y[5*i+2] += sum3;
742     y[5*i+3] += sum4;
743     y[5*i+4] += sum5;
744   }
745 
746   ierr = PetscLogFlops(10*a->nz);CHKERRQ(ierr);
747   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
748   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
749   PetscFunctionReturn(0);
750 }
751 
752 #undef __FUNCT__
753 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_5"
754 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
755 {
756   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
757   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
758   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5;
759   PetscErrorCode ierr;
760   PetscInt       m = b->AIJ->m,n,i,*idx;
761 
762   PetscFunctionBegin;
763   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
764   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
765   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
766 
767   for (i=0; i<m; i++) {
768     idx    = a->j + a->i[i] ;
769     v      = a->a + a->i[i] ;
770     n      = a->i[i+1] - a->i[i];
771     alpha1 = x[5*i];
772     alpha2 = x[5*i+1];
773     alpha3 = x[5*i+2];
774     alpha4 = x[5*i+3];
775     alpha5 = x[5*i+4];
776     while (n-->0) {
777       y[5*(*idx)]   += alpha1*(*v);
778       y[5*(*idx)+1] += alpha2*(*v);
779       y[5*(*idx)+2] += alpha3*(*v);
780       y[5*(*idx)+3] += alpha4*(*v);
781       y[5*(*idx)+4] += alpha5*(*v);
782       idx++; v++;
783     }
784   }
785   ierr = PetscLogFlops(10*a->nz);CHKERRQ(ierr);
786   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
787   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
788   PetscFunctionReturn(0);
789 }
790 
791 /* ------------------------------------------------------------------------------*/
792 #undef __FUNCT__
793 #define __FUNCT__ "MatMult_SeqMAIJ_6"
794 PetscErrorCode MatMult_SeqMAIJ_6(Mat A,Vec xx,Vec yy)
795 {
796   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
797   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
798   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6;
799   PetscErrorCode ierr;
800   PetscInt       m = b->AIJ->m,*idx,*ii;
801   PetscInt       n,i,jrow,j;
802 
803   PetscFunctionBegin;
804   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
805   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
806   idx  = a->j;
807   v    = a->a;
808   ii   = a->i;
809 
810   for (i=0; i<m; i++) {
811     jrow = ii[i];
812     n    = ii[i+1] - jrow;
813     sum1  = 0.0;
814     sum2  = 0.0;
815     sum3  = 0.0;
816     sum4  = 0.0;
817     sum5  = 0.0;
818     sum6  = 0.0;
819     for (j=0; j<n; j++) {
820       sum1 += v[jrow]*x[6*idx[jrow]];
821       sum2 += v[jrow]*x[6*idx[jrow]+1];
822       sum3 += v[jrow]*x[6*idx[jrow]+2];
823       sum4 += v[jrow]*x[6*idx[jrow]+3];
824       sum5 += v[jrow]*x[6*idx[jrow]+4];
825       sum6 += v[jrow]*x[6*idx[jrow]+5];
826       jrow++;
827      }
828     y[6*i]   = sum1;
829     y[6*i+1] = sum2;
830     y[6*i+2] = sum3;
831     y[6*i+3] = sum4;
832     y[6*i+4] = sum5;
833     y[6*i+5] = sum6;
834   }
835 
836   ierr = PetscLogFlops(12*a->nz - 6*m);CHKERRQ(ierr);
837   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
838   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
839   PetscFunctionReturn(0);
840 }
841 
842 #undef __FUNCT__
843 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_6"
844 PetscErrorCode MatMultTranspose_SeqMAIJ_6(Mat A,Vec xx,Vec yy)
845 {
846   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
847   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
848   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,zero = 0.0;
849   PetscErrorCode ierr;
850   PetscInt       m = b->AIJ->m,n,i,*idx;
851 
852   PetscFunctionBegin;
853   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
854   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
855   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
856 
857   for (i=0; i<m; i++) {
858     idx    = a->j + a->i[i] ;
859     v      = a->a + a->i[i] ;
860     n      = a->i[i+1] - a->i[i];
861     alpha1 = x[6*i];
862     alpha2 = x[6*i+1];
863     alpha3 = x[6*i+2];
864     alpha4 = x[6*i+3];
865     alpha5 = x[6*i+4];
866     alpha6 = x[6*i+5];
867     while (n-->0) {
868       y[6*(*idx)]   += alpha1*(*v);
869       y[6*(*idx)+1] += alpha2*(*v);
870       y[6*(*idx)+2] += alpha3*(*v);
871       y[6*(*idx)+3] += alpha4*(*v);
872       y[6*(*idx)+4] += alpha5*(*v);
873       y[6*(*idx)+5] += alpha6*(*v);
874       idx++; v++;
875     }
876   }
877   ierr = PetscLogFlops(12*a->nz - 6*b->AIJ->n);CHKERRQ(ierr);
878   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
879   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
880   PetscFunctionReturn(0);
881 }
882 
883 #undef __FUNCT__
884 #define __FUNCT__ "MatMultAdd_SeqMAIJ_6"
885 PetscErrorCode MatMultAdd_SeqMAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
886 {
887   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
888   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
889   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6;
890   PetscErrorCode ierr;
891   PetscInt       m = b->AIJ->m,*idx,*ii;
892   PetscInt       n,i,jrow,j;
893 
894   PetscFunctionBegin;
895   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
896   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
897   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
898   idx  = a->j;
899   v    = a->a;
900   ii   = a->i;
901 
902   for (i=0; i<m; i++) {
903     jrow = ii[i];
904     n    = ii[i+1] - jrow;
905     sum1  = 0.0;
906     sum2  = 0.0;
907     sum3  = 0.0;
908     sum4  = 0.0;
909     sum5  = 0.0;
910     sum6  = 0.0;
911     for (j=0; j<n; j++) {
912       sum1 += v[jrow]*x[6*idx[jrow]];
913       sum2 += v[jrow]*x[6*idx[jrow]+1];
914       sum3 += v[jrow]*x[6*idx[jrow]+2];
915       sum4 += v[jrow]*x[6*idx[jrow]+3];
916       sum5 += v[jrow]*x[6*idx[jrow]+4];
917       sum6 += v[jrow]*x[6*idx[jrow]+5];
918       jrow++;
919      }
920     y[6*i]   += sum1;
921     y[6*i+1] += sum2;
922     y[6*i+2] += sum3;
923     y[6*i+3] += sum4;
924     y[6*i+4] += sum5;
925     y[6*i+5] += sum6;
926   }
927 
928   ierr = PetscLogFlops(12*a->nz);CHKERRQ(ierr);
929   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
930   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
931   PetscFunctionReturn(0);
932 }
933 
934 #undef __FUNCT__
935 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_6"
936 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
937 {
938   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
939   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
940   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6;
941   PetscErrorCode ierr;
942   PetscInt       m = b->AIJ->m,n,i,*idx;
943 
944   PetscFunctionBegin;
945   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
946   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
947   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
948 
949   for (i=0; i<m; i++) {
950     idx    = a->j + a->i[i] ;
951     v      = a->a + a->i[i] ;
952     n      = a->i[i+1] - a->i[i];
953     alpha1 = x[6*i];
954     alpha2 = x[6*i+1];
955     alpha3 = x[6*i+2];
956     alpha4 = x[6*i+3];
957     alpha5 = x[6*i+4];
958     alpha6 = x[6*i+5];
959     while (n-->0) {
960       y[6*(*idx)]   += alpha1*(*v);
961       y[6*(*idx)+1] += alpha2*(*v);
962       y[6*(*idx)+2] += alpha3*(*v);
963       y[6*(*idx)+3] += alpha4*(*v);
964       y[6*(*idx)+4] += alpha5*(*v);
965       y[6*(*idx)+5] += alpha6*(*v);
966       idx++; v++;
967     }
968   }
969   ierr = PetscLogFlops(12*a->nz);CHKERRQ(ierr);
970   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
971   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
972   PetscFunctionReturn(0);
973 }
974 
975 /* ------------------------------------------------------------------------------*/
976 #undef __FUNCT__
977 #define __FUNCT__ "MatMult_SeqMAIJ_7"
978 PetscErrorCode MatMult_SeqMAIJ_7(Mat A,Vec xx,Vec yy)
979 {
980   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
981   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
982   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7;
983   PetscErrorCode ierr;
984   PetscInt       m = b->AIJ->m,*idx,*ii;
985   PetscInt       n,i,jrow,j;
986 
987   PetscFunctionBegin;
988   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
989   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
990   idx  = a->j;
991   v    = a->a;
992   ii   = a->i;
993 
994   for (i=0; i<m; i++) {
995     jrow = ii[i];
996     n    = ii[i+1] - jrow;
997     sum1  = 0.0;
998     sum2  = 0.0;
999     sum3  = 0.0;
1000     sum4  = 0.0;
1001     sum5  = 0.0;
1002     sum6  = 0.0;
1003     sum7  = 0.0;
1004     for (j=0; j<n; j++) {
1005       sum1 += v[jrow]*x[7*idx[jrow]];
1006       sum2 += v[jrow]*x[7*idx[jrow]+1];
1007       sum3 += v[jrow]*x[7*idx[jrow]+2];
1008       sum4 += v[jrow]*x[7*idx[jrow]+3];
1009       sum5 += v[jrow]*x[7*idx[jrow]+4];
1010       sum6 += v[jrow]*x[7*idx[jrow]+5];
1011       sum7 += v[jrow]*x[7*idx[jrow]+6];
1012       jrow++;
1013      }
1014     y[7*i]   = sum1;
1015     y[7*i+1] = sum2;
1016     y[7*i+2] = sum3;
1017     y[7*i+3] = sum4;
1018     y[7*i+4] = sum5;
1019     y[7*i+5] = sum6;
1020     y[7*i+6] = sum7;
1021   }
1022 
1023   ierr = PetscLogFlops(14*a->nz - 7*m);CHKERRQ(ierr);
1024   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1025   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 #undef __FUNCT__
1030 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_7"
1031 PetscErrorCode MatMultTranspose_SeqMAIJ_7(Mat A,Vec xx,Vec yy)
1032 {
1033   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1034   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1035   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,zero = 0.0;
1036   PetscErrorCode ierr;
1037   PetscInt       m = b->AIJ->m,n,i,*idx;
1038 
1039   PetscFunctionBegin;
1040   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
1041   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1042   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
1043 
1044   for (i=0; i<m; i++) {
1045     idx    = a->j + a->i[i] ;
1046     v      = a->a + a->i[i] ;
1047     n      = a->i[i+1] - a->i[i];
1048     alpha1 = x[7*i];
1049     alpha2 = x[7*i+1];
1050     alpha3 = x[7*i+2];
1051     alpha4 = x[7*i+3];
1052     alpha5 = x[7*i+4];
1053     alpha6 = x[7*i+5];
1054     alpha7 = x[7*i+6];
1055     while (n-->0) {
1056       y[7*(*idx)]   += alpha1*(*v);
1057       y[7*(*idx)+1] += alpha2*(*v);
1058       y[7*(*idx)+2] += alpha3*(*v);
1059       y[7*(*idx)+3] += alpha4*(*v);
1060       y[7*(*idx)+4] += alpha5*(*v);
1061       y[7*(*idx)+5] += alpha6*(*v);
1062       y[7*(*idx)+6] += alpha7*(*v);
1063       idx++; v++;
1064     }
1065   }
1066   ierr = PetscLogFlops(14*a->nz - 7*b->AIJ->n);CHKERRQ(ierr);
1067   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1068   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1069   PetscFunctionReturn(0);
1070 }
1071 
1072 #undef __FUNCT__
1073 #define __FUNCT__ "MatMultAdd_SeqMAIJ_7"
1074 PetscErrorCode MatMultAdd_SeqMAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
1075 {
1076   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1077   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1078   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7;
1079   PetscErrorCode ierr;
1080   PetscInt       m = b->AIJ->m,*idx,*ii;
1081   PetscInt       n,i,jrow,j;
1082 
1083   PetscFunctionBegin;
1084   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
1085   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1086   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
1087   idx  = a->j;
1088   v    = a->a;
1089   ii   = a->i;
1090 
1091   for (i=0; i<m; i++) {
1092     jrow = ii[i];
1093     n    = ii[i+1] - jrow;
1094     sum1  = 0.0;
1095     sum2  = 0.0;
1096     sum3  = 0.0;
1097     sum4  = 0.0;
1098     sum5  = 0.0;
1099     sum6  = 0.0;
1100     sum7  = 0.0;
1101     for (j=0; j<n; j++) {
1102       sum1 += v[jrow]*x[7*idx[jrow]];
1103       sum2 += v[jrow]*x[7*idx[jrow]+1];
1104       sum3 += v[jrow]*x[7*idx[jrow]+2];
1105       sum4 += v[jrow]*x[7*idx[jrow]+3];
1106       sum5 += v[jrow]*x[7*idx[jrow]+4];
1107       sum6 += v[jrow]*x[7*idx[jrow]+5];
1108       sum7 += v[jrow]*x[7*idx[jrow]+6];
1109       jrow++;
1110      }
1111     y[7*i]   += sum1;
1112     y[7*i+1] += sum2;
1113     y[7*i+2] += sum3;
1114     y[7*i+3] += sum4;
1115     y[7*i+4] += sum5;
1116     y[7*i+5] += sum6;
1117     y[7*i+6] += sum7;
1118   }
1119 
1120   ierr = PetscLogFlops(14*a->nz);CHKERRQ(ierr);
1121   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1122   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 #undef __FUNCT__
1127 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_7"
1128 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
1129 {
1130   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1131   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1132   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7;
1133   PetscErrorCode ierr;
1134   PetscInt       m = b->AIJ->m,n,i,*idx;
1135 
1136   PetscFunctionBegin;
1137   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
1138   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1139   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
1140   for (i=0; i<m; i++) {
1141     idx    = a->j + a->i[i] ;
1142     v      = a->a + a->i[i] ;
1143     n      = a->i[i+1] - a->i[i];
1144     alpha1 = x[7*i];
1145     alpha2 = x[7*i+1];
1146     alpha3 = x[7*i+2];
1147     alpha4 = x[7*i+3];
1148     alpha5 = x[7*i+4];
1149     alpha6 = x[7*i+5];
1150     alpha7 = x[7*i+6];
1151     while (n-->0) {
1152       y[7*(*idx)]   += alpha1*(*v);
1153       y[7*(*idx)+1] += alpha2*(*v);
1154       y[7*(*idx)+2] += alpha3*(*v);
1155       y[7*(*idx)+3] += alpha4*(*v);
1156       y[7*(*idx)+4] += alpha5*(*v);
1157       y[7*(*idx)+5] += alpha6*(*v);
1158       y[7*(*idx)+6] += alpha7*(*v);
1159       idx++; v++;
1160     }
1161   }
1162   ierr = PetscLogFlops(14*a->nz);CHKERRQ(ierr);
1163   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1164   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
1165   PetscFunctionReturn(0);
1166 }
1167 
1168 #undef __FUNCT__
1169 #define __FUNCT__ "MatMult_SeqMAIJ_8"
1170 PetscErrorCode MatMult_SeqMAIJ_8(Mat A,Vec xx,Vec yy)
1171 {
1172   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1173   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1174   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7, sum8;
1175   PetscErrorCode ierr;
1176   PetscInt       m = b->AIJ->m,*idx,*ii;
1177   PetscInt       n,i,jrow,j;
1178 
1179   PetscFunctionBegin;
1180   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1181   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
1182   idx  = a->j;
1183   v    = a->a;
1184   ii   = a->i;
1185 
1186   for (i=0; i<m; i++) {
1187     jrow = ii[i];
1188     n    = ii[i+1] - jrow;
1189     sum1  = 0.0;
1190     sum2  = 0.0;
1191     sum3  = 0.0;
1192     sum4  = 0.0;
1193     sum5  = 0.0;
1194     sum6  = 0.0;
1195     sum7  = 0.0;
1196     sum8  = 0.0;
1197     for (j=0; j<n; j++) {
1198       sum1 += v[jrow]*x[8*idx[jrow]];
1199       sum2 += v[jrow]*x[8*idx[jrow]+1];
1200       sum3 += v[jrow]*x[8*idx[jrow]+2];
1201       sum4 += v[jrow]*x[8*idx[jrow]+3];
1202       sum5 += v[jrow]*x[8*idx[jrow]+4];
1203       sum6 += v[jrow]*x[8*idx[jrow]+5];
1204       sum7 += v[jrow]*x[8*idx[jrow]+6];
1205       sum8 += v[jrow]*x[8*idx[jrow]+7];
1206       jrow++;
1207      }
1208     y[8*i]   = sum1;
1209     y[8*i+1] = sum2;
1210     y[8*i+2] = sum3;
1211     y[8*i+3] = sum4;
1212     y[8*i+4] = sum5;
1213     y[8*i+5] = sum6;
1214     y[8*i+6] = sum7;
1215     y[8*i+7] = sum8;
1216   }
1217 
1218   ierr = PetscLogFlops(16*a->nz - 8*m);CHKERRQ(ierr);
1219   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1220   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1221   PetscFunctionReturn(0);
1222 }
1223 
1224 #undef __FUNCT__
1225 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_8"
1226 PetscErrorCode MatMultTranspose_SeqMAIJ_8(Mat A,Vec xx,Vec yy)
1227 {
1228   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1229   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1230   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,alpha8,zero = 0.0;
1231   PetscErrorCode ierr;
1232   PetscInt       m = b->AIJ->m,n,i,*idx;
1233 
1234   PetscFunctionBegin;
1235   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
1236   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1237   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
1238 
1239   for (i=0; i<m; i++) {
1240     idx    = a->j + a->i[i] ;
1241     v      = a->a + a->i[i] ;
1242     n      = a->i[i+1] - a->i[i];
1243     alpha1 = x[8*i];
1244     alpha2 = x[8*i+1];
1245     alpha3 = x[8*i+2];
1246     alpha4 = x[8*i+3];
1247     alpha5 = x[8*i+4];
1248     alpha6 = x[8*i+5];
1249     alpha7 = x[8*i+6];
1250     alpha8 = x[8*i+7];
1251     while (n-->0) {
1252       y[8*(*idx)]   += alpha1*(*v);
1253       y[8*(*idx)+1] += alpha2*(*v);
1254       y[8*(*idx)+2] += alpha3*(*v);
1255       y[8*(*idx)+3] += alpha4*(*v);
1256       y[8*(*idx)+4] += alpha5*(*v);
1257       y[8*(*idx)+5] += alpha6*(*v);
1258       y[8*(*idx)+6] += alpha7*(*v);
1259       y[8*(*idx)+7] += alpha8*(*v);
1260       idx++; v++;
1261     }
1262   }
1263   ierr = PetscLogFlops(16*a->nz - 8*b->AIJ->n);CHKERRQ(ierr);
1264   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1265   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1266   PetscFunctionReturn(0);
1267 }
1268 
1269 #undef __FUNCT__
1270 #define __FUNCT__ "MatMultAdd_SeqMAIJ_8"
1271 PetscErrorCode MatMultAdd_SeqMAIJ_8(Mat A,Vec xx,Vec yy,Vec zz)
1272 {
1273   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1274   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1275   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7, sum8;
1276   PetscErrorCode ierr;
1277   PetscInt       m = b->AIJ->m,*idx,*ii;
1278   PetscInt       n,i,jrow,j;
1279 
1280   PetscFunctionBegin;
1281   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
1282   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1283   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
1284   idx  = a->j;
1285   v    = a->a;
1286   ii   = a->i;
1287 
1288   for (i=0; i<m; i++) {
1289     jrow = ii[i];
1290     n    = ii[i+1] - jrow;
1291     sum1  = 0.0;
1292     sum2  = 0.0;
1293     sum3  = 0.0;
1294     sum4  = 0.0;
1295     sum5  = 0.0;
1296     sum6  = 0.0;
1297     sum7  = 0.0;
1298     sum8  = 0.0;
1299     for (j=0; j<n; j++) {
1300       sum1 += v[jrow]*x[8*idx[jrow]];
1301       sum2 += v[jrow]*x[8*idx[jrow]+1];
1302       sum3 += v[jrow]*x[8*idx[jrow]+2];
1303       sum4 += v[jrow]*x[8*idx[jrow]+3];
1304       sum5 += v[jrow]*x[8*idx[jrow]+4];
1305       sum6 += v[jrow]*x[8*idx[jrow]+5];
1306       sum7 += v[jrow]*x[8*idx[jrow]+6];
1307       sum8 += v[jrow]*x[8*idx[jrow]+7];
1308       jrow++;
1309      }
1310     y[8*i]   += sum1;
1311     y[8*i+1] += sum2;
1312     y[8*i+2] += sum3;
1313     y[8*i+3] += sum4;
1314     y[8*i+4] += sum5;
1315     y[8*i+5] += sum6;
1316     y[8*i+6] += sum7;
1317     y[8*i+7] += sum8;
1318   }
1319 
1320   ierr = PetscLogFlops(16*a->nz);CHKERRQ(ierr);
1321   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1322   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
1323   PetscFunctionReturn(0);
1324 }
1325 
1326 #undef __FUNCT__
1327 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_8"
1328 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_8(Mat A,Vec xx,Vec yy,Vec zz)
1329 {
1330   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1331   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1332   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,alpha8;
1333   PetscErrorCode ierr;
1334   PetscInt       m = b->AIJ->m,n,i,*idx;
1335 
1336   PetscFunctionBegin;
1337   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
1338   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1339   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
1340   for (i=0; i<m; i++) {
1341     idx    = a->j + a->i[i] ;
1342     v      = a->a + a->i[i] ;
1343     n      = a->i[i+1] - a->i[i];
1344     alpha1 = x[8*i];
1345     alpha2 = x[8*i+1];
1346     alpha3 = x[8*i+2];
1347     alpha4 = x[8*i+3];
1348     alpha5 = x[8*i+4];
1349     alpha6 = x[8*i+5];
1350     alpha7 = x[8*i+6];
1351     alpha8 = x[8*i+7];
1352     while (n-->0) {
1353       y[8*(*idx)]   += alpha1*(*v);
1354       y[8*(*idx)+1] += alpha2*(*v);
1355       y[8*(*idx)+2] += alpha3*(*v);
1356       y[8*(*idx)+3] += alpha4*(*v);
1357       y[8*(*idx)+4] += alpha5*(*v);
1358       y[8*(*idx)+5] += alpha6*(*v);
1359       y[8*(*idx)+6] += alpha7*(*v);
1360       y[8*(*idx)+7] += alpha8*(*v);
1361       idx++; v++;
1362     }
1363   }
1364   ierr = PetscLogFlops(16*a->nz);CHKERRQ(ierr);
1365   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1366   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
1367   PetscFunctionReturn(0);
1368 }
1369 
1370 
1371 /* ------------------------------------------------------------------------------*/
1372 #undef __FUNCT__
1373 #define __FUNCT__ "MatMult_SeqMAIJ_9"
1374 PetscErrorCode MatMult_SeqMAIJ_9(Mat A,Vec xx,Vec yy)
1375 {
1376   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1377   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1378   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7, sum8, sum9;
1379   PetscErrorCode ierr;
1380   PetscInt       m = b->AIJ->m,*idx,*ii;
1381   PetscInt       n,i,jrow,j;
1382 
1383   PetscFunctionBegin;
1384   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1385   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
1386   idx  = a->j;
1387   v    = a->a;
1388   ii   = a->i;
1389 
1390   for (i=0; i<m; i++) {
1391     jrow = ii[i];
1392     n    = ii[i+1] - jrow;
1393     sum1  = 0.0;
1394     sum2  = 0.0;
1395     sum3  = 0.0;
1396     sum4  = 0.0;
1397     sum5  = 0.0;
1398     sum6  = 0.0;
1399     sum7  = 0.0;
1400     sum8  = 0.0;
1401     sum9  = 0.0;
1402     for (j=0; j<n; j++) {
1403       sum1 += v[jrow]*x[9*idx[jrow]];
1404       sum2 += v[jrow]*x[9*idx[jrow]+1];
1405       sum3 += v[jrow]*x[9*idx[jrow]+2];
1406       sum4 += v[jrow]*x[9*idx[jrow]+3];
1407       sum5 += v[jrow]*x[9*idx[jrow]+4];
1408       sum6 += v[jrow]*x[9*idx[jrow]+5];
1409       sum7 += v[jrow]*x[9*idx[jrow]+6];
1410       sum8 += v[jrow]*x[9*idx[jrow]+7];
1411       sum9 += v[jrow]*x[9*idx[jrow]+8];
1412       jrow++;
1413      }
1414     y[9*i]   = sum1;
1415     y[9*i+1] = sum2;
1416     y[9*i+2] = sum3;
1417     y[9*i+3] = sum4;
1418     y[9*i+4] = sum5;
1419     y[9*i+5] = sum6;
1420     y[9*i+6] = sum7;
1421     y[9*i+7] = sum8;
1422     y[9*i+8] = sum9;
1423   }
1424 
1425   ierr = PetscLogFlops(18*a->nz - 9*m);CHKERRQ(ierr);
1426   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1427   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1428   PetscFunctionReturn(0);
1429 }
1430 
1431 #undef __FUNCT__
1432 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_9"
1433 PetscErrorCode MatMultTranspose_SeqMAIJ_9(Mat A,Vec xx,Vec yy)
1434 {
1435   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1436   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1437   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,alpha8,alpha9,zero = 0.0;
1438   PetscErrorCode ierr;
1439   PetscInt       m = b->AIJ->m,n,i,*idx;
1440 
1441   PetscFunctionBegin;
1442   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
1443   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1444   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
1445 
1446   for (i=0; i<m; i++) {
1447     idx    = a->j + a->i[i] ;
1448     v      = a->a + a->i[i] ;
1449     n      = a->i[i+1] - a->i[i];
1450     alpha1 = x[9*i];
1451     alpha2 = x[9*i+1];
1452     alpha3 = x[9*i+2];
1453     alpha4 = x[9*i+3];
1454     alpha5 = x[9*i+4];
1455     alpha6 = x[9*i+5];
1456     alpha7 = x[9*i+6];
1457     alpha8 = x[9*i+7];
1458     alpha9 = x[9*i+8];
1459     while (n-->0) {
1460       y[9*(*idx)]   += alpha1*(*v);
1461       y[9*(*idx)+1] += alpha2*(*v);
1462       y[9*(*idx)+2] += alpha3*(*v);
1463       y[9*(*idx)+3] += alpha4*(*v);
1464       y[9*(*idx)+4] += alpha5*(*v);
1465       y[9*(*idx)+5] += alpha6*(*v);
1466       y[9*(*idx)+6] += alpha7*(*v);
1467       y[9*(*idx)+7] += alpha8*(*v);
1468       y[9*(*idx)+8] += alpha9*(*v);
1469       idx++; v++;
1470     }
1471   }
1472   ierr = PetscLogFlops(18*a->nz - 9*b->AIJ->n);CHKERRQ(ierr);
1473   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1474   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1475   PetscFunctionReturn(0);
1476 }
1477 
1478 #undef __FUNCT__
1479 #define __FUNCT__ "MatMultAdd_SeqMAIJ_9"
1480 PetscErrorCode MatMultAdd_SeqMAIJ_9(Mat A,Vec xx,Vec yy,Vec zz)
1481 {
1482   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1483   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1484   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7, sum8, sum9;
1485   PetscErrorCode ierr;
1486   PetscInt       m = b->AIJ->m,*idx,*ii;
1487   PetscInt       n,i,jrow,j;
1488 
1489   PetscFunctionBegin;
1490   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
1491   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1492   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
1493   idx  = a->j;
1494   v    = a->a;
1495   ii   = a->i;
1496 
1497   for (i=0; i<m; i++) {
1498     jrow = ii[i];
1499     n    = ii[i+1] - jrow;
1500     sum1  = 0.0;
1501     sum2  = 0.0;
1502     sum3  = 0.0;
1503     sum4  = 0.0;
1504     sum5  = 0.0;
1505     sum6  = 0.0;
1506     sum7  = 0.0;
1507     sum8  = 0.0;
1508     sum9  = 0.0;
1509     for (j=0; j<n; j++) {
1510       sum1 += v[jrow]*x[9*idx[jrow]];
1511       sum2 += v[jrow]*x[9*idx[jrow]+1];
1512       sum3 += v[jrow]*x[9*idx[jrow]+2];
1513       sum4 += v[jrow]*x[9*idx[jrow]+3];
1514       sum5 += v[jrow]*x[9*idx[jrow]+4];
1515       sum6 += v[jrow]*x[9*idx[jrow]+5];
1516       sum7 += v[jrow]*x[9*idx[jrow]+6];
1517       sum8 += v[jrow]*x[9*idx[jrow]+7];
1518       sum9 += v[jrow]*x[9*idx[jrow]+8];
1519       jrow++;
1520      }
1521     y[9*i]   += sum1;
1522     y[9*i+1] += sum2;
1523     y[9*i+2] += sum3;
1524     y[9*i+3] += sum4;
1525     y[9*i+4] += sum5;
1526     y[9*i+5] += sum6;
1527     y[9*i+6] += sum7;
1528     y[9*i+7] += sum8;
1529     y[9*i+8] += sum9;
1530   }
1531 
1532   ierr = PetscLogFlops(18*a->nz);CHKERRQ(ierr);
1533   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1534   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
1535   PetscFunctionReturn(0);
1536 }
1537 
1538 #undef __FUNCT__
1539 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_9"
1540 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_9(Mat A,Vec xx,Vec yy,Vec zz)
1541 {
1542   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1543   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1544   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,alpha8,alpha9;
1545   PetscErrorCode ierr;
1546   PetscInt       m = b->AIJ->m,n,i,*idx;
1547 
1548   PetscFunctionBegin;
1549   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
1550   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1551   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
1552   for (i=0; i<m; i++) {
1553     idx    = a->j + a->i[i] ;
1554     v      = a->a + a->i[i] ;
1555     n      = a->i[i+1] - a->i[i];
1556     alpha1 = x[9*i];
1557     alpha2 = x[9*i+1];
1558     alpha3 = x[9*i+2];
1559     alpha4 = x[9*i+3];
1560     alpha5 = x[9*i+4];
1561     alpha6 = x[9*i+5];
1562     alpha7 = x[9*i+6];
1563     alpha8 = x[9*i+7];
1564     alpha9 = x[9*i+8];
1565     while (n-->0) {
1566       y[9*(*idx)]   += alpha1*(*v);
1567       y[9*(*idx)+1] += alpha2*(*v);
1568       y[9*(*idx)+2] += alpha3*(*v);
1569       y[9*(*idx)+3] += alpha4*(*v);
1570       y[9*(*idx)+4] += alpha5*(*v);
1571       y[9*(*idx)+5] += alpha6*(*v);
1572       y[9*(*idx)+6] += alpha7*(*v);
1573       y[9*(*idx)+7] += alpha8*(*v);
1574       y[9*(*idx)+8] += alpha9*(*v);
1575       idx++; v++;
1576     }
1577   }
1578   ierr = PetscLogFlops(18*a->nz);CHKERRQ(ierr);
1579   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1580   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
1581   PetscFunctionReturn(0);
1582 }
1583 
1584 /*--------------------------------------------------------------------------------------------*/
1585 #undef __FUNCT__
1586 #define __FUNCT__ "MatMult_SeqMAIJ_16"
1587 PetscErrorCode MatMult_SeqMAIJ_16(Mat A,Vec xx,Vec yy)
1588 {
1589   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1590   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1591   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7, sum8;
1592   PetscScalar    sum9, sum10, sum11, sum12, sum13, sum14, sum15, sum16;
1593   PetscErrorCode ierr;
1594   PetscInt       m = b->AIJ->m,*idx,*ii;
1595   PetscInt       n,i,jrow,j;
1596 
1597   PetscFunctionBegin;
1598   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1599   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
1600   idx  = a->j;
1601   v    = a->a;
1602   ii   = a->i;
1603 
1604   for (i=0; i<m; i++) {
1605     jrow = ii[i];
1606     n    = ii[i+1] - jrow;
1607     sum1  = 0.0;
1608     sum2  = 0.0;
1609     sum3  = 0.0;
1610     sum4  = 0.0;
1611     sum5  = 0.0;
1612     sum6  = 0.0;
1613     sum7  = 0.0;
1614     sum8  = 0.0;
1615     sum9  = 0.0;
1616     sum10 = 0.0;
1617     sum11 = 0.0;
1618     sum12 = 0.0;
1619     sum13 = 0.0;
1620     sum14 = 0.0;
1621     sum15 = 0.0;
1622     sum16 = 0.0;
1623     for (j=0; j<n; j++) {
1624       sum1  += v[jrow]*x[16*idx[jrow]];
1625       sum2  += v[jrow]*x[16*idx[jrow]+1];
1626       sum3  += v[jrow]*x[16*idx[jrow]+2];
1627       sum4  += v[jrow]*x[16*idx[jrow]+3];
1628       sum5  += v[jrow]*x[16*idx[jrow]+4];
1629       sum6  += v[jrow]*x[16*idx[jrow]+5];
1630       sum7  += v[jrow]*x[16*idx[jrow]+6];
1631       sum8  += v[jrow]*x[16*idx[jrow]+7];
1632       sum9  += v[jrow]*x[16*idx[jrow]+8];
1633       sum10 += v[jrow]*x[16*idx[jrow]+9];
1634       sum11 += v[jrow]*x[16*idx[jrow]+10];
1635       sum12 += v[jrow]*x[16*idx[jrow]+11];
1636       sum13 += v[jrow]*x[16*idx[jrow]+12];
1637       sum14 += v[jrow]*x[16*idx[jrow]+13];
1638       sum15 += v[jrow]*x[16*idx[jrow]+14];
1639       sum16 += v[jrow]*x[16*idx[jrow]+15];
1640       jrow++;
1641      }
1642     y[16*i]    = sum1;
1643     y[16*i+1]  = sum2;
1644     y[16*i+2]  = sum3;
1645     y[16*i+3]  = sum4;
1646     y[16*i+4]  = sum5;
1647     y[16*i+5]  = sum6;
1648     y[16*i+6]  = sum7;
1649     y[16*i+7]  = sum8;
1650     y[16*i+8]  = sum9;
1651     y[16*i+9]  = sum10;
1652     y[16*i+10] = sum11;
1653     y[16*i+11] = sum12;
1654     y[16*i+12] = sum13;
1655     y[16*i+13] = sum14;
1656     y[16*i+14] = sum15;
1657     y[16*i+15] = sum16;
1658   }
1659 
1660   ierr = PetscLogFlops(32*a->nz - 16*m);CHKERRQ(ierr);
1661   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1662   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1663   PetscFunctionReturn(0);
1664 }
1665 
1666 #undef __FUNCT__
1667 #define __FUNCT__ "MatMultTranspose_SeqMAIJ_16"
1668 PetscErrorCode MatMultTranspose_SeqMAIJ_16(Mat A,Vec xx,Vec yy)
1669 {
1670   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1671   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1672   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,alpha8,zero = 0.0;
1673   PetscScalar    alpha9,alpha10,alpha11,alpha12,alpha13,alpha14,alpha15,alpha16;
1674   PetscErrorCode ierr;
1675   PetscInt       m = b->AIJ->m,n,i,*idx;
1676 
1677   PetscFunctionBegin;
1678   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
1679   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1680   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
1681 
1682   for (i=0; i<m; i++) {
1683     idx    = a->j + a->i[i] ;
1684     v      = a->a + a->i[i] ;
1685     n      = a->i[i+1] - a->i[i];
1686     alpha1  = x[16*i];
1687     alpha2  = x[16*i+1];
1688     alpha3  = x[16*i+2];
1689     alpha4  = x[16*i+3];
1690     alpha5  = x[16*i+4];
1691     alpha6  = x[16*i+5];
1692     alpha7  = x[16*i+6];
1693     alpha8  = x[16*i+7];
1694     alpha9  = x[16*i+8];
1695     alpha10 = x[16*i+9];
1696     alpha11 = x[16*i+10];
1697     alpha12 = x[16*i+11];
1698     alpha13 = x[16*i+12];
1699     alpha14 = x[16*i+13];
1700     alpha15 = x[16*i+14];
1701     alpha16 = x[16*i+15];
1702     while (n-->0) {
1703       y[16*(*idx)]    += alpha1*(*v);
1704       y[16*(*idx)+1]  += alpha2*(*v);
1705       y[16*(*idx)+2]  += alpha3*(*v);
1706       y[16*(*idx)+3]  += alpha4*(*v);
1707       y[16*(*idx)+4]  += alpha5*(*v);
1708       y[16*(*idx)+5]  += alpha6*(*v);
1709       y[16*(*idx)+6]  += alpha7*(*v);
1710       y[16*(*idx)+7]  += alpha8*(*v);
1711       y[16*(*idx)+8]  += alpha9*(*v);
1712       y[16*(*idx)+9]  += alpha10*(*v);
1713       y[16*(*idx)+10] += alpha11*(*v);
1714       y[16*(*idx)+11] += alpha12*(*v);
1715       y[16*(*idx)+12] += alpha13*(*v);
1716       y[16*(*idx)+13] += alpha14*(*v);
1717       y[16*(*idx)+14] += alpha15*(*v);
1718       y[16*(*idx)+15] += alpha16*(*v);
1719       idx++; v++;
1720     }
1721   }
1722   ierr = PetscLogFlops(32*a->nz - 16*b->AIJ->n);CHKERRQ(ierr);
1723   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1724   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1725   PetscFunctionReturn(0);
1726 }
1727 
1728 #undef __FUNCT__
1729 #define __FUNCT__ "MatMultAdd_SeqMAIJ_16"
1730 PetscErrorCode MatMultAdd_SeqMAIJ_16(Mat A,Vec xx,Vec yy,Vec zz)
1731 {
1732   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1733   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1734   PetscScalar    *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7, sum8;
1735   PetscScalar    sum9, sum10, sum11, sum12, sum13, sum14, sum15, sum16;
1736   PetscErrorCode ierr;
1737   PetscInt       m = b->AIJ->m,*idx,*ii;
1738   PetscInt       n,i,jrow,j;
1739 
1740   PetscFunctionBegin;
1741   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
1742   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1743   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
1744   idx  = a->j;
1745   v    = a->a;
1746   ii   = a->i;
1747 
1748   for (i=0; i<m; i++) {
1749     jrow = ii[i];
1750     n    = ii[i+1] - jrow;
1751     sum1  = 0.0;
1752     sum2  = 0.0;
1753     sum3  = 0.0;
1754     sum4  = 0.0;
1755     sum5  = 0.0;
1756     sum6  = 0.0;
1757     sum7  = 0.0;
1758     sum8  = 0.0;
1759     sum9  = 0.0;
1760     sum10 = 0.0;
1761     sum11 = 0.0;
1762     sum12 = 0.0;
1763     sum13 = 0.0;
1764     sum14 = 0.0;
1765     sum15 = 0.0;
1766     sum16 = 0.0;
1767     for (j=0; j<n; j++) {
1768       sum1  += v[jrow]*x[16*idx[jrow]];
1769       sum2  += v[jrow]*x[16*idx[jrow]+1];
1770       sum3  += v[jrow]*x[16*idx[jrow]+2];
1771       sum4  += v[jrow]*x[16*idx[jrow]+3];
1772       sum5  += v[jrow]*x[16*idx[jrow]+4];
1773       sum6  += v[jrow]*x[16*idx[jrow]+5];
1774       sum7  += v[jrow]*x[16*idx[jrow]+6];
1775       sum8  += v[jrow]*x[16*idx[jrow]+7];
1776       sum9  += v[jrow]*x[16*idx[jrow]+8];
1777       sum10 += v[jrow]*x[16*idx[jrow]+9];
1778       sum11 += v[jrow]*x[16*idx[jrow]+10];
1779       sum12 += v[jrow]*x[16*idx[jrow]+11];
1780       sum13 += v[jrow]*x[16*idx[jrow]+12];
1781       sum14 += v[jrow]*x[16*idx[jrow]+13];
1782       sum15 += v[jrow]*x[16*idx[jrow]+14];
1783       sum16 += v[jrow]*x[16*idx[jrow]+15];
1784       jrow++;
1785      }
1786     y[16*i]    += sum1;
1787     y[16*i+1]  += sum2;
1788     y[16*i+2]  += sum3;
1789     y[16*i+3]  += sum4;
1790     y[16*i+4]  += sum5;
1791     y[16*i+5]  += sum6;
1792     y[16*i+6]  += sum7;
1793     y[16*i+7]  += sum8;
1794     y[16*i+8]  += sum9;
1795     y[16*i+9]  += sum10;
1796     y[16*i+10] += sum11;
1797     y[16*i+11] += sum12;
1798     y[16*i+12] += sum13;
1799     y[16*i+13] += sum14;
1800     y[16*i+14] += sum15;
1801     y[16*i+15] += sum16;
1802   }
1803 
1804   ierr = PetscLogFlops(32*a->nz);CHKERRQ(ierr);
1805   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1806   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
1807   PetscFunctionReturn(0);
1808 }
1809 
1810 #undef __FUNCT__
1811 #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_16"
1812 PetscErrorCode MatMultTransposeAdd_SeqMAIJ_16(Mat A,Vec xx,Vec yy,Vec zz)
1813 {
1814   Mat_SeqMAIJ    *b = (Mat_SeqMAIJ*)A->data;
1815   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)b->AIJ->data;
1816   PetscScalar    *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,alpha8;
1817   PetscScalar    alpha9,alpha10,alpha11,alpha12,alpha13,alpha14,alpha15,alpha16;
1818   PetscErrorCode ierr;
1819   PetscInt       m = b->AIJ->m,n,i,*idx;
1820 
1821   PetscFunctionBegin;
1822   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
1823   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1824   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
1825   for (i=0; i<m; i++) {
1826     idx    = a->j + a->i[i] ;
1827     v      = a->a + a->i[i] ;
1828     n      = a->i[i+1] - a->i[i];
1829     alpha1 = x[16*i];
1830     alpha2 = x[16*i+1];
1831     alpha3 = x[16*i+2];
1832     alpha4 = x[16*i+3];
1833     alpha5 = x[16*i+4];
1834     alpha6 = x[16*i+5];
1835     alpha7 = x[16*i+6];
1836     alpha8 = x[16*i+7];
1837     alpha9  = x[16*i+8];
1838     alpha10 = x[16*i+9];
1839     alpha11 = x[16*i+10];
1840     alpha12 = x[16*i+11];
1841     alpha13 = x[16*i+12];
1842     alpha14 = x[16*i+13];
1843     alpha15 = x[16*i+14];
1844     alpha16 = x[16*i+15];
1845     while (n-->0) {
1846       y[16*(*idx)]   += alpha1*(*v);
1847       y[16*(*idx)+1] += alpha2*(*v);
1848       y[16*(*idx)+2] += alpha3*(*v);
1849       y[16*(*idx)+3] += alpha4*(*v);
1850       y[16*(*idx)+4] += alpha5*(*v);
1851       y[16*(*idx)+5] += alpha6*(*v);
1852       y[16*(*idx)+6] += alpha7*(*v);
1853       y[16*(*idx)+7] += alpha8*(*v);
1854       y[16*(*idx)+8]  += alpha9*(*v);
1855       y[16*(*idx)+9]  += alpha10*(*v);
1856       y[16*(*idx)+10] += alpha11*(*v);
1857       y[16*(*idx)+11] += alpha12*(*v);
1858       y[16*(*idx)+12] += alpha13*(*v);
1859       y[16*(*idx)+13] += alpha14*(*v);
1860       y[16*(*idx)+14] += alpha15*(*v);
1861       y[16*(*idx)+15] += alpha16*(*v);
1862       idx++; v++;
1863     }
1864   }
1865   ierr = PetscLogFlops(32*a->nz);CHKERRQ(ierr);
1866   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1867   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
1868   PetscFunctionReturn(0);
1869 }
1870 
1871 /*===================================================================================*/
1872 #undef __FUNCT__
1873 #define __FUNCT__ "MatMult_MPIMAIJ_dof"
1874 PetscErrorCode MatMult_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
1875 {
1876   Mat_MPIMAIJ    *b = (Mat_MPIMAIJ*)A->data;
1877   PetscErrorCode ierr;
1878 
1879   PetscFunctionBegin;
1880   /* start the scatter */
1881   ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
1882   ierr = (*b->AIJ->ops->mult)(b->AIJ,xx,yy);CHKERRQ(ierr);
1883   ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
1884   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr);
1885   PetscFunctionReturn(0);
1886 }
1887 
1888 #undef __FUNCT__
1889 #define __FUNCT__ "MatMultTranspose_MPIMAIJ_dof"
1890 PetscErrorCode MatMultTranspose_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
1891 {
1892   Mat_MPIMAIJ    *b = (Mat_MPIMAIJ*)A->data;
1893   PetscErrorCode ierr;
1894 
1895   PetscFunctionBegin;
1896   ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr);
1897   ierr = (*b->AIJ->ops->multtranspose)(b->AIJ,xx,yy);CHKERRQ(ierr);
1898   ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
1899   ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
1900   PetscFunctionReturn(0);
1901 }
1902 
1903 #undef __FUNCT__
1904 #define __FUNCT__ "MatMultAdd_MPIMAIJ_dof"
1905 PetscErrorCode MatMultAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz)
1906 {
1907   Mat_MPIMAIJ    *b = (Mat_MPIMAIJ*)A->data;
1908   PetscErrorCode ierr;
1909 
1910   PetscFunctionBegin;
1911   /* start the scatter */
1912   ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
1913   ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr);
1914   ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
1915   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,zz,zz);CHKERRQ(ierr);
1916   PetscFunctionReturn(0);
1917 }
1918 
1919 #undef __FUNCT__
1920 #define __FUNCT__ "MatMultTransposeAdd_MPIMAIJ_dof"
1921 PetscErrorCode MatMultTransposeAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz)
1922 {
1923   Mat_MPIMAIJ    *b = (Mat_MPIMAIJ*)A->data;
1924   PetscErrorCode ierr;
1925 
1926   PetscFunctionBegin;
1927   ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr);
1928   ierr = VecScatterBegin(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
1929   ierr = (*b->AIJ->ops->multtransposeadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr);
1930   ierr = VecScatterEnd(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
1931   PetscFunctionReturn(0);
1932 }
1933 
1934 #undef __FUNCT__
1935 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ"
1936 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,PetscReal fill,Mat *C)
1937 {
1938   /* This routine requires testing -- but it's getting better. */
1939   PetscErrorCode ierr;
1940   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
1941   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
1942   Mat            P=pp->AIJ;
1943   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
1944   PetscInt       *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
1945   PetscInt       *ci,*cj,*ptadenserow,*ptasparserow,*denserow,*sparserow,*ptaj;
1946   PetscInt       an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof,cn;
1947   PetscInt       i,j,k,dof,pshift,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
1948   MatScalar      *ca;
1949 
1950   PetscFunctionBegin;
1951   /* Start timer */
1952   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
1953 
1954   /* Get ij structure of P^T */
1955   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
1956 
1957   cn = pn*ppdof;
1958   /* Allocate ci array, arrays for fill computation and */
1959   /* free space for accumulating nonzero column info */
1960   ierr = PetscMalloc((cn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
1961   ci[0] = 0;
1962 
1963   /* Work arrays for rows of P^T*A */
1964   ierr = PetscMalloc((2*cn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
1965   ierr = PetscMemzero(ptadenserow,(2*cn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
1966   ptasparserow = ptadenserow  + an;
1967   denserow     = ptasparserow + an;
1968   sparserow    = denserow     + cn;
1969 
1970   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
1971   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
1972   /* Note, aspect ratio of P is the same as the aspect ratio of SeqAIJ inside P */
1973   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
1974   current_space = free_space;
1975 
1976   /* Determine symbolic info for each row of C: */
1977   for (i=0;i<pn;i++) {
1978     ptnzi  = pti[i+1] - pti[i];
1979     ptJ    = ptj + pti[i];
1980     for (dof=0;dof<ppdof;dof++) {
1981       ptanzi = 0;
1982       /* Determine symbolic row of PtA: */
1983       for (j=0;j<ptnzi;j++) {
1984         /* Expand ptJ[j] by block size and shift by dof to get the right row of A */
1985         arow = ptJ[j]*ppdof + dof;
1986         /* Nonzeros of P^T*A will be in same locations as any element of A in that row */
1987         anzj = ai[arow+1] - ai[arow];
1988         ajj  = aj + ai[arow];
1989         for (k=0;k<anzj;k++) {
1990           if (!ptadenserow[ajj[k]]) {
1991             ptadenserow[ajj[k]]    = -1;
1992             ptasparserow[ptanzi++] = ajj[k];
1993           }
1994         }
1995       }
1996       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
1997       ptaj = ptasparserow;
1998       cnzi   = 0;
1999       for (j=0;j<ptanzi;j++) {
2000         /* Get offset within block of P */
2001         pshift = *ptaj%ppdof;
2002         /* Get block row of P */
2003         prow = (*ptaj++)/ppdof; /* integer division */
2004         /* P has same number of nonzeros per row as the compressed form */
2005         pnzj = pi[prow+1] - pi[prow];
2006         pjj  = pj + pi[prow];
2007         for (k=0;k<pnzj;k++) {
2008           /* Locations in C are shifted by the offset within the block */
2009           /* Note: we cannot use PetscLLAdd here because of the additional offset for the write location */
2010           if (!denserow[pjj[k]*ppdof+pshift]) {
2011             denserow[pjj[k]*ppdof+pshift] = -1;
2012             sparserow[cnzi++]             = pjj[k]*ppdof+pshift;
2013           }
2014         }
2015       }
2016 
2017       /* sort sparserow */
2018       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
2019 
2020       /* If free space is not available, make more free space */
2021       /* Double the amount of total space in the list */
2022       if (current_space->local_remaining<cnzi) {
2023         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
2024       }
2025 
2026       /* Copy data into free space, and zero out denserows */
2027       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr);
2028       current_space->array           += cnzi;
2029       current_space->local_used      += cnzi;
2030       current_space->local_remaining -= cnzi;
2031 
2032       for (j=0;j<ptanzi;j++) {
2033         ptadenserow[ptasparserow[j]] = 0;
2034       }
2035       for (j=0;j<cnzi;j++) {
2036         denserow[sparserow[j]] = 0;
2037       }
2038       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
2039       /*        For now, we will recompute what is needed. */
2040       ci[i*ppdof+1+dof] = ci[i*ppdof+dof] + cnzi;
2041     }
2042   }
2043   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
2044   /* Allocate space for cj, initialize cj, and */
2045   /* destroy list of free space and other temporary array(s) */
2046   ierr = PetscMalloc((ci[cn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
2047   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
2048   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
2049 
2050   /* Allocate space for ca */
2051   ierr = PetscMalloc((ci[cn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
2052   ierr = PetscMemzero(ca,(ci[cn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
2053 
2054   /* put together the new matrix */
2055   ierr = MatCreateSeqAIJWithArrays(A->comm,cn,cn,ci,cj,ca,C);CHKERRQ(ierr);
2056 
2057   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2058   /* Since these are PETSc arrays, change flags to free them as necessary. */
2059   c = (Mat_SeqAIJ *)((*C)->data);
2060   c->freedata = PETSC_TRUE;
2061   c->nonew    = 0;
2062 
2063   /* Clean up. */
2064   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
2065 
2066   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
2067   PetscFunctionReturn(0);
2068 }
2069 
2070 #undef __FUNCT__
2071 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqMAIJ"
2072 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat C)
2073 {
2074   /* This routine requires testing -- first draft only */
2075   PetscErrorCode ierr;
2076   PetscInt       flops=0;
2077   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
2078   Mat            P=pp->AIJ;
2079   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
2080   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
2081   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
2082   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
2083   PetscInt       *ci=c->i,*cj=c->j,*cjj;
2084   PetscInt       am=A->M,cn=C->N,cm=C->M,ppdof=pp->dof;
2085   PetscInt       i,j,k,pshift,poffset,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
2086   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
2087 
2088   PetscFunctionBegin;
2089   /* Allocate temporary array for storage of one row of A*P */
2090   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
2091   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
2092 
2093   apj      = (PetscInt *)(apa + cn);
2094   apjdense = apj + cn;
2095 
2096   /* Clear old values in C */
2097   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
2098 
2099   for (i=0;i<am;i++) {
2100     /* Form sparse row of A*P */
2101     anzi  = ai[i+1] - ai[i];
2102     apnzj = 0;
2103     for (j=0;j<anzi;j++) {
2104       /* Get offset within block of P */
2105       pshift = *aj%ppdof;
2106       /* Get block row of P */
2107       prow   = *aj++/ppdof; /* integer division */
2108       pnzj = pi[prow+1] - pi[prow];
2109       pjj  = pj + pi[prow];
2110       paj  = pa + pi[prow];
2111       for (k=0;k<pnzj;k++) {
2112         poffset = pjj[k]*ppdof+pshift;
2113         if (!apjdense[poffset]) {
2114           apjdense[poffset] = -1;
2115           apj[apnzj++]      = poffset;
2116         }
2117         apa[poffset] += (*aa)*paj[k];
2118       }
2119       flops += 2*pnzj;
2120       aa++;
2121     }
2122 
2123     /* Sort the j index array for quick sparse axpy. */
2124     /* Note: a array does not need sorting as it is in dense storage locations. */
2125     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
2126 
2127     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
2128     prow    = i/ppdof; /* integer division */
2129     pshift  = i%ppdof;
2130     poffset = pi[prow];
2131     pnzi = pi[prow+1] - poffset;
2132     /* Reset pJ and pA so we can traverse the same row of P 'dof' times. */
2133     pJ   = pj+poffset;
2134     pA   = pa+poffset;
2135     for (j=0;j<pnzi;j++) {
2136       crow   = (*pJ)*ppdof+pshift;
2137       cjj    = cj + ci[crow];
2138       caj    = ca + ci[crow];
2139       pJ++;
2140       /* Perform sparse axpy operation.  Note cjj includes apj. */
2141       for (k=0,nextap=0;nextap<apnzj;k++) {
2142         if (cjj[k]==apj[nextap]) {
2143           caj[k] += (*pA)*apa[apj[nextap++]];
2144         }
2145       }
2146       flops += 2*apnzj;
2147       pA++;
2148     }
2149 
2150     /* Zero the current row info for A*P */
2151     for (j=0;j<apnzj;j++) {
2152       apa[apj[j]]      = 0.;
2153       apjdense[apj[j]] = 0;
2154     }
2155   }
2156 
2157   /* Assemble the final matrix and clean up */
2158   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2159   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2160   ierr = PetscFree(apa);CHKERRQ(ierr);
2161   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
2162 
2163   PetscFunctionReturn(0);
2164 }
2165 
2166 #undef __FUNCT__
2167 #define __FUNCT__ "MatConvert_SeqMAIJ_SeqAIJ"
2168 PetscErrorCode MatConvert_SeqMAIJ_SeqAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat)
2169 {
2170   Mat_SeqMAIJ       *b = (Mat_SeqMAIJ*)A->data;
2171   Mat               a = b->AIJ,B;
2172   Mat_SeqAIJ        *aij = (Mat_SeqAIJ*)a->data;
2173   PetscErrorCode    ierr;
2174   PetscInt          m,n,i,ncols,*ilen,nmax = 0,*icols,j,k,ii,dof = b->dof;
2175   PetscInt          *cols;
2176   PetscScalar       *vals;
2177 
2178   PetscFunctionBegin;
2179   ierr = MatGetSize(a,&m,&n);CHKERRQ(ierr);
2180   ierr = PetscMalloc(dof*m*sizeof(int),&ilen);CHKERRQ(ierr);
2181   for (i=0; i<m; i++) {
2182     nmax = PetscMax(nmax,aij->ilen[i]);
2183     for (j=0; j<dof; j++) {
2184       ilen[dof*i+j] = aij->ilen[i];
2185     }
2186   }
2187   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,dof*m,dof*n,0,ilen,&B);CHKERRQ(ierr);
2188   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
2189   ierr = PetscFree(ilen);CHKERRQ(ierr);
2190   ierr = PetscMalloc(nmax*sizeof(PetscInt),&icols);CHKERRQ(ierr);
2191   ii   = 0;
2192   for (i=0; i<m; i++) {
2193     ierr = MatGetRow_SeqAIJ(a,i,&ncols,&cols,&vals);CHKERRQ(ierr);
2194     for (j=0; j<dof; j++) {
2195       for (k=0; k<ncols; k++) {
2196         icols[k] = dof*cols[k]+j;
2197       }
2198       ierr = MatSetValues_SeqAIJ(B,1,&ii,ncols,icols,vals,INSERT_VALUES);CHKERRQ(ierr);
2199       ii++;
2200     }
2201     ierr = MatRestoreRow_SeqAIJ(a,i,&ncols,&cols,&vals);CHKERRQ(ierr);
2202   }
2203   ierr = PetscFree(icols);CHKERRQ(ierr);
2204   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2205   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2206 
2207   if (reuse == MAT_REUSE_MATRIX) {
2208     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
2209   } else {
2210     *newmat = B;
2211   }
2212   PetscFunctionReturn(0);
2213 }
2214 
2215 #include "src/mat/impls/aij/mpi/mpiaij.h"
2216 #undef __FUNCT__
2217 #define __FUNCT__ "MatConvert_MPIMAIJ_MPIAIJ"
2218 PetscErrorCode MatConvert_MPIMAIJ_MPIAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat)
2219 {
2220   Mat_MPIMAIJ       *maij = (Mat_MPIMAIJ*)A->data;
2221   Mat               MatAIJ  = ((Mat_SeqMAIJ*)maij->AIJ->data)->AIJ,B;
2222   Mat               MatOAIJ = ((Mat_SeqMAIJ*)maij->OAIJ->data)->AIJ;
2223   Mat_SeqAIJ        *AIJ = (Mat_SeqAIJ*) MatAIJ->data;
2224   Mat_SeqAIJ        *OAIJ =(Mat_SeqAIJ*) MatOAIJ->data;
2225   Mat_MPIAIJ        *mpiaij = (Mat_MPIAIJ*) maij->A->data;
2226   PetscInt          dof = maij->dof,i,j,*dnz,*onz,nmax = 0,onmax = 0,*oicols,*icols,ncols,*cols,oncols,*ocols;
2227   PetscInt          rstart,cstart,*garray,ii,k;
2228   PetscErrorCode    ierr;
2229   PetscScalar       *vals,*ovals;
2230 
2231   PetscFunctionBegin;
2232   ierr = PetscMalloc2(A->m,PetscInt,&dnz,A->m,PetscInt,&onz);CHKERRQ(ierr);
2233   for (i=0; i<A->m/dof; i++) {
2234     nmax  = PetscMax(nmax,AIJ->ilen[i]);
2235     onmax = PetscMax(onmax,OAIJ->ilen[i]);
2236     for (j=0; j<dof; j++) {
2237       dnz[dof*i+j] = AIJ->ilen[i];
2238       onz[dof*i+j] = OAIJ->ilen[i];
2239     }
2240   }
2241   ierr = MatCreateMPIAIJ(A->comm,A->m,A->n,A->M,A->N,0,dnz,0,onz,&B);CHKERRQ(ierr);
2242   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
2243   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
2244 
2245   ierr   = PetscMalloc2(nmax,PetscInt,&icols,onmax,PetscInt,&oicols);CHKERRQ(ierr);
2246   rstart = dof*mpiaij->rstart;
2247   cstart = dof*mpiaij->cstart;
2248   garray = mpiaij->garray;
2249 
2250   ii = rstart;
2251   for (i=0; i<A->m/dof; i++) {
2252     ierr = MatGetRow_SeqAIJ(MatAIJ,i,&ncols,&cols,&vals);CHKERRQ(ierr);
2253     ierr = MatGetRow_SeqAIJ(MatOAIJ,i,&oncols,&ocols,&ovals);CHKERRQ(ierr);
2254     for (j=0; j<dof; j++) {
2255       for (k=0; k<ncols; k++) {
2256         icols[k] = cstart + dof*cols[k]+j;
2257       }
2258       for (k=0; k<oncols; k++) {
2259         oicols[k] = dof*garray[ocols[k]]+j;
2260       }
2261       ierr = MatSetValues_MPIAIJ(B,1,&ii,ncols,icols,vals,INSERT_VALUES);CHKERRQ(ierr);
2262       ierr = MatSetValues_MPIAIJ(B,1,&ii,oncols,oicols,ovals,INSERT_VALUES);CHKERRQ(ierr);
2263       ii++;
2264     }
2265     ierr = MatRestoreRow_SeqAIJ(MatAIJ,i,&ncols,&cols,&vals);CHKERRQ(ierr);
2266     ierr = MatRestoreRow_SeqAIJ(MatOAIJ,i,&oncols,&ocols,&ovals);CHKERRQ(ierr);
2267   }
2268   ierr = PetscFree2(icols,oicols);CHKERRQ(ierr);
2269 
2270   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2271   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2272 
2273   if (reuse == MAT_REUSE_MATRIX) {
2274     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
2275   } else {
2276     *newmat = B;
2277   }
2278   PetscFunctionReturn(0);
2279 }
2280 
2281 
2282 
2283 /* ---------------------------------------------------------------------------------- */
2284 /*MC
2285   MatCreateMAIJ - Creates a matrix type providing restriction and interpolation
2286   operations for multicomponent problems.  It interpolates each component the same
2287   way independently.  The matrix type is based on MATSEQAIJ for sequential matrices,
2288   and MATMPIAIJ for distributed matrices.
2289 
2290   Operations provided:
2291 + MatMult
2292 . MatMultTranspose
2293 . MatMultAdd
2294 . MatMultTransposeAdd
2295 - MatView
2296 
2297   Level: advanced
2298 
2299 M*/
2300 #undef __FUNCT__
2301 #define __FUNCT__ "MatCreateMAIJ"
2302 PetscErrorCode MatCreateMAIJ(Mat A,PetscInt dof,Mat *maij)
2303 {
2304   PetscErrorCode ierr;
2305   PetscMPIInt    size;
2306   PetscInt       n;
2307   Mat_MPIMAIJ    *b;
2308   Mat            B;
2309 
2310   PetscFunctionBegin;
2311   ierr   = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
2312 
2313   if (dof == 1) {
2314     *maij = A;
2315   } else {
2316     ierr = MatCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr);
2317     B->assembled    = PETSC_TRUE;
2318 
2319     ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2320     if (size == 1) {
2321       ierr = MatSetType(B,MATSEQMAIJ);CHKERRQ(ierr);
2322       B->ops->destroy = MatDestroy_SeqMAIJ;
2323       B->ops->view    = MatView_SeqMAIJ;
2324       b      = (Mat_MPIMAIJ*)B->data;
2325       b->dof = dof;
2326       b->AIJ = A;
2327       if (dof == 2) {
2328         B->ops->mult             = MatMult_SeqMAIJ_2;
2329         B->ops->multadd          = MatMultAdd_SeqMAIJ_2;
2330         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_2;
2331         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2;
2332       } else if (dof == 3) {
2333         B->ops->mult             = MatMult_SeqMAIJ_3;
2334         B->ops->multadd          = MatMultAdd_SeqMAIJ_3;
2335         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_3;
2336         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3;
2337       } else if (dof == 4) {
2338         B->ops->mult             = MatMult_SeqMAIJ_4;
2339         B->ops->multadd          = MatMultAdd_SeqMAIJ_4;
2340         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_4;
2341         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4;
2342       } else if (dof == 5) {
2343         B->ops->mult             = MatMult_SeqMAIJ_5;
2344         B->ops->multadd          = MatMultAdd_SeqMAIJ_5;
2345         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_5;
2346         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5;
2347       } else if (dof == 6) {
2348         B->ops->mult             = MatMult_SeqMAIJ_6;
2349         B->ops->multadd          = MatMultAdd_SeqMAIJ_6;
2350         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_6;
2351         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_6;
2352       } else if (dof == 7) {
2353         B->ops->mult             = MatMult_SeqMAIJ_7;
2354         B->ops->multadd          = MatMultAdd_SeqMAIJ_7;
2355         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_7;
2356         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_7;
2357       } else if (dof == 8) {
2358         B->ops->mult             = MatMult_SeqMAIJ_8;
2359         B->ops->multadd          = MatMultAdd_SeqMAIJ_8;
2360         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_8;
2361         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_8;
2362       } else if (dof == 9) {
2363         B->ops->mult             = MatMult_SeqMAIJ_9;
2364         B->ops->multadd          = MatMultAdd_SeqMAIJ_9;
2365         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_9;
2366         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_9;
2367       } else if (dof == 16) {
2368         B->ops->mult             = MatMult_SeqMAIJ_16;
2369         B->ops->multadd          = MatMultAdd_SeqMAIJ_16;
2370         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_16;
2371         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_16;
2372       } else {
2373         SETERRQ1(PETSC_ERR_SUP,"Cannot handle a dof of %D. Send request for code to petsc-maint@mcs.anl.gov\n",dof);
2374       }
2375       B->ops->ptapsymbolic_seqaij = MatPtAPSymbolic_SeqAIJ_SeqMAIJ;
2376       B->ops->ptapnumeric_seqaij  = MatPtAPNumeric_SeqAIJ_SeqMAIJ;
2377       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqmaij_seqaij_C","MatConvert_SeqMAIJ_SeqAIJ",MatConvert_SeqMAIJ_SeqAIJ);CHKERRQ(ierr);
2378     } else {
2379       Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data;
2380       IS         from,to;
2381       Vec        gvec;
2382       PetscInt   *garray,i;
2383 
2384       ierr = MatSetType(B,MATMPIMAIJ);CHKERRQ(ierr);
2385       B->ops->destroy = MatDestroy_MPIMAIJ;
2386       B->ops->view    = MatView_MPIMAIJ;
2387       b      = (Mat_MPIMAIJ*)B->data;
2388       b->dof = dof;
2389       b->A   = A;
2390       ierr = MatCreateMAIJ(mpiaij->A,dof,&b->AIJ);CHKERRQ(ierr);
2391       ierr = MatCreateMAIJ(mpiaij->B,dof,&b->OAIJ);CHKERRQ(ierr);
2392 
2393       ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr);
2394       ierr = VecCreateSeq(PETSC_COMM_SELF,n*dof,&b->w);CHKERRQ(ierr);
2395       ierr = VecSetBlockSize(b->w,dof);CHKERRQ(ierr);
2396 
2397       /* create two temporary Index sets for build scatter gather */
2398       ierr = PetscMalloc((n+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
2399       for (i=0; i<n; i++) garray[i] = dof*mpiaij->garray[i];
2400       ierr = ISCreateBlock(A->comm,dof,n,garray,&from);CHKERRQ(ierr);
2401       ierr = PetscFree(garray);CHKERRQ(ierr);
2402       ierr = ISCreateStride(PETSC_COMM_SELF,n*dof,0,1,&to);CHKERRQ(ierr);
2403 
2404       /* create temporary global vector to generate scatter context */
2405       ierr = VecCreateMPI(A->comm,dof*A->n,dof*A->N,&gvec);CHKERRQ(ierr);
2406       ierr = VecSetBlockSize(gvec,dof);CHKERRQ(ierr);
2407 
2408       /* generate the scatter context */
2409       ierr = VecScatterCreate(gvec,from,b->w,to,&b->ctx);CHKERRQ(ierr);
2410 
2411       ierr = ISDestroy(from);CHKERRQ(ierr);
2412       ierr = ISDestroy(to);CHKERRQ(ierr);
2413       ierr = VecDestroy(gvec);CHKERRQ(ierr);
2414 
2415       B->ops->mult             = MatMult_MPIMAIJ_dof;
2416       B->ops->multtranspose    = MatMultTranspose_MPIMAIJ_dof;
2417       B->ops->multadd          = MatMultAdd_MPIMAIJ_dof;
2418       B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof;
2419       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpimaij_mpiaij_C","MatConvert_MPIMAIJ_MPIAIJ",MatConvert_MPIMAIJ_MPIAIJ);CHKERRQ(ierr);
2420     }
2421     *maij = B;
2422     ierr = MatView_Private(B);CHKERRQ(ierr);
2423   }
2424   PetscFunctionReturn(0);
2425 }
2426 
2427 
2428 
2429 
2430 
2431 
2432 
2433 
2434 
2435 
2436 
2437 
2438