xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision 904d1e7050d8278cb5f67df70e5edfa4eaab048b)
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;
30cf97cf3cSHong Zhang   PetscInt    algType;         /* implementation algorithm */
31f8487c73SHong Zhang 
32f8487c73SHong Zhang   Mat_Merge_SeqsToMPI *merge;
33f8487c73SHong Zhang   PetscErrorCode (*destroy)(Mat);
344ae0847bSHong Zhang   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
35cf97cf3cSHong Zhang   PetscErrorCode (*view)(Mat,PetscViewer);
36f8487c73SHong Zhang } Mat_PtAPMPI;
37f8487c73SHong Zhang 
38f8487c73SHong Zhang typedef struct {
39f8487c73SHong Zhang   Mat A,B;                             /* local submatrices: A (diag part),
40f8487c73SHong Zhang                                            B (off-diag part) */
41f8487c73SHong Zhang   PetscMPIInt size;                     /* size of communicator */
42f8487c73SHong Zhang   PetscMPIInt rank;                     /* rank of proc in communicator */
43f8487c73SHong Zhang 
44f8487c73SHong Zhang   /* The following variables are used for matrix assembly */
45f8487c73SHong Zhang   PetscBool   donotstash;               /* PETSC_TRUE if off processor entries dropped */
46f8487c73SHong Zhang   MPI_Request *send_waits;              /* array of send requests */
47f8487c73SHong Zhang   MPI_Request *recv_waits;              /* array of receive requests */
48f8487c73SHong Zhang   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
49f8487c73SHong Zhang   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
50f8487c73SHong Zhang   PetscInt    rmax;                     /* maximum message length */
51f8487c73SHong Zhang #if defined(PETSC_USE_CTABLE)
52f8487c73SHong Zhang   PetscTable colmap;
53f8487c73SHong Zhang #else
54f8487c73SHong Zhang   PetscInt *colmap;                     /* local col number of off-diag col */
55f8487c73SHong Zhang #endif
56f8487c73SHong Zhang   PetscInt *garray;                     /* global index of all off-processor columns */
57f8487c73SHong Zhang 
58f8487c73SHong Zhang   /* The following variables are used for matrix-vector products */
59f8487c73SHong Zhang   Vec        lvec;                 /* local vector */
60f8487c73SHong Zhang   Vec        diag;
614b8d542aSHong Zhang   VecScatter Mvctx,Mvctx_mpi1;     /* scatter context for vector */
62fa110396SHong Zhang   PetscBool  Mvctx_mpi1_flg;       /* if true, additional Mvctx_mpi1 is requested for mat-mat ops, default false */
63f8487c73SHong Zhang   PetscBool  roworiented;          /* if true, row-oriented input, default true */
64f8487c73SHong Zhang 
65f8487c73SHong Zhang   /* The following variables are for MatGetRow() */
66f8487c73SHong Zhang   PetscInt    *rowindices;         /* column indices for row */
67f8487c73SHong Zhang   PetscScalar *rowvalues;          /* nonzero values in row */
68f8487c73SHong Zhang   PetscBool   getrowactive;        /* indicates MatGetRow(), not restored */
69f8487c73SHong Zhang 
70f8487c73SHong Zhang   /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation  and nonzero pattern */
71f8487c73SHong Zhang   PetscInt *ld;                    /* number of entries per row left of diagona block */
72f8487c73SHong Zhang 
73a1a4e84aSHong Zhang   /* Used by MatMatMult() and MatPtAP() */
74f8487c73SHong Zhang   Mat_PtAPMPI *ptap;
75ca45077fSPaul Mullowney 
76f996eeb8SHong Zhang   /* used by MatMatMatMult() */
77f996eeb8SHong Zhang   Mat_MatMatMatMult *matmatmatmult;
78f996eeb8SHong Zhang 
79bbf3fe20SPaul Mullowney   /* Used by MPICUSP and MPICUSPARSE classes */
80bbf3fe20SPaul Mullowney   void * spptr;
815cc03489SHong Zhang 
82f8487c73SHong Zhang } Mat_MPIAIJ;
83f8487c73SHong Zhang 
848cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
852327d61dSSatish Balay 
8693283c78SAlejandro Lamas Daviña PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat,MatAssemblyType);
8793283c78SAlejandro Lamas Daviña 
885a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
895a576424SJed Brown PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
905a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*);
915a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
92b1b1104fSBarry Smith PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat,PetscInt,IS [],PetscInt);
9393dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
94f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
957dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
967dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
977dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat,MatCreateSubMatrixOption,MatReuse,Mat *[]);
9853dd7562SDmitry Karpeev 
995494a064SHong Zhang 
1007dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*);
101618cbb4aSHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_nonscalable(Mat,IS,IS,PetscInt,MatReuse,Mat*);
1021358a193SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowDist(Mat,IS,IS,IS,MatReuse,Mat*);
1033b00a383SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowColDist(Mat,IS,IS,MatReuse,Mat*);
1043233f0ddSHong Zhang PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*);
105d6037b41SHong Zhang 
1065a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer);
107e3225e9dSHong Zhang PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1095a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1100fc8cf34SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
111a9d7e0ccSandi selinger PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat,Mat,PetscReal,Mat*);
112b2405163SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1135a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
114f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
115f8487c73SHong Zhang 
1165a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
1175a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*);
1185a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);
119f996eeb8SHong Zhang 
1205a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1215a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1225a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
1230d3441aeSHong Zhang 
124a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat,Mat,PetscReal,Mat*);
125a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat,Mat,Mat);
126a9d7e0ccSandi selinger 
127a9d7e0ccSandi selinger 
128a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE)
129a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
130a4ffb656SHong Zhang #endif
1310d3441aeSHong Zhang 
1325a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat);
1335a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
134c5af039cSHong Zhang 
1358d45306eSHong Zhang PETSC_INTERN PetscErrorCode MatRARt_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1368d45306eSHong Zhang 
1375a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
1385a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
139*904d1e70Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat(Mat,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscInt[],const PetscInt[]);
140*904d1e70Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Mat,const PetscInt[],const PetscInt[]);
1415a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat);
1425a576424SJed Brown PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*);
1435a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);
144187b3c17SHong Zhang 
1455a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1462bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
1476da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1486da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
1492bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
150c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
1518949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
1528949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
1538949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
1545a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);
155c0d3702cSSatish Balay 
1564416b707SBarry Smith PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat);
15710e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1583a7fca6bSBarry Smith 
159570b7f6dSBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
1605a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*);
16197304618SKris Buschelman #endif
1625a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec);
1635a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*);
16497304618SKris Buschelman 
165001ddc4fSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
166001ddc4fSHong Zhang 
16711bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*);
1687087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);
1693a7fca6bSBarry Smith 
17053dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*);
17153dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat);
17253dd7562SDmitry Karpeev 
173cf742fccSHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
174cf742fccSHong Zhang #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \
175cf742fccSHong Zhang {\
176cf742fccSHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ;      \
177cf742fccSHong Zhang   PetscScalar *_aa,_valtmp,*_pa;                             \
178cf742fccSHong Zhang   _apJ = apj + api[i];\
179cf742fccSHong Zhang   /* diagonal portion of A */\
180cf742fccSHong Zhang   _ai  = ad->i;\
181cf742fccSHong Zhang   _anz = _ai[i+1] - _ai[i];\
182cf742fccSHong Zhang   _aj  = ad->j + _ai[i];\
183cf742fccSHong Zhang   _aa  = ad->a + _ai[i];\
184cf742fccSHong Zhang   for (_j=0; _j<_anz; _j++) {\
185cf742fccSHong Zhang     _row = _aj[_j]; \
186cf742fccSHong Zhang     _pi  = p_loc->i;                                 \
187cf742fccSHong Zhang     _pnz = _pi[_row+1] - _pi[_row];         \
188cf742fccSHong Zhang     _pj  = p_loc->j + _pi[_row];                 \
189cf742fccSHong Zhang     _pa  = p_loc->a + _pi[_row];                 \
190cf742fccSHong Zhang     /* perform sparse axpy */                    \
191cf742fccSHong Zhang     _valtmp = _aa[_j];                           \
192cf742fccSHong Zhang     _nextp  = 0; \
193cf742fccSHong Zhang     for (_k=0; _nextp<_pnz; _k++) {                    \
194cf742fccSHong Zhang       if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */   \
195cf742fccSHong Zhang         apa[_k] += _valtmp*_pa[_nextp++];                                \
196cf742fccSHong Zhang       } \
197cf742fccSHong Zhang     }                                           \
19846eb3cd7SSatish Balay     (void)PetscLogFlops(2.0*_pnz);              \
199cf742fccSHong Zhang   }                                             \
200cf742fccSHong Zhang   /* off-diagonal portion of A */               \
201cf742fccSHong Zhang   _ai  = ao->i;\
202cf742fccSHong Zhang   _anz = _ai[i+1] - _ai[i];                     \
203cf742fccSHong Zhang   _aj  = ao->j + _ai[i];                         \
204cf742fccSHong Zhang   _aa  = ao->a + _ai[i];                         \
205cf742fccSHong Zhang   for (_j=0; _j<_anz; _j++) {                      \
206cf742fccSHong Zhang     _row = _aj[_j];    \
207cf742fccSHong Zhang     _pi  = p_oth->i;                         \
208cf742fccSHong Zhang     _pnz = _pi[_row+1] - _pi[_row];          \
209cf742fccSHong Zhang     _pj  = p_oth->j + _pi[_row];                  \
210cf742fccSHong Zhang     _pa  = p_oth->a + _pi[_row];                  \
211cf742fccSHong Zhang     /* perform sparse axpy */                     \
212cf742fccSHong Zhang     _valtmp = _aa[_j];                             \
213cf742fccSHong Zhang     _nextp  = 0; \
214cf742fccSHong Zhang     for (_k=0; _nextp<_pnz; _k++) {                     \
215cf742fccSHong Zhang       if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
216cf742fccSHong Zhang         apa[_k] += _valtmp*_pa[_nextp++];                       \
217cf742fccSHong Zhang       }                                                     \
218cf742fccSHong Zhang     }                                            \
21946eb3cd7SSatish Balay     (void)PetscLogFlops(2.0*_pnz);               \
220cf742fccSHong Zhang   } \
221cf742fccSHong Zhang }
222cf742fccSHong Zhang 
223e497d3c8SHong Zhang #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
224e497d3c8SHong Zhang {\
225e497d3c8SHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;      \
226e497d3c8SHong Zhang   PetscScalar *_aa,_valtmp,*_pa;                             \
227e497d3c8SHong Zhang   /* diagonal portion of A */\
228e497d3c8SHong Zhang   _ai  = ad->i;\
229e497d3c8SHong Zhang   _anz = _ai[i+1] - _ai[i];\
230e497d3c8SHong Zhang   _aj  = ad->j + _ai[i];\
231e497d3c8SHong Zhang   _aa  = ad->a + _ai[i];\
232e497d3c8SHong Zhang   for (_j=0; _j<_anz; _j++) {\
233e497d3c8SHong Zhang     _row = _aj[_j]; \
234e497d3c8SHong Zhang     _pi  = p_loc->i;                                 \
235e497d3c8SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];         \
236e497d3c8SHong Zhang     _pj  = p_loc->j + _pi[_row];                 \
237e497d3c8SHong Zhang     _pa  = p_loc->a + _pi[_row];                 \
238e497d3c8SHong Zhang     /* perform dense axpy */                    \
239e497d3c8SHong Zhang     _valtmp = _aa[_j];                           \
240e497d3c8SHong Zhang     for (_k=0; _k<_pnz; _k++) {                    \
241e497d3c8SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];               \
242e497d3c8SHong Zhang     }                                           \
24346eb3cd7SSatish Balay     (void)PetscLogFlops(2.0*_pnz);              \
244e497d3c8SHong Zhang   }                                             \
245e497d3c8SHong Zhang   /* off-diagonal portion of A */               \
246e497d3c8SHong Zhang   _ai  = ao->i;\
247e497d3c8SHong Zhang   _anz = _ai[i+1] - _ai[i];                     \
248e497d3c8SHong Zhang   _aj  = ao->j + _ai[i];                         \
249e497d3c8SHong Zhang   _aa  = ao->a + _ai[i];                         \
250e497d3c8SHong Zhang   for (_j=0; _j<_anz; _j++) {                      \
251e497d3c8SHong Zhang     _row = _aj[_j];    \
252e497d3c8SHong Zhang     _pi  = p_oth->i;                         \
253e497d3c8SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];          \
254e497d3c8SHong Zhang     _pj  = p_oth->j + _pi[_row];                  \
255e497d3c8SHong Zhang     _pa  = p_oth->a + _pi[_row];                  \
256e497d3c8SHong Zhang     /* perform dense axpy */                     \
257e497d3c8SHong Zhang     _valtmp = _aa[_j];                             \
258e497d3c8SHong Zhang     for (_k=0; _k<_pnz; _k++) {                     \
259e497d3c8SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];                \
260e497d3c8SHong Zhang     }                                            \
26146eb3cd7SSatish Balay     (void)PetscLogFlops(2.0*_pnz);               \
262e497d3c8SHong Zhang   } \
263e497d3c8SHong Zhang }
264cf742fccSHong Zhang 
2653369ce9aSBarry Smith #endif
266