xref: /petsc/src/mat/impls/maij/maij.c (revision 87828ca270d8140797fd4271705413c4ecfcb535)
1*87828ca2SBarry Smith /*$Id: maij.c,v 1.17 2001/05/29 03:29:32 bsmith Exp bsmith $*/
282b1193eSBarry Smith /*
3cd3bbe55SBarry Smith     Defines the basic matrix operations for the MAIJ  matrix storage format.
4cd3bbe55SBarry Smith   This format is used for restriction and interpolation operations for
5cd3bbe55SBarry Smith   multicomponent problems. It interpolates each component the same way
6cd3bbe55SBarry Smith   independently.
7cd3bbe55SBarry Smith 
8cd3bbe55SBarry Smith      We provide:
9cd3bbe55SBarry Smith          MatMult()
10cd3bbe55SBarry Smith          MatMultTranspose()
11cd3bbe55SBarry Smith          MatMultTransposeAdd()
12cd3bbe55SBarry Smith          MatMultAdd()
13cd3bbe55SBarry Smith           and
14cd3bbe55SBarry Smith          MatCreateMAIJ(Mat,dof,Mat*)
15f4a53059SBarry Smith 
16f4a53059SBarry Smith      This single directory handles both the sequential and parallel codes
1782b1193eSBarry Smith */
1882b1193eSBarry Smith 
19f4a53059SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h"
2082b1193eSBarry Smith 
21cd3bbe55SBarry Smith typedef struct {
22cd3bbe55SBarry Smith   int        dof;         /* number of components */
234cb9d9b8SBarry Smith   Mat        AIJ;        /* representation of interpolation for one component */
244cb9d9b8SBarry Smith } Mat_SeqMAIJ;
254cb9d9b8SBarry Smith 
264cb9d9b8SBarry Smith typedef struct {
274cb9d9b8SBarry Smith   int        dof;         /* number of components */
28f4a53059SBarry Smith   Mat        AIJ,OAIJ;    /* representation of interpolation for one component */
294cb9d9b8SBarry Smith   Mat        A;
30f4a53059SBarry Smith   VecScatter ctx;         /* update ghost points for parallel case */
31f4a53059SBarry Smith   Vec        w;           /* work space for ghost values for parallel case */
324cb9d9b8SBarry Smith } Mat_MPIMAIJ;
3382b1193eSBarry Smith 
344a2ae208SSatish Balay #undef __FUNCT__
354a2ae208SSatish Balay #define __FUNCT__ "MatMAIJGetAIJ"
36b9b97703SBarry Smith int MatMAIJGetAIJ(Mat A,Mat *B)
37b9b97703SBarry Smith {
38b9b97703SBarry Smith   int         ierr;
39b9b97703SBarry Smith   PetscTruth  ismpimaij,isseqmaij;
40b9b97703SBarry Smith 
41b9b97703SBarry Smith   PetscFunctionBegin;
42b9b97703SBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATMPIMAIJ,&ismpimaij);CHKERRQ(ierr);
43b9b97703SBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATSEQMAIJ,&isseqmaij);CHKERRQ(ierr);
44b9b97703SBarry Smith   if (ismpimaij) {
45b9b97703SBarry Smith     Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
46b9b97703SBarry Smith 
47b9b97703SBarry Smith     *B = b->A;
48b9b97703SBarry Smith   } else if (isseqmaij) {
49b9b97703SBarry Smith     Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
50b9b97703SBarry Smith 
51b9b97703SBarry Smith     *B = b->AIJ;
52b9b97703SBarry Smith   } else {
53b9b97703SBarry Smith     *B = A;
54b9b97703SBarry Smith   }
55b9b97703SBarry Smith   PetscFunctionReturn(0);
56b9b97703SBarry Smith }
57b9b97703SBarry Smith 
584a2ae208SSatish Balay #undef __FUNCT__
594a2ae208SSatish Balay #define __FUNCT__ "MatMAIJRedimension"
60b9b97703SBarry Smith int MatMAIJRedimension(Mat A,int dof,Mat *B)
61b9b97703SBarry Smith {
62b9b97703SBarry Smith   int ierr;
63b9b97703SBarry Smith   Mat Aij;
64b9b97703SBarry Smith 
65b9b97703SBarry Smith   PetscFunctionBegin;
66b9b97703SBarry Smith   ierr = MatMAIJGetAIJ(A,&Aij);CHKERRQ(ierr);
67b9b97703SBarry Smith   ierr = MatCreateMAIJ(Aij,dof,B);CHKERRQ(ierr);
68b9b97703SBarry Smith   PetscFunctionReturn(0);
69b9b97703SBarry Smith }
70b9b97703SBarry Smith 
714a2ae208SSatish Balay #undef __FUNCT__
724a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqMAIJ"
734cb9d9b8SBarry Smith int MatDestroy_SeqMAIJ(Mat A)
7482b1193eSBarry Smith {
7582b1193eSBarry Smith   int         ierr;
764cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
7782b1193eSBarry Smith 
7882b1193eSBarry Smith   PetscFunctionBegin;
79cd3bbe55SBarry Smith   if (b->AIJ) {
80cd3bbe55SBarry Smith     ierr = MatDestroy(b->AIJ);CHKERRQ(ierr);
8182b1193eSBarry Smith   }
824cb9d9b8SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
834cb9d9b8SBarry Smith   PetscFunctionReturn(0);
844cb9d9b8SBarry Smith }
854cb9d9b8SBarry Smith 
864a2ae208SSatish Balay #undef __FUNCT__
874a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIMAIJ"
884cb9d9b8SBarry Smith int MatDestroy_MPIMAIJ(Mat A)
894cb9d9b8SBarry Smith {
904cb9d9b8SBarry Smith   int         ierr;
914cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
924cb9d9b8SBarry Smith 
934cb9d9b8SBarry Smith   PetscFunctionBegin;
944cb9d9b8SBarry Smith   if (b->AIJ) {
954cb9d9b8SBarry Smith     ierr = MatDestroy(b->AIJ);CHKERRQ(ierr);
964cb9d9b8SBarry Smith   }
97f4a53059SBarry Smith   if (b->OAIJ) {
98f4a53059SBarry Smith     ierr = MatDestroy(b->OAIJ);CHKERRQ(ierr);
99f4a53059SBarry Smith   }
100b9b97703SBarry Smith   if (b->A) {
101b9b97703SBarry Smith     ierr = MatDestroy(b->A);CHKERRQ(ierr);
102b9b97703SBarry Smith   }
103f4a53059SBarry Smith   if (b->ctx) {
104f4a53059SBarry Smith     ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr);
105f4a53059SBarry Smith   }
106f4a53059SBarry Smith   if (b->w) {
107f4a53059SBarry Smith     ierr = VecDestroy(b->w);CHKERRQ(ierr);
108f4a53059SBarry Smith   }
109cd3bbe55SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
110b9b97703SBarry Smith   PetscFunctionReturn(0);
111b9b97703SBarry Smith }
112b9b97703SBarry Smith 
11382b1193eSBarry Smith EXTERN_C_BEGIN
1144a2ae208SSatish Balay #undef __FUNCT__
1154a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MAIJ"
116f4a53059SBarry Smith int MatCreate_MAIJ(Mat A)
11782b1193eSBarry Smith {
118cd3bbe55SBarry Smith   int         ierr;
1194cb9d9b8SBarry Smith   Mat_MPIMAIJ *b;
12082b1193eSBarry Smith 
12182b1193eSBarry Smith   PetscFunctionBegin;
122b0a32e0cSBarry Smith   ierr     = PetscNew(Mat_MPIMAIJ,&b);CHKERRQ(ierr);
123b0a32e0cSBarry Smith   A->data  = (void*)b;
1244cb9d9b8SBarry Smith   ierr = PetscMemzero(b,sizeof(Mat_MPIMAIJ));CHKERRQ(ierr);
125cd3bbe55SBarry Smith   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
126cd3bbe55SBarry Smith   A->factor           = 0;
127cd3bbe55SBarry Smith   A->mapping          = 0;
128f4a53059SBarry Smith 
129cd3bbe55SBarry Smith   b->AIJ  = 0;
130cd3bbe55SBarry Smith   b->dof  = 0;
131f4a53059SBarry Smith   b->OAIJ = 0;
132f4a53059SBarry Smith   b->ctx  = 0;
133f4a53059SBarry Smith   b->w    = 0;
13482b1193eSBarry Smith   PetscFunctionReturn(0);
13582b1193eSBarry Smith }
13682b1193eSBarry Smith EXTERN_C_END
13782b1193eSBarry Smith 
138cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/
1394a2ae208SSatish Balay #undef __FUNCT__
140b9617806SBarry Smith #define __FUNCT__ "MatMult_SeqMAIJ_2"
141cd3bbe55SBarry Smith int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
14282b1193eSBarry Smith {
1434cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
144bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
145*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2;
146273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
147bcc973b7SBarry Smith   int         n,i,jrow,j;
14882b1193eSBarry Smith 
149bcc973b7SBarry Smith   PetscFunctionBegin;
150bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
151bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
152bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
153bcc973b7SBarry Smith   idx  = a->j;
154bcc973b7SBarry Smith   v    = a->a;
155bcc973b7SBarry Smith   ii   = a->i;
156bcc973b7SBarry Smith 
157bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
158bcc973b7SBarry Smith   idx  += shift;
159bcc973b7SBarry Smith   for (i=0; i<m; i++) {
160bcc973b7SBarry Smith     jrow = ii[i];
161bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
162bcc973b7SBarry Smith     sum1  = 0.0;
163bcc973b7SBarry Smith     sum2  = 0.0;
164bcc973b7SBarry Smith     for (j=0; j<n; j++) {
165bcc973b7SBarry Smith       sum1 += v[jrow]*x[2*idx[jrow]];
166bcc973b7SBarry Smith       sum2 += v[jrow]*x[2*idx[jrow]+1];
167bcc973b7SBarry Smith       jrow++;
168bcc973b7SBarry Smith      }
169bcc973b7SBarry Smith     y[2*i]   = sum1;
170bcc973b7SBarry Smith     y[2*i+1] = sum2;
171bcc973b7SBarry Smith   }
172bcc973b7SBarry Smith 
173b0a32e0cSBarry Smith   PetscLogFlops(4*a->nz - 2*m);
174bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
175bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
17682b1193eSBarry Smith   PetscFunctionReturn(0);
17782b1193eSBarry Smith }
178bcc973b7SBarry Smith 
1794a2ae208SSatish Balay #undef __FUNCT__
180b9617806SBarry Smith #define __FUNCT__ "MatMultTranspose_SeqMAIJ_2"
181cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
18282b1193eSBarry Smith {
1834cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
184bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
185*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,zero = 0.0;
186273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
18782b1193eSBarry Smith 
188bcc973b7SBarry Smith   PetscFunctionBegin;
189bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
190bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
191bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
192bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
193bcc973b7SBarry Smith   for (i=0; i<m; i++) {
194bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
195bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
196bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
197bcc973b7SBarry Smith     alpha1 = x[2*i];
198bcc973b7SBarry Smith     alpha2 = x[2*i+1];
199bcc973b7SBarry Smith     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
200bcc973b7SBarry Smith   }
201b0a32e0cSBarry Smith   PetscLogFlops(4*a->nz - 2*b->AIJ->n);
202bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
203bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
20482b1193eSBarry Smith   PetscFunctionReturn(0);
20582b1193eSBarry Smith }
206bcc973b7SBarry Smith 
2074a2ae208SSatish Balay #undef __FUNCT__
208b9617806SBarry Smith #define __FUNCT__ "MatMultAdd_SeqMAIJ_2"
209cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
21082b1193eSBarry Smith {
2114cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
212bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
213*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2;
214273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
215bcc973b7SBarry Smith   int         n,i,jrow,j;
21682b1193eSBarry Smith 
217bcc973b7SBarry Smith   PetscFunctionBegin;
218f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
219bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
220f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
221bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
222bcc973b7SBarry Smith   idx  = a->j;
223bcc973b7SBarry Smith   v    = a->a;
224bcc973b7SBarry Smith   ii   = a->i;
225bcc973b7SBarry Smith 
226bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
227bcc973b7SBarry Smith   idx  += shift;
228bcc973b7SBarry Smith   for (i=0; i<m; i++) {
229bcc973b7SBarry Smith     jrow = ii[i];
230bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
231bcc973b7SBarry Smith     sum1  = 0.0;
232bcc973b7SBarry Smith     sum2  = 0.0;
233bcc973b7SBarry Smith     for (j=0; j<n; j++) {
234bcc973b7SBarry Smith       sum1 += v[jrow]*x[2*idx[jrow]];
235bcc973b7SBarry Smith       sum2 += v[jrow]*x[2*idx[jrow]+1];
236bcc973b7SBarry Smith       jrow++;
237bcc973b7SBarry Smith      }
238bcc973b7SBarry Smith     y[2*i]   += sum1;
239bcc973b7SBarry Smith     y[2*i+1] += sum2;
240bcc973b7SBarry Smith   }
241bcc973b7SBarry Smith 
242b0a32e0cSBarry Smith   PetscLogFlops(4*a->nz - 2*m);
243bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
244f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
245bcc973b7SBarry Smith   PetscFunctionReturn(0);
24682b1193eSBarry Smith }
2474a2ae208SSatish Balay #undef __FUNCT__
248b9617806SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_2"
249cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
25082b1193eSBarry Smith {
2514cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
252bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
253*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2;
254273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
25582b1193eSBarry Smith 
256bcc973b7SBarry Smith   PetscFunctionBegin;
257f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
258bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
259f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
260bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
261bcc973b7SBarry Smith   for (i=0; i<m; i++) {
262bcc973b7SBarry Smith     idx   = a->j + a->i[i] + shift;
263bcc973b7SBarry Smith     v     = a->a + a->i[i] + shift;
264bcc973b7SBarry Smith     n     = a->i[i+1] - a->i[i];
265bcc973b7SBarry Smith     alpha1 = x[2*i];
266bcc973b7SBarry Smith     alpha2 = x[2*i+1];
267bcc973b7SBarry Smith     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
268bcc973b7SBarry Smith   }
269b0a32e0cSBarry Smith   PetscLogFlops(4*a->nz - 2*b->AIJ->n);
270bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
271f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
272bcc973b7SBarry Smith   PetscFunctionReturn(0);
27382b1193eSBarry Smith }
274cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/
2754a2ae208SSatish Balay #undef __FUNCT__
276b9617806SBarry Smith #define __FUNCT__ "MatMult_SeqMAIJ_3"
277bcc973b7SBarry Smith int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
278bcc973b7SBarry Smith {
2794cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
280bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
281*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3;
282273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
283bcc973b7SBarry Smith   int         n,i,jrow,j;
28482b1193eSBarry Smith 
285bcc973b7SBarry Smith   PetscFunctionBegin;
286bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
287bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
288bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
289bcc973b7SBarry Smith   idx  = a->j;
290bcc973b7SBarry Smith   v    = a->a;
291bcc973b7SBarry Smith   ii   = a->i;
292bcc973b7SBarry Smith 
293bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
294bcc973b7SBarry Smith   idx  += shift;
295bcc973b7SBarry Smith   for (i=0; i<m; i++) {
296bcc973b7SBarry Smith     jrow = ii[i];
297bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
298bcc973b7SBarry Smith     sum1  = 0.0;
299bcc973b7SBarry Smith     sum2  = 0.0;
300bcc973b7SBarry Smith     sum3  = 0.0;
301bcc973b7SBarry Smith     for (j=0; j<n; j++) {
302bcc973b7SBarry Smith       sum1 += v[jrow]*x[3*idx[jrow]];
303bcc973b7SBarry Smith       sum2 += v[jrow]*x[3*idx[jrow]+1];
304bcc973b7SBarry Smith       sum3 += v[jrow]*x[3*idx[jrow]+2];
305bcc973b7SBarry Smith       jrow++;
306bcc973b7SBarry Smith      }
307bcc973b7SBarry Smith     y[3*i]   = sum1;
308bcc973b7SBarry Smith     y[3*i+1] = sum2;
309bcc973b7SBarry Smith     y[3*i+2] = sum3;
310bcc973b7SBarry Smith   }
311bcc973b7SBarry Smith 
312b0a32e0cSBarry Smith   PetscLogFlops(6*a->nz - 3*m);
313bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
314bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
315bcc973b7SBarry Smith   PetscFunctionReturn(0);
316bcc973b7SBarry Smith }
317bcc973b7SBarry Smith 
3184a2ae208SSatish Balay #undef __FUNCT__
319b9617806SBarry Smith #define __FUNCT__ "MatMultTranspose_SeqMAIJ_3"
320bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
321bcc973b7SBarry Smith {
3224cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
323bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
324*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0;
325273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
326bcc973b7SBarry Smith 
327bcc973b7SBarry Smith   PetscFunctionBegin;
328bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
329bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
330bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
331bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
332bcc973b7SBarry Smith   for (i=0; i<m; i++) {
333bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
334bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
335bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
336bcc973b7SBarry Smith     alpha1 = x[3*i];
337bcc973b7SBarry Smith     alpha2 = x[3*i+1];
338bcc973b7SBarry Smith     alpha3 = x[3*i+2];
339bcc973b7SBarry Smith     while (n-->0) {
340bcc973b7SBarry Smith       y[3*(*idx)]   += alpha1*(*v);
341bcc973b7SBarry Smith       y[3*(*idx)+1] += alpha2*(*v);
342bcc973b7SBarry Smith       y[3*(*idx)+2] += alpha3*(*v);
343bcc973b7SBarry Smith       idx++; v++;
344bcc973b7SBarry Smith     }
345bcc973b7SBarry Smith   }
346b0a32e0cSBarry Smith   PetscLogFlops(6*a->nz - 3*b->AIJ->n);
347bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
348bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
349bcc973b7SBarry Smith   PetscFunctionReturn(0);
350bcc973b7SBarry Smith }
351bcc973b7SBarry Smith 
3524a2ae208SSatish Balay #undef __FUNCT__
353b9617806SBarry Smith #define __FUNCT__ "MatMultAdd_SeqMAIJ_3"
354bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
355bcc973b7SBarry Smith {
3564cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
357bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
358*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3;
359273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
360bcc973b7SBarry Smith   int         n,i,jrow,j;
361bcc973b7SBarry Smith 
362bcc973b7SBarry Smith   PetscFunctionBegin;
363f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
364bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
365f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
366bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
367bcc973b7SBarry Smith   idx  = a->j;
368bcc973b7SBarry Smith   v    = a->a;
369bcc973b7SBarry Smith   ii   = a->i;
370bcc973b7SBarry Smith 
371bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
372bcc973b7SBarry Smith   idx  += shift;
373bcc973b7SBarry Smith   for (i=0; i<m; i++) {
374bcc973b7SBarry Smith     jrow = ii[i];
375bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
376bcc973b7SBarry Smith     sum1  = 0.0;
377bcc973b7SBarry Smith     sum2  = 0.0;
378bcc973b7SBarry Smith     sum3  = 0.0;
379bcc973b7SBarry Smith     for (j=0; j<n; j++) {
380bcc973b7SBarry Smith       sum1 += v[jrow]*x[3*idx[jrow]];
381bcc973b7SBarry Smith       sum2 += v[jrow]*x[3*idx[jrow]+1];
382bcc973b7SBarry Smith       sum3 += v[jrow]*x[3*idx[jrow]+2];
383bcc973b7SBarry Smith       jrow++;
384bcc973b7SBarry Smith      }
385bcc973b7SBarry Smith     y[3*i]   += sum1;
386bcc973b7SBarry Smith     y[3*i+1] += sum2;
387bcc973b7SBarry Smith     y[3*i+2] += sum3;
388bcc973b7SBarry Smith   }
389bcc973b7SBarry Smith 
390b0a32e0cSBarry Smith   PetscLogFlops(6*a->nz);
391bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
392f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
393bcc973b7SBarry Smith   PetscFunctionReturn(0);
394bcc973b7SBarry Smith }
3954a2ae208SSatish Balay #undef __FUNCT__
396b9617806SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_3"
397bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
398bcc973b7SBarry Smith {
3994cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
400bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
401*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3;
402273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
403bcc973b7SBarry Smith 
404bcc973b7SBarry Smith   PetscFunctionBegin;
405f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
406bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
407f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
408bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
409bcc973b7SBarry Smith   for (i=0; i<m; i++) {
410bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
411bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
412bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
413bcc973b7SBarry Smith     alpha1 = x[3*i];
414bcc973b7SBarry Smith     alpha2 = x[3*i+1];
415bcc973b7SBarry Smith     alpha3 = x[3*i+2];
416bcc973b7SBarry Smith     while (n-->0) {
417bcc973b7SBarry Smith       y[3*(*idx)]   += alpha1*(*v);
418bcc973b7SBarry Smith       y[3*(*idx)+1] += alpha2*(*v);
419bcc973b7SBarry Smith       y[3*(*idx)+2] += alpha3*(*v);
420bcc973b7SBarry Smith       idx++; v++;
421bcc973b7SBarry Smith     }
422bcc973b7SBarry Smith   }
423b0a32e0cSBarry Smith   PetscLogFlops(6*a->nz);
424bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
425f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
426bcc973b7SBarry Smith   PetscFunctionReturn(0);
427bcc973b7SBarry Smith }
428bcc973b7SBarry Smith 
429bcc973b7SBarry Smith /* ------------------------------------------------------------------------------*/
4304a2ae208SSatish Balay #undef __FUNCT__
431b9617806SBarry Smith #define __FUNCT__ "MatMult_SeqMAIJ_4"
432bcc973b7SBarry Smith int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
433bcc973b7SBarry Smith {
4344cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
435bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
436*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3, sum4;
437273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
438bcc973b7SBarry Smith   int         n,i,jrow,j;
439bcc973b7SBarry Smith 
440bcc973b7SBarry Smith   PetscFunctionBegin;
441bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
442bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
443bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
444bcc973b7SBarry Smith   idx  = a->j;
445bcc973b7SBarry Smith   v    = a->a;
446bcc973b7SBarry Smith   ii   = a->i;
447bcc973b7SBarry Smith 
448bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
449bcc973b7SBarry Smith   idx  += shift;
450bcc973b7SBarry Smith   for (i=0; i<m; i++) {
451bcc973b7SBarry Smith     jrow = ii[i];
452bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
453bcc973b7SBarry Smith     sum1  = 0.0;
454bcc973b7SBarry Smith     sum2  = 0.0;
455bcc973b7SBarry Smith     sum3  = 0.0;
456bcc973b7SBarry Smith     sum4  = 0.0;
457bcc973b7SBarry Smith     for (j=0; j<n; j++) {
458bcc973b7SBarry Smith       sum1 += v[jrow]*x[4*idx[jrow]];
459bcc973b7SBarry Smith       sum2 += v[jrow]*x[4*idx[jrow]+1];
460bcc973b7SBarry Smith       sum3 += v[jrow]*x[4*idx[jrow]+2];
461bcc973b7SBarry Smith       sum4 += v[jrow]*x[4*idx[jrow]+3];
462bcc973b7SBarry Smith       jrow++;
463bcc973b7SBarry Smith      }
464bcc973b7SBarry Smith     y[4*i]   = sum1;
465bcc973b7SBarry Smith     y[4*i+1] = sum2;
466bcc973b7SBarry Smith     y[4*i+2] = sum3;
467bcc973b7SBarry Smith     y[4*i+3] = sum4;
468bcc973b7SBarry Smith   }
469bcc973b7SBarry Smith 
470b0a32e0cSBarry Smith   PetscLogFlops(8*a->nz - 4*m);
471bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
472bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
473bcc973b7SBarry Smith   PetscFunctionReturn(0);
474bcc973b7SBarry Smith }
475bcc973b7SBarry Smith 
4764a2ae208SSatish Balay #undef __FUNCT__
477b9617806SBarry Smith #define __FUNCT__ "MatMultTranspose_SeqMAIJ_4"
478bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
479bcc973b7SBarry Smith {
4804cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
481bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
482*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0;
483273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
484bcc973b7SBarry Smith 
485bcc973b7SBarry Smith   PetscFunctionBegin;
486bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
487bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
488bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
489bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
490bcc973b7SBarry Smith   for (i=0; i<m; i++) {
491bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
492bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
493bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
494bcc973b7SBarry Smith     alpha1 = x[4*i];
495bcc973b7SBarry Smith     alpha2 = x[4*i+1];
496bcc973b7SBarry Smith     alpha3 = x[4*i+2];
497bcc973b7SBarry Smith     alpha4 = x[4*i+3];
498bcc973b7SBarry Smith     while (n-->0) {
499bcc973b7SBarry Smith       y[4*(*idx)]   += alpha1*(*v);
500bcc973b7SBarry Smith       y[4*(*idx)+1] += alpha2*(*v);
501bcc973b7SBarry Smith       y[4*(*idx)+2] += alpha3*(*v);
502bcc973b7SBarry Smith       y[4*(*idx)+3] += alpha4*(*v);
503bcc973b7SBarry Smith       idx++; v++;
504bcc973b7SBarry Smith     }
505bcc973b7SBarry Smith   }
506b0a32e0cSBarry Smith   PetscLogFlops(8*a->nz - 4*b->AIJ->n);
507bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
508bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
509bcc973b7SBarry Smith   PetscFunctionReturn(0);
510bcc973b7SBarry Smith }
511bcc973b7SBarry Smith 
5124a2ae208SSatish Balay #undef __FUNCT__
513b9617806SBarry Smith #define __FUNCT__ "MatMultAdd_SeqMAIJ_4"
514bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
515bcc973b7SBarry Smith {
5164cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
517f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
518*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3, sum4;
519273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
520f9fae5adSBarry Smith   int         n,i,jrow,j;
521f9fae5adSBarry Smith 
522f9fae5adSBarry Smith   PetscFunctionBegin;
523f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
524f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
525f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
526f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
527f9fae5adSBarry Smith   idx  = a->j;
528f9fae5adSBarry Smith   v    = a->a;
529f9fae5adSBarry Smith   ii   = a->i;
530f9fae5adSBarry Smith 
531f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
532f9fae5adSBarry Smith   idx  += shift;
533f9fae5adSBarry Smith   for (i=0; i<m; i++) {
534f9fae5adSBarry Smith     jrow = ii[i];
535f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
536f9fae5adSBarry Smith     sum1  = 0.0;
537f9fae5adSBarry Smith     sum2  = 0.0;
538f9fae5adSBarry Smith     sum3  = 0.0;
539f9fae5adSBarry Smith     sum4  = 0.0;
540f9fae5adSBarry Smith     for (j=0; j<n; j++) {
541f9fae5adSBarry Smith       sum1 += v[jrow]*x[4*idx[jrow]];
542f9fae5adSBarry Smith       sum2 += v[jrow]*x[4*idx[jrow]+1];
543f9fae5adSBarry Smith       sum3 += v[jrow]*x[4*idx[jrow]+2];
544f9fae5adSBarry Smith       sum4 += v[jrow]*x[4*idx[jrow]+3];
545f9fae5adSBarry Smith       jrow++;
546f9fae5adSBarry Smith      }
547f9fae5adSBarry Smith     y[4*i]   += sum1;
548f9fae5adSBarry Smith     y[4*i+1] += sum2;
549f9fae5adSBarry Smith     y[4*i+2] += sum3;
550f9fae5adSBarry Smith     y[4*i+3] += sum4;
551f9fae5adSBarry Smith   }
552f9fae5adSBarry Smith 
553b0a32e0cSBarry Smith   PetscLogFlops(8*a->nz - 4*m);
554f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
555f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
556f9fae5adSBarry Smith   PetscFunctionReturn(0);
557bcc973b7SBarry Smith }
5584a2ae208SSatish Balay #undef __FUNCT__
559b9617806SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_4"
560bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
561bcc973b7SBarry Smith {
5624cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
563f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
564*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4;
565273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
566f9fae5adSBarry Smith 
567f9fae5adSBarry Smith   PetscFunctionBegin;
568f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
569f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
570f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
571f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
572f9fae5adSBarry Smith   for (i=0; i<m; i++) {
573f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
574f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
575f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
576f9fae5adSBarry Smith     alpha1 = x[4*i];
577f9fae5adSBarry Smith     alpha2 = x[4*i+1];
578f9fae5adSBarry Smith     alpha3 = x[4*i+2];
579f9fae5adSBarry Smith     alpha4 = x[4*i+3];
580f9fae5adSBarry Smith     while (n-->0) {
581f9fae5adSBarry Smith       y[4*(*idx)]   += alpha1*(*v);
582f9fae5adSBarry Smith       y[4*(*idx)+1] += alpha2*(*v);
583f9fae5adSBarry Smith       y[4*(*idx)+2] += alpha3*(*v);
584f9fae5adSBarry Smith       y[4*(*idx)+3] += alpha4*(*v);
585f9fae5adSBarry Smith       idx++; v++;
586f9fae5adSBarry Smith     }
587f9fae5adSBarry Smith   }
588b0a32e0cSBarry Smith   PetscLogFlops(8*a->nz - 4*b->AIJ->n);
589f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
590f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
591f9fae5adSBarry Smith   PetscFunctionReturn(0);
592f9fae5adSBarry Smith }
593f9fae5adSBarry Smith /* ------------------------------------------------------------------------------*/
5946cd98798SBarry Smith 
5954a2ae208SSatish Balay #undef __FUNCT__
596b9617806SBarry Smith #define __FUNCT__ "MatMult_SeqMAIJ_5"
597f9fae5adSBarry Smith int MatMult_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
598f9fae5adSBarry Smith {
5994cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
600f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
601*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
602273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
603f9fae5adSBarry Smith   int         n,i,jrow,j;
604f9fae5adSBarry Smith 
605f9fae5adSBarry Smith   PetscFunctionBegin;
606f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
607f9fae5adSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
608f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
609f9fae5adSBarry Smith   idx  = a->j;
610f9fae5adSBarry Smith   v    = a->a;
611f9fae5adSBarry Smith   ii   = a->i;
612f9fae5adSBarry Smith 
613f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
614f9fae5adSBarry Smith   idx  += shift;
615f9fae5adSBarry Smith   for (i=0; i<m; i++) {
616f9fae5adSBarry Smith     jrow = ii[i];
617f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
618f9fae5adSBarry Smith     sum1  = 0.0;
619f9fae5adSBarry Smith     sum2  = 0.0;
620f9fae5adSBarry Smith     sum3  = 0.0;
621f9fae5adSBarry Smith     sum4  = 0.0;
622f9fae5adSBarry Smith     sum5  = 0.0;
623f9fae5adSBarry Smith     for (j=0; j<n; j++) {
624f9fae5adSBarry Smith       sum1 += v[jrow]*x[5*idx[jrow]];
625f9fae5adSBarry Smith       sum2 += v[jrow]*x[5*idx[jrow]+1];
626f9fae5adSBarry Smith       sum3 += v[jrow]*x[5*idx[jrow]+2];
627f9fae5adSBarry Smith       sum4 += v[jrow]*x[5*idx[jrow]+3];
628f9fae5adSBarry Smith       sum5 += v[jrow]*x[5*idx[jrow]+4];
629f9fae5adSBarry Smith       jrow++;
630f9fae5adSBarry Smith      }
631f9fae5adSBarry Smith     y[5*i]   = sum1;
632f9fae5adSBarry Smith     y[5*i+1] = sum2;
633f9fae5adSBarry Smith     y[5*i+2] = sum3;
634f9fae5adSBarry Smith     y[5*i+3] = sum4;
635f9fae5adSBarry Smith     y[5*i+4] = sum5;
636f9fae5adSBarry Smith   }
637f9fae5adSBarry Smith 
638b0a32e0cSBarry Smith   PetscLogFlops(10*a->nz - 5*m);
639f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
640f9fae5adSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
641f9fae5adSBarry Smith   PetscFunctionReturn(0);
642f9fae5adSBarry Smith }
643f9fae5adSBarry Smith 
6444a2ae208SSatish Balay #undef __FUNCT__
645b9617806SBarry Smith #define __FUNCT__ "MatMultTranspose_SeqMAIJ_5"
646f9fae5adSBarry Smith int MatMultTranspose_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
647f9fae5adSBarry Smith {
6484cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
649f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
650*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,zero = 0.0;
651273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
652f9fae5adSBarry Smith 
653f9fae5adSBarry Smith   PetscFunctionBegin;
654f9fae5adSBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
655f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
656f9fae5adSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
657f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
658f9fae5adSBarry Smith   for (i=0; i<m; i++) {
659f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
660f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
661f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
662f9fae5adSBarry Smith     alpha1 = x[5*i];
663f9fae5adSBarry Smith     alpha2 = x[5*i+1];
664f9fae5adSBarry Smith     alpha3 = x[5*i+2];
665f9fae5adSBarry Smith     alpha4 = x[5*i+3];
666f9fae5adSBarry Smith     alpha5 = x[5*i+4];
667f9fae5adSBarry Smith     while (n-->0) {
668f9fae5adSBarry Smith       y[5*(*idx)]   += alpha1*(*v);
669f9fae5adSBarry Smith       y[5*(*idx)+1] += alpha2*(*v);
670f9fae5adSBarry Smith       y[5*(*idx)+2] += alpha3*(*v);
671f9fae5adSBarry Smith       y[5*(*idx)+3] += alpha4*(*v);
672f9fae5adSBarry Smith       y[5*(*idx)+4] += alpha5*(*v);
673f9fae5adSBarry Smith       idx++; v++;
674f9fae5adSBarry Smith     }
675f9fae5adSBarry Smith   }
676b0a32e0cSBarry Smith   PetscLogFlops(10*a->nz - 5*b->AIJ->n);
677f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
678f9fae5adSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
679f9fae5adSBarry Smith   PetscFunctionReturn(0);
680f9fae5adSBarry Smith }
681f9fae5adSBarry Smith 
6824a2ae208SSatish Balay #undef __FUNCT__
683b9617806SBarry Smith #define __FUNCT__ "MatMultAdd_SeqMAIJ_5"
684f9fae5adSBarry Smith int MatMultAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
685f9fae5adSBarry Smith {
6864cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
687f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
688*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
689273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
690f9fae5adSBarry Smith   int         n,i,jrow,j;
691f9fae5adSBarry Smith 
692f9fae5adSBarry Smith   PetscFunctionBegin;
693f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
694f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
695f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
696f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
697f9fae5adSBarry Smith   idx  = a->j;
698f9fae5adSBarry Smith   v    = a->a;
699f9fae5adSBarry Smith   ii   = a->i;
700f9fae5adSBarry Smith 
701f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
702f9fae5adSBarry Smith   idx  += shift;
703f9fae5adSBarry Smith   for (i=0; i<m; i++) {
704f9fae5adSBarry Smith     jrow = ii[i];
705f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
706f9fae5adSBarry Smith     sum1  = 0.0;
707f9fae5adSBarry Smith     sum2  = 0.0;
708f9fae5adSBarry Smith     sum3  = 0.0;
709f9fae5adSBarry Smith     sum4  = 0.0;
710f9fae5adSBarry Smith     sum5  = 0.0;
711f9fae5adSBarry Smith     for (j=0; j<n; j++) {
712f9fae5adSBarry Smith       sum1 += v[jrow]*x[5*idx[jrow]];
713f9fae5adSBarry Smith       sum2 += v[jrow]*x[5*idx[jrow]+1];
714f9fae5adSBarry Smith       sum3 += v[jrow]*x[5*idx[jrow]+2];
715f9fae5adSBarry Smith       sum4 += v[jrow]*x[5*idx[jrow]+3];
716f9fae5adSBarry Smith       sum5 += v[jrow]*x[5*idx[jrow]+4];
717f9fae5adSBarry Smith       jrow++;
718f9fae5adSBarry Smith      }
719f9fae5adSBarry Smith     y[5*i]   += sum1;
720f9fae5adSBarry Smith     y[5*i+1] += sum2;
721f9fae5adSBarry Smith     y[5*i+2] += sum3;
722f9fae5adSBarry Smith     y[5*i+3] += sum4;
723f9fae5adSBarry Smith     y[5*i+4] += sum5;
724f9fae5adSBarry Smith   }
725f9fae5adSBarry Smith 
726b0a32e0cSBarry Smith   PetscLogFlops(10*a->nz);
727f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
728f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
729f9fae5adSBarry Smith   PetscFunctionReturn(0);
730f9fae5adSBarry Smith }
731f9fae5adSBarry Smith 
7324a2ae208SSatish Balay #undef __FUNCT__
733b9617806SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_5"
734f9fae5adSBarry Smith int MatMultTransposeAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
735f9fae5adSBarry Smith {
7364cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
737f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
738*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5;
739273d9f13SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
740f9fae5adSBarry Smith 
741f9fae5adSBarry Smith   PetscFunctionBegin;
742f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
743f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
744f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
745f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
746f9fae5adSBarry Smith   for (i=0; i<m; i++) {
747f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
748f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
749f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
750f9fae5adSBarry Smith     alpha1 = x[5*i];
751f9fae5adSBarry Smith     alpha2 = x[5*i+1];
752f9fae5adSBarry Smith     alpha3 = x[5*i+2];
753f9fae5adSBarry Smith     alpha4 = x[5*i+3];
754f9fae5adSBarry Smith     alpha5 = x[5*i+4];
755f9fae5adSBarry Smith     while (n-->0) {
756f9fae5adSBarry Smith       y[5*(*idx)]   += alpha1*(*v);
757f9fae5adSBarry Smith       y[5*(*idx)+1] += alpha2*(*v);
758f9fae5adSBarry Smith       y[5*(*idx)+2] += alpha3*(*v);
759f9fae5adSBarry Smith       y[5*(*idx)+3] += alpha4*(*v);
760f9fae5adSBarry Smith       y[5*(*idx)+4] += alpha5*(*v);
761f9fae5adSBarry Smith       idx++; v++;
762f9fae5adSBarry Smith     }
763f9fae5adSBarry Smith   }
764b0a32e0cSBarry Smith   PetscLogFlops(10*a->nz);
765f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
766f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
767f9fae5adSBarry Smith   PetscFunctionReturn(0);
768bcc973b7SBarry Smith }
769bcc973b7SBarry Smith 
7706cd98798SBarry Smith /* ------------------------------------------------------------------------------*/
7714a2ae208SSatish Balay #undef __FUNCT__
772b9617806SBarry Smith #define __FUNCT__ "MatMult_SeqMAIJ_6"
7736cd98798SBarry Smith int MatMult_SeqMAIJ_6(Mat A,Vec xx,Vec yy)
7746cd98798SBarry Smith {
7756cd98798SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
7766cd98798SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
777*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6;
7786cd98798SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
7796cd98798SBarry Smith   int         n,i,jrow,j;
7806cd98798SBarry Smith 
7816cd98798SBarry Smith   PetscFunctionBegin;
7826cd98798SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
7836cd98798SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
7846cd98798SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
7856cd98798SBarry Smith   idx  = a->j;
7866cd98798SBarry Smith   v    = a->a;
7876cd98798SBarry Smith   ii   = a->i;
7886cd98798SBarry Smith 
7896cd98798SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
7906cd98798SBarry Smith   idx  += shift;
7916cd98798SBarry Smith   for (i=0; i<m; i++) {
7926cd98798SBarry Smith     jrow = ii[i];
7936cd98798SBarry Smith     n    = ii[i+1] - jrow;
7946cd98798SBarry Smith     sum1  = 0.0;
7956cd98798SBarry Smith     sum2  = 0.0;
7966cd98798SBarry Smith     sum3  = 0.0;
7976cd98798SBarry Smith     sum4  = 0.0;
7986cd98798SBarry Smith     sum5  = 0.0;
7996cd98798SBarry Smith     sum6  = 0.0;
8006cd98798SBarry Smith     for (j=0; j<n; j++) {
8016cd98798SBarry Smith       sum1 += v[jrow]*x[6*idx[jrow]];
8026cd98798SBarry Smith       sum2 += v[jrow]*x[6*idx[jrow]+1];
8036cd98798SBarry Smith       sum3 += v[jrow]*x[6*idx[jrow]+2];
8046cd98798SBarry Smith       sum4 += v[jrow]*x[6*idx[jrow]+3];
8056cd98798SBarry Smith       sum5 += v[jrow]*x[6*idx[jrow]+4];
8066cd98798SBarry Smith       sum6 += v[jrow]*x[6*idx[jrow]+5];
8076cd98798SBarry Smith       jrow++;
8086cd98798SBarry Smith      }
8096cd98798SBarry Smith     y[6*i]   = sum1;
8106cd98798SBarry Smith     y[6*i+1] = sum2;
8116cd98798SBarry Smith     y[6*i+2] = sum3;
8126cd98798SBarry Smith     y[6*i+3] = sum4;
8136cd98798SBarry Smith     y[6*i+4] = sum5;
8146cd98798SBarry Smith     y[6*i+5] = sum6;
8156cd98798SBarry Smith   }
8166cd98798SBarry Smith 
817b0a32e0cSBarry Smith   PetscLogFlops(12*a->nz - 6*m);
8186cd98798SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8196cd98798SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8206cd98798SBarry Smith   PetscFunctionReturn(0);
8216cd98798SBarry Smith }
8226cd98798SBarry Smith 
8234a2ae208SSatish Balay #undef __FUNCT__
824b9617806SBarry Smith #define __FUNCT__ "MatMultTranspose_SeqMAIJ_6"
8256cd98798SBarry Smith int MatMultTranspose_SeqMAIJ_6(Mat A,Vec xx,Vec yy)
8266cd98798SBarry Smith {
8276cd98798SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
8286cd98798SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
829*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,zero = 0.0;
8306cd98798SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
8316cd98798SBarry Smith 
8326cd98798SBarry Smith   PetscFunctionBegin;
8336cd98798SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
8346cd98798SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
8356cd98798SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
8366cd98798SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
8376cd98798SBarry Smith   for (i=0; i<m; i++) {
8386cd98798SBarry Smith     idx    = a->j + a->i[i] + shift;
8396cd98798SBarry Smith     v      = a->a + a->i[i] + shift;
8406cd98798SBarry Smith     n      = a->i[i+1] - a->i[i];
8416cd98798SBarry Smith     alpha1 = x[6*i];
8426cd98798SBarry Smith     alpha2 = x[6*i+1];
8436cd98798SBarry Smith     alpha3 = x[6*i+2];
8446cd98798SBarry Smith     alpha4 = x[6*i+3];
8456cd98798SBarry Smith     alpha5 = x[6*i+4];
8466cd98798SBarry Smith     alpha6 = x[6*i+5];
8476cd98798SBarry Smith     while (n-->0) {
8486cd98798SBarry Smith       y[6*(*idx)]   += alpha1*(*v);
8496cd98798SBarry Smith       y[6*(*idx)+1] += alpha2*(*v);
8506cd98798SBarry Smith       y[6*(*idx)+2] += alpha3*(*v);
8516cd98798SBarry Smith       y[6*(*idx)+3] += alpha4*(*v);
8526cd98798SBarry Smith       y[6*(*idx)+4] += alpha5*(*v);
8536cd98798SBarry Smith       y[6*(*idx)+5] += alpha6*(*v);
8546cd98798SBarry Smith       idx++; v++;
8556cd98798SBarry Smith     }
8566cd98798SBarry Smith   }
857b0a32e0cSBarry Smith   PetscLogFlops(12*a->nz - 6*b->AIJ->n);
8586cd98798SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8596cd98798SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8606cd98798SBarry Smith   PetscFunctionReturn(0);
8616cd98798SBarry Smith }
8626cd98798SBarry Smith 
8634a2ae208SSatish Balay #undef __FUNCT__
864b9617806SBarry Smith #define __FUNCT__ "MatMultAdd_SeqMAIJ_6"
8656cd98798SBarry Smith int MatMultAdd_SeqMAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
8666cd98798SBarry Smith {
8676cd98798SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
8686cd98798SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
869*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6;
8706cd98798SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
8716cd98798SBarry Smith   int         n,i,jrow,j;
8726cd98798SBarry Smith 
8736cd98798SBarry Smith   PetscFunctionBegin;
8746cd98798SBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
8756cd98798SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
8766cd98798SBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
8776cd98798SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
8786cd98798SBarry Smith   idx  = a->j;
8796cd98798SBarry Smith   v    = a->a;
8806cd98798SBarry Smith   ii   = a->i;
8816cd98798SBarry Smith 
8826cd98798SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
8836cd98798SBarry Smith   idx  += shift;
8846cd98798SBarry Smith   for (i=0; i<m; i++) {
8856cd98798SBarry Smith     jrow = ii[i];
8866cd98798SBarry Smith     n    = ii[i+1] - jrow;
8876cd98798SBarry Smith     sum1  = 0.0;
8886cd98798SBarry Smith     sum2  = 0.0;
8896cd98798SBarry Smith     sum3  = 0.0;
8906cd98798SBarry Smith     sum4  = 0.0;
8916cd98798SBarry Smith     sum5  = 0.0;
8926cd98798SBarry Smith     sum6  = 0.0;
8936cd98798SBarry Smith     for (j=0; j<n; j++) {
8946cd98798SBarry Smith       sum1 += v[jrow]*x[6*idx[jrow]];
8956cd98798SBarry Smith       sum2 += v[jrow]*x[6*idx[jrow]+1];
8966cd98798SBarry Smith       sum3 += v[jrow]*x[6*idx[jrow]+2];
8976cd98798SBarry Smith       sum4 += v[jrow]*x[6*idx[jrow]+3];
8986cd98798SBarry Smith       sum5 += v[jrow]*x[6*idx[jrow]+4];
8996cd98798SBarry Smith       sum6 += v[jrow]*x[6*idx[jrow]+5];
9006cd98798SBarry Smith       jrow++;
9016cd98798SBarry Smith      }
9026cd98798SBarry Smith     y[6*i]   += sum1;
9036cd98798SBarry Smith     y[6*i+1] += sum2;
9046cd98798SBarry Smith     y[6*i+2] += sum3;
9056cd98798SBarry Smith     y[6*i+3] += sum4;
9066cd98798SBarry Smith     y[6*i+4] += sum5;
9076cd98798SBarry Smith     y[6*i+5] += sum6;
9086cd98798SBarry Smith   }
9096cd98798SBarry Smith 
910b0a32e0cSBarry Smith   PetscLogFlops(12*a->nz);
9116cd98798SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9126cd98798SBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
9136cd98798SBarry Smith   PetscFunctionReturn(0);
9146cd98798SBarry Smith }
9156cd98798SBarry Smith 
9164a2ae208SSatish Balay #undef __FUNCT__
917b9617806SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_6"
9186cd98798SBarry Smith int MatMultTransposeAdd_SeqMAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
9196cd98798SBarry Smith {
9206cd98798SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
9216cd98798SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
922*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6;
9236cd98798SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
9246cd98798SBarry Smith 
9256cd98798SBarry Smith   PetscFunctionBegin;
9266cd98798SBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
9276cd98798SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
9286cd98798SBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
9296cd98798SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
9306cd98798SBarry Smith   for (i=0; i<m; i++) {
9316cd98798SBarry Smith     idx    = a->j + a->i[i] + shift;
9326cd98798SBarry Smith     v      = a->a + a->i[i] + shift;
9336cd98798SBarry Smith     n      = a->i[i+1] - a->i[i];
9346cd98798SBarry Smith     alpha1 = x[6*i];
9356cd98798SBarry Smith     alpha2 = x[6*i+1];
9366cd98798SBarry Smith     alpha3 = x[6*i+2];
9376cd98798SBarry Smith     alpha4 = x[6*i+3];
9386cd98798SBarry Smith     alpha5 = x[6*i+4];
9396cd98798SBarry Smith     alpha6 = x[6*i+5];
9406cd98798SBarry Smith     while (n-->0) {
9416cd98798SBarry Smith       y[6*(*idx)]   += alpha1*(*v);
9426cd98798SBarry Smith       y[6*(*idx)+1] += alpha2*(*v);
9436cd98798SBarry Smith       y[6*(*idx)+2] += alpha3*(*v);
9446cd98798SBarry Smith       y[6*(*idx)+3] += alpha4*(*v);
9456cd98798SBarry Smith       y[6*(*idx)+4] += alpha5*(*v);
9466cd98798SBarry Smith       y[6*(*idx)+5] += alpha6*(*v);
9476cd98798SBarry Smith       idx++; v++;
9486cd98798SBarry Smith     }
9496cd98798SBarry Smith   }
950b0a32e0cSBarry Smith   PetscLogFlops(12*a->nz);
9516cd98798SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9526cd98798SBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
9536cd98798SBarry Smith   PetscFunctionReturn(0);
9546cd98798SBarry Smith }
9556cd98798SBarry Smith 
95666ed3db0SBarry Smith /* ------------------------------------------------------------------------------*/
95766ed3db0SBarry Smith #undef __FUNCT__
95866ed3db0SBarry Smith #define __FUNCT__ "MatMult_SeqMAIJ_8"
95966ed3db0SBarry Smith int MatMult_SeqMAIJ_8(Mat A,Vec xx,Vec yy)
96066ed3db0SBarry Smith {
96166ed3db0SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
96266ed3db0SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
963*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7, sum8;
96466ed3db0SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
96566ed3db0SBarry Smith   int         n,i,jrow,j;
96666ed3db0SBarry Smith 
96766ed3db0SBarry Smith   PetscFunctionBegin;
96866ed3db0SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
96966ed3db0SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
97066ed3db0SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
97166ed3db0SBarry Smith   idx  = a->j;
97266ed3db0SBarry Smith   v    = a->a;
97366ed3db0SBarry Smith   ii   = a->i;
97466ed3db0SBarry Smith 
97566ed3db0SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
97666ed3db0SBarry Smith   idx  += shift;
97766ed3db0SBarry Smith   for (i=0; i<m; i++) {
97866ed3db0SBarry Smith     jrow = ii[i];
97966ed3db0SBarry Smith     n    = ii[i+1] - jrow;
98066ed3db0SBarry Smith     sum1  = 0.0;
98166ed3db0SBarry Smith     sum2  = 0.0;
98266ed3db0SBarry Smith     sum3  = 0.0;
98366ed3db0SBarry Smith     sum4  = 0.0;
98466ed3db0SBarry Smith     sum5  = 0.0;
98566ed3db0SBarry Smith     sum6  = 0.0;
98666ed3db0SBarry Smith     sum7  = 0.0;
98766ed3db0SBarry Smith     sum8  = 0.0;
98866ed3db0SBarry Smith     for (j=0; j<n; j++) {
98966ed3db0SBarry Smith       sum1 += v[jrow]*x[8*idx[jrow]];
99066ed3db0SBarry Smith       sum2 += v[jrow]*x[8*idx[jrow]+1];
99166ed3db0SBarry Smith       sum3 += v[jrow]*x[8*idx[jrow]+2];
99266ed3db0SBarry Smith       sum4 += v[jrow]*x[8*idx[jrow]+3];
99366ed3db0SBarry Smith       sum5 += v[jrow]*x[8*idx[jrow]+4];
99466ed3db0SBarry Smith       sum6 += v[jrow]*x[8*idx[jrow]+5];
99566ed3db0SBarry Smith       sum7 += v[jrow]*x[8*idx[jrow]+6];
99666ed3db0SBarry Smith       sum8 += v[jrow]*x[8*idx[jrow]+7];
99766ed3db0SBarry Smith       jrow++;
99866ed3db0SBarry Smith      }
99966ed3db0SBarry Smith     y[8*i]   = sum1;
100066ed3db0SBarry Smith     y[8*i+1] = sum2;
100166ed3db0SBarry Smith     y[8*i+2] = sum3;
100266ed3db0SBarry Smith     y[8*i+3] = sum4;
100366ed3db0SBarry Smith     y[8*i+4] = sum5;
100466ed3db0SBarry Smith     y[8*i+5] = sum6;
100566ed3db0SBarry Smith     y[8*i+6] = sum7;
100666ed3db0SBarry Smith     y[8*i+7] = sum8;
100766ed3db0SBarry Smith   }
100866ed3db0SBarry Smith 
100966ed3db0SBarry Smith   PetscLogFlops(16*a->nz - 8*m);
101066ed3db0SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
101166ed3db0SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
101266ed3db0SBarry Smith   PetscFunctionReturn(0);
101366ed3db0SBarry Smith }
101466ed3db0SBarry Smith 
101566ed3db0SBarry Smith #undef __FUNCT__
101666ed3db0SBarry Smith #define __FUNCT__ "MatMultTranspose_SeqMAIJ_8"
101766ed3db0SBarry Smith int MatMultTranspose_SeqMAIJ_8(Mat A,Vec xx,Vec yy)
101866ed3db0SBarry Smith {
101966ed3db0SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
102066ed3db0SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
1021*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,alpha8,zero = 0.0;
102266ed3db0SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
102366ed3db0SBarry Smith 
102466ed3db0SBarry Smith   PetscFunctionBegin;
102566ed3db0SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
102666ed3db0SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
102766ed3db0SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
102866ed3db0SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
102966ed3db0SBarry Smith   for (i=0; i<m; i++) {
103066ed3db0SBarry Smith     idx    = a->j + a->i[i] + shift;
103166ed3db0SBarry Smith     v      = a->a + a->i[i] + shift;
103266ed3db0SBarry Smith     n      = a->i[i+1] - a->i[i];
103366ed3db0SBarry Smith     alpha1 = x[8*i];
103466ed3db0SBarry Smith     alpha2 = x[8*i+1];
103566ed3db0SBarry Smith     alpha3 = x[8*i+2];
103666ed3db0SBarry Smith     alpha4 = x[8*i+3];
103766ed3db0SBarry Smith     alpha5 = x[8*i+4];
103866ed3db0SBarry Smith     alpha6 = x[8*i+5];
103966ed3db0SBarry Smith     alpha7 = x[8*i+6];
104066ed3db0SBarry Smith     alpha8 = x[8*i+7];
104166ed3db0SBarry Smith     while (n-->0) {
104266ed3db0SBarry Smith       y[8*(*idx)]   += alpha1*(*v);
104366ed3db0SBarry Smith       y[8*(*idx)+1] += alpha2*(*v);
104466ed3db0SBarry Smith       y[8*(*idx)+2] += alpha3*(*v);
104566ed3db0SBarry Smith       y[8*(*idx)+3] += alpha4*(*v);
104666ed3db0SBarry Smith       y[8*(*idx)+4] += alpha5*(*v);
104766ed3db0SBarry Smith       y[8*(*idx)+5] += alpha6*(*v);
104866ed3db0SBarry Smith       y[8*(*idx)+6] += alpha7*(*v);
104966ed3db0SBarry Smith       y[8*(*idx)+7] += alpha8*(*v);
105066ed3db0SBarry Smith       idx++; v++;
105166ed3db0SBarry Smith     }
105266ed3db0SBarry Smith   }
105366ed3db0SBarry Smith   PetscLogFlops(16*a->nz - 8*b->AIJ->n);
105466ed3db0SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
105566ed3db0SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
105666ed3db0SBarry Smith   PetscFunctionReturn(0);
105766ed3db0SBarry Smith }
105866ed3db0SBarry Smith 
105966ed3db0SBarry Smith #undef __FUNCT__
106066ed3db0SBarry Smith #define __FUNCT__ "MatMultAdd_SeqMAIJ_8"
106166ed3db0SBarry Smith int MatMultAdd_SeqMAIJ_8(Mat A,Vec xx,Vec yy,Vec zz)
106266ed3db0SBarry Smith {
106366ed3db0SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
106466ed3db0SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
1065*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6, sum7, sum8;
106666ed3db0SBarry Smith   int         ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii;
106766ed3db0SBarry Smith   int         n,i,jrow,j;
106866ed3db0SBarry Smith 
106966ed3db0SBarry Smith   PetscFunctionBegin;
107066ed3db0SBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
107166ed3db0SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
107266ed3db0SBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
107366ed3db0SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
107466ed3db0SBarry Smith   idx  = a->j;
107566ed3db0SBarry Smith   v    = a->a;
107666ed3db0SBarry Smith   ii   = a->i;
107766ed3db0SBarry Smith 
107866ed3db0SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
107966ed3db0SBarry Smith   idx  += shift;
108066ed3db0SBarry Smith   for (i=0; i<m; i++) {
108166ed3db0SBarry Smith     jrow = ii[i];
108266ed3db0SBarry Smith     n    = ii[i+1] - jrow;
108366ed3db0SBarry Smith     sum1  = 0.0;
108466ed3db0SBarry Smith     sum2  = 0.0;
108566ed3db0SBarry Smith     sum3  = 0.0;
108666ed3db0SBarry Smith     sum4  = 0.0;
108766ed3db0SBarry Smith     sum5  = 0.0;
108866ed3db0SBarry Smith     sum6  = 0.0;
108966ed3db0SBarry Smith     sum7  = 0.0;
109066ed3db0SBarry Smith     sum8  = 0.0;
109166ed3db0SBarry Smith     for (j=0; j<n; j++) {
109266ed3db0SBarry Smith       sum1 += v[jrow]*x[8*idx[jrow]];
109366ed3db0SBarry Smith       sum2 += v[jrow]*x[8*idx[jrow]+1];
109466ed3db0SBarry Smith       sum3 += v[jrow]*x[8*idx[jrow]+2];
109566ed3db0SBarry Smith       sum4 += v[jrow]*x[8*idx[jrow]+3];
109666ed3db0SBarry Smith       sum5 += v[jrow]*x[8*idx[jrow]+4];
109766ed3db0SBarry Smith       sum6 += v[jrow]*x[8*idx[jrow]+5];
109866ed3db0SBarry Smith       sum7 += v[jrow]*x[8*idx[jrow]+6];
109966ed3db0SBarry Smith       sum8 += v[jrow]*x[8*idx[jrow]+7];
110066ed3db0SBarry Smith       jrow++;
110166ed3db0SBarry Smith      }
110266ed3db0SBarry Smith     y[8*i]   += sum1;
110366ed3db0SBarry Smith     y[8*i+1] += sum2;
110466ed3db0SBarry Smith     y[8*i+2] += sum3;
110566ed3db0SBarry Smith     y[8*i+3] += sum4;
110666ed3db0SBarry Smith     y[8*i+4] += sum5;
110766ed3db0SBarry Smith     y[8*i+5] += sum6;
110866ed3db0SBarry Smith     y[8*i+6] += sum7;
110966ed3db0SBarry Smith     y[8*i+7] += sum8;
111066ed3db0SBarry Smith   }
111166ed3db0SBarry Smith 
111266ed3db0SBarry Smith   PetscLogFlops(16*a->nz);
111366ed3db0SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
111466ed3db0SBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
111566ed3db0SBarry Smith   PetscFunctionReturn(0);
111666ed3db0SBarry Smith }
111766ed3db0SBarry Smith 
111866ed3db0SBarry Smith #undef __FUNCT__
111966ed3db0SBarry Smith #define __FUNCT__ "MatMultTransposeAdd_SeqMAIJ_8"
112066ed3db0SBarry Smith int MatMultTransposeAdd_SeqMAIJ_8(Mat A,Vec xx,Vec yy,Vec zz)
112166ed3db0SBarry Smith {
112266ed3db0SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
112366ed3db0SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
1124*87828ca2SBarry Smith   PetscScalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,alpha7,alpha8;
112566ed3db0SBarry Smith   int         ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift;
112666ed3db0SBarry Smith 
112766ed3db0SBarry Smith   PetscFunctionBegin;
112866ed3db0SBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
112966ed3db0SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
113066ed3db0SBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
113166ed3db0SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
113266ed3db0SBarry Smith   for (i=0; i<m; i++) {
113366ed3db0SBarry Smith     idx    = a->j + a->i[i] + shift;
113466ed3db0SBarry Smith     v      = a->a + a->i[i] + shift;
113566ed3db0SBarry Smith     n      = a->i[i+1] - a->i[i];
113666ed3db0SBarry Smith     alpha1 = x[8*i];
113766ed3db0SBarry Smith     alpha2 = x[8*i+1];
113866ed3db0SBarry Smith     alpha3 = x[8*i+2];
113966ed3db0SBarry Smith     alpha4 = x[8*i+3];
114066ed3db0SBarry Smith     alpha5 = x[8*i+4];
114166ed3db0SBarry Smith     alpha6 = x[8*i+5];
114266ed3db0SBarry Smith     alpha7 = x[8*i+6];
114366ed3db0SBarry Smith     alpha8 = x[8*i+7];
114466ed3db0SBarry Smith     while (n-->0) {
114566ed3db0SBarry Smith       y[8*(*idx)]   += alpha1*(*v);
114666ed3db0SBarry Smith       y[8*(*idx)+1] += alpha2*(*v);
114766ed3db0SBarry Smith       y[8*(*idx)+2] += alpha3*(*v);
114866ed3db0SBarry Smith       y[8*(*idx)+3] += alpha4*(*v);
114966ed3db0SBarry Smith       y[8*(*idx)+4] += alpha5*(*v);
115066ed3db0SBarry Smith       y[8*(*idx)+5] += alpha6*(*v);
115166ed3db0SBarry Smith       y[8*(*idx)+6] += alpha7*(*v);
115266ed3db0SBarry Smith       y[8*(*idx)+7] += alpha8*(*v);
115366ed3db0SBarry Smith       idx++; v++;
115466ed3db0SBarry Smith     }
115566ed3db0SBarry Smith   }
115666ed3db0SBarry Smith   PetscLogFlops(16*a->nz);
115766ed3db0SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
115866ed3db0SBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
115966ed3db0SBarry Smith   PetscFunctionReturn(0);
116066ed3db0SBarry Smith }
116166ed3db0SBarry Smith 
1162f4a53059SBarry Smith /*===================================================================================*/
11634a2ae208SSatish Balay #undef __FUNCT__
11644a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIMAIJ_dof"
1165f4a53059SBarry Smith int MatMult_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
1166f4a53059SBarry Smith {
11674cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
1168f4a53059SBarry Smith   int         ierr;
1169f4a53059SBarry Smith   PetscFunctionBegin;
1170f4a53059SBarry Smith 
1171f4a53059SBarry Smith   /* start the scatter */
1172f4a53059SBarry Smith   ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
11734cb9d9b8SBarry Smith   ierr = (*b->AIJ->ops->mult)(b->AIJ,xx,yy);CHKERRQ(ierr);
1174f4a53059SBarry Smith   ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
11754cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr);
1176f4a53059SBarry Smith   PetscFunctionReturn(0);
1177f4a53059SBarry Smith }
1178f4a53059SBarry Smith 
11794a2ae208SSatish Balay #undef __FUNCT__
11804a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIMAIJ_dof"
11814cb9d9b8SBarry Smith int MatMultTranspose_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
11824cb9d9b8SBarry Smith {
11834cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
11844cb9d9b8SBarry Smith   int         ierr;
11854cb9d9b8SBarry Smith   PetscFunctionBegin;
11864cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr);
11874cb9d9b8SBarry Smith   ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
11884cb9d9b8SBarry Smith   ierr = (*b->AIJ->ops->multtranspose)(b->AIJ,xx,yy);CHKERRQ(ierr);
11894cb9d9b8SBarry Smith   ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
11904cb9d9b8SBarry Smith   PetscFunctionReturn(0);
11914cb9d9b8SBarry Smith }
11924cb9d9b8SBarry Smith 
11934a2ae208SSatish Balay #undef __FUNCT__
11944a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIMAIJ_dof"
1195d72c5c08SBarry Smith int MatMultAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz)
11964cb9d9b8SBarry Smith {
11974cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
11984cb9d9b8SBarry Smith   int         ierr;
11994cb9d9b8SBarry Smith   PetscFunctionBegin;
12004cb9d9b8SBarry Smith 
12014cb9d9b8SBarry Smith   /* start the scatter */
12024cb9d9b8SBarry Smith   ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
1203d72c5c08SBarry Smith   ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr);
12044cb9d9b8SBarry Smith   ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
1205d72c5c08SBarry Smith   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,zz);CHKERRQ(ierr);
12064cb9d9b8SBarry Smith   PetscFunctionReturn(0);
12074cb9d9b8SBarry Smith }
12084cb9d9b8SBarry Smith 
12094a2ae208SSatish Balay #undef __FUNCT__
12104a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIMAIJ_dof"
1211d72c5c08SBarry Smith int MatMultTransposeAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz)
12124cb9d9b8SBarry Smith {
12134cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
12144cb9d9b8SBarry Smith   int         ierr;
12154cb9d9b8SBarry Smith   PetscFunctionBegin;
12164cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr);
1217d72c5c08SBarry Smith   ierr = VecScatterBegin(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
1218d72c5c08SBarry Smith   ierr = (*b->AIJ->ops->multtransposeadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr);
1219d72c5c08SBarry Smith   ierr = VecScatterEnd(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
12204cb9d9b8SBarry Smith   PetscFunctionReturn(0);
12214cb9d9b8SBarry Smith }
12224cb9d9b8SBarry Smith 
1223bcc973b7SBarry Smith /* ---------------------------------------------------------------------------------- */
12244a2ae208SSatish Balay #undef __FUNCT__
12254a2ae208SSatish Balay #define __FUNCT__ "MatCreateMAIJ"
1226cd3bbe55SBarry Smith int MatCreateMAIJ(Mat A,int dof,Mat *maij)
122782b1193eSBarry Smith {
1228f4a53059SBarry Smith   int         ierr,size,n;
12294cb9d9b8SBarry Smith   Mat_MPIMAIJ *b;
123082b1193eSBarry Smith   Mat         B;
123182b1193eSBarry Smith 
123282b1193eSBarry Smith   PetscFunctionBegin;
1233d72c5c08SBarry Smith   ierr   = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1234d72c5c08SBarry Smith 
1235d72c5c08SBarry Smith   if (dof == 1) {
1236d72c5c08SBarry Smith     *maij = A;
1237d72c5c08SBarry Smith   } else {
123883903688SBarry Smith     ierr = MatCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr);
1239cd3bbe55SBarry Smith     B->assembled    = PETSC_TRUE;
1240d72c5c08SBarry Smith 
1241f4a53059SBarry Smith     ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1242f4a53059SBarry Smith     if (size == 1) {
1243b9b97703SBarry Smith       ierr = MatSetType(B,MATSEQMAIJ);CHKERRQ(ierr);
12444cb9d9b8SBarry Smith       B->ops->destroy = MatDestroy_SeqMAIJ;
1245b9b97703SBarry Smith       b      = (Mat_MPIMAIJ*)B->data;
1246b9b97703SBarry Smith       b->dof = dof;
12474cb9d9b8SBarry Smith       b->AIJ = A;
1248d72c5c08SBarry Smith       if (dof == 2) {
1249cd3bbe55SBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_2;
1250cd3bbe55SBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_2;
1251cd3bbe55SBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_2;
1252cd3bbe55SBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2;
1253bcc973b7SBarry Smith       } else if (dof == 3) {
1254bcc973b7SBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_3;
1255bcc973b7SBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_3;
1256bcc973b7SBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_3;
1257bcc973b7SBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3;
1258bcc973b7SBarry Smith       } else if (dof == 4) {
1259bcc973b7SBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_4;
1260bcc973b7SBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_4;
1261bcc973b7SBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_4;
1262bcc973b7SBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4;
1263f9fae5adSBarry Smith       } else if (dof == 5) {
1264f9fae5adSBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_5;
1265f9fae5adSBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_5;
1266f9fae5adSBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_5;
1267f9fae5adSBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5;
12686cd98798SBarry Smith       } else if (dof == 6) {
12696cd98798SBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_6;
12706cd98798SBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_6;
12716cd98798SBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_6;
12726cd98798SBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_6;
127366ed3db0SBarry Smith       } else if (dof == 8) {
127466ed3db0SBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_8;
127566ed3db0SBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_8;
127666ed3db0SBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_8;
127766ed3db0SBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_8;
127882b1193eSBarry Smith       } else {
127966ed3db0SBarry Smith         SETERRQ1(1,"Cannot handle a dof of %d. Send request for code to petsc-maint@mcs.anl.gov\n",dof);
128082b1193eSBarry Smith       }
1281f4a53059SBarry Smith     } else {
1282f4a53059SBarry Smith       Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data;
1283f4a53059SBarry Smith       IS         from,to;
1284f4a53059SBarry Smith       Vec        gvec;
1285f4a53059SBarry Smith       int        *garray,i;
1286f4a53059SBarry Smith 
1287b9b97703SBarry Smith       ierr = MatSetType(B,MATMPIMAIJ);CHKERRQ(ierr);
1288d72c5c08SBarry Smith       B->ops->destroy = MatDestroy_MPIMAIJ;
1289b9b97703SBarry Smith       b      = (Mat_MPIMAIJ*)B->data;
1290b9b97703SBarry Smith       b->dof = dof;
1291b9b97703SBarry Smith       b->A   = A;
12924cb9d9b8SBarry Smith       ierr = MatCreateMAIJ(mpiaij->A,dof,&b->AIJ);CHKERRQ(ierr);
12934cb9d9b8SBarry Smith       ierr = MatCreateMAIJ(mpiaij->B,dof,&b->OAIJ);CHKERRQ(ierr);
12944cb9d9b8SBarry Smith 
1295f4a53059SBarry Smith       ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr);
1296f4a53059SBarry Smith       ierr = VecCreateSeq(PETSC_COMM_SELF,n*dof,&b->w);CHKERRQ(ierr);
1297f4a53059SBarry Smith 
1298f4a53059SBarry Smith       /* create two temporary Index sets for build scatter gather */
1299b0a32e0cSBarry Smith       ierr = PetscMalloc((n+1)*sizeof(int),&garray);CHKERRQ(ierr);
1300f4a53059SBarry Smith       for (i=0; i<n; i++) garray[i] = dof*mpiaij->garray[i];
1301f4a53059SBarry Smith       ierr = ISCreateBlock(A->comm,dof,n,garray,&from);CHKERRQ(ierr);
1302f4a53059SBarry Smith       ierr = PetscFree(garray);CHKERRQ(ierr);
1303f4a53059SBarry Smith       ierr = ISCreateStride(PETSC_COMM_SELF,n*dof,0,1,&to);CHKERRQ(ierr);
1304f4a53059SBarry Smith 
1305f4a53059SBarry Smith       /* create temporary global vector to generate scatter context */
1306f4a53059SBarry Smith       ierr = VecCreateMPI(A->comm,dof*A->n,dof*A->N,&gvec);CHKERRQ(ierr);
1307f4a53059SBarry Smith 
1308f4a53059SBarry Smith       /* generate the scatter context */
1309f4a53059SBarry Smith       ierr = VecScatterCreate(gvec,from,b->w,to,&b->ctx);CHKERRQ(ierr);
1310f4a53059SBarry Smith 
1311f4a53059SBarry Smith       ierr = ISDestroy(from);CHKERRQ(ierr);
1312f4a53059SBarry Smith       ierr = ISDestroy(to);CHKERRQ(ierr);
1313f4a53059SBarry Smith       ierr = VecDestroy(gvec);CHKERRQ(ierr);
1314f4a53059SBarry Smith 
1315f4a53059SBarry Smith       B->ops->mult             = MatMult_MPIMAIJ_dof;
13164cb9d9b8SBarry Smith       B->ops->multtranspose    = MatMultTranspose_MPIMAIJ_dof;
13174cb9d9b8SBarry Smith       B->ops->multadd          = MatMultAdd_MPIMAIJ_dof;
13184cb9d9b8SBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof;
1319f4a53059SBarry Smith     }
1320cd3bbe55SBarry Smith     *maij = B;
1321d72c5c08SBarry Smith   }
132282b1193eSBarry Smith   PetscFunctionReturn(0);
132382b1193eSBarry Smith }
132482b1193eSBarry Smith 
132582b1193eSBarry Smith 
132682b1193eSBarry Smith 
132782b1193eSBarry Smith 
132882b1193eSBarry Smith 
132982b1193eSBarry Smith 
133082b1193eSBarry Smith 
133182b1193eSBarry Smith 
133282b1193eSBarry Smith 
133382b1193eSBarry Smith 
133482b1193eSBarry Smith 
133582b1193eSBarry Smith 
1336