xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision cf742fcc470053b8746f003f5618f6b40d3ea504)
18a729477SBarry Smith 
23369ce9aSBarry Smith #if !defined(__MPIAIJ_H)
33369ce9aSBarry Smith #define __MPIAIJ_H
43369ce9aSBarry Smith 
5c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
68a729477SBarry Smith 
790431a8fSHong Zhang typedef struct { /* used by MatCreateMPIAIJSumSeqAIJ for reusing the merged matrix */
826283091SBarry Smith   PetscLayout rowmap;
9de0260b3SHong Zhang   PetscInt    **buf_ri,**buf_rj;
100febcb4bSHong Zhang   PetscMPIInt *len_s,*len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
110febcb4bSHong Zhang   PetscMPIInt nsend,nrecv;
12fc08c53fSHong Zhang   PetscInt    *bi,*bj;    /* i and j array of the local portion of mpi C (matrix product) - rename to ci, cj! */
13de0260b3SHong Zhang   PetscInt    *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
14dce485f0SBarry Smith   PetscErrorCode (*destroy)(Mat);
15dce485f0SBarry Smith   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
16b90dcfe3SHong Zhang } Mat_Merge_SeqsToMPI;
17b90dcfe3SHong Zhang 
18a1a4e84aSHong Zhang typedef struct { /* used by MatPtAP_MPIAIJ_MPIAIJ() and MatMatMult_MPIAIJ_MPIAIJ() */
19b7f45c76SHong Zhang   PetscInt    *startsj_s,*startsj_r;    /* used by MatGetBrowsOfAoCols_MPIAIJ */
20b7f45c76SHong Zhang   PetscScalar *bufa;                    /* used by MatGetBrowsOfAoCols_MPIAIJ */
21a1a4e84aSHong Zhang   Mat         P_loc,P_oth;     /* partial B_seq -- intend to replace B_seq */
22a1a4e84aSHong Zhang   PetscInt    *api,*apj;       /* symbolic i and j arrays of the local product A_loc*B_seq */
23445158ffSHong Zhang   PetscScalar *apv;
24b7f45c76SHong Zhang   MatReuse    reuse;           /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */
25598bc09dSHong Zhang   PetscScalar *apa;            /* tmp array for store a row of A*P used in MatMatMult() */
26c5af039cSHong Zhang   Mat         A_loc;           /* used by MatTransposeMatMult(), contains api and apj */
27ab1b013aSHong Zhang   Mat         Pt;              /* used by MatTransposeMatMult(), Pt = P^T */
28414894bdSHong Zhang   PetscBool   scalable;        /* flag determines scalable or non-scalable implementation */
299ce11a7cSHong Zhang   Mat         Rd,Ro,AP_loc,C_loc,C_oth;
30f8487c73SHong Zhang 
31f8487c73SHong Zhang   Mat_Merge_SeqsToMPI *merge;
32f8487c73SHong Zhang   PetscErrorCode (*destroy)(Mat);
334ae0847bSHong Zhang   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
34f8487c73SHong Zhang } Mat_PtAPMPI;
35f8487c73SHong Zhang 
36f8487c73SHong Zhang typedef struct {
37f8487c73SHong Zhang   Mat A,B;                             /* local submatrices: A (diag part),
38f8487c73SHong Zhang                                            B (off-diag part) */
39f8487c73SHong Zhang   PetscMPIInt size;                     /* size of communicator */
40f8487c73SHong Zhang   PetscMPIInt rank;                     /* rank of proc in communicator */
41f8487c73SHong Zhang 
42f8487c73SHong Zhang   /* The following variables are used for matrix assembly */
43f8487c73SHong Zhang   PetscBool   donotstash;               /* PETSC_TRUE if off processor entries dropped */
44f8487c73SHong Zhang   MPI_Request *send_waits;              /* array of send requests */
45f8487c73SHong Zhang   MPI_Request *recv_waits;              /* array of receive requests */
46f8487c73SHong Zhang   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
47f8487c73SHong Zhang   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
48f8487c73SHong Zhang   PetscInt    rmax;                     /* maximum message length */
49f8487c73SHong Zhang #if defined(PETSC_USE_CTABLE)
50f8487c73SHong Zhang   PetscTable colmap;
51f8487c73SHong Zhang #else
52f8487c73SHong Zhang   PetscInt *colmap;                     /* local col number of off-diag col */
53f8487c73SHong Zhang #endif
54f8487c73SHong Zhang   PetscInt *garray;                     /* global index of all off-processor columns */
55f8487c73SHong Zhang 
56f8487c73SHong Zhang   /* The following variables are used for matrix-vector products */
57f8487c73SHong Zhang   Vec        lvec;                 /* local vector */
58f8487c73SHong Zhang   Vec        diag;
59f8487c73SHong Zhang   VecScatter Mvctx;                /* scatter context for vector */
60f8487c73SHong Zhang   PetscBool  roworiented;          /* if true, row-oriented input, default true */
61f8487c73SHong Zhang 
62f8487c73SHong Zhang   /* The following variables are for MatGetRow() */
63f8487c73SHong Zhang   PetscInt    *rowindices;         /* column indices for row */
64f8487c73SHong Zhang   PetscScalar *rowvalues;          /* nonzero values in row */
65f8487c73SHong Zhang   PetscBool   getrowactive;        /* indicates MatGetRow(), not restored */
66f8487c73SHong Zhang 
67f8487c73SHong Zhang   /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation  and nonzero pattern */
68f8487c73SHong Zhang   PetscInt *ld;                    /* number of entries per row left of diagona block */
69f8487c73SHong Zhang 
70a1a4e84aSHong Zhang   /* Used by MatMatMult() and MatPtAP() */
71f8487c73SHong Zhang   Mat_PtAPMPI *ptap;
72ca45077fSPaul Mullowney 
73f996eeb8SHong Zhang   /* used by MatMatMatMult() */
74f996eeb8SHong Zhang   Mat_MatMatMatMult *matmatmatmult;
75f996eeb8SHong Zhang 
76bbf3fe20SPaul Mullowney   /* Used by MPICUSP and MPICUSPARSE classes */
77bbf3fe20SPaul Mullowney   void * spptr;
785cc03489SHong Zhang 
79f8487c73SHong Zhang } Mat_MPIAIJ;
80f8487c73SHong Zhang 
818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
822327d61dSSatish Balay 
835a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
845a576424SJed Brown PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
855a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*);
865a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
87b1b1104fSBarry Smith PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat,PetscInt,IS [],PetscInt);
8893dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
89f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
905a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
9153dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
925a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_All(Mat,MatGetSubMatrixOption,MatReuse,Mat *[]);
9353dd7562SDmitry Karpeev 
945494a064SHong Zhang 
955a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*);
965a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_Private (Mat,IS,IS,PetscInt,MatReuse,Mat*);
973233f0ddSHong Zhang PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*);
98d6037b41SHong Zhang 
995a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer);
100e3225e9dSHong Zhang PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1025a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1030fc8cf34SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
104b2405163SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1055a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
106f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
107f8487c73SHong Zhang 
1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
1095a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*);
1105a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);
111f996eeb8SHong Zhang 
1125a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1135a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1145a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
1150d3441aeSHong Zhang 
1168403a639SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat,Mat,PetscReal,Mat*);
1178403a639SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat,Mat,Mat);
1180d3441aeSHong Zhang 
1195a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat);
1205a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
121c5af039cSHong Zhang 
1225a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
1235a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
1245a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat);
1255a576424SJed Brown PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*);
1265a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);
127187b3c17SHong Zhang 
1285a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1292bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
1306da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1316da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
1322bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
133c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
1348949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
1358949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
1368949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
1375a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);
138c0d3702cSSatish Balay 
1394416b707SBarry Smith PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat);
14010e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1413a7fca6bSBarry Smith 
142570b7f6dSBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
1435a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*);
14497304618SKris Buschelman #endif
1455a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec);
1465a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*);
14797304618SKris Buschelman 
148001ddc4fSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
149001ddc4fSHong Zhang 
15011bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*);
1517087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);
1523a7fca6bSBarry Smith 
15353dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*);
15453dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat);
15553dd7562SDmitry Karpeev 
156*cf742fccSHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
157*cf742fccSHong Zhang #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \
158*cf742fccSHong Zhang {\
159*cf742fccSHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ;      \
160*cf742fccSHong Zhang   PetscScalar *_aa,_valtmp,*_pa;                             \
161*cf742fccSHong Zhang   _apJ = apj + api[i];\
162*cf742fccSHong Zhang   /* diagonal portion of A */\
163*cf742fccSHong Zhang   _ai  = ad->i;\
164*cf742fccSHong Zhang   _anz = _ai[i+1] - _ai[i];\
165*cf742fccSHong Zhang   _aj  = ad->j + _ai[i];\
166*cf742fccSHong Zhang   _aa  = ad->a + _ai[i];\
167*cf742fccSHong Zhang   for (_j=0; _j<_anz; _j++) {\
168*cf742fccSHong Zhang     _row = _aj[_j]; \
169*cf742fccSHong Zhang     _pi  = p_loc->i;                                 \
170*cf742fccSHong Zhang     _pnz = _pi[_row+1] - _pi[_row];         \
171*cf742fccSHong Zhang     _pj  = p_loc->j + _pi[_row];                 \
172*cf742fccSHong Zhang     _pa  = p_loc->a + _pi[_row];                 \
173*cf742fccSHong Zhang     /* perform sparse axpy */                    \
174*cf742fccSHong Zhang     _valtmp = _aa[_j];                           \
175*cf742fccSHong Zhang     _nextp  = 0; \
176*cf742fccSHong Zhang     for (_k=0; _nextp<_pnz; _k++) {                    \
177*cf742fccSHong Zhang       if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */   \
178*cf742fccSHong Zhang         apa[_k] += _valtmp*_pa[_nextp++];                                \
179*cf742fccSHong Zhang       } \
180*cf742fccSHong Zhang     }                                           \
181*cf742fccSHong Zhang     PetscLogFlops(2.0*_pnz);                    \
182*cf742fccSHong Zhang   }                                             \
183*cf742fccSHong Zhang   /* off-diagonal portion of A */               \
184*cf742fccSHong Zhang   _ai  = ao->i;\
185*cf742fccSHong Zhang   _anz = _ai[i+1] - _ai[i];                     \
186*cf742fccSHong Zhang   _aj  = ao->j + _ai[i];                         \
187*cf742fccSHong Zhang   _aa  = ao->a + _ai[i];                         \
188*cf742fccSHong Zhang   for (_j=0; _j<_anz; _j++) {                      \
189*cf742fccSHong Zhang     _row = _aj[_j];    \
190*cf742fccSHong Zhang     _pi  = p_oth->i;                         \
191*cf742fccSHong Zhang     _pnz = _pi[_row+1] - _pi[_row];          \
192*cf742fccSHong Zhang     _pj  = p_oth->j + _pi[_row];                  \
193*cf742fccSHong Zhang     _pa  = p_oth->a + _pi[_row];                  \
194*cf742fccSHong Zhang     /* perform sparse axpy */                     \
195*cf742fccSHong Zhang     _valtmp = _aa[_j];                             \
196*cf742fccSHong Zhang     _nextp  = 0; \
197*cf742fccSHong Zhang     for (_k=0; _nextp<_pnz; _k++) {                     \
198*cf742fccSHong Zhang       if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
199*cf742fccSHong Zhang         apa[_k] += _valtmp*_pa[_nextp++];                       \
200*cf742fccSHong Zhang       }                                                     \
201*cf742fccSHong Zhang     }                                            \
202*cf742fccSHong Zhang     PetscLogFlops(2.0*_pnz);                     \
203*cf742fccSHong Zhang   } \
204*cf742fccSHong Zhang }
205*cf742fccSHong Zhang 
206e497d3c8SHong Zhang #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
207e497d3c8SHong Zhang {\
208e497d3c8SHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;      \
209e497d3c8SHong Zhang   PetscScalar *_aa,_valtmp,*_pa;                             \
210e497d3c8SHong Zhang   /* diagonal portion of A */\
211e497d3c8SHong Zhang   _ai  = ad->i;\
212e497d3c8SHong Zhang   _anz = _ai[i+1] - _ai[i];\
213e497d3c8SHong Zhang   _aj  = ad->j + _ai[i];\
214e497d3c8SHong Zhang   _aa  = ad->a + _ai[i];\
215e497d3c8SHong Zhang   for (_j=0; _j<_anz; _j++) {\
216e497d3c8SHong Zhang     _row = _aj[_j]; \
217e497d3c8SHong Zhang     _pi  = p_loc->i;                                 \
218e497d3c8SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];         \
219e497d3c8SHong Zhang     _pj  = p_loc->j + _pi[_row];                 \
220e497d3c8SHong Zhang     _pa  = p_loc->a + _pi[_row];                 \
221e497d3c8SHong Zhang     /* perform dense axpy */                    \
222e497d3c8SHong Zhang     _valtmp = _aa[_j];                           \
223e497d3c8SHong Zhang     for (_k=0; _k<_pnz; _k++) {                    \
224e497d3c8SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];               \
225e497d3c8SHong Zhang     }                                           \
226e497d3c8SHong Zhang     PetscLogFlops(2.0*_pnz);                    \
227e497d3c8SHong Zhang   }                                             \
228e497d3c8SHong Zhang   /* off-diagonal portion of A */               \
229e497d3c8SHong Zhang   _ai  = ao->i;\
230e497d3c8SHong Zhang   _anz = _ai[i+1] - _ai[i];                     \
231e497d3c8SHong Zhang   _aj  = ao->j + _ai[i];                         \
232e497d3c8SHong Zhang   _aa  = ao->a + _ai[i];                         \
233e497d3c8SHong Zhang   for (_j=0; _j<_anz; _j++) {                      \
234e497d3c8SHong Zhang     _row = _aj[_j];    \
235e497d3c8SHong Zhang     _pi  = p_oth->i;                         \
236e497d3c8SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];          \
237e497d3c8SHong Zhang     _pj  = p_oth->j + _pi[_row];                  \
238e497d3c8SHong Zhang     _pa  = p_oth->a + _pi[_row];                  \
239e497d3c8SHong Zhang     /* perform dense axpy */                     \
240e497d3c8SHong Zhang     _valtmp = _aa[_j];                             \
241e497d3c8SHong Zhang     for (_k=0; _k<_pnz; _k++) {                     \
242e497d3c8SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];                \
243e497d3c8SHong Zhang     }                                            \
244e497d3c8SHong Zhang     PetscLogFlops(2.0*_pnz);                     \
245e497d3c8SHong Zhang   } \
246e497d3c8SHong Zhang }
247*cf742fccSHong Zhang 
2483369ce9aSBarry Smith #endif
249