xref: /petsc/src/mat/impls/maij/maij.c (revision f9fae5ad43ef59e3b226cdede57a43592e8b35bd)
1*f9fae5adSBarry Smith /*$Id: maij.c,v 1.3 2000/06/01 19:56:50 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*)
1582b1193eSBarry Smith */
1682b1193eSBarry Smith 
1782b1193eSBarry Smith #include "src/mat/impls/aij/seq/aij.h"
1882b1193eSBarry Smith 
19cd3bbe55SBarry Smith typedef struct {
20cd3bbe55SBarry Smith   int dof;               /* number of components */
21cd3bbe55SBarry Smith   Mat AIJ;               /* representation of interpolation for one component */
22cd3bbe55SBarry Smith } Mat_SeqMAIJ;
2382b1193eSBarry Smith 
2482b1193eSBarry Smith #undef __FUNC__
25cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatDestroy_SeqMAIJ"></a>*/"MatDestroy_SeqMAIJ"
26cd3bbe55SBarry Smith int MatDestroy_SeqMAIJ(Mat A)
2782b1193eSBarry Smith {
2882b1193eSBarry Smith   int         ierr;
29cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
3082b1193eSBarry Smith 
3182b1193eSBarry Smith   PetscFunctionBegin;
32cd3bbe55SBarry Smith   if (b->AIJ) {
33cd3bbe55SBarry Smith     ierr = MatDestroy(b->AIJ);CHKERRQ(ierr);
3482b1193eSBarry Smith   }
35cd3bbe55SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
3682b1193eSBarry Smith   PetscHeaderDestroy(A);
3782b1193eSBarry Smith   PetscFunctionReturn(0);
3882b1193eSBarry Smith }
3982b1193eSBarry Smith 
4082b1193eSBarry Smith EXTERN_C_BEGIN
4182b1193eSBarry Smith #undef __FUNC__
42cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatCreate_SeqMAIJ"></a>*/"MatCreate_SeqMAIJ"
43cd3bbe55SBarry Smith int MatCreate_SeqMAIJ(Mat A)
4482b1193eSBarry Smith {
45cd3bbe55SBarry Smith   int         ierr;
46cd3bbe55SBarry Smith   Mat_SeqMAIJ *b;
4782b1193eSBarry Smith 
4882b1193eSBarry Smith   PetscFunctionBegin;
49cd3bbe55SBarry Smith   A->data             = (void*)(b = PetscNew(Mat_SeqMAIJ));CHKPTRQ(b);
50cd3bbe55SBarry Smith   ierr = PetscMemzero(b,sizeof(Mat_SeqMAIJ));CHKERRQ(ierr);
51cd3bbe55SBarry Smith   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
52cd3bbe55SBarry Smith   A->factor           = 0;
53cd3bbe55SBarry Smith   A->mapping          = 0;
54cd3bbe55SBarry Smith   b->AIJ = 0;
55cd3bbe55SBarry Smith   b->dof = 0;
5682b1193eSBarry Smith   PetscFunctionReturn(0);
5782b1193eSBarry Smith }
5882b1193eSBarry Smith EXTERN_C_END
5982b1193eSBarry Smith 
60bcc973b7SBarry Smith /* ----------------------------------------------------------------------------------*/
61cd3bbe55SBarry Smith EXTERN int MatMult_SeqAIJ(Mat,Vec,Vec);
62cd3bbe55SBarry Smith EXTERN int MatMultTranspose_SeqAIJ(Mat,Vec,Vec);
63cd3bbe55SBarry Smith EXTERN int MatMultTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
64cd3bbe55SBarry Smith EXTERN int matmulttransposeadd_seqaijMatMultAdd_SeqAIJ(Mat,Vec,Vec,Vec);
65cd3bbe55SBarry Smith 
6682b1193eSBarry Smith #undef __FUNC__
67cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_1"></a>*/"MatMult_SeqMAIJ_1"
68cd3bbe55SBarry Smith int MatMult_SeqMAIJ_1(Mat A,Vec xx,Vec yy)
6982b1193eSBarry Smith {
70cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
71cd3bbe55SBarry Smith   int         ierr;
7282b1193eSBarry Smith   PetscFunctionBegin;
73cd3bbe55SBarry Smith   ierr = MatMult_SeqAIJ(b->AIJ,xx,yy);
74cd3bbe55SBarry Smith   PetscFunctionReturn(0);
7582b1193eSBarry Smith }
76cd3bbe55SBarry Smith #undef __FUNC__
77cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_1"></a>*/"MatMultTranspose_SeqMAIJ_1"
78cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_1(Mat A,Vec xx,Vec yy)
79cd3bbe55SBarry Smith {
80cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
81cd3bbe55SBarry Smith   int         ierr;
82cd3bbe55SBarry Smith   PetscFunctionBegin;
83cd3bbe55SBarry Smith   ierr = MatMultTranspose_SeqAIJ(b->AIJ,xx,yy);
84cd3bbe55SBarry Smith   PetscFunctionReturn(0);
85cd3bbe55SBarry Smith }
86cd3bbe55SBarry Smith #undef __FUNC__
87cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_1"></a>*/"MatMultAdd_SeqMAIJ_1"
88cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
89cd3bbe55SBarry Smith {
90cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
91cd3bbe55SBarry Smith   int         ierr;
92cd3bbe55SBarry Smith   PetscFunctionBegin;
93cd3bbe55SBarry Smith   ierr = MatMultAdd_SeqAIJ(b->AIJ,xx,yy,zz);
94cd3bbe55SBarry Smith   PetscFunctionReturn(0);
95cd3bbe55SBarry Smith }
96cd3bbe55SBarry Smith #undef __FUNC__
97cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_1"></a>*/"MatMultTransposeAdd_SeqMAIJ_1"
98cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
99cd3bbe55SBarry Smith {
100cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
101cd3bbe55SBarry Smith   int         ierr;
102cd3bbe55SBarry Smith   PetscFunctionBegin;
103cd3bbe55SBarry Smith   ierr = MatMultTransposeAdd_SeqAIJ(b->AIJ,xx,yy,zz);
10482b1193eSBarry Smith   PetscFunctionReturn(0);
10582b1193eSBarry Smith }
10682b1193eSBarry Smith 
107cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/
10882b1193eSBarry Smith #undef __FUNC__
109cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_2"></a>*/"MatMult_SeqMAIJ_2"
110cd3bbe55SBarry Smith int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
11182b1193eSBarry Smith {
112cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
113bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
114bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2;
115bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
116bcc973b7SBarry Smith   int         n,i,jrow,j;
11782b1193eSBarry Smith 
118bcc973b7SBarry Smith   PetscFunctionBegin;
119bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
120bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
121bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
122bcc973b7SBarry Smith   idx  = a->j;
123bcc973b7SBarry Smith   v    = a->a;
124bcc973b7SBarry Smith   ii   = a->i;
125bcc973b7SBarry Smith 
126bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
127bcc973b7SBarry Smith   idx  += shift;
128bcc973b7SBarry Smith   for (i=0; i<m; i++) {
129bcc973b7SBarry Smith     jrow = ii[i];
130bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
131bcc973b7SBarry Smith     sum1  = 0.0;
132bcc973b7SBarry Smith     sum2  = 0.0;
133bcc973b7SBarry Smith     for (j=0; j<n; j++) {
134bcc973b7SBarry Smith       sum1 += v[jrow]*x[2*idx[jrow]];
135bcc973b7SBarry Smith       sum2 += v[jrow]*x[2*idx[jrow]+1];
136bcc973b7SBarry Smith       jrow++;
137bcc973b7SBarry Smith      }
138bcc973b7SBarry Smith     y[2*i]   = sum1;
139bcc973b7SBarry Smith     y[2*i+1] = sum2;
140bcc973b7SBarry Smith   }
141bcc973b7SBarry Smith 
142bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*m);
143bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
144bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
14582b1193eSBarry Smith   PetscFunctionReturn(0);
14682b1193eSBarry Smith }
147bcc973b7SBarry Smith 
14882b1193eSBarry Smith #undef __FUNC__
149cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_2"></a>*/"MatMultTranspose_SeqMAIJ_2"
150cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
15182b1193eSBarry Smith {
152cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
153bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
154bcc973b7SBarry Smith   Scalar     *x,*y,*v,alpha1,alpha2,zero = 0.0;
155bcc973b7SBarry Smith   int        ierr,m = a->m,n,i,*idx,shift = a->indexshift;
15682b1193eSBarry Smith 
157bcc973b7SBarry Smith   PetscFunctionBegin;
158bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
159bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
160bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
161bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
162bcc973b7SBarry Smith   for (i=0; i<m; i++) {
163bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
164bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
165bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
166bcc973b7SBarry Smith     alpha1 = x[2*i];
167bcc973b7SBarry Smith     alpha2 = x[2*i+1];
168bcc973b7SBarry Smith     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
169bcc973b7SBarry Smith   }
170bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*a->n);
171bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
172bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
17382b1193eSBarry Smith   PetscFunctionReturn(0);
17482b1193eSBarry Smith }
175bcc973b7SBarry Smith 
17682b1193eSBarry Smith #undef __FUNC__
177cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_2"></a>*/"MatMultAdd_SeqMAIJ_2"
178cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
17982b1193eSBarry Smith {
180cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
181bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
182bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2;
183bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
184bcc973b7SBarry Smith   int         n,i,jrow,j;
18582b1193eSBarry Smith 
186bcc973b7SBarry Smith   PetscFunctionBegin;
187*f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
188bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
189*f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
190bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
191bcc973b7SBarry Smith   idx  = a->j;
192bcc973b7SBarry Smith   v    = a->a;
193bcc973b7SBarry Smith   ii   = a->i;
194bcc973b7SBarry Smith 
195bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
196bcc973b7SBarry Smith   idx  += shift;
197bcc973b7SBarry Smith   for (i=0; i<m; i++) {
198bcc973b7SBarry Smith     jrow = ii[i];
199bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
200bcc973b7SBarry Smith     sum1  = 0.0;
201bcc973b7SBarry Smith     sum2  = 0.0;
202bcc973b7SBarry Smith     for (j=0; j<n; j++) {
203bcc973b7SBarry Smith       sum1 += v[jrow]*x[2*idx[jrow]];
204bcc973b7SBarry Smith       sum2 += v[jrow]*x[2*idx[jrow]+1];
205bcc973b7SBarry Smith       jrow++;
206bcc973b7SBarry Smith      }
207bcc973b7SBarry Smith     y[2*i]   += sum1;
208bcc973b7SBarry Smith     y[2*i+1] += sum2;
209bcc973b7SBarry Smith   }
210bcc973b7SBarry Smith 
211bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*m);
212bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
213*f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
214bcc973b7SBarry Smith   PetscFunctionReturn(0);
21582b1193eSBarry Smith   PetscFunctionReturn(0);
21682b1193eSBarry Smith }
21782b1193eSBarry Smith #undef __FUNC__
218cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_2"></a>*/"MatMultTransposeAdd_SeqMAIJ_2"
219cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
22082b1193eSBarry Smith {
221cd3bbe55SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
222bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
223bcc973b7SBarry Smith   Scalar     *x,*y,*v,alpha1,alpha2;
224bcc973b7SBarry Smith   int        ierr,m = a->m,n,i,*idx,shift = a->indexshift;
22582b1193eSBarry Smith 
226bcc973b7SBarry Smith   PetscFunctionBegin;
227*f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
228bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
229*f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
230bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
231bcc973b7SBarry Smith   for (i=0; i<m; i++) {
232bcc973b7SBarry Smith     idx   = a->j + a->i[i] + shift;
233bcc973b7SBarry Smith     v     = a->a + a->i[i] + shift;
234bcc973b7SBarry Smith     n     = a->i[i+1] - a->i[i];
235bcc973b7SBarry Smith     alpha1 = x[2*i];
236bcc973b7SBarry Smith     alpha2 = x[2*i+1];
237bcc973b7SBarry Smith     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
238bcc973b7SBarry Smith   }
239bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*a->n);
240bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
241*f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
242bcc973b7SBarry Smith   PetscFunctionReturn(0);
24382b1193eSBarry Smith   PetscFunctionReturn(0);
24482b1193eSBarry Smith }
245cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/
246bcc973b7SBarry Smith #undef __FUNC__
247bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_3"></a>*/"MatMult_SeqMAIJ_3"
248bcc973b7SBarry Smith int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
249bcc973b7SBarry Smith {
250bcc973b7SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
251bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
252bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3;
253bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
254bcc973b7SBarry Smith   int         n,i,jrow,j;
25582b1193eSBarry Smith 
256bcc973b7SBarry Smith   PetscFunctionBegin;
257bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
258bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
259bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
260bcc973b7SBarry Smith   idx  = a->j;
261bcc973b7SBarry Smith   v    = a->a;
262bcc973b7SBarry Smith   ii   = a->i;
263bcc973b7SBarry Smith 
264bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
265bcc973b7SBarry Smith   idx  += shift;
266bcc973b7SBarry Smith   for (i=0; i<m; i++) {
267bcc973b7SBarry Smith     jrow = ii[i];
268bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
269bcc973b7SBarry Smith     sum1  = 0.0;
270bcc973b7SBarry Smith     sum2  = 0.0;
271bcc973b7SBarry Smith     sum3  = 0.0;
272bcc973b7SBarry Smith     for (j=0; j<n; j++) {
273bcc973b7SBarry Smith       sum1 += v[jrow]*x[3*idx[jrow]];
274bcc973b7SBarry Smith       sum2 += v[jrow]*x[3*idx[jrow]+1];
275bcc973b7SBarry Smith       sum3 += v[jrow]*x[3*idx[jrow]+2];
276bcc973b7SBarry Smith       jrow++;
277bcc973b7SBarry Smith      }
278bcc973b7SBarry Smith     y[3*i]   = sum1;
279bcc973b7SBarry Smith     y[3*i+1] = sum2;
280bcc973b7SBarry Smith     y[3*i+2] = sum3;
281bcc973b7SBarry Smith   }
282bcc973b7SBarry Smith 
283bcc973b7SBarry Smith   PLogFlops(6*a->nz - 3*m);
284bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
285bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
286bcc973b7SBarry Smith   PetscFunctionReturn(0);
287bcc973b7SBarry Smith }
288bcc973b7SBarry Smith 
289bcc973b7SBarry Smith #undef __FUNC__
290bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_3"></a>*/"MatMultTranspose_SeqMAIJ_3"
291bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
292bcc973b7SBarry Smith {
293bcc973b7SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
294bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
295bcc973b7SBarry Smith   Scalar     *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0;
296bcc973b7SBarry Smith   int        ierr,m = a->m,n,i,*idx,shift = a->indexshift;
297bcc973b7SBarry Smith 
298bcc973b7SBarry Smith   PetscFunctionBegin;
299bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
300bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
301bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
302bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
303bcc973b7SBarry Smith   for (i=0; i<m; i++) {
304bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
305bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
306bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
307bcc973b7SBarry Smith     alpha1 = x[3*i];
308bcc973b7SBarry Smith     alpha2 = x[3*i+1];
309bcc973b7SBarry Smith     alpha3 = x[3*i+2];
310bcc973b7SBarry Smith     while (n-->0) {
311bcc973b7SBarry Smith       y[3*(*idx)]   += alpha1*(*v);
312bcc973b7SBarry Smith       y[3*(*idx)+1] += alpha2*(*v);
313bcc973b7SBarry Smith       y[3*(*idx)+2] += alpha3*(*v);
314bcc973b7SBarry Smith       idx++; v++;
315bcc973b7SBarry Smith     }
316bcc973b7SBarry Smith   }
317bcc973b7SBarry Smith   PLogFlops(6*a->nz - 3*a->n);
318bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
319bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
320bcc973b7SBarry Smith   PetscFunctionReturn(0);
321bcc973b7SBarry Smith }
322bcc973b7SBarry Smith 
323bcc973b7SBarry Smith #undef __FUNC__
324bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_3"></a>*/"MatMultAdd_SeqMAIJ_3"
325bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
326bcc973b7SBarry Smith {
327bcc973b7SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
328bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
329bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3;
330bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
331bcc973b7SBarry Smith   int         n,i,jrow,j;
332bcc973b7SBarry Smith 
333bcc973b7SBarry Smith   PetscFunctionBegin;
334*f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
335bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
336*f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
337bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
338bcc973b7SBarry Smith   idx  = a->j;
339bcc973b7SBarry Smith   v    = a->a;
340bcc973b7SBarry Smith   ii   = a->i;
341bcc973b7SBarry Smith 
342bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
343bcc973b7SBarry Smith   idx  += shift;
344bcc973b7SBarry Smith   for (i=0; i<m; i++) {
345bcc973b7SBarry Smith     jrow = ii[i];
346bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
347bcc973b7SBarry Smith     sum1  = 0.0;
348bcc973b7SBarry Smith     sum2  = 0.0;
349bcc973b7SBarry Smith     sum3  = 0.0;
350bcc973b7SBarry Smith     for (j=0; j<n; j++) {
351bcc973b7SBarry Smith       sum1 += v[jrow]*x[3*idx[jrow]];
352bcc973b7SBarry Smith       sum2 += v[jrow]*x[3*idx[jrow]+1];
353bcc973b7SBarry Smith       sum3 += v[jrow]*x[3*idx[jrow]+2];
354bcc973b7SBarry Smith       jrow++;
355bcc973b7SBarry Smith      }
356bcc973b7SBarry Smith     y[3*i]   += sum1;
357bcc973b7SBarry Smith     y[3*i+1] += sum2;
358bcc973b7SBarry Smith     y[3*i+2] += sum3;
359bcc973b7SBarry Smith   }
360bcc973b7SBarry Smith 
361bcc973b7SBarry Smith   PLogFlops(6*a->nz);
362bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
363*f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
364bcc973b7SBarry Smith   PetscFunctionReturn(0);
365bcc973b7SBarry Smith }
366bcc973b7SBarry Smith #undef __FUNC__
367bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_3"></a>*/"MatMultTransposeAdd_SeqMAIJ_3"
368bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
369bcc973b7SBarry Smith {
370bcc973b7SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
371bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
372bcc973b7SBarry Smith   Scalar     *x,*y,*v,alpha1,alpha2,alpha3;
373bcc973b7SBarry Smith   int        ierr,m = a->m,n,i,*idx,shift = a->indexshift;
374bcc973b7SBarry Smith 
375bcc973b7SBarry Smith   PetscFunctionBegin;
376*f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
377bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
378*f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
379bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
380bcc973b7SBarry Smith   for (i=0; i<m; i++) {
381bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
382bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
383bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
384bcc973b7SBarry Smith     alpha1 = x[3*i];
385bcc973b7SBarry Smith     alpha2 = x[3*i+1];
386bcc973b7SBarry Smith     alpha3 = x[3*i+2];
387bcc973b7SBarry Smith     while (n-->0) {
388bcc973b7SBarry Smith       y[3*(*idx)]   += alpha1*(*v);
389bcc973b7SBarry Smith       y[3*(*idx)+1] += alpha2*(*v);
390bcc973b7SBarry Smith       y[3*(*idx)+2] += alpha3*(*v);
391bcc973b7SBarry Smith       idx++; v++;
392bcc973b7SBarry Smith     }
393bcc973b7SBarry Smith   }
394bcc973b7SBarry Smith   PLogFlops(6*a->nz);
395bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
396*f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
397bcc973b7SBarry Smith   PetscFunctionReturn(0);
398bcc973b7SBarry Smith }
399bcc973b7SBarry Smith 
400bcc973b7SBarry Smith /* ------------------------------------------------------------------------------*/
401bcc973b7SBarry Smith #undef __FUNC__
402bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_4"></a>*/"MatMult_SeqMAIJ_4"
403bcc973b7SBarry Smith int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
404bcc973b7SBarry Smith {
405bcc973b7SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
406bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
407bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4;
408bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
409bcc973b7SBarry Smith   int         n,i,jrow,j;
410bcc973b7SBarry Smith 
411bcc973b7SBarry Smith   PetscFunctionBegin;
412bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
413bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
414bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
415bcc973b7SBarry Smith   idx  = a->j;
416bcc973b7SBarry Smith   v    = a->a;
417bcc973b7SBarry Smith   ii   = a->i;
418bcc973b7SBarry Smith 
419bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
420bcc973b7SBarry Smith   idx  += shift;
421bcc973b7SBarry Smith   for (i=0; i<m; i++) {
422bcc973b7SBarry Smith     jrow = ii[i];
423bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
424bcc973b7SBarry Smith     sum1  = 0.0;
425bcc973b7SBarry Smith     sum2  = 0.0;
426bcc973b7SBarry Smith     sum3  = 0.0;
427bcc973b7SBarry Smith     sum4  = 0.0;
428bcc973b7SBarry Smith     for (j=0; j<n; j++) {
429bcc973b7SBarry Smith       sum1 += v[jrow]*x[4*idx[jrow]];
430bcc973b7SBarry Smith       sum2 += v[jrow]*x[4*idx[jrow]+1];
431bcc973b7SBarry Smith       sum3 += v[jrow]*x[4*idx[jrow]+2];
432bcc973b7SBarry Smith       sum4 += v[jrow]*x[4*idx[jrow]+3];
433bcc973b7SBarry Smith       jrow++;
434bcc973b7SBarry Smith      }
435bcc973b7SBarry Smith     y[4*i]   = sum1;
436bcc973b7SBarry Smith     y[4*i+1] = sum2;
437bcc973b7SBarry Smith     y[4*i+2] = sum3;
438bcc973b7SBarry Smith     y[4*i+3] = sum4;
439bcc973b7SBarry Smith   }
440bcc973b7SBarry Smith 
441bcc973b7SBarry Smith   PLogFlops(8*a->nz - 4*m);
442bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
443bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
444bcc973b7SBarry Smith   PetscFunctionReturn(0);
445bcc973b7SBarry Smith }
446bcc973b7SBarry Smith 
447bcc973b7SBarry Smith #undef __FUNC__
448bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_4"></a>*/"MatMultTranspose_SeqMAIJ_4"
449bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
450bcc973b7SBarry Smith {
451bcc973b7SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
452bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
453bcc973b7SBarry Smith   Scalar     *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0;
454bcc973b7SBarry Smith   int        ierr,m = a->m,n,i,*idx,shift = a->indexshift;
455bcc973b7SBarry Smith 
456bcc973b7SBarry Smith   PetscFunctionBegin;
457bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
458bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
459bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
460bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
461bcc973b7SBarry Smith   for (i=0; i<m; i++) {
462bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
463bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
464bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
465bcc973b7SBarry Smith     alpha1 = x[4*i];
466bcc973b7SBarry Smith     alpha2 = x[4*i+1];
467bcc973b7SBarry Smith     alpha3 = x[4*i+2];
468bcc973b7SBarry Smith     alpha4 = x[4*i+3];
469bcc973b7SBarry Smith     while (n-->0) {
470bcc973b7SBarry Smith       y[4*(*idx)]   += alpha1*(*v);
471bcc973b7SBarry Smith       y[4*(*idx)+1] += alpha2*(*v);
472bcc973b7SBarry Smith       y[4*(*idx)+2] += alpha3*(*v);
473bcc973b7SBarry Smith       y[4*(*idx)+3] += alpha4*(*v);
474bcc973b7SBarry Smith       idx++; v++;
475bcc973b7SBarry Smith     }
476bcc973b7SBarry Smith   }
477bcc973b7SBarry Smith   PLogFlops(8*a->nz - 4*a->n);
478bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
479bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
480bcc973b7SBarry Smith   PetscFunctionReturn(0);
481bcc973b7SBarry Smith }
482bcc973b7SBarry Smith 
483bcc973b7SBarry Smith #undef __FUNC__
484bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_4"></a>*/"MatMultAdd_SeqMAIJ_4"
485bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
486bcc973b7SBarry Smith {
487*f9fae5adSBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
488*f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
489*f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4;
490*f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
491*f9fae5adSBarry Smith   int         n,i,jrow,j;
492*f9fae5adSBarry Smith 
493*f9fae5adSBarry Smith   PetscFunctionBegin;
494*f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
495*f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
496*f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
497*f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
498*f9fae5adSBarry Smith   idx  = a->j;
499*f9fae5adSBarry Smith   v    = a->a;
500*f9fae5adSBarry Smith   ii   = a->i;
501*f9fae5adSBarry Smith 
502*f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
503*f9fae5adSBarry Smith   idx  += shift;
504*f9fae5adSBarry Smith   for (i=0; i<m; i++) {
505*f9fae5adSBarry Smith     jrow = ii[i];
506*f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
507*f9fae5adSBarry Smith     sum1  = 0.0;
508*f9fae5adSBarry Smith     sum2  = 0.0;
509*f9fae5adSBarry Smith     sum3  = 0.0;
510*f9fae5adSBarry Smith     sum4  = 0.0;
511*f9fae5adSBarry Smith     for (j=0; j<n; j++) {
512*f9fae5adSBarry Smith       sum1 += v[jrow]*x[4*idx[jrow]];
513*f9fae5adSBarry Smith       sum2 += v[jrow]*x[4*idx[jrow]+1];
514*f9fae5adSBarry Smith       sum3 += v[jrow]*x[4*idx[jrow]+2];
515*f9fae5adSBarry Smith       sum4 += v[jrow]*x[4*idx[jrow]+3];
516*f9fae5adSBarry Smith       jrow++;
517*f9fae5adSBarry Smith      }
518*f9fae5adSBarry Smith     y[4*i]   += sum1;
519*f9fae5adSBarry Smith     y[4*i+1] += sum2;
520*f9fae5adSBarry Smith     y[4*i+2] += sum3;
521*f9fae5adSBarry Smith     y[4*i+3] += sum4;
522*f9fae5adSBarry Smith   }
523*f9fae5adSBarry Smith 
524*f9fae5adSBarry Smith   PLogFlops(8*a->nz - 4*m);
525*f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
526*f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
527*f9fae5adSBarry Smith   PetscFunctionReturn(0);
528bcc973b7SBarry Smith }
529bcc973b7SBarry Smith #undef __FUNC__
530bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_4"></a>*/"MatMultTransposeAdd_SeqMAIJ_4"
531bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
532bcc973b7SBarry Smith {
533*f9fae5adSBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
534*f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
535*f9fae5adSBarry Smith   Scalar     *x,*y,*v,alpha1,alpha2,alpha3,alpha4;
536*f9fae5adSBarry Smith   int        ierr,m = a->m,n,i,*idx,shift = a->indexshift;
537*f9fae5adSBarry Smith 
538*f9fae5adSBarry Smith   PetscFunctionBegin;
539*f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
540*f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
541*f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
542*f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
543*f9fae5adSBarry Smith   for (i=0; i<m; i++) {
544*f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
545*f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
546*f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
547*f9fae5adSBarry Smith     alpha1 = x[4*i];
548*f9fae5adSBarry Smith     alpha2 = x[4*i+1];
549*f9fae5adSBarry Smith     alpha3 = x[4*i+2];
550*f9fae5adSBarry Smith     alpha4 = x[4*i+3];
551*f9fae5adSBarry Smith     while (n-->0) {
552*f9fae5adSBarry Smith       y[4*(*idx)]   += alpha1*(*v);
553*f9fae5adSBarry Smith       y[4*(*idx)+1] += alpha2*(*v);
554*f9fae5adSBarry Smith       y[4*(*idx)+2] += alpha3*(*v);
555*f9fae5adSBarry Smith       y[4*(*idx)+3] += alpha4*(*v);
556*f9fae5adSBarry Smith       idx++; v++;
557*f9fae5adSBarry Smith     }
558*f9fae5adSBarry Smith   }
559*f9fae5adSBarry Smith   PLogFlops(8*a->nz - 4*a->n);
560*f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
561*f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
562*f9fae5adSBarry Smith   PetscFunctionReturn(0);
563*f9fae5adSBarry Smith 
564*f9fae5adSBarry Smith }
565*f9fae5adSBarry Smith 
566*f9fae5adSBarry Smith /* ------------------------------------------------------------------------------*/
567*f9fae5adSBarry Smith #undef __FUNC__
568*f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_5"></a>*/"MatMult_SeqMAIJ_5"
569*f9fae5adSBarry Smith int MatMult_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
570*f9fae5adSBarry Smith {
571*f9fae5adSBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
572*f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
573*f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
574*f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
575*f9fae5adSBarry Smith   int         n,i,jrow,j;
576*f9fae5adSBarry Smith 
577*f9fae5adSBarry Smith   PetscFunctionBegin;
578*f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
579*f9fae5adSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
580*f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
581*f9fae5adSBarry Smith   idx  = a->j;
582*f9fae5adSBarry Smith   v    = a->a;
583*f9fae5adSBarry Smith   ii   = a->i;
584*f9fae5adSBarry Smith 
585*f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
586*f9fae5adSBarry Smith   idx  += shift;
587*f9fae5adSBarry Smith   for (i=0; i<m; i++) {
588*f9fae5adSBarry Smith     jrow = ii[i];
589*f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
590*f9fae5adSBarry Smith     sum1  = 0.0;
591*f9fae5adSBarry Smith     sum2  = 0.0;
592*f9fae5adSBarry Smith     sum3  = 0.0;
593*f9fae5adSBarry Smith     sum4  = 0.0;
594*f9fae5adSBarry Smith     sum5  = 0.0;
595*f9fae5adSBarry Smith     for (j=0; j<n; j++) {
596*f9fae5adSBarry Smith       sum1 += v[jrow]*x[5*idx[jrow]];
597*f9fae5adSBarry Smith       sum2 += v[jrow]*x[5*idx[jrow]+1];
598*f9fae5adSBarry Smith       sum3 += v[jrow]*x[5*idx[jrow]+2];
599*f9fae5adSBarry Smith       sum4 += v[jrow]*x[5*idx[jrow]+3];
600*f9fae5adSBarry Smith       sum5 += v[jrow]*x[5*idx[jrow]+4];
601*f9fae5adSBarry Smith       jrow++;
602*f9fae5adSBarry Smith      }
603*f9fae5adSBarry Smith     y[5*i]   = sum1;
604*f9fae5adSBarry Smith     y[5*i+1] = sum2;
605*f9fae5adSBarry Smith     y[5*i+2] = sum3;
606*f9fae5adSBarry Smith     y[5*i+3] = sum4;
607*f9fae5adSBarry Smith     y[5*i+4] = sum5;
608*f9fae5adSBarry Smith   }
609*f9fae5adSBarry Smith 
610*f9fae5adSBarry Smith   PLogFlops(10*a->nz - 5*m);
611*f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
612*f9fae5adSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
613*f9fae5adSBarry Smith   PetscFunctionReturn(0);
614*f9fae5adSBarry Smith }
615*f9fae5adSBarry Smith 
616*f9fae5adSBarry Smith #undef __FUNC__
617*f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_5"></a>*/"MatMultTranspose_SeqMAIJ_5"
618*f9fae5adSBarry Smith int MatMultTranspose_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
619*f9fae5adSBarry Smith {
620*f9fae5adSBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
621*f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
622*f9fae5adSBarry Smith   Scalar     *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,zero = 0.0;
623*f9fae5adSBarry Smith   int        ierr,m = a->m,n,i,*idx,shift = a->indexshift;
624*f9fae5adSBarry Smith 
625*f9fae5adSBarry Smith   PetscFunctionBegin;
626*f9fae5adSBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
627*f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
628*f9fae5adSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
629*f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
630*f9fae5adSBarry Smith   for (i=0; i<m; i++) {
631*f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
632*f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
633*f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
634*f9fae5adSBarry Smith     alpha1 = x[5*i];
635*f9fae5adSBarry Smith     alpha2 = x[5*i+1];
636*f9fae5adSBarry Smith     alpha3 = x[5*i+2];
637*f9fae5adSBarry Smith     alpha4 = x[5*i+3];
638*f9fae5adSBarry Smith     alpha5 = x[5*i+4];
639*f9fae5adSBarry Smith     while (n-->0) {
640*f9fae5adSBarry Smith       y[5*(*idx)]   += alpha1*(*v);
641*f9fae5adSBarry Smith       y[5*(*idx)+1] += alpha2*(*v);
642*f9fae5adSBarry Smith       y[5*(*idx)+2] += alpha3*(*v);
643*f9fae5adSBarry Smith       y[5*(*idx)+3] += alpha4*(*v);
644*f9fae5adSBarry Smith       y[5*(*idx)+4] += alpha5*(*v);
645*f9fae5adSBarry Smith       idx++; v++;
646*f9fae5adSBarry Smith     }
647*f9fae5adSBarry Smith   }
648*f9fae5adSBarry Smith   PLogFlops(10*a->nz - 5*a->n);
649*f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
650*f9fae5adSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
651*f9fae5adSBarry Smith   PetscFunctionReturn(0);
652*f9fae5adSBarry Smith }
653*f9fae5adSBarry Smith 
654*f9fae5adSBarry Smith #undef __FUNC__
655*f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_5"></a>*/"MatMultAdd_SeqMAIJ_5"
656*f9fae5adSBarry Smith int MatMultAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
657*f9fae5adSBarry Smith {
658*f9fae5adSBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
659*f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
660*f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
661*f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
662*f9fae5adSBarry Smith   int         n,i,jrow,j;
663*f9fae5adSBarry Smith 
664*f9fae5adSBarry Smith   PetscFunctionBegin;
665*f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
666*f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
667*f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
668*f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
669*f9fae5adSBarry Smith   idx  = a->j;
670*f9fae5adSBarry Smith   v    = a->a;
671*f9fae5adSBarry Smith   ii   = a->i;
672*f9fae5adSBarry Smith 
673*f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
674*f9fae5adSBarry Smith   idx  += shift;
675*f9fae5adSBarry Smith   for (i=0; i<m; i++) {
676*f9fae5adSBarry Smith     jrow = ii[i];
677*f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
678*f9fae5adSBarry Smith     sum1  = 0.0;
679*f9fae5adSBarry Smith     sum2  = 0.0;
680*f9fae5adSBarry Smith     sum3  = 0.0;
681*f9fae5adSBarry Smith     sum4  = 0.0;
682*f9fae5adSBarry Smith     sum5  = 0.0;
683*f9fae5adSBarry Smith     for (j=0; j<n; j++) {
684*f9fae5adSBarry Smith       sum1 += v[jrow]*x[5*idx[jrow]];
685*f9fae5adSBarry Smith       sum2 += v[jrow]*x[5*idx[jrow]+1];
686*f9fae5adSBarry Smith       sum3 += v[jrow]*x[5*idx[jrow]+2];
687*f9fae5adSBarry Smith       sum4 += v[jrow]*x[5*idx[jrow]+3];
688*f9fae5adSBarry Smith       sum5 += v[jrow]*x[5*idx[jrow]+4];
689*f9fae5adSBarry Smith       jrow++;
690*f9fae5adSBarry Smith      }
691*f9fae5adSBarry Smith     y[5*i]   += sum1;
692*f9fae5adSBarry Smith     y[5*i+1] += sum2;
693*f9fae5adSBarry Smith     y[5*i+2] += sum3;
694*f9fae5adSBarry Smith     y[5*i+3] += sum4;
695*f9fae5adSBarry Smith     y[5*i+4] += sum5;
696*f9fae5adSBarry Smith   }
697*f9fae5adSBarry Smith 
698*f9fae5adSBarry Smith   PLogFlops(10*a->nz);
699*f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
700*f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
701*f9fae5adSBarry Smith   PetscFunctionReturn(0);
702*f9fae5adSBarry Smith }
703*f9fae5adSBarry Smith 
704*f9fae5adSBarry Smith #undef __FUNC__
705*f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_5"></a>*/"MatMultTransposeAdd_SeqMAIJ_5"
706*f9fae5adSBarry Smith int MatMultTransposeAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
707*f9fae5adSBarry Smith {
708*f9fae5adSBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
709*f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
710*f9fae5adSBarry Smith   Scalar     *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5;
711*f9fae5adSBarry Smith   int        ierr,m = a->m,n,i,*idx,shift = a->indexshift;
712*f9fae5adSBarry Smith 
713*f9fae5adSBarry Smith   PetscFunctionBegin;
714*f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
715*f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
716*f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
717*f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
718*f9fae5adSBarry Smith   for (i=0; i<m; i++) {
719*f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
720*f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
721*f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
722*f9fae5adSBarry Smith     alpha1 = x[5*i];
723*f9fae5adSBarry Smith     alpha2 = x[5*i+1];
724*f9fae5adSBarry Smith     alpha3 = x[5*i+2];
725*f9fae5adSBarry Smith     alpha4 = x[5*i+3];
726*f9fae5adSBarry Smith     alpha5 = x[5*i+4];
727*f9fae5adSBarry Smith     while (n-->0) {
728*f9fae5adSBarry Smith       y[5*(*idx)]   += alpha1*(*v);
729*f9fae5adSBarry Smith       y[5*(*idx)+1] += alpha2*(*v);
730*f9fae5adSBarry Smith       y[5*(*idx)+2] += alpha3*(*v);
731*f9fae5adSBarry Smith       y[5*(*idx)+3] += alpha4*(*v);
732*f9fae5adSBarry Smith       y[5*(*idx)+4] += alpha5*(*v);
733*f9fae5adSBarry Smith       idx++; v++;
734*f9fae5adSBarry Smith     }
735*f9fae5adSBarry Smith   }
736*f9fae5adSBarry Smith   PLogFlops(10*a->nz);
737*f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
738*f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
739*f9fae5adSBarry Smith   PetscFunctionReturn(0);
740bcc973b7SBarry Smith }
741bcc973b7SBarry Smith 
742bcc973b7SBarry Smith /* ---------------------------------------------------------------------------------- */
74382b1193eSBarry Smith #undef __FUNC__
744cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatCreateMAIJ"></a>*/"MatCreateMAIJ"
745cd3bbe55SBarry Smith int MatCreateMAIJ(Mat A,int dof,Mat *maij)
74682b1193eSBarry Smith {
747cd3bbe55SBarry Smith   int         ierr;
748cd3bbe55SBarry Smith   Mat_SeqMAIJ *b;
74982b1193eSBarry Smith   Mat         B;
75082b1193eSBarry Smith 
75182b1193eSBarry Smith   PetscFunctionBegin;
752cd3bbe55SBarry Smith   ierr = MATCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr);
753cd3bbe55SBarry Smith   ierr = MatSetType(B,MATSEQMAIJ);CHKERRQ(ierr);
75482b1193eSBarry Smith 
755cd3bbe55SBarry Smith   B->assembled    = PETSC_TRUE;
756cd3bbe55SBarry Smith   B->ops->destroy = MatDestroy_SeqMAIJ;
757cd3bbe55SBarry Smith   b = (Mat_SeqMAIJ*)B->data;
75882b1193eSBarry Smith 
759cd3bbe55SBarry Smith   b->AIJ = A;
760cd3bbe55SBarry Smith   b->dof = dof;
761cd3bbe55SBarry Smith   ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
762cd3bbe55SBarry Smith   if (dof == 1) {
763cd3bbe55SBarry Smith     B->ops->mult             = MatMult_SeqMAIJ_1;
764cd3bbe55SBarry Smith     B->ops->multadd          = MatMultAdd_SeqMAIJ_1;
765cd3bbe55SBarry Smith     B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_1;
766cd3bbe55SBarry Smith     B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_1;
767cd3bbe55SBarry Smith   } else if (dof == 2) {
768cd3bbe55SBarry Smith     B->ops->mult             = MatMult_SeqMAIJ_2;
769cd3bbe55SBarry Smith     B->ops->multadd          = MatMultAdd_SeqMAIJ_2;
770cd3bbe55SBarry Smith     B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_2;
771cd3bbe55SBarry Smith     B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2;
772bcc973b7SBarry Smith   } else if (dof == 3) {
773bcc973b7SBarry Smith     B->ops->mult             = MatMult_SeqMAIJ_3;
774bcc973b7SBarry Smith     B->ops->multadd          = MatMultAdd_SeqMAIJ_3;
775bcc973b7SBarry Smith     B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_3;
776bcc973b7SBarry Smith     B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3;
777bcc973b7SBarry Smith   } else if (dof == 4) {
778bcc973b7SBarry Smith     B->ops->mult             = MatMult_SeqMAIJ_4;
779bcc973b7SBarry Smith     B->ops->multadd          = MatMultAdd_SeqMAIJ_4;
780bcc973b7SBarry Smith     B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_4;
781bcc973b7SBarry Smith     B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4;
782*f9fae5adSBarry Smith   } else if (dof == 5) {
783*f9fae5adSBarry Smith     B->ops->mult             = MatMult_SeqMAIJ_5;
784*f9fae5adSBarry Smith     B->ops->multadd          = MatMultAdd_SeqMAIJ_5;
785*f9fae5adSBarry Smith     B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_5;
786*f9fae5adSBarry Smith     B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5;
78782b1193eSBarry Smith   } else {
788cd3bbe55SBarry Smith     SETERRQ1(1,1,"Cannot handle a dof of %d\n",dof);
78982b1193eSBarry Smith   }
790cd3bbe55SBarry Smith   *maij = B;
79182b1193eSBarry Smith   PetscFunctionReturn(0);
79282b1193eSBarry Smith }
79382b1193eSBarry Smith 
79482b1193eSBarry Smith 
79582b1193eSBarry Smith 
79682b1193eSBarry Smith 
79782b1193eSBarry Smith 
79882b1193eSBarry Smith 
79982b1193eSBarry Smith 
80082b1193eSBarry Smith 
80182b1193eSBarry Smith 
80282b1193eSBarry Smith 
80382b1193eSBarry Smith 
80482b1193eSBarry Smith 
805