xref: /petsc/src/mat/impls/maij/maij.c (revision 4cb9d9b824e32e79ec1cb9b3d471027df5d623c0)
1*4cb9d9b8SBarry Smith /*$Id: maij.c,v 1.5 2000/06/05 19:38:31 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 */
23*4cb9d9b8SBarry Smith   Mat        AIJ;        /* representation of interpolation for one component */
24*4cb9d9b8SBarry Smith } Mat_SeqMAIJ;
25*4cb9d9b8SBarry Smith 
26*4cb9d9b8SBarry Smith typedef struct {
27*4cb9d9b8SBarry Smith   int        dof;         /* number of components */
28f4a53059SBarry Smith   Mat        AIJ,OAIJ;    /* representation of interpolation for one component */
29*4cb9d9b8SBarry 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 */
32*4cb9d9b8SBarry Smith } Mat_MPIMAIJ;
3382b1193eSBarry Smith 
3482b1193eSBarry Smith #undef __FUNC__
35*4cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatDestroy_SeqMAIJ"></a>*/"MatDestroy_SeqMAIJ"
36*4cb9d9b8SBarry Smith int MatDestroy_SeqMAIJ(Mat A)
3782b1193eSBarry Smith {
3882b1193eSBarry Smith   int         ierr;
39*4cb9d9b8SBarry 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   }
45*4cb9d9b8SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
46*4cb9d9b8SBarry Smith   PetscHeaderDestroy(A);
47*4cb9d9b8SBarry Smith   PetscFunctionReturn(0);
48*4cb9d9b8SBarry Smith }
49*4cb9d9b8SBarry Smith 
50*4cb9d9b8SBarry Smith #undef __FUNC__
51*4cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatDestroy_MPIMAIJ"></a>*/"MatDestroy_MPIMAIJ"
52*4cb9d9b8SBarry Smith int MatDestroy_MPIMAIJ(Mat A)
53*4cb9d9b8SBarry Smith {
54*4cb9d9b8SBarry Smith   int         ierr;
55*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
56*4cb9d9b8SBarry Smith 
57*4cb9d9b8SBarry Smith   PetscFunctionBegin;
58*4cb9d9b8SBarry Smith   if (b->A) {
59*4cb9d9b8SBarry Smith     ierr = MatDestroy(b->A);CHKERRQ(ierr);
60*4cb9d9b8SBarry Smith   }
61*4cb9d9b8SBarry Smith   if (b->AIJ) {
62*4cb9d9b8SBarry Smith     ierr = MatDestroy(b->AIJ);CHKERRQ(ierr);
63*4cb9d9b8SBarry 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;
84*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b;
8582b1193eSBarry Smith 
8682b1193eSBarry Smith   PetscFunctionBegin;
87*4cb9d9b8SBarry Smith   A->data             = (void*)(b = PetscNew(Mat_MPIMAIJ));CHKPTRQ(b);
88*4cb9d9b8SBarry 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 
102bcc973b7SBarry Smith /* ----------------------------------------------------------------------------------*/
103cd3bbe55SBarry Smith EXTERN int MatMult_SeqAIJ(Mat,Vec,Vec);
104f4a53059SBarry Smith EXTERN int MatMultAdd_SeqAIJ(Mat,Vec,Vec,Vec);
105cd3bbe55SBarry Smith EXTERN int MatMultTranspose_SeqAIJ(Mat,Vec,Vec);
106cd3bbe55SBarry Smith EXTERN int MatMultTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
107cd3bbe55SBarry Smith 
10882b1193eSBarry Smith #undef __FUNC__
109cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_1"></a>*/"MatMult_SeqMAIJ_1"
110cd3bbe55SBarry Smith int MatMult_SeqMAIJ_1(Mat A,Vec xx,Vec yy)
11182b1193eSBarry Smith {
112*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
113cd3bbe55SBarry Smith   int         ierr;
11482b1193eSBarry Smith   PetscFunctionBegin;
115cd3bbe55SBarry Smith   ierr = MatMult_SeqAIJ(b->AIJ,xx,yy);
116cd3bbe55SBarry Smith   PetscFunctionReturn(0);
11782b1193eSBarry Smith }
118cd3bbe55SBarry Smith #undef __FUNC__
119cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_1"></a>*/"MatMultTranspose_SeqMAIJ_1"
120cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_1(Mat A,Vec xx,Vec yy)
121cd3bbe55SBarry Smith {
122*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
123cd3bbe55SBarry Smith   int         ierr;
124cd3bbe55SBarry Smith   PetscFunctionBegin;
125cd3bbe55SBarry Smith   ierr = MatMultTranspose_SeqAIJ(b->AIJ,xx,yy);
126cd3bbe55SBarry Smith   PetscFunctionReturn(0);
127cd3bbe55SBarry Smith }
128cd3bbe55SBarry Smith #undef __FUNC__
129cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_1"></a>*/"MatMultAdd_SeqMAIJ_1"
130cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
131cd3bbe55SBarry Smith {
132*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
133cd3bbe55SBarry Smith   int         ierr;
134cd3bbe55SBarry Smith   PetscFunctionBegin;
135cd3bbe55SBarry Smith   ierr = MatMultAdd_SeqAIJ(b->AIJ,xx,yy,zz);
136cd3bbe55SBarry Smith   PetscFunctionReturn(0);
137cd3bbe55SBarry Smith }
138cd3bbe55SBarry Smith #undef __FUNC__
139cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_1"></a>*/"MatMultTransposeAdd_SeqMAIJ_1"
140cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
141cd3bbe55SBarry Smith {
142*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
143cd3bbe55SBarry Smith   int         ierr;
144cd3bbe55SBarry Smith   PetscFunctionBegin;
145cd3bbe55SBarry Smith   ierr = MatMultTransposeAdd_SeqAIJ(b->AIJ,xx,yy,zz);
14682b1193eSBarry Smith   PetscFunctionReturn(0);
14782b1193eSBarry Smith }
14882b1193eSBarry Smith 
149cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/
15082b1193eSBarry Smith #undef __FUNC__
151cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_2"></a>*/"MatMult_SeqMAIJ_2"
152cd3bbe55SBarry Smith int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
15382b1193eSBarry Smith {
154*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
155bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
156bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2;
157bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
158bcc973b7SBarry Smith   int         n,i,jrow,j;
15982b1193eSBarry Smith 
160bcc973b7SBarry Smith   PetscFunctionBegin;
161bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
162bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
163bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
164bcc973b7SBarry Smith   idx  = a->j;
165bcc973b7SBarry Smith   v    = a->a;
166bcc973b7SBarry Smith   ii   = a->i;
167bcc973b7SBarry Smith 
168bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
169bcc973b7SBarry Smith   idx  += shift;
170bcc973b7SBarry Smith   for (i=0; i<m; i++) {
171bcc973b7SBarry Smith     jrow = ii[i];
172bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
173bcc973b7SBarry Smith     sum1  = 0.0;
174bcc973b7SBarry Smith     sum2  = 0.0;
175bcc973b7SBarry Smith     for (j=0; j<n; j++) {
176bcc973b7SBarry Smith       sum1 += v[jrow]*x[2*idx[jrow]];
177bcc973b7SBarry Smith       sum2 += v[jrow]*x[2*idx[jrow]+1];
178bcc973b7SBarry Smith       jrow++;
179bcc973b7SBarry Smith      }
180bcc973b7SBarry Smith     y[2*i]   = sum1;
181bcc973b7SBarry Smith     y[2*i+1] = sum2;
182bcc973b7SBarry Smith   }
183bcc973b7SBarry Smith 
184bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*m);
185bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
186bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
18782b1193eSBarry Smith   PetscFunctionReturn(0);
18882b1193eSBarry Smith }
189bcc973b7SBarry Smith 
19082b1193eSBarry Smith #undef __FUNC__
191cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_2"></a>*/"MatMultTranspose_SeqMAIJ_2"
192cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy)
19382b1193eSBarry Smith {
194*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
195bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
196bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,zero = 0.0;
197bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
19882b1193eSBarry Smith 
199bcc973b7SBarry Smith   PetscFunctionBegin;
200bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
201bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
202bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
203bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
204bcc973b7SBarry Smith   for (i=0; i<m; i++) {
205bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
206bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
207bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
208bcc973b7SBarry Smith     alpha1 = x[2*i];
209bcc973b7SBarry Smith     alpha2 = x[2*i+1];
210bcc973b7SBarry Smith     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
211bcc973b7SBarry Smith   }
212bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*a->n);
213bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
214bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
21582b1193eSBarry Smith   PetscFunctionReturn(0);
21682b1193eSBarry Smith }
217bcc973b7SBarry Smith 
21882b1193eSBarry Smith #undef __FUNC__
219cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_2"></a>*/"MatMultAdd_SeqMAIJ_2"
220cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
22182b1193eSBarry Smith {
222*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
223bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
224bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2;
225bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
226bcc973b7SBarry Smith   int         n,i,jrow,j;
22782b1193eSBarry Smith 
228bcc973b7SBarry Smith   PetscFunctionBegin;
229f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
230bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
231f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
232bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
233bcc973b7SBarry Smith   idx  = a->j;
234bcc973b7SBarry Smith   v    = a->a;
235bcc973b7SBarry Smith   ii   = a->i;
236bcc973b7SBarry Smith 
237bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
238bcc973b7SBarry Smith   idx  += shift;
239bcc973b7SBarry Smith   for (i=0; i<m; i++) {
240bcc973b7SBarry Smith     jrow = ii[i];
241bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
242bcc973b7SBarry Smith     sum1  = 0.0;
243bcc973b7SBarry Smith     sum2  = 0.0;
244bcc973b7SBarry Smith     for (j=0; j<n; j++) {
245bcc973b7SBarry Smith       sum1 += v[jrow]*x[2*idx[jrow]];
246bcc973b7SBarry Smith       sum2 += v[jrow]*x[2*idx[jrow]+1];
247bcc973b7SBarry Smith       jrow++;
248bcc973b7SBarry Smith      }
249bcc973b7SBarry Smith     y[2*i]   += sum1;
250bcc973b7SBarry Smith     y[2*i+1] += sum2;
251bcc973b7SBarry Smith   }
252bcc973b7SBarry Smith 
253bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*m);
254bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
255f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
256bcc973b7SBarry Smith   PetscFunctionReturn(0);
25782b1193eSBarry Smith }
25882b1193eSBarry Smith #undef __FUNC__
259cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_2"></a>*/"MatMultTransposeAdd_SeqMAIJ_2"
260cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
26182b1193eSBarry Smith {
262*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
263bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
264bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2;
265bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
26682b1193eSBarry Smith 
267bcc973b7SBarry Smith   PetscFunctionBegin;
268f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
269bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
270f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
271bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
272bcc973b7SBarry Smith   for (i=0; i<m; i++) {
273bcc973b7SBarry Smith     idx   = a->j + a->i[i] + shift;
274bcc973b7SBarry Smith     v     = a->a + a->i[i] + shift;
275bcc973b7SBarry Smith     n     = a->i[i+1] - a->i[i];
276bcc973b7SBarry Smith     alpha1 = x[2*i];
277bcc973b7SBarry Smith     alpha2 = x[2*i+1];
278bcc973b7SBarry Smith     while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;}
279bcc973b7SBarry Smith   }
280bcc973b7SBarry Smith   PLogFlops(4*a->nz - 2*a->n);
281bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
282f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
283bcc973b7SBarry Smith   PetscFunctionReturn(0);
28482b1193eSBarry Smith }
285cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/
286bcc973b7SBarry Smith #undef __FUNC__
287bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_3"></a>*/"MatMult_SeqMAIJ_3"
288bcc973b7SBarry Smith int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
289bcc973b7SBarry Smith {
290*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
291bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
292bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3;
293bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
294bcc973b7SBarry Smith   int         n,i,jrow,j;
29582b1193eSBarry Smith 
296bcc973b7SBarry Smith   PetscFunctionBegin;
297bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
298bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
299bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
300bcc973b7SBarry Smith   idx  = a->j;
301bcc973b7SBarry Smith   v    = a->a;
302bcc973b7SBarry Smith   ii   = a->i;
303bcc973b7SBarry Smith 
304bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
305bcc973b7SBarry Smith   idx  += shift;
306bcc973b7SBarry Smith   for (i=0; i<m; i++) {
307bcc973b7SBarry Smith     jrow = ii[i];
308bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
309bcc973b7SBarry Smith     sum1  = 0.0;
310bcc973b7SBarry Smith     sum2  = 0.0;
311bcc973b7SBarry Smith     sum3  = 0.0;
312bcc973b7SBarry Smith     for (j=0; j<n; j++) {
313bcc973b7SBarry Smith       sum1 += v[jrow]*x[3*idx[jrow]];
314bcc973b7SBarry Smith       sum2 += v[jrow]*x[3*idx[jrow]+1];
315bcc973b7SBarry Smith       sum3 += v[jrow]*x[3*idx[jrow]+2];
316bcc973b7SBarry Smith       jrow++;
317bcc973b7SBarry Smith      }
318bcc973b7SBarry Smith     y[3*i]   = sum1;
319bcc973b7SBarry Smith     y[3*i+1] = sum2;
320bcc973b7SBarry Smith     y[3*i+2] = sum3;
321bcc973b7SBarry Smith   }
322bcc973b7SBarry Smith 
323bcc973b7SBarry Smith   PLogFlops(6*a->nz - 3*m);
324bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
325bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
326bcc973b7SBarry Smith   PetscFunctionReturn(0);
327bcc973b7SBarry Smith }
328bcc973b7SBarry Smith 
329bcc973b7SBarry Smith #undef __FUNC__
330bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_3"></a>*/"MatMultTranspose_SeqMAIJ_3"
331bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy)
332bcc973b7SBarry Smith {
333*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
334bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
335bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0;
336bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
337bcc973b7SBarry Smith 
338bcc973b7SBarry Smith   PetscFunctionBegin;
339bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
340bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
341bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
342bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
343bcc973b7SBarry Smith   for (i=0; i<m; i++) {
344bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
345bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
346bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
347bcc973b7SBarry Smith     alpha1 = x[3*i];
348bcc973b7SBarry Smith     alpha2 = x[3*i+1];
349bcc973b7SBarry Smith     alpha3 = x[3*i+2];
350bcc973b7SBarry Smith     while (n-->0) {
351bcc973b7SBarry Smith       y[3*(*idx)]   += alpha1*(*v);
352bcc973b7SBarry Smith       y[3*(*idx)+1] += alpha2*(*v);
353bcc973b7SBarry Smith       y[3*(*idx)+2] += alpha3*(*v);
354bcc973b7SBarry Smith       idx++; v++;
355bcc973b7SBarry Smith     }
356bcc973b7SBarry Smith   }
357bcc973b7SBarry Smith   PLogFlops(6*a->nz - 3*a->n);
358bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
359bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
360bcc973b7SBarry Smith   PetscFunctionReturn(0);
361bcc973b7SBarry Smith }
362bcc973b7SBarry Smith 
363bcc973b7SBarry Smith #undef __FUNC__
364bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_3"></a>*/"MatMultAdd_SeqMAIJ_3"
365bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
366bcc973b7SBarry Smith {
367*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
368bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
369bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3;
370bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
371bcc973b7SBarry Smith   int         n,i,jrow,j;
372bcc973b7SBarry Smith 
373bcc973b7SBarry Smith   PetscFunctionBegin;
374f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
375bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
376f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
377bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
378bcc973b7SBarry Smith   idx  = a->j;
379bcc973b7SBarry Smith   v    = a->a;
380bcc973b7SBarry Smith   ii   = a->i;
381bcc973b7SBarry Smith 
382bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
383bcc973b7SBarry Smith   idx  += shift;
384bcc973b7SBarry Smith   for (i=0; i<m; i++) {
385bcc973b7SBarry Smith     jrow = ii[i];
386bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
387bcc973b7SBarry Smith     sum1  = 0.0;
388bcc973b7SBarry Smith     sum2  = 0.0;
389bcc973b7SBarry Smith     sum3  = 0.0;
390bcc973b7SBarry Smith     for (j=0; j<n; j++) {
391bcc973b7SBarry Smith       sum1 += v[jrow]*x[3*idx[jrow]];
392bcc973b7SBarry Smith       sum2 += v[jrow]*x[3*idx[jrow]+1];
393bcc973b7SBarry Smith       sum3 += v[jrow]*x[3*idx[jrow]+2];
394bcc973b7SBarry Smith       jrow++;
395bcc973b7SBarry Smith      }
396bcc973b7SBarry Smith     y[3*i]   += sum1;
397bcc973b7SBarry Smith     y[3*i+1] += sum2;
398bcc973b7SBarry Smith     y[3*i+2] += sum3;
399bcc973b7SBarry Smith   }
400bcc973b7SBarry Smith 
401bcc973b7SBarry Smith   PLogFlops(6*a->nz);
402bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
403f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
404bcc973b7SBarry Smith   PetscFunctionReturn(0);
405bcc973b7SBarry Smith }
406bcc973b7SBarry Smith #undef __FUNC__
407bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_3"></a>*/"MatMultTransposeAdd_SeqMAIJ_3"
408bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
409bcc973b7SBarry Smith {
410*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
411bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
412bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3;
413bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
414bcc973b7SBarry Smith 
415bcc973b7SBarry Smith   PetscFunctionBegin;
416f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
417bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
418f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
419bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
420bcc973b7SBarry Smith   for (i=0; i<m; i++) {
421bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
422bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
423bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
424bcc973b7SBarry Smith     alpha1 = x[3*i];
425bcc973b7SBarry Smith     alpha2 = x[3*i+1];
426bcc973b7SBarry Smith     alpha3 = x[3*i+2];
427bcc973b7SBarry Smith     while (n-->0) {
428bcc973b7SBarry Smith       y[3*(*idx)]   += alpha1*(*v);
429bcc973b7SBarry Smith       y[3*(*idx)+1] += alpha2*(*v);
430bcc973b7SBarry Smith       y[3*(*idx)+2] += alpha3*(*v);
431bcc973b7SBarry Smith       idx++; v++;
432bcc973b7SBarry Smith     }
433bcc973b7SBarry Smith   }
434bcc973b7SBarry Smith   PLogFlops(6*a->nz);
435bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
436f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
437bcc973b7SBarry Smith   PetscFunctionReturn(0);
438bcc973b7SBarry Smith }
439bcc973b7SBarry Smith 
440bcc973b7SBarry Smith /* ------------------------------------------------------------------------------*/
441bcc973b7SBarry Smith #undef __FUNC__
442bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_4"></a>*/"MatMult_SeqMAIJ_4"
443bcc973b7SBarry Smith int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
444bcc973b7SBarry Smith {
445*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
446bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
447bcc973b7SBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4;
448bcc973b7SBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
449bcc973b7SBarry Smith   int         n,i,jrow,j;
450bcc973b7SBarry Smith 
451bcc973b7SBarry Smith   PetscFunctionBegin;
452bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
453bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
454bcc973b7SBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
455bcc973b7SBarry Smith   idx  = a->j;
456bcc973b7SBarry Smith   v    = a->a;
457bcc973b7SBarry Smith   ii   = a->i;
458bcc973b7SBarry Smith 
459bcc973b7SBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
460bcc973b7SBarry Smith   idx  += shift;
461bcc973b7SBarry Smith   for (i=0; i<m; i++) {
462bcc973b7SBarry Smith     jrow = ii[i];
463bcc973b7SBarry Smith     n    = ii[i+1] - jrow;
464bcc973b7SBarry Smith     sum1  = 0.0;
465bcc973b7SBarry Smith     sum2  = 0.0;
466bcc973b7SBarry Smith     sum3  = 0.0;
467bcc973b7SBarry Smith     sum4  = 0.0;
468bcc973b7SBarry Smith     for (j=0; j<n; j++) {
469bcc973b7SBarry Smith       sum1 += v[jrow]*x[4*idx[jrow]];
470bcc973b7SBarry Smith       sum2 += v[jrow]*x[4*idx[jrow]+1];
471bcc973b7SBarry Smith       sum3 += v[jrow]*x[4*idx[jrow]+2];
472bcc973b7SBarry Smith       sum4 += v[jrow]*x[4*idx[jrow]+3];
473bcc973b7SBarry Smith       jrow++;
474bcc973b7SBarry Smith      }
475bcc973b7SBarry Smith     y[4*i]   = sum1;
476bcc973b7SBarry Smith     y[4*i+1] = sum2;
477bcc973b7SBarry Smith     y[4*i+2] = sum3;
478bcc973b7SBarry Smith     y[4*i+3] = sum4;
479bcc973b7SBarry Smith   }
480bcc973b7SBarry Smith 
481bcc973b7SBarry Smith   PLogFlops(8*a->nz - 4*m);
482bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
483bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
484bcc973b7SBarry Smith   PetscFunctionReturn(0);
485bcc973b7SBarry Smith }
486bcc973b7SBarry Smith 
487bcc973b7SBarry Smith #undef __FUNC__
488bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_4"></a>*/"MatMultTranspose_SeqMAIJ_4"
489bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy)
490bcc973b7SBarry Smith {
491*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
492bcc973b7SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
493bcc973b7SBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0;
494bcc973b7SBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
495bcc973b7SBarry Smith 
496bcc973b7SBarry Smith   PetscFunctionBegin;
497bcc973b7SBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
498bcc973b7SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
499bcc973b7SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
500bcc973b7SBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
501bcc973b7SBarry Smith   for (i=0; i<m; i++) {
502bcc973b7SBarry Smith     idx    = a->j + a->i[i] + shift;
503bcc973b7SBarry Smith     v      = a->a + a->i[i] + shift;
504bcc973b7SBarry Smith     n      = a->i[i+1] - a->i[i];
505bcc973b7SBarry Smith     alpha1 = x[4*i];
506bcc973b7SBarry Smith     alpha2 = x[4*i+1];
507bcc973b7SBarry Smith     alpha3 = x[4*i+2];
508bcc973b7SBarry Smith     alpha4 = x[4*i+3];
509bcc973b7SBarry Smith     while (n-->0) {
510bcc973b7SBarry Smith       y[4*(*idx)]   += alpha1*(*v);
511bcc973b7SBarry Smith       y[4*(*idx)+1] += alpha2*(*v);
512bcc973b7SBarry Smith       y[4*(*idx)+2] += alpha3*(*v);
513bcc973b7SBarry Smith       y[4*(*idx)+3] += alpha4*(*v);
514bcc973b7SBarry Smith       idx++; v++;
515bcc973b7SBarry Smith     }
516bcc973b7SBarry Smith   }
517bcc973b7SBarry Smith   PLogFlops(8*a->nz - 4*a->n);
518bcc973b7SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
519bcc973b7SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
520bcc973b7SBarry Smith   PetscFunctionReturn(0);
521bcc973b7SBarry Smith }
522bcc973b7SBarry Smith 
523bcc973b7SBarry Smith #undef __FUNC__
524bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_4"></a>*/"MatMultAdd_SeqMAIJ_4"
525bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
526bcc973b7SBarry Smith {
527*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
528f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
529f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4;
530f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
531f9fae5adSBarry Smith   int         n,i,jrow,j;
532f9fae5adSBarry Smith 
533f9fae5adSBarry Smith   PetscFunctionBegin;
534f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
535f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
536f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
537f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
538f9fae5adSBarry Smith   idx  = a->j;
539f9fae5adSBarry Smith   v    = a->a;
540f9fae5adSBarry Smith   ii   = a->i;
541f9fae5adSBarry Smith 
542f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
543f9fae5adSBarry Smith   idx  += shift;
544f9fae5adSBarry Smith   for (i=0; i<m; i++) {
545f9fae5adSBarry Smith     jrow = ii[i];
546f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
547f9fae5adSBarry Smith     sum1  = 0.0;
548f9fae5adSBarry Smith     sum2  = 0.0;
549f9fae5adSBarry Smith     sum3  = 0.0;
550f9fae5adSBarry Smith     sum4  = 0.0;
551f9fae5adSBarry Smith     for (j=0; j<n; j++) {
552f9fae5adSBarry Smith       sum1 += v[jrow]*x[4*idx[jrow]];
553f9fae5adSBarry Smith       sum2 += v[jrow]*x[4*idx[jrow]+1];
554f9fae5adSBarry Smith       sum3 += v[jrow]*x[4*idx[jrow]+2];
555f9fae5adSBarry Smith       sum4 += v[jrow]*x[4*idx[jrow]+3];
556f9fae5adSBarry Smith       jrow++;
557f9fae5adSBarry Smith      }
558f9fae5adSBarry Smith     y[4*i]   += sum1;
559f9fae5adSBarry Smith     y[4*i+1] += sum2;
560f9fae5adSBarry Smith     y[4*i+2] += sum3;
561f9fae5adSBarry Smith     y[4*i+3] += sum4;
562f9fae5adSBarry Smith   }
563f9fae5adSBarry Smith 
564f9fae5adSBarry Smith   PLogFlops(8*a->nz - 4*m);
565f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
566f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
567f9fae5adSBarry Smith   PetscFunctionReturn(0);
568bcc973b7SBarry Smith }
569bcc973b7SBarry Smith #undef __FUNC__
570bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_4"></a>*/"MatMultTransposeAdd_SeqMAIJ_4"
571bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
572bcc973b7SBarry Smith {
573*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
574f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
575f9fae5adSBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4;
576f9fae5adSBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
577f9fae5adSBarry Smith 
578f9fae5adSBarry Smith   PetscFunctionBegin;
579f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
580f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
581f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
582f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
583f9fae5adSBarry Smith   for (i=0; i<m; i++) {
584f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
585f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
586f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
587f9fae5adSBarry Smith     alpha1 = x[4*i];
588f9fae5adSBarry Smith     alpha2 = x[4*i+1];
589f9fae5adSBarry Smith     alpha3 = x[4*i+2];
590f9fae5adSBarry Smith     alpha4 = x[4*i+3];
591f9fae5adSBarry Smith     while (n-->0) {
592f9fae5adSBarry Smith       y[4*(*idx)]   += alpha1*(*v);
593f9fae5adSBarry Smith       y[4*(*idx)+1] += alpha2*(*v);
594f9fae5adSBarry Smith       y[4*(*idx)+2] += alpha3*(*v);
595f9fae5adSBarry Smith       y[4*(*idx)+3] += alpha4*(*v);
596f9fae5adSBarry Smith       idx++; v++;
597f9fae5adSBarry Smith     }
598f9fae5adSBarry Smith   }
599f9fae5adSBarry Smith   PLogFlops(8*a->nz - 4*a->n);
600f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
601f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
602f9fae5adSBarry Smith   PetscFunctionReturn(0);
603f9fae5adSBarry Smith 
604f9fae5adSBarry Smith }
605f9fae5adSBarry Smith 
606f9fae5adSBarry Smith /* ------------------------------------------------------------------------------*/
607f9fae5adSBarry Smith #undef __FUNC__
608f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_5"></a>*/"MatMult_SeqMAIJ_5"
609f9fae5adSBarry Smith int MatMult_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
610f9fae5adSBarry Smith {
611*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
612f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
613f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
614f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
615f9fae5adSBarry Smith   int         n,i,jrow,j;
616f9fae5adSBarry Smith 
617f9fae5adSBarry Smith   PetscFunctionBegin;
618f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
619f9fae5adSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
620f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
621f9fae5adSBarry Smith   idx  = a->j;
622f9fae5adSBarry Smith   v    = a->a;
623f9fae5adSBarry Smith   ii   = a->i;
624f9fae5adSBarry Smith 
625f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
626f9fae5adSBarry Smith   idx  += shift;
627f9fae5adSBarry Smith   for (i=0; i<m; i++) {
628f9fae5adSBarry Smith     jrow = ii[i];
629f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
630f9fae5adSBarry Smith     sum1  = 0.0;
631f9fae5adSBarry Smith     sum2  = 0.0;
632f9fae5adSBarry Smith     sum3  = 0.0;
633f9fae5adSBarry Smith     sum4  = 0.0;
634f9fae5adSBarry Smith     sum5  = 0.0;
635f9fae5adSBarry Smith     for (j=0; j<n; j++) {
636f9fae5adSBarry Smith       sum1 += v[jrow]*x[5*idx[jrow]];
637f9fae5adSBarry Smith       sum2 += v[jrow]*x[5*idx[jrow]+1];
638f9fae5adSBarry Smith       sum3 += v[jrow]*x[5*idx[jrow]+2];
639f9fae5adSBarry Smith       sum4 += v[jrow]*x[5*idx[jrow]+3];
640f9fae5adSBarry Smith       sum5 += v[jrow]*x[5*idx[jrow]+4];
641f9fae5adSBarry Smith       jrow++;
642f9fae5adSBarry Smith      }
643f9fae5adSBarry Smith     y[5*i]   = sum1;
644f9fae5adSBarry Smith     y[5*i+1] = sum2;
645f9fae5adSBarry Smith     y[5*i+2] = sum3;
646f9fae5adSBarry Smith     y[5*i+3] = sum4;
647f9fae5adSBarry Smith     y[5*i+4] = sum5;
648f9fae5adSBarry Smith   }
649f9fae5adSBarry Smith 
650f9fae5adSBarry Smith   PLogFlops(10*a->nz - 5*m);
651f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
652f9fae5adSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
653f9fae5adSBarry Smith   PetscFunctionReturn(0);
654f9fae5adSBarry Smith }
655f9fae5adSBarry Smith 
656f9fae5adSBarry Smith #undef __FUNC__
657f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_5"></a>*/"MatMultTranspose_SeqMAIJ_5"
658f9fae5adSBarry Smith int MatMultTranspose_SeqMAIJ_5(Mat A,Vec xx,Vec yy)
659f9fae5adSBarry Smith {
660*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
661f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
662f9fae5adSBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,zero = 0.0;
663f9fae5adSBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
664f9fae5adSBarry Smith 
665f9fae5adSBarry Smith   PetscFunctionBegin;
666f9fae5adSBarry Smith   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
667f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
668f9fae5adSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
669f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
670f9fae5adSBarry Smith   for (i=0; i<m; i++) {
671f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
672f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
673f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
674f9fae5adSBarry Smith     alpha1 = x[5*i];
675f9fae5adSBarry Smith     alpha2 = x[5*i+1];
676f9fae5adSBarry Smith     alpha3 = x[5*i+2];
677f9fae5adSBarry Smith     alpha4 = x[5*i+3];
678f9fae5adSBarry Smith     alpha5 = x[5*i+4];
679f9fae5adSBarry Smith     while (n-->0) {
680f9fae5adSBarry Smith       y[5*(*idx)]   += alpha1*(*v);
681f9fae5adSBarry Smith       y[5*(*idx)+1] += alpha2*(*v);
682f9fae5adSBarry Smith       y[5*(*idx)+2] += alpha3*(*v);
683f9fae5adSBarry Smith       y[5*(*idx)+3] += alpha4*(*v);
684f9fae5adSBarry Smith       y[5*(*idx)+4] += alpha5*(*v);
685f9fae5adSBarry Smith       idx++; v++;
686f9fae5adSBarry Smith     }
687f9fae5adSBarry Smith   }
688f9fae5adSBarry Smith   PLogFlops(10*a->nz - 5*a->n);
689f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
690f9fae5adSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
691f9fae5adSBarry Smith   PetscFunctionReturn(0);
692f9fae5adSBarry Smith }
693f9fae5adSBarry Smith 
694f9fae5adSBarry Smith #undef __FUNC__
695f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_5"></a>*/"MatMultAdd_SeqMAIJ_5"
696f9fae5adSBarry Smith int MatMultAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
697f9fae5adSBarry Smith {
698*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
699f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
700f9fae5adSBarry Smith   Scalar      *x,*y,*v,sum1, sum2, sum3, sum4, sum5;
701f9fae5adSBarry Smith   int         ierr,m = a->m,*idx,shift = a->indexshift,*ii;
702f9fae5adSBarry Smith   int         n,i,jrow,j;
703f9fae5adSBarry Smith 
704f9fae5adSBarry Smith   PetscFunctionBegin;
705f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
706f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
707f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
708f9fae5adSBarry Smith   x    = x + shift;    /* shift for Fortran start by 1 indexing */
709f9fae5adSBarry Smith   idx  = a->j;
710f9fae5adSBarry Smith   v    = a->a;
711f9fae5adSBarry Smith   ii   = a->i;
712f9fae5adSBarry Smith 
713f9fae5adSBarry Smith   v    += shift; /* shift for Fortran start by 1 indexing */
714f9fae5adSBarry Smith   idx  += shift;
715f9fae5adSBarry Smith   for (i=0; i<m; i++) {
716f9fae5adSBarry Smith     jrow = ii[i];
717f9fae5adSBarry Smith     n    = ii[i+1] - jrow;
718f9fae5adSBarry Smith     sum1  = 0.0;
719f9fae5adSBarry Smith     sum2  = 0.0;
720f9fae5adSBarry Smith     sum3  = 0.0;
721f9fae5adSBarry Smith     sum4  = 0.0;
722f9fae5adSBarry Smith     sum5  = 0.0;
723f9fae5adSBarry Smith     for (j=0; j<n; j++) {
724f9fae5adSBarry Smith       sum1 += v[jrow]*x[5*idx[jrow]];
725f9fae5adSBarry Smith       sum2 += v[jrow]*x[5*idx[jrow]+1];
726f9fae5adSBarry Smith       sum3 += v[jrow]*x[5*idx[jrow]+2];
727f9fae5adSBarry Smith       sum4 += v[jrow]*x[5*idx[jrow]+3];
728f9fae5adSBarry Smith       sum5 += v[jrow]*x[5*idx[jrow]+4];
729f9fae5adSBarry Smith       jrow++;
730f9fae5adSBarry Smith      }
731f9fae5adSBarry Smith     y[5*i]   += sum1;
732f9fae5adSBarry Smith     y[5*i+1] += sum2;
733f9fae5adSBarry Smith     y[5*i+2] += sum3;
734f9fae5adSBarry Smith     y[5*i+3] += sum4;
735f9fae5adSBarry Smith     y[5*i+4] += sum5;
736f9fae5adSBarry Smith   }
737f9fae5adSBarry Smith 
738f9fae5adSBarry Smith   PLogFlops(10*a->nz);
739f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
740f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
741f9fae5adSBarry Smith   PetscFunctionReturn(0);
742f9fae5adSBarry Smith }
743f9fae5adSBarry Smith 
744f9fae5adSBarry Smith #undef __FUNC__
745f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_5"></a>*/"MatMultTransposeAdd_SeqMAIJ_5"
746f9fae5adSBarry Smith int MatMultTransposeAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
747f9fae5adSBarry Smith {
748*4cb9d9b8SBarry Smith   Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data;
749f9fae5adSBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)b->AIJ->data;
750f9fae5adSBarry Smith   Scalar      *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5;
751f9fae5adSBarry Smith   int         ierr,m = a->m,n,i,*idx,shift = a->indexshift;
752f9fae5adSBarry Smith 
753f9fae5adSBarry Smith   PetscFunctionBegin;
754f9fae5adSBarry Smith   if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);}
755f9fae5adSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
756f9fae5adSBarry Smith   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
757f9fae5adSBarry Smith   y = y + shift; /* shift for Fortran start by 1 indexing */
758f9fae5adSBarry Smith   for (i=0; i<m; i++) {
759f9fae5adSBarry Smith     idx    = a->j + a->i[i] + shift;
760f9fae5adSBarry Smith     v      = a->a + a->i[i] + shift;
761f9fae5adSBarry Smith     n      = a->i[i+1] - a->i[i];
762f9fae5adSBarry Smith     alpha1 = x[5*i];
763f9fae5adSBarry Smith     alpha2 = x[5*i+1];
764f9fae5adSBarry Smith     alpha3 = x[5*i+2];
765f9fae5adSBarry Smith     alpha4 = x[5*i+3];
766f9fae5adSBarry Smith     alpha5 = x[5*i+4];
767f9fae5adSBarry Smith     while (n-->0) {
768f9fae5adSBarry Smith       y[5*(*idx)]   += alpha1*(*v);
769f9fae5adSBarry Smith       y[5*(*idx)+1] += alpha2*(*v);
770f9fae5adSBarry Smith       y[5*(*idx)+2] += alpha3*(*v);
771f9fae5adSBarry Smith       y[5*(*idx)+3] += alpha4*(*v);
772f9fae5adSBarry Smith       y[5*(*idx)+4] += alpha5*(*v);
773f9fae5adSBarry Smith       idx++; v++;
774f9fae5adSBarry Smith     }
775f9fae5adSBarry Smith   }
776f9fae5adSBarry Smith   PLogFlops(10*a->nz);
777f9fae5adSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
778f9fae5adSBarry Smith   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
779f9fae5adSBarry Smith   PetscFunctionReturn(0);
780bcc973b7SBarry Smith }
781bcc973b7SBarry Smith 
782f4a53059SBarry Smith /*===================================================================================*/
783f4a53059SBarry Smith EXTERN int MatMult_MPIAIJ(Mat,Vec,Vec);
784f4a53059SBarry Smith EXTERN int MatMultAdd_MPIAIJ(Mat,Vec,Vec,Vec);
785f4a53059SBarry Smith EXTERN int MatMultTranspose_MPIAIJ(Mat,Vec,Vec);
786f4a53059SBarry Smith EXTERN int MatMultTransposeAdd_MPIAIJ(Mat,Vec,Vec,Vec);
787f4a53059SBarry Smith 
788f4a53059SBarry Smith #undef __FUNC__
789f4a53059SBarry Smith #define __FUNC__ /*<a name="MatMult_MPIMAIJ_1"></a>*/"MatMult_MPIMAIJ_1"
790f4a53059SBarry Smith int MatMult_MPIMAIJ_1(Mat A,Vec xx,Vec yy)
791f4a53059SBarry Smith {
792*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
793f4a53059SBarry Smith   int         ierr;
794f4a53059SBarry Smith   PetscFunctionBegin;
795f4a53059SBarry Smith   ierr = MatMult_MPIAIJ(b->AIJ,xx,yy);
796f4a53059SBarry Smith   PetscFunctionReturn(0);
797f4a53059SBarry Smith }
798f4a53059SBarry Smith #undef __FUNC__
799f4a53059SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_MPIMAIJ_1"></a>*/"MatMultTranspose_MPIMAIJ_1"
800f4a53059SBarry Smith int MatMultTranspose_MPIMAIJ_1(Mat A,Vec xx,Vec yy)
801f4a53059SBarry Smith {
802*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
803f4a53059SBarry Smith   int         ierr;
804f4a53059SBarry Smith   PetscFunctionBegin;
805f4a53059SBarry Smith   ierr = MatMultTranspose_MPIAIJ(b->AIJ,xx,yy);
806f4a53059SBarry Smith   PetscFunctionReturn(0);
807f4a53059SBarry Smith }
808f4a53059SBarry Smith #undef __FUNC__
809f4a53059SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_MPIMAIJ_1"></a>*/"MatMultAdd_MPIMAIJ_1"
810f4a53059SBarry Smith int MatMultAdd_MPIMAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
811f4a53059SBarry Smith {
812*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
813f4a53059SBarry Smith   int         ierr;
814f4a53059SBarry Smith   PetscFunctionBegin;
815f4a53059SBarry Smith   ierr = MatMultAdd_MPIAIJ(b->AIJ,xx,yy,zz);
816f4a53059SBarry Smith   PetscFunctionReturn(0);
817f4a53059SBarry Smith }
818f4a53059SBarry Smith #undef __FUNC__
819f4a53059SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_MPIMAIJ_1"></a>*/"MatMultTransposeAdd_MPIMAIJ_1"
820f4a53059SBarry Smith int MatMultTransposeAdd_MPIMAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
821f4a53059SBarry Smith {
822*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
823f4a53059SBarry Smith   int         ierr;
824f4a53059SBarry Smith   PetscFunctionBegin;
825f4a53059SBarry Smith   ierr = MatMultTransposeAdd_MPIAIJ(b->AIJ,xx,yy,zz);
826f4a53059SBarry Smith   PetscFunctionReturn(0);
827f4a53059SBarry Smith }
828f4a53059SBarry Smith 
829f4a53059SBarry Smith /* ---------------------------------------------------------------------------------- */
830f4a53059SBarry Smith #undef __FUNC__
831f4a53059SBarry Smith #define __FUNC__ /*<a name="MatMult_MPIMAIJ_dof"></a>*/"MatMult_MPIMAIJ_dof"
832f4a53059SBarry Smith int MatMult_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
833f4a53059SBarry Smith {
834*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
835f4a53059SBarry Smith   int         ierr;
836f4a53059SBarry Smith   PetscFunctionBegin;
837f4a53059SBarry Smith 
838f4a53059SBarry Smith   /* start the scatter */
839f4a53059SBarry Smith   ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
840*4cb9d9b8SBarry Smith   ierr = (*b->AIJ->ops->mult)(b->AIJ,xx,yy);CHKERRQ(ierr);
841f4a53059SBarry Smith   ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
842*4cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr);
843f4a53059SBarry Smith   PetscFunctionReturn(0);
844f4a53059SBarry Smith }
845f4a53059SBarry Smith 
846*4cb9d9b8SBarry Smith #undef __FUNC__
847*4cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_MPIMAIJ_dof"></a>*/"MatMultTranspose_MPIMAIJ_dof"
848*4cb9d9b8SBarry Smith int MatMultTranspose_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
849*4cb9d9b8SBarry Smith {
850*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
851*4cb9d9b8SBarry Smith   int         ierr;
852*4cb9d9b8SBarry Smith   PetscFunctionBegin;
853*4cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr);
854*4cb9d9b8SBarry Smith   ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
855*4cb9d9b8SBarry Smith   ierr = (*b->AIJ->ops->multtranspose)(b->AIJ,xx,yy);CHKERRQ(ierr);
856*4cb9d9b8SBarry Smith   ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
857*4cb9d9b8SBarry Smith   PetscFunctionReturn(0);
858*4cb9d9b8SBarry Smith }
859*4cb9d9b8SBarry Smith 
860*4cb9d9b8SBarry Smith #undef __FUNC__
861*4cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_MPIMAIJ_dof"></a>*/"MatMultAdd_MPIMAIJ_dof"
862*4cb9d9b8SBarry Smith int MatMultAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
863*4cb9d9b8SBarry Smith {
864*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
865*4cb9d9b8SBarry Smith   int         ierr;
866*4cb9d9b8SBarry Smith   PetscFunctionBegin;
867*4cb9d9b8SBarry Smith 
868*4cb9d9b8SBarry Smith   /* start the scatter */
869*4cb9d9b8SBarry Smith   ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
870*4cb9d9b8SBarry Smith   ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,yy,yy);CHKERRQ(ierr);
871*4cb9d9b8SBarry Smith   ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr);
872*4cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr);
873*4cb9d9b8SBarry Smith   PetscFunctionReturn(0);
874*4cb9d9b8SBarry Smith }
875*4cb9d9b8SBarry Smith 
876*4cb9d9b8SBarry Smith #undef __FUNC__
877*4cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_MPIMAIJ_dof"></a>*/"MatMultTransposeAdd_MPIMAIJ_dof"
878*4cb9d9b8SBarry Smith int MatMultTransposeAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy)
879*4cb9d9b8SBarry Smith {
880*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data;
881*4cb9d9b8SBarry Smith   int         ierr;
882*4cb9d9b8SBarry Smith   PetscFunctionBegin;
883*4cb9d9b8SBarry Smith   ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr);
884*4cb9d9b8SBarry Smith   ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
885*4cb9d9b8SBarry Smith   ierr = (*b->AIJ->ops->multtransposeadd)(b->AIJ,xx,yy,yy);CHKERRQ(ierr);
886*4cb9d9b8SBarry Smith   ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr);
887*4cb9d9b8SBarry Smith   PetscFunctionReturn(0);
888*4cb9d9b8SBarry Smith }
889*4cb9d9b8SBarry Smith 
890*4cb9d9b8SBarry Smith 
891bcc973b7SBarry Smith /* ---------------------------------------------------------------------------------- */
89282b1193eSBarry Smith #undef __FUNC__
893cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatCreateMAIJ"></a>*/"MatCreateMAIJ"
894cd3bbe55SBarry Smith int MatCreateMAIJ(Mat A,int dof,Mat *maij)
89582b1193eSBarry Smith {
896f4a53059SBarry Smith   int         ierr,size,n;
897*4cb9d9b8SBarry Smith   Mat_MPIMAIJ *b;
89882b1193eSBarry Smith   Mat         B;
89982b1193eSBarry Smith 
90082b1193eSBarry Smith   PetscFunctionBegin;
901cd3bbe55SBarry Smith   ierr = MATCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr);
902f4a53059SBarry Smith   ierr = MatSetType(B,MATMAIJ);CHKERRQ(ierr);
90382b1193eSBarry Smith 
904cd3bbe55SBarry Smith   B->assembled    = PETSC_TRUE;
905*4cb9d9b8SBarry Smith   b = (Mat_MPIMAIJ*)B->data;
90682b1193eSBarry Smith 
907cd3bbe55SBarry Smith   b->dof = dof;
908cd3bbe55SBarry Smith   ierr   = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
909f4a53059SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
910f4a53059SBarry Smith   if (size == 1) {
911*4cb9d9b8SBarry Smith     B->ops->destroy = MatDestroy_SeqMAIJ;
912*4cb9d9b8SBarry Smith     b->AIJ = A;
913cd3bbe55SBarry Smith     if (dof == 1) {
914cd3bbe55SBarry Smith       B->ops->mult             = MatMult_SeqMAIJ_1;
915cd3bbe55SBarry Smith       B->ops->multadd          = MatMultAdd_SeqMAIJ_1;
916cd3bbe55SBarry Smith       B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_1;
917cd3bbe55SBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_1;
918cd3bbe55SBarry Smith     } else if (dof == 2) {
919cd3bbe55SBarry Smith       B->ops->mult             = MatMult_SeqMAIJ_2;
920cd3bbe55SBarry Smith       B->ops->multadd          = MatMultAdd_SeqMAIJ_2;
921cd3bbe55SBarry Smith       B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_2;
922cd3bbe55SBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2;
923bcc973b7SBarry Smith     } else if (dof == 3) {
924bcc973b7SBarry Smith       B->ops->mult             = MatMult_SeqMAIJ_3;
925bcc973b7SBarry Smith       B->ops->multadd          = MatMultAdd_SeqMAIJ_3;
926bcc973b7SBarry Smith       B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_3;
927bcc973b7SBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3;
928bcc973b7SBarry Smith     } else if (dof == 4) {
929bcc973b7SBarry Smith       B->ops->mult             = MatMult_SeqMAIJ_4;
930bcc973b7SBarry Smith       B->ops->multadd          = MatMultAdd_SeqMAIJ_4;
931bcc973b7SBarry Smith       B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_4;
932bcc973b7SBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4;
933f9fae5adSBarry Smith     } else if (dof == 5) {
934f9fae5adSBarry Smith       B->ops->mult             = MatMult_SeqMAIJ_5;
935f9fae5adSBarry Smith       B->ops->multadd          = MatMultAdd_SeqMAIJ_5;
936f9fae5adSBarry Smith       B->ops->multtranspose    = MatMultTranspose_SeqMAIJ_5;
937f9fae5adSBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5;
93882b1193eSBarry Smith     } else {
939cd3bbe55SBarry Smith       SETERRQ1(1,1,"Cannot handle a dof of %d\n",dof);
94082b1193eSBarry Smith     }
941f4a53059SBarry Smith   } else {
942*4cb9d9b8SBarry Smith     b->A = A;
943*4cb9d9b8SBarry Smith     B->ops->destroy = MatDestroy_MPIMAIJ;
944f4a53059SBarry Smith     if (dof == 1) {
945*4cb9d9b8SBarry Smith       b->AIJ = A;
946f4a53059SBarry Smith       B->ops->mult             = MatMult_MPIMAIJ_1;
947f4a53059SBarry Smith       B->ops->multadd          = MatMultAdd_MPIMAIJ_1;
948f4a53059SBarry Smith       B->ops->multtranspose    = MatMultTranspose_MPIMAIJ_1;
949f4a53059SBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_1;
950f4a53059SBarry Smith     } else {
951f4a53059SBarry Smith       Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data;
952f4a53059SBarry Smith       IS         from,to;
953f4a53059SBarry Smith       Vec        gvec;
954f4a53059SBarry Smith       int        *garray,i;
955f4a53059SBarry Smith 
956*4cb9d9b8SBarry Smith       ierr = MatCreateMAIJ(mpiaij->A,dof,&b->AIJ);CHKERRQ(ierr);
957*4cb9d9b8SBarry Smith       ierr = MatCreateMAIJ(mpiaij->B,dof,&b->OAIJ);CHKERRQ(ierr);
958*4cb9d9b8SBarry Smith 
959f4a53059SBarry Smith       ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr);
960f4a53059SBarry Smith       ierr = VecCreateSeq(PETSC_COMM_SELF,n*dof,&b->w);CHKERRQ(ierr);
961f4a53059SBarry Smith 
962f4a53059SBarry Smith       /* create two temporary Index sets for build scatter gather */
963f4a53059SBarry Smith       garray = (int*)PetscMalloc((n+1)*sizeof(int));CHKPTRQ(garray);
964f4a53059SBarry Smith       for (i=0; i<n; i++) garray[i] = dof*mpiaij->garray[i];
965f4a53059SBarry Smith       ierr = ISCreateBlock(A->comm,dof,n,garray,&from);CHKERRQ(ierr);
966f4a53059SBarry Smith       ierr = PetscFree(garray);CHKERRQ(ierr);
967f4a53059SBarry Smith       ierr = ISCreateStride(PETSC_COMM_SELF,n*dof,0,1,&to);CHKERRQ(ierr);
968f4a53059SBarry Smith 
969f4a53059SBarry Smith       /* create temporary global vector to generate scatter context */
970f4a53059SBarry Smith       ierr = VecCreateMPI(A->comm,dof*A->n,dof*A->N,&gvec);CHKERRQ(ierr);
971f4a53059SBarry Smith 
972f4a53059SBarry Smith       /* generate the scatter context */
973f4a53059SBarry Smith       ierr = VecScatterCreate(gvec,from,b->w,to,&b->ctx);CHKERRQ(ierr);
974f4a53059SBarry Smith 
975f4a53059SBarry Smith       ierr = ISDestroy(from);CHKERRQ(ierr);
976f4a53059SBarry Smith       ierr = ISDestroy(to);CHKERRQ(ierr);
977f4a53059SBarry Smith       ierr = VecDestroy(gvec);CHKERRQ(ierr);
978f4a53059SBarry Smith 
979f4a53059SBarry Smith       B->ops->mult             = MatMult_MPIMAIJ_dof;
980*4cb9d9b8SBarry Smith       B->ops->multtranspose    = MatMultTranspose_MPIMAIJ_dof;
981*4cb9d9b8SBarry Smith       B->ops->multadd          = MatMultAdd_MPIMAIJ_dof;
982*4cb9d9b8SBarry Smith       B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof;
983f4a53059SBarry Smith     }
984f4a53059SBarry Smith   }
985cd3bbe55SBarry Smith   *maij = B;
98682b1193eSBarry Smith   PetscFunctionReturn(0);
98782b1193eSBarry Smith }
98882b1193eSBarry Smith 
98982b1193eSBarry Smith 
99082b1193eSBarry Smith 
99182b1193eSBarry Smith 
99282b1193eSBarry Smith 
99382b1193eSBarry Smith 
99482b1193eSBarry Smith 
99582b1193eSBarry Smith 
99682b1193eSBarry Smith 
99782b1193eSBarry Smith 
99882b1193eSBarry Smith 
99982b1193eSBarry Smith 
1000