xref: /petsc/src/mat/impls/maij/maij.c (revision d72c5c086a9ad2ca28c7153cd0684f41320be008)
1*d72c5c08SBarry Smith /*$Id: maij.c,v 1.6 2000/06/27 17:46:45 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 
3482b1193eSBarry Smith #undef __FUNC__
354cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatDestroy_SeqMAIJ"></a>*/"MatDestroy_SeqMAIJ"
364cb9d9b8SBarry Smith int MatDestroy_SeqMAIJ(Mat A)
3782b1193eSBarry Smith {
3882b1193eSBarry Smith   int         ierr;
394cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
4082b1193eSBarry Smith 
4182b1193eSBarry Smith   PetscFunctionBegin;
42cd3bbe55SBarry Smith   if (b->AIJ) {
43cd3bbe55SBarry Smith     ierr = MatDestroy(b->AIJ);CHKERRQ(ierr);
4482b1193eSBarry Smith   }
454cb9d9b8SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
464cb9d9b8SBarry Smith   PetscHeaderDestroy(A);
474cb9d9b8SBarry Smith   PetscFunctionReturn(0);
484cb9d9b8SBarry Smith }
494cb9d9b8SBarry Smith 
504cb9d9b8SBarry Smith #undef __FUNC__
514cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatDestroy_MPIMAIJ"></a>*/"MatDestroy_MPIMAIJ"
524cb9d9b8SBarry Smith int MatDestroy_MPIMAIJ(Mat A)
534cb9d9b8SBarry Smith {
544cb9d9b8SBarry Smith   int         ierr;
554cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
564cb9d9b8SBarry Smith 
574cb9d9b8SBarry Smith   PetscFunctionBegin;
584cb9d9b8SBarry Smith   if (b->A) {
594cb9d9b8SBarry Smith     ierr = MatDestroy(b->A);CHKERRQ(ierr);
604cb9d9b8SBarry Smith   }
614cb9d9b8SBarry Smith   if (b->AIJ) {
624cb9d9b8SBarry Smith     ierr = MatDestroy(b->AIJ);CHKERRQ(ierr);
634cb9d9b8SBarry Smith   }
64f4a53059SBarry Smith   if (b->OAIJ) {
65f4a53059SBarry Smith     ierr = MatDestroy(b->OAIJ);CHKERRQ(ierr);
66f4a53059SBarry Smith   }
67f4a53059SBarry Smith   if (b->ctx) {
68f4a53059SBarry Smith     ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr);
69f4a53059SBarry Smith   }
70f4a53059SBarry Smith   if (b->w) {
71f4a53059SBarry Smith     ierr = VecDestroy(b->w);CHKERRQ(ierr);
72f4a53059SBarry Smith   }
73cd3bbe55SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
7482b1193eSBarry Smith   PetscHeaderDestroy(A);
7582b1193eSBarry Smith   PetscFunctionReturn(0);
7682b1193eSBarry Smith }
7782b1193eSBarry Smith 
7882b1193eSBarry Smith EXTERN_C_BEGIN
7982b1193eSBarry Smith #undef __FUNC__
80f4a53059SBarry Smith #define __FUNC__ /*<a name="MatCreate_MAIJ"></a>*/"MatCreate_MAIJ"
81f4a53059SBarry Smith int MatCreate_MAIJ(Mat A)
8282b1193eSBarry Smith {
83cd3bbe55SBarry Smith   int         ierr;
844cb9d9b8SBarry Smith   Mat_MPIMAIJ *b;
8582b1193eSBarry Smith 
8682b1193eSBarry Smith   PetscFunctionBegin;
874cb9d9b8SBarry Smith   A->data             = (void*)(b = PetscNew(Mat_MPIMAIJ));CHKPTRQ(b);
884cb9d9b8SBarry Smith   ierr = PetscMemzero(b,sizeof(Mat_MPIMAIJ));CHKERRQ(ierr);
89cd3bbe55SBarry Smith   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
90cd3bbe55SBarry Smith   A->factor           = 0;
91cd3bbe55SBarry Smith   A->mapping          = 0;
92f4a53059SBarry Smith 
93cd3bbe55SBarry Smith   b->AIJ  = 0;
94cd3bbe55SBarry Smith   b->dof  = 0;
95f4a53059SBarry Smith   b->OAIJ = 0;
96f4a53059SBarry Smith   b->ctx  = 0;
97f4a53059SBarry Smith   b->w    = 0;
9882b1193eSBarry Smith   PetscFunctionReturn(0);
9982b1193eSBarry Smith }
10082b1193eSBarry Smith EXTERN_C_END
10182b1193eSBarry Smith 
102cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/
10382b1193eSBarry Smith #undef __FUNC__
104cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_2"></a>*/"MatMult_SeqMAIJ_2"
105cd3bbe55SBarry Smith int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
10682b1193eSBarry Smith {
1074cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
108bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
109bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2;
110bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
111bcc973b7SBarry Smith   int         n,i,jrow,j;
11282b1193eSBarry Smith 
113bcc973b7SBarry Smith   PetscFunctionBegin;
114bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
115bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
116bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
117bcc973b7SBarry Smith   idx  = a->j;
118bcc973b7SBarry Smith   v    = a->a;
119bcc973b7SBarry Smith   ii   = a->i;
120bcc973b7SBarry Smith 
121bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
122bcc973b7SBarry Smith   idx  += shift;
123bcc973b7SBarry Smith   for (i=0; i<m; i++) {
124bcc973b7SBarry Smith     jrow = ii[i];
125bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
126bcc973b7SBarry Smith     sum1  = 0.0;
127bcc973b7SBarry Smith     sum2  = 0.0;
128bcc973b7SBarry Smith     for (j=0; j<n; j++) {
129bcc973b7SBarry Smith       sum1 += v[jrow]*x[2*idx[jrow]];
130bcc973b7SBarry Smith       sum2 += v[jrow]*x[2*idx[jrow]+1];
131bcc973b7SBarry Smith       jrow++;
132bcc973b7SBarry Smith      }
133bcc973b7SBarry Smith     y[2*i]   = sum1;
134bcc973b7SBarry Smith     y[2*i+1] = sum2;
135bcc973b7SBarry Smith   }
136bcc973b7SBarry Smith 
137bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*m);
138bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
139bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
14082b1193eSBarry Smith   PetscFunctionReturn(0);
14182b1193eSBarry Smith }
142bcc973b7SBarry Smith 
14382b1193eSBarry Smith #undef __FUNC__
144cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_2"></a>*/"MatMultTranspose_SeqMAIJ_2"
145cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
14682b1193eSBarry Smith {
1474cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
148bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
149bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,zero = 0.0;
150bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
15182b1193eSBarry Smith 
152bcc973b7SBarry Smith   PetscFunctionBegin;
153bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
154bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
155bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
156bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
157bcc973b7SBarry Smith   for (i=0; i<m; i++) {
158bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
159bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
160bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
161bcc973b7SBarry Smith     alpha1 = x[2*i];
162bcc973b7SBarry Smith     alpha2 = x[2*i+1];
163bcc973b7SBarry Smith     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
164bcc973b7SBarry Smith   }
165bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*a->n);
166bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
167bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
16882b1193eSBarry Smith   PetscFunctionReturn(0);
16982b1193eSBarry Smith }
170bcc973b7SBarry Smith 
17182b1193eSBarry Smith #undef __FUNC__
172cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_2"></a>*/"MatMultAdd_SeqMAIJ_2"
173cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
17482b1193eSBarry Smith {
1754cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
176bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
177bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2;
178bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
179bcc973b7SBarry Smith   int         n,i,jrow,j;
18082b1193eSBarry Smith 
181bcc973b7SBarry Smith   PetscFunctionBegin;
182f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
183bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
184f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
185bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
186bcc973b7SBarry Smith   idx  = a->j;
187bcc973b7SBarry Smith   v    = a->a;
188bcc973b7SBarry Smith   ii   = a->i;
189bcc973b7SBarry Smith 
190bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
191bcc973b7SBarry Smith   idx  += shift;
192bcc973b7SBarry Smith   for (i=0; i<m; i++) {
193bcc973b7SBarry Smith     jrow = ii[i];
194bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
195bcc973b7SBarry Smith     sum1  = 0.0;
196bcc973b7SBarry Smith     sum2  = 0.0;
197bcc973b7SBarry Smith     for (j=0; j<n; j++) {
198bcc973b7SBarry Smith       sum1 += v[jrow]*x[2*idx[jrow]];
199bcc973b7SBarry Smith       sum2 += v[jrow]*x[2*idx[jrow]+1];
200bcc973b7SBarry Smith       jrow++;
201bcc973b7SBarry Smith      }
202bcc973b7SBarry Smith     y[2*i]   += sum1;
203bcc973b7SBarry Smith     y[2*i+1] += sum2;
204bcc973b7SBarry Smith   }
205bcc973b7SBarry Smith 
206bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*m);
207bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
208f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
209bcc973b7SBarry Smith   PetscFunctionReturn(0);
21082b1193eSBarry Smith }
21182b1193eSBarry Smith #undef __FUNC__
212cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_2"></a>*/"MatMultTransposeAdd_SeqMAIJ_2"
213cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
21482b1193eSBarry Smith {
2154cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
216bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
217bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2;
218bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
21982b1193eSBarry Smith 
220bcc973b7SBarry Smith   PetscFunctionBegin;
221f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
222bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
223f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
224bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
225bcc973b7SBarry Smith   for (i=0; i<m; i++) {
226bcc973b7SBarry Smith     idx   = a->j + a->i[i] + shift;
227bcc973b7SBarry Smith     v     = a->a + a->i[i] + shift;
228bcc973b7SBarry Smith     n     = a->i[i+1] - a->i[i];
229bcc973b7SBarry Smith     alpha1 = x[2*i];
230bcc973b7SBarry Smith     alpha2 = x[2*i+1];
231bcc973b7SBarry Smith     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
232bcc973b7SBarry Smith   }
233bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*a->n);
234bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
235f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
236bcc973b7SBarry Smith   PetscFunctionReturn(0);
23782b1193eSBarry Smith }
238cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/
239bcc973b7SBarry Smith #undef __FUNC__
240bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_3"></a>*/"MatMult_SeqMAIJ_3"
241bcc973b7SBarry Smith int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
242bcc973b7SBarry Smith {
2434cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
244bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
245bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3;
246bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
247bcc973b7SBarry Smith   int         n,i,jrow,j;
24882b1193eSBarry Smith 
249bcc973b7SBarry Smith   PetscFunctionBegin;
250bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
251bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
252bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
253bcc973b7SBarry Smith   idx  = a->j;
254bcc973b7SBarry Smith   v    = a->a;
255bcc973b7SBarry Smith   ii   = a->i;
256bcc973b7SBarry Smith 
257bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
258bcc973b7SBarry Smith   idx  += shift;
259bcc973b7SBarry Smith   for (i=0; i<m; i++) {
260bcc973b7SBarry Smith     jrow = ii[i];
261bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
262bcc973b7SBarry Smith     sum1  = 0.0;
263bcc973b7SBarry Smith     sum2  = 0.0;
264bcc973b7SBarry Smith     sum3  = 0.0;
265bcc973b7SBarry Smith     for (j=0; j<n; j++) {
266bcc973b7SBarry Smith       sum1 += v[jrow]*x[3*idx[jrow]];
267bcc973b7SBarry Smith       sum2 += v[jrow]*x[3*idx[jrow]+1];
268bcc973b7SBarry Smith       sum3 += v[jrow]*x[3*idx[jrow]+2];
269bcc973b7SBarry Smith       jrow++;
270bcc973b7SBarry Smith      }
271bcc973b7SBarry Smith     y[3*i]   = sum1;
272bcc973b7SBarry Smith     y[3*i+1] = sum2;
273bcc973b7SBarry Smith     y[3*i+2] = sum3;
274bcc973b7SBarry Smith   }
275bcc973b7SBarry Smith 
276bcc973b7SBarry Smith   PLogFlops(6*a->nz - 3*m);
277bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
278bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
279bcc973b7SBarry Smith   PetscFunctionReturn(0);
280bcc973b7SBarry Smith }
281bcc973b7SBarry Smith 
282bcc973b7SBarry Smith #undef __FUNC__
283bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_3"></a>*/"MatMultTranspose_SeqMAIJ_3"
284bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
285bcc973b7SBarry Smith {
2864cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
287bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
288bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0;
289bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
290bcc973b7SBarry Smith 
291bcc973b7SBarry Smith   PetscFunctionBegin;
292bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
293bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
294bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
295bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
296bcc973b7SBarry Smith   for (i=0; i<m; i++) {
297bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
298bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
299bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
300bcc973b7SBarry Smith     alpha1 = x[3*i];
301bcc973b7SBarry Smith     alpha2 = x[3*i+1];
302bcc973b7SBarry Smith     alpha3 = x[3*i+2];
303bcc973b7SBarry Smith     while (n-->0) {
304bcc973b7SBarry Smith       y[3*(*idx)]   += alpha1*(*v);
305bcc973b7SBarry Smith       y[3*(*idx)+1] += alpha2*(*v);
306bcc973b7SBarry Smith       y[3*(*idx)+2] += alpha3*(*v);
307bcc973b7SBarry Smith       idx++; v++;
308bcc973b7SBarry Smith     }
309bcc973b7SBarry Smith   }
310bcc973b7SBarry Smith   PLogFlops(6*a->nz - 3*a->n);
311bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
312bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
313bcc973b7SBarry Smith   PetscFunctionReturn(0);
314bcc973b7SBarry Smith }
315bcc973b7SBarry Smith 
316bcc973b7SBarry Smith #undef __FUNC__
317bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_3"></a>*/"MatMultAdd_SeqMAIJ_3"
318bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
319bcc973b7SBarry Smith {
3204cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
321bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
322bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3;
323bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
324bcc973b7SBarry Smith   int         n,i,jrow,j;
325bcc973b7SBarry Smith 
326bcc973b7SBarry Smith   PetscFunctionBegin;
327f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
328bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
329f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
330bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
331bcc973b7SBarry Smith   idx  = a->j;
332bcc973b7SBarry Smith   v    = a->a;
333bcc973b7SBarry Smith   ii   = a->i;
334bcc973b7SBarry Smith 
335bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
336bcc973b7SBarry Smith   idx  += shift;
337bcc973b7SBarry Smith   for (i=0; i<m; i++) {
338bcc973b7SBarry Smith     jrow = ii[i];
339bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
340bcc973b7SBarry Smith     sum1  = 0.0;
341bcc973b7SBarry Smith     sum2  = 0.0;
342bcc973b7SBarry Smith     sum3  = 0.0;
343bcc973b7SBarry Smith     for (j=0; j<n; j++) {
344bcc973b7SBarry Smith       sum1 += v[jrow]*x[3*idx[jrow]];
345bcc973b7SBarry Smith       sum2 += v[jrow]*x[3*idx[jrow]+1];
346bcc973b7SBarry Smith       sum3 += v[jrow]*x[3*idx[jrow]+2];
347bcc973b7SBarry Smith       jrow++;
348bcc973b7SBarry Smith      }
349bcc973b7SBarry Smith     y[3*i]   += sum1;
350bcc973b7SBarry Smith     y[3*i+1] += sum2;
351bcc973b7SBarry Smith     y[3*i+2] += sum3;
352bcc973b7SBarry Smith   }
353bcc973b7SBarry Smith 
354bcc973b7SBarry Smith   PLogFlops(6*a->nz);
355bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
356f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
357bcc973b7SBarry Smith   PetscFunctionReturn(0);
358bcc973b7SBarry Smith }
359bcc973b7SBarry Smith #undef __FUNC__
360bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_3"></a>*/"MatMultTransposeAdd_SeqMAIJ_3"
361bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
362bcc973b7SBarry Smith {
3634cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
364bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
365bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3;
366bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
367bcc973b7SBarry Smith 
368bcc973b7SBarry Smith   PetscFunctionBegin;
369f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
370bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
371f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
372bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
373bcc973b7SBarry Smith   for (i=0; i<m; i++) {
374bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
375bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
376bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
377bcc973b7SBarry Smith     alpha1 = x[3*i];
378bcc973b7SBarry Smith     alpha2 = x[3*i+1];
379bcc973b7SBarry Smith     alpha3 = x[3*i+2];
380bcc973b7SBarry Smith     while (n-->0) {
381bcc973b7SBarry Smith       y[3*(*idx)]   += alpha1*(*v);
382bcc973b7SBarry Smith       y[3*(*idx)+1] += alpha2*(*v);
383bcc973b7SBarry Smith       y[3*(*idx)+2] += alpha3*(*v);
384bcc973b7SBarry Smith       idx++; v++;
385bcc973b7SBarry Smith     }
386bcc973b7SBarry Smith   }
387bcc973b7SBarry Smith   PLogFlops(6*a->nz);
388bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
389f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
390bcc973b7SBarry Smith   PetscFunctionReturn(0);
391bcc973b7SBarry Smith }
392bcc973b7SBarry Smith 
393bcc973b7SBarry Smith /* ------------------------------------------------------------------------------*/
394bcc973b7SBarry Smith #undef __FUNC__
395bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_4"></a>*/"MatMult_SeqMAIJ_4"
396bcc973b7SBarry Smith int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
397bcc973b7SBarry Smith {
3984cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
399bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
400bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4;
401bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
402bcc973b7SBarry Smith   int         n,i,jrow,j;
403bcc973b7SBarry Smith 
404bcc973b7SBarry Smith   PetscFunctionBegin;
405bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
406bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
407bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
408bcc973b7SBarry Smith   idx  = a->j;
409bcc973b7SBarry Smith   v    = a->a;
410bcc973b7SBarry Smith   ii   = a->i;
411bcc973b7SBarry Smith 
412bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
413bcc973b7SBarry Smith   idx  += shift;
414bcc973b7SBarry Smith   for (i=0; i<m; i++) {
415bcc973b7SBarry Smith     jrow = ii[i];
416bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
417bcc973b7SBarry Smith     sum1  = 0.0;
418bcc973b7SBarry Smith     sum2  = 0.0;
419bcc973b7SBarry Smith     sum3  = 0.0;
420bcc973b7SBarry Smith     sum4  = 0.0;
421bcc973b7SBarry Smith     for (j=0; j<n; j++) {
422bcc973b7SBarry Smith       sum1 += v[jrow]*x[4*idx[jrow]];
423bcc973b7SBarry Smith       sum2 += v[jrow]*x[4*idx[jrow]+1];
424bcc973b7SBarry Smith       sum3 += v[jrow]*x[4*idx[jrow]+2];
425bcc973b7SBarry Smith       sum4 += v[jrow]*x[4*idx[jrow]+3];
426bcc973b7SBarry Smith       jrow++;
427bcc973b7SBarry Smith      }
428bcc973b7SBarry Smith     y[4*i]   = sum1;
429bcc973b7SBarry Smith     y[4*i+1] = sum2;
430bcc973b7SBarry Smith     y[4*i+2] = sum3;
431bcc973b7SBarry Smith     y[4*i+3] = sum4;
432bcc973b7SBarry Smith   }
433bcc973b7SBarry Smith 
434bcc973b7SBarry Smith   PLogFlops(8*a->nz - 4*m);
435bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
436bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
437bcc973b7SBarry Smith   PetscFunctionReturn(0);
438bcc973b7SBarry Smith }
439bcc973b7SBarry Smith 
440bcc973b7SBarry Smith #undef __FUNC__
441bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_4"></a>*/"MatMultTranspose_SeqMAIJ_4"
442bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
443bcc973b7SBarry Smith {
4444cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
445bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
446bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0;
447bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
448bcc973b7SBarry Smith 
449bcc973b7SBarry Smith   PetscFunctionBegin;
450bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
451bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
452bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
453bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
454bcc973b7SBarry Smith   for (i=0; i<m; i++) {
455bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
456bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
457bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
458bcc973b7SBarry Smith     alpha1 = x[4*i];
459bcc973b7SBarry Smith     alpha2 = x[4*i+1];
460bcc973b7SBarry Smith     alpha3 = x[4*i+2];
461bcc973b7SBarry Smith     alpha4 = x[4*i+3];
462bcc973b7SBarry Smith     while (n-->0) {
463bcc973b7SBarry Smith       y[4*(*idx)]   += alpha1*(*v);
464bcc973b7SBarry Smith       y[4*(*idx)+1] += alpha2*(*v);
465bcc973b7SBarry Smith       y[4*(*idx)+2] += alpha3*(*v);
466bcc973b7SBarry Smith       y[4*(*idx)+3] += alpha4*(*v);
467bcc973b7SBarry Smith       idx++; v++;
468bcc973b7SBarry Smith     }
469bcc973b7SBarry Smith   }
470bcc973b7SBarry Smith   PLogFlops(8*a->nz - 4*a->n);
471bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
472bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
473bcc973b7SBarry Smith   PetscFunctionReturn(0);
474bcc973b7SBarry Smith }
475bcc973b7SBarry Smith 
476bcc973b7SBarry Smith #undef __FUNC__
477bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_4"></a>*/"MatMultAdd_SeqMAIJ_4"
478bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
479bcc973b7SBarry Smith {
4804cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
481f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
482f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4;
483f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
484f9fae5adSBarry Smith   int         n,i,jrow,j;
485f9fae5adSBarry Smith 
486f9fae5adSBarry Smith   PetscFunctionBegin;
487f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
488f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
489f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
490f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
491f9fae5adSBarry Smith   idx  = a->j;
492f9fae5adSBarry Smith   v    = a->a;
493f9fae5adSBarry Smith   ii   = a->i;
494f9fae5adSBarry Smith 
495f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
496f9fae5adSBarry Smith   idx  += shift;
497f9fae5adSBarry Smith   for (i=0; i<m; i++) {
498f9fae5adSBarry Smith     jrow = ii[i];
499f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
500f9fae5adSBarry Smith     sum1  = 0.0;
501f9fae5adSBarry Smith     sum2  = 0.0;
502f9fae5adSBarry Smith     sum3  = 0.0;
503f9fae5adSBarry Smith     sum4  = 0.0;
504f9fae5adSBarry Smith     for (j=0; j<n; j++) {
505f9fae5adSBarry Smith       sum1 += v[jrow]*x[4*idx[jrow]];
506f9fae5adSBarry Smith       sum2 += v[jrow]*x[4*idx[jrow]+1];
507f9fae5adSBarry Smith       sum3 += v[jrow]*x[4*idx[jrow]+2];
508f9fae5adSBarry Smith       sum4 += v[jrow]*x[4*idx[jrow]+3];
509f9fae5adSBarry Smith       jrow++;
510f9fae5adSBarry Smith      }
511f9fae5adSBarry Smith     y[4*i]   += sum1;
512f9fae5adSBarry Smith     y[4*i+1] += sum2;
513f9fae5adSBarry Smith     y[4*i+2] += sum3;
514f9fae5adSBarry Smith     y[4*i+3] += sum4;
515f9fae5adSBarry Smith   }
516f9fae5adSBarry Smith 
517f9fae5adSBarry Smith   PLogFlops(8*a->nz - 4*m);
518f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
519f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
520f9fae5adSBarry Smith   PetscFunctionReturn(0);
521bcc973b7SBarry Smith }
522bcc973b7SBarry Smith #undef __FUNC__
523bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_4"></a>*/"MatMultTransposeAdd_SeqMAIJ_4"
524bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
525bcc973b7SBarry Smith {
5264cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
527f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
528f9fae5adSBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4;
529f9fae5adSBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
530f9fae5adSBarry Smith 
531f9fae5adSBarry Smith   PetscFunctionBegin;
532f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
533f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
534f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
535f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
536f9fae5adSBarry Smith   for (i=0; i<m; i++) {
537f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
538f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
539f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
540f9fae5adSBarry Smith     alpha1 = x[4*i];
541f9fae5adSBarry Smith     alpha2 = x[4*i+1];
542f9fae5adSBarry Smith     alpha3 = x[4*i+2];
543f9fae5adSBarry Smith     alpha4 = x[4*i+3];
544f9fae5adSBarry Smith     while (n-->0) {
545f9fae5adSBarry Smith       y[4*(*idx)]   += alpha1*(*v);
546f9fae5adSBarry Smith       y[4*(*idx)+1] += alpha2*(*v);
547f9fae5adSBarry Smith       y[4*(*idx)+2] += alpha3*(*v);
548f9fae5adSBarry Smith       y[4*(*idx)+3] += alpha4*(*v);
549f9fae5adSBarry Smith       idx++; v++;
550f9fae5adSBarry Smith     }
551f9fae5adSBarry Smith   }
552f9fae5adSBarry Smith   PLogFlops(8*a->nz - 4*a->n);
553f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
554f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
555f9fae5adSBarry Smith   PetscFunctionReturn(0);
556f9fae5adSBarry Smith 
557f9fae5adSBarry Smith }
558f9fae5adSBarry Smith 
559f9fae5adSBarry Smith /* ------------------------------------------------------------------------------*/
560f9fae5adSBarry Smith #undef __FUNC__
561f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_5"></a>*/"MatMult_SeqMAIJ_5"
562f9fae5adSBarry Smith int MatMult_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
563f9fae5adSBarry Smith {
5644cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
565f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
566f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
567f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
568f9fae5adSBarry Smith   int         n,i,jrow,j;
569f9fae5adSBarry Smith 
570f9fae5adSBarry Smith   PetscFunctionBegin;
571f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
572f9fae5adSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
573f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
574f9fae5adSBarry Smith   idx  = a->j;
575f9fae5adSBarry Smith   v    = a->a;
576f9fae5adSBarry Smith   ii   = a->i;
577f9fae5adSBarry Smith 
578f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
579f9fae5adSBarry Smith   idx  += shift;
580f9fae5adSBarry Smith   for (i=0; i<m; i++) {
581f9fae5adSBarry Smith     jrow = ii[i];
582f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
583f9fae5adSBarry Smith     sum1  = 0.0;
584f9fae5adSBarry Smith     sum2  = 0.0;
585f9fae5adSBarry Smith     sum3  = 0.0;
586f9fae5adSBarry Smith     sum4  = 0.0;
587f9fae5adSBarry Smith     sum5  = 0.0;
588f9fae5adSBarry Smith     for (j=0; j<n; j++) {
589f9fae5adSBarry Smith       sum1 += v[jrow]*x[5*idx[jrow]];
590f9fae5adSBarry Smith       sum2 += v[jrow]*x[5*idx[jrow]+1];
591f9fae5adSBarry Smith       sum3 += v[jrow]*x[5*idx[jrow]+2];
592f9fae5adSBarry Smith       sum4 += v[jrow]*x[5*idx[jrow]+3];
593f9fae5adSBarry Smith       sum5 += v[jrow]*x[5*idx[jrow]+4];
594f9fae5adSBarry Smith       jrow++;
595f9fae5adSBarry Smith      }
596f9fae5adSBarry Smith     y[5*i]   = sum1;
597f9fae5adSBarry Smith     y[5*i+1] = sum2;
598f9fae5adSBarry Smith     y[5*i+2] = sum3;
599f9fae5adSBarry Smith     y[5*i+3] = sum4;
600f9fae5adSBarry Smith     y[5*i+4] = sum5;
601f9fae5adSBarry Smith   }
602f9fae5adSBarry Smith 
603f9fae5adSBarry Smith   PLogFlops(10*a->nz - 5*m);
604f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
605f9fae5adSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
606f9fae5adSBarry Smith   PetscFunctionReturn(0);
607f9fae5adSBarry Smith }
608f9fae5adSBarry Smith 
609f9fae5adSBarry Smith #undef __FUNC__
610f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_5"></a>*/"MatMultTranspose_SeqMAIJ_5"
611f9fae5adSBarry Smith int MatMultTranspose_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
612f9fae5adSBarry Smith {
6134cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
614f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
615f9fae5adSBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,zero = 0.0;
616f9fae5adSBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
617f9fae5adSBarry Smith 
618f9fae5adSBarry Smith   PetscFunctionBegin;
619f9fae5adSBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
620f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
621f9fae5adSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
622f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
623f9fae5adSBarry Smith   for (i=0; i<m; i++) {
624f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
625f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
626f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
627f9fae5adSBarry Smith     alpha1 = x[5*i];
628f9fae5adSBarry Smith     alpha2 = x[5*i+1];
629f9fae5adSBarry Smith     alpha3 = x[5*i+2];
630f9fae5adSBarry Smith     alpha4 = x[5*i+3];
631f9fae5adSBarry Smith     alpha5 = x[5*i+4];
632f9fae5adSBarry Smith     while (n-->0) {
633f9fae5adSBarry Smith       y[5*(*idx)]   += alpha1*(*v);
634f9fae5adSBarry Smith       y[5*(*idx)+1] += alpha2*(*v);
635f9fae5adSBarry Smith       y[5*(*idx)+2] += alpha3*(*v);
636f9fae5adSBarry Smith       y[5*(*idx)+3] += alpha4*(*v);
637f9fae5adSBarry Smith       y[5*(*idx)+4] += alpha5*(*v);
638f9fae5adSBarry Smith       idx++; v++;
639f9fae5adSBarry Smith     }
640f9fae5adSBarry Smith   }
641f9fae5adSBarry Smith   PLogFlops(10*a->nz - 5*a->n);
642f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
643f9fae5adSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
644f9fae5adSBarry Smith   PetscFunctionReturn(0);
645f9fae5adSBarry Smith }
646f9fae5adSBarry Smith 
647f9fae5adSBarry Smith #undef __FUNC__
648f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_5"></a>*/"MatMultAdd_SeqMAIJ_5"
649f9fae5adSBarry Smith int MatMultAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
650f9fae5adSBarry Smith {
6514cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
652f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
653f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
654f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
655f9fae5adSBarry Smith   int         n,i,jrow,j;
656f9fae5adSBarry Smith 
657f9fae5adSBarry Smith   PetscFunctionBegin;
658f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
659f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
660f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
661f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
662f9fae5adSBarry Smith   idx  = a->j;
663f9fae5adSBarry Smith   v    = a->a;
664f9fae5adSBarry Smith   ii   = a->i;
665f9fae5adSBarry Smith 
666f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
667f9fae5adSBarry Smith   idx  += shift;
668f9fae5adSBarry Smith   for (i=0; i<m; i++) {
669f9fae5adSBarry Smith     jrow = ii[i];
670f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
671f9fae5adSBarry Smith     sum1  = 0.0;
672f9fae5adSBarry Smith     sum2  = 0.0;
673f9fae5adSBarry Smith     sum3  = 0.0;
674f9fae5adSBarry Smith     sum4  = 0.0;
675f9fae5adSBarry Smith     sum5  = 0.0;
676f9fae5adSBarry Smith     for (j=0; j<n; j++) {
677f9fae5adSBarry Smith       sum1 += v[jrow]*x[5*idx[jrow]];
678f9fae5adSBarry Smith       sum2 += v[jrow]*x[5*idx[jrow]+1];
679f9fae5adSBarry Smith       sum3 += v[jrow]*x[5*idx[jrow]+2];
680f9fae5adSBarry Smith       sum4 += v[jrow]*x[5*idx[jrow]+3];
681f9fae5adSBarry Smith       sum5 += v[jrow]*x[5*idx[jrow]+4];
682f9fae5adSBarry Smith       jrow++;
683f9fae5adSBarry Smith      }
684f9fae5adSBarry Smith     y[5*i]   += sum1;
685f9fae5adSBarry Smith     y[5*i+1] += sum2;
686f9fae5adSBarry Smith     y[5*i+2] += sum3;
687f9fae5adSBarry Smith     y[5*i+3] += sum4;
688f9fae5adSBarry Smith     y[5*i+4] += sum5;
689f9fae5adSBarry Smith   }
690f9fae5adSBarry Smith 
691f9fae5adSBarry Smith   PLogFlops(10*a->nz);
692f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
693f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
694f9fae5adSBarry Smith   PetscFunctionReturn(0);
695f9fae5adSBarry Smith }
696f9fae5adSBarry Smith 
697f9fae5adSBarry Smith #undef __FUNC__
698f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_5"></a>*/"MatMultTransposeAdd_SeqMAIJ_5"
699f9fae5adSBarry Smith int MatMultTransposeAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
700f9fae5adSBarry Smith {
7014cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
702f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
703f9fae5adSBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5;
704f9fae5adSBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
705f9fae5adSBarry Smith 
706f9fae5adSBarry Smith   PetscFunctionBegin;
707f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
708f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
709f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
710f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
711f9fae5adSBarry Smith   for (i=0; i<m; i++) {
712f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
713f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
714f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
715f9fae5adSBarry Smith     alpha1 = x[5*i];
716f9fae5adSBarry Smith     alpha2 = x[5*i+1];
717f9fae5adSBarry Smith     alpha3 = x[5*i+2];
718f9fae5adSBarry Smith     alpha4 = x[5*i+3];
719f9fae5adSBarry Smith     alpha5 = x[5*i+4];
720f9fae5adSBarry Smith     while (n-->0) {
721f9fae5adSBarry Smith       y[5*(*idx)]   += alpha1*(*v);
722f9fae5adSBarry Smith       y[5*(*idx)+1] += alpha2*(*v);
723f9fae5adSBarry Smith       y[5*(*idx)+2] += alpha3*(*v);
724f9fae5adSBarry Smith       y[5*(*idx)+3] += alpha4*(*v);
725f9fae5adSBarry Smith       y[5*(*idx)+4] += alpha5*(*v);
726f9fae5adSBarry Smith       idx++; v++;
727f9fae5adSBarry Smith     }
728f9fae5adSBarry Smith   }
729f9fae5adSBarry Smith   PLogFlops(10*a->nz);
730f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
731f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
732f9fae5adSBarry Smith   PetscFunctionReturn(0);
733bcc973b7SBarry Smith }
734bcc973b7SBarry Smith 
735f4a53059SBarry Smith /*===================================================================================*/
736f4a53059SBarry Smith #undef __FUNC__
737f4a53059SBarry Smith #define __FUNC__ /*<a name="MatMult_MPIMAIJ_dof"></a>*/"MatMult_MPIMAIJ_dof"
738f4a53059SBarry Smith int MatMult_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
739f4a53059SBarry Smith {
7404cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
741f4a53059SBarry Smith   int         ierr;
742f4a53059SBarry Smith   PetscFunctionBegin;
743f4a53059SBarry Smith 
744f4a53059SBarry Smith   /* start the scatter */
745f4a53059SBarry Smith   ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
7464cb9d9b8SBarry Smith   ierr = (*b->AIJ->ops->mult)(b->AIJ,xx,yy);CHKERRQ(ierr);
747f4a53059SBarry Smith   ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
7484cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr);
749f4a53059SBarry Smith   PetscFunctionReturn(0);
750f4a53059SBarry Smith }
751f4a53059SBarry Smith 
7524cb9d9b8SBarry Smith #undef __FUNC__
7534cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_MPIMAIJ_dof"></a>*/"MatMultTranspose_MPIMAIJ_dof"
7544cb9d9b8SBarry Smith int MatMultTranspose_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
7554cb9d9b8SBarry Smith {
7564cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
7574cb9d9b8SBarry Smith   int         ierr;
7584cb9d9b8SBarry Smith   PetscFunctionBegin;
7594cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr);
7604cb9d9b8SBarry Smith   ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
7614cb9d9b8SBarry Smith   ierr = (*b->AIJ->ops->multtranspose)(b->AIJ,xx,yy);CHKERRQ(ierr);
7624cb9d9b8SBarry Smith   ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
7634cb9d9b8SBarry Smith   PetscFunctionReturn(0);
7644cb9d9b8SBarry Smith }
7654cb9d9b8SBarry Smith 
7664cb9d9b8SBarry Smith #undef __FUNC__
7674cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_MPIMAIJ_dof"></a>*/"MatMultAdd_MPIMAIJ_dof"
768*d72c5c08SBarry Smith int MatMultAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz)
7694cb9d9b8SBarry Smith {
7704cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
7714cb9d9b8SBarry Smith   int         ierr;
7724cb9d9b8SBarry Smith   PetscFunctionBegin;
7734cb9d9b8SBarry Smith 
7744cb9d9b8SBarry Smith   /* start the scatter */
7754cb9d9b8SBarry Smith   ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
776*d72c5c08SBarry Smith   ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr);
7774cb9d9b8SBarry Smith   ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
778*d72c5c08SBarry Smith   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,zz);CHKERRQ(ierr);
7794cb9d9b8SBarry Smith   PetscFunctionReturn(0);
7804cb9d9b8SBarry Smith }
7814cb9d9b8SBarry Smith 
7824cb9d9b8SBarry Smith #undef __FUNC__
7834cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_MPIMAIJ_dof"></a>*/"MatMultTransposeAdd_MPIMAIJ_dof"
784*d72c5c08SBarry Smith int MatMultTransposeAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz)
7854cb9d9b8SBarry Smith {
7864cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
7874cb9d9b8SBarry Smith   int         ierr;
7884cb9d9b8SBarry Smith   PetscFunctionBegin;
7894cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr);
790*d72c5c08SBarry Smith   ierr = VecScatterBegin(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
791*d72c5c08SBarry Smith   ierr = (*b->AIJ->ops->multtransposeadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr);
792*d72c5c08SBarry Smith   ierr = VecScatterEnd(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
7934cb9d9b8SBarry Smith   PetscFunctionReturn(0);
7944cb9d9b8SBarry Smith }
7954cb9d9b8SBarry Smith 
796bcc973b7SBarry Smith /* ---------------------------------------------------------------------------------- */
79782b1193eSBarry Smith #undef __FUNC__
798cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatCreateMAIJ"></a>*/"MatCreateMAIJ"
799cd3bbe55SBarry Smith int MatCreateMAIJ(Mat A,int dof,Mat *maij)
80082b1193eSBarry Smith {
801f4a53059SBarry Smith   int         ierr,size,n;
8024cb9d9b8SBarry Smith   Mat_MPIMAIJ *b;
80382b1193eSBarry Smith   Mat         B;
80482b1193eSBarry Smith 
80582b1193eSBarry Smith   PetscFunctionBegin;
806*d72c5c08SBarry Smith   ierr   = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
807*d72c5c08SBarry Smith 
808*d72c5c08SBarry Smith   if (dof == 1) {
809*d72c5c08SBarry Smith     *maij = A;
810*d72c5c08SBarry Smith   } else {
811cd3bbe55SBarry Smith     ierr = MATCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr);
812f4a53059SBarry Smith     ierr = MatSetType(B,MATMAIJ);CHKERRQ(ierr);
813cd3bbe55SBarry Smith     B->assembled    = PETSC_TRUE;
8144cb9d9b8SBarry Smith     b = (Mat_MPIMAIJ*)B->data;
815cd3bbe55SBarry Smith     b->dof = dof;
816*d72c5c08SBarry Smith 
817f4a53059SBarry Smith     ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
818f4a53059SBarry Smith     if (size == 1) {
8194cb9d9b8SBarry Smith       B->ops->destroy = MatDestroy_SeqMAIJ;
8204cb9d9b8SBarry Smith       b->AIJ = A;
821*d72c5c08SBarry Smith       if (dof == 2) {
822cd3bbe55SBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_2;
823cd3bbe55SBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_2;
824cd3bbe55SBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_2;
825cd3bbe55SBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2;
826bcc973b7SBarry Smith       } else if (dof == 3) {
827bcc973b7SBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_3;
828bcc973b7SBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_3;
829bcc973b7SBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_3;
830bcc973b7SBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3;
831bcc973b7SBarry Smith       } else if (dof == 4) {
832bcc973b7SBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_4;
833bcc973b7SBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_4;
834bcc973b7SBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_4;
835bcc973b7SBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4;
836f9fae5adSBarry Smith       } else if (dof == 5) {
837f9fae5adSBarry Smith         B->ops->mult             = MatMult_SeqMAIJ_5;
838f9fae5adSBarry Smith         B->ops->multadd          = MatMultAdd_SeqMAIJ_5;
839f9fae5adSBarry Smith         B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_5;
840f9fae5adSBarry Smith         B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5;
84182b1193eSBarry Smith       } else {
842cd3bbe55SBarry Smith         SETERRQ1(1,1,"Cannot handle a dof of %d\n",dof);
84382b1193eSBarry Smith       }
844f4a53059SBarry Smith     } else {
845f4a53059SBarry Smith       Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data;
846f4a53059SBarry Smith       IS         from,to;
847f4a53059SBarry Smith       Vec        gvec;
848f4a53059SBarry Smith       int        *garray,i;
849f4a53059SBarry Smith 
850*d72c5c08SBarry Smith       b->A = A;
851*d72c5c08SBarry Smith       B->ops->destroy = MatDestroy_MPIMAIJ;
8524cb9d9b8SBarry Smith       ierr = MatCreateMAIJ(mpiaij->A,dof,&b->AIJ);CHKERRQ(ierr);
8534cb9d9b8SBarry Smith       ierr = MatCreateMAIJ(mpiaij->B,dof,&b->OAIJ);CHKERRQ(ierr);
8544cb9d9b8SBarry Smith 
855f4a53059SBarry Smith       ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr);
856f4a53059SBarry Smith       ierr = VecCreateSeq(PETSC_COMM_SELF,n*dof,&b->w);CHKERRQ(ierr);
857f4a53059SBarry Smith 
858f4a53059SBarry Smith       /* create two temporary Index sets for build scatter gather */
859f4a53059SBarry Smith       garray = (int*)PetscMalloc((n+1)*sizeof(int));CHKPTRQ(garray);
860f4a53059SBarry Smith       for (i=0; i<n; i++) garray[i] = dof*mpiaij->garray[i];
861f4a53059SBarry Smith       ierr = ISCreateBlock(A->comm,dof,n,garray,&from);CHKERRQ(ierr);
862f4a53059SBarry Smith       ierr = PetscFree(garray);CHKERRQ(ierr);
863f4a53059SBarry Smith       ierr = ISCreateStride(PETSC_COMM_SELF,n*dof,0,1,&to);CHKERRQ(ierr);
864f4a53059SBarry Smith 
865f4a53059SBarry Smith       /* create temporary global vector to generate scatter context */
866f4a53059SBarry Smith       ierr = VecCreateMPI(A->comm,dof*A->n,dof*A->N,&gvec);CHKERRQ(ierr);
867f4a53059SBarry Smith 
868f4a53059SBarry Smith       /* generate the scatter context */
869f4a53059SBarry Smith       ierr = VecScatterCreate(gvec,from,b->w,to,&b->ctx);CHKERRQ(ierr);
870f4a53059SBarry Smith 
871f4a53059SBarry Smith       ierr = ISDestroy(from);CHKERRQ(ierr);
872f4a53059SBarry Smith       ierr = ISDestroy(to);CHKERRQ(ierr);
873f4a53059SBarry Smith       ierr = VecDestroy(gvec);CHKERRQ(ierr);
874f4a53059SBarry Smith 
875f4a53059SBarry Smith       B->ops->mult             = MatMult_MPIMAIJ_dof;
8764cb9d9b8SBarry Smith       B->ops->multtranspose    = MatMultTranspose_MPIMAIJ_dof;
8774cb9d9b8SBarry Smith       B->ops->multadd          = MatMultAdd_MPIMAIJ_dof;
8784cb9d9b8SBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof;
879f4a53059SBarry Smith     }
880cd3bbe55SBarry Smith     *maij = B;
881*d72c5c08SBarry Smith   }
88282b1193eSBarry Smith   PetscFunctionReturn(0);
88382b1193eSBarry Smith }
88482b1193eSBarry Smith 
88582b1193eSBarry Smith 
88682b1193eSBarry Smith 
88782b1193eSBarry Smith 
88882b1193eSBarry Smith 
88982b1193eSBarry Smith 
89082b1193eSBarry Smith 
89182b1193eSBarry Smith 
89282b1193eSBarry Smith 
89382b1193eSBarry Smith 
89482b1193eSBarry Smith 
89582b1193eSBarry Smith 
896