xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision eec179cf895b6fbcd6a0b58694319392b06e361b)
16524c165SJacob Faibussowitsch #ifndef __MPIAIJ_H
23369ce9aSBarry Smith #define __MPIAIJ_H
33369ce9aSBarry Smith 
4c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
58a729477SBarry Smith 
690431a8fSHong Zhang typedef struct { /* used by MatCreateMPIAIJSumSeqAIJ for reusing the merged matrix */
726283091SBarry Smith   PetscLayout  rowmap;
8de0260b3SHong Zhang   PetscInt   **buf_ri, **buf_rj;
90febcb4bSHong Zhang   PetscMPIInt *len_s, *len_r, *id_r; /* array of length of comm->size, store send/recv matrix values */
100febcb4bSHong Zhang   PetscMPIInt  nsend, nrecv;
11fc08c53fSHong Zhang   PetscInt    *bi, *bj;               /* i and j array of the local portion of mpi C (matrix product) - rename to ci, cj! */
12de0260b3SHong Zhang   PetscInt    *owners_co, *coi, *coj; /* i and j array of (p->B)^T*A*P - used in the communication */
13b90dcfe3SHong Zhang } Mat_Merge_SeqsToMPI;
14b90dcfe3SHong Zhang 
154222ddf1SHong Zhang typedef struct {                                /* used by MatPtAPXXX_MPIAIJ_MPIAIJ() and MatMatMultXXX_MPIAIJ_MPIAIJ() */
16b7f45c76SHong Zhang   PetscInt              *startsj_s, *startsj_r; /* used by MatGetBrowsOfAoCols_MPIAIJ */
17b7f45c76SHong Zhang   PetscScalar           *bufa;                  /* used by MatGetBrowsOfAoCols_MPIAIJ */
18a1a4e84aSHong Zhang   Mat                    P_loc, P_oth;          /* partial B_seq -- intend to replace B_seq */
19a1a4e84aSHong Zhang   PetscInt              *api, *apj;             /* symbolic i and j arrays of the local product A_loc*B_seq */
20445158ffSHong Zhang   PetscScalar           *apv;
21b7f45c76SHong Zhang   MatReuse               reuse; /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */
22598bc09dSHong Zhang   PetscScalar           *apa;   /* tmp array for store a row of A*P used in MatMatMult() */
23c5af039cSHong Zhang   Mat                    A_loc; /* used by MatTransposeMatMult(), contains api and apj */
24a3bb6f32SFande Kong   ISLocalToGlobalMapping ltog;  /* mapping from local column indices to global column indices for A_loc */
25ab1b013aSHong Zhang   Mat                    Pt;    /* used by MatTransposeMatMult(), Pt = P^T */
269ce11a7cSHong Zhang   Mat                    Rd, Ro, AP_loc, C_loc, C_oth;
27cf97cf3cSHong Zhang   PetscInt               algType; /* implementation algorithm */
284a56b808SFande Kong   PetscSF                sf;      /* use it to communicate remote part of C */
294a56b808SFande Kong   PetscInt              *c_othi, *c_rmti;
30f8487c73SHong Zhang 
31f8487c73SHong Zhang   Mat_Merge_SeqsToMPI *merge;
323cdca5ebSHong Zhang } Mat_APMPI;
33f8487c73SHong Zhang 
34f8487c73SHong Zhang typedef struct {
35f8487c73SHong Zhang   Mat         A, B; /* local submatrices: A (diag part),
36f8487c73SHong Zhang                                            B (off-diag part) */
37f8487c73SHong Zhang   PetscMPIInt size; /* size of communicator */
38f8487c73SHong Zhang   PetscMPIInt rank; /* rank of proc in communicator */
39f8487c73SHong Zhang 
40f8487c73SHong Zhang   /* The following variables are used for matrix assembly */
41f8487c73SHong Zhang   PetscBool    donotstash;        /* PETSC_TRUE if off processor entries dropped */
42f8487c73SHong Zhang   MPI_Request *send_waits;        /* array of send requests */
43f8487c73SHong Zhang   MPI_Request *recv_waits;        /* array of receive requests */
44f8487c73SHong Zhang   PetscInt     nsends, nrecvs;    /* numbers of sends and receives */
45f8487c73SHong Zhang   PetscScalar *svalues, *rvalues; /* sending and receiving data */
46f8487c73SHong Zhang   PetscInt     rmax;              /* maximum message length */
47f8487c73SHong Zhang #if defined(PETSC_USE_CTABLE)
48*eec179cfSJacob Faibussowitsch   PetscHMapI colmap;
49f8487c73SHong Zhang #else
50f8487c73SHong Zhang   PetscInt *colmap; /* local col number of off-diag col */
51f8487c73SHong Zhang #endif
52f8487c73SHong Zhang   PetscInt *garray; /* global index of all off-processor columns */
53f8487c73SHong Zhang 
54f8487c73SHong Zhang   /* The following variables are used for matrix-vector products */
55f8487c73SHong Zhang   Vec        lvec; /* local vector */
56f8487c73SHong Zhang   Vec        diag;
5797929ea7SJunchao Zhang   VecScatter Mvctx;       /* scatter context for vector */
58f8487c73SHong Zhang   PetscBool  roworiented; /* if true, row-oriented input, default true */
59f8487c73SHong Zhang 
60f8487c73SHong Zhang   /* The following variables are for MatGetRow() */
61f8487c73SHong Zhang   PetscInt    *rowindices;   /* column indices for row */
62f8487c73SHong Zhang   PetscScalar *rowvalues;    /* nonzero values in row */
63f8487c73SHong Zhang   PetscBool    getrowactive; /* indicates MatGetRow(), not restored */
64f8487c73SHong Zhang 
65ce496241SStefano Zampini   PetscInt *ld; /* number of entries per row left of diagonal block */
66f8487c73SHong Zhang 
67ce496241SStefano Zampini   /* Used by device classes */
68bbf3fe20SPaul Mullowney   void *spptr;
695cc03489SHong Zhang 
70394ed5ebSJunchao Zhang   /* MatSetValuesCOO() related stuff */
71394ed5ebSJunchao Zhang   PetscCount   coo_n;                      /* Number of COOs passed to MatSetPreallocationCOO)() */
72394ed5ebSJunchao Zhang   PetscSF      coo_sf;                     /* SF to send/recv remote values in MatSetValuesCOO() */
73158ec288SJunchao Zhang   PetscCount   Annz, Bnnz;                 /* Number of entries in diagonal A and off-diagonal B */
74158ec288SJunchao Zhang   PetscCount   Annz2, Bnnz2;               /* Number of unique remote entries belonging to A and B */
75158ec288SJunchao Zhang   PetscCount   Atot1, Atot2, Btot1, Btot2; /* Total local (tot1) and remote (tot2) entries (which might contain repeats) belonging to A and B */
76158ec288SJunchao Zhang   PetscCount  *Ajmap1, *Aperm1;            /* Lengths: [Annz+1], [Atot1]. Local entries to diag */
77158ec288SJunchao Zhang   PetscCount  *Bjmap1, *Bperm1;            /* Lengths: [Bnnz+1], [Btot1]. Local entries to offdiag */
78394ed5ebSJunchao Zhang   PetscCount  *Aimap2, *Ajmap2, *Aperm2;   /* Lengths: [Annz2], [Annz2+1], [Atot2]. Remote entries to diag */
79394ed5ebSJunchao Zhang   PetscCount  *Bimap2, *Bjmap2, *Bperm2;   /* Lengths: [Bnnz2], [Bnnz2+1], [Btot2]. Remote entries to offdiag */
80394ed5ebSJunchao Zhang   PetscCount  *Cperm1;                     /* [sendlen] Permutation to fill MPI send buffer. 'C' for communication */
81394ed5ebSJunchao Zhang   PetscScalar *sendbuf, *recvbuf;          /* Buffers for remote values in MatSetValuesCOO() */
82394ed5ebSJunchao Zhang   PetscInt     sendlen, recvlen;           /* Lengths (in unit of PetscScalar) of send/recvbuf */
83f8487c73SHong Zhang } Mat_MPIAIJ;
84f8487c73SHong Zhang 
858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
862327d61dSSatish Balay 
8793283c78SAlejandro Lamas Daviña PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat, MatAssemblyType);
8893283c78SAlejandro Lamas Daviña 
895a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
905a576424SJed Brown PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
915a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat, MatDuplicateOption, Mat *);
925a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat, PetscInt, IS[], PetscInt);
93b1b1104fSBarry Smith PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat, PetscInt, IS[], PetscInt);
9493dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat, ISColoring, MatFDColoring);
95f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat, ISColoring, MatFDColoring);
967dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIAIJ(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
977dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
987dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat, MatCreateSubMatrixOption, MatReuse, Mat *[]);
99baa1ba78SHong Zhang PETSC_INTERN PetscErrorCode MatView_MPIAIJ(Mat, PetscViewer);
1005494a064SHong Zhang 
1017dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ(Mat, IS, IS, MatReuse, Mat *);
102618cbb4aSHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_nonscalable(Mat, IS, IS, PetscInt, MatReuse, Mat *);
1031358a193SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowDist(Mat, IS, IS, IS, MatReuse, Mat *);
1043b00a383SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowColDist(Mat, IS, IS, MatReuse, Mat *);
1053233f0ddSHong Zhang PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat, MPI_Comm, MatReuse, Mat *);
106d6037b41SHong Zhang 
1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat, PetscViewer);
10852f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ_Binary(Mat, PetscViewer);
109e3225e9dSHong Zhang PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
1104222ddf1SHong Zhang 
1114222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat);
1124e84afc0SStefano Zampini PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJBACKEND(Mat);
1134e84afc0SStefano Zampini PETSC_INTERN PetscErrorCode MatProductSymbolic_MPIAIJBACKEND(Mat);
1144222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat);
1154222ddf1SHong Zhang 
1164222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ(Mat);
1174222ddf1SHong Zhang 
1184222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_MPIAIJ_MPIAIJ(Mat);
1194222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductNumeric_RARt_MPIAIJ_MPIAIJ(Mat);
1204222ddf1SHong Zhang 
1214222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat, Mat, PetscReal, Mat);
1224222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat, Mat, PetscReal, Mat);
1234222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat, Mat, PetscReal, Mat);
1245a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat, Mat, Mat);
125f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat, Mat, Mat);
126f8487c73SHong Zhang 
1274222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat, Mat, Mat, PetscReal, Mat);
1285a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat, Mat, Mat, Mat);
129f996eeb8SHong Zhang 
1304222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat, Mat, PetscReal, Mat);
1315a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat, Mat, Mat);
1320d3441aeSHong Zhang 
1334222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat, Mat, PetscReal, Mat);
1344222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat, Mat, PetscReal, Mat);
1354222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat, Mat, PetscReal, Mat);
136a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat, Mat, Mat);
1374a56b808SFande Kong PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat, Mat, Mat);
1384a56b808SFande Kong PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat, Mat, Mat);
139a9d7e0ccSandi selinger 
140a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE)
1414222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
142a4ffb656SHong Zhang #endif
1434e29119aSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPIDense(Mat, MatType, MatReuse, Mat *);
144d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
145d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
146d24d4204SJose E. Roman #endif
1470d3441aeSHong Zhang 
1485a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
1496718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(void *);
1506718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void *);
151c5af039cSHong Zhang 
1525a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat, Mat, MatReuse, PetscInt **, PetscInt **, MatScalar **, Mat *);
1535a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
154e9ede7d0Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat(Mat, const PetscInt[], const PetscInt[], const PetscScalar[]);
155904d1e70Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Mat, const PetscInt[], const PetscInt[]);
1565a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat, MatOption, PetscBool);
157187b3c17SHong Zhang 
1584222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat, Mat, PetscReal, Mat);
1594222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat, Mat, PetscReal, Mat);
1606da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat, Mat, Mat);
1612bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat, Mat, Mat);
162c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat, Mat, Mat);
1634222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat, Mat, PetscReal, Mat);
1645a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat, Mat *);
165c0d3702cSSatish Balay 
166dbbe0bcdSBarry Smith PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(Mat, PetscOptionItems *);
16710e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]);
1683a7fca6bSBarry Smith 
169570b7f6dSBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
1705a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat, IS, IS, const MatFactorInfo *, Mat *);
17197304618SKris Buschelman #endif
1725a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat, Vec, Vec);
1735a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat, IS, IS, const MatFactorInfo *);
17497304618SKris Buschelman 
175001ddc4fSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *);
176001ddc4fSHong Zhang 
17711bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat, Mat *);
1787087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat, Vec);
1793a7fca6bSBarry Smith 
18053dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat, Mat *, Mat *);
18153dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat, IS, IS, IS, MatStructure, Mat, Mat);
18253dd7562SDmitry Karpeev 
183e8729f6fSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_MPIAIJ(Mat, PetscCount, PetscInt[], PetscInt[]);
184cbc6b225SStefano Zampini PETSC_INTERN PetscErrorCode MatResetPreallocationCOO_MPIAIJ(Mat);
185394ed5ebSJunchao Zhang 
186cf742fccSHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
187cf742fccSHong Zhang #define AProw_scalable(i, ad, ao, p_loc, p_oth, api, apj, apa) \
188cf742fccSHong Zhang   { \
189cf742fccSHong Zhang     PetscInt     _anz, _pnz, _j, _k, *_ai, *_aj, _row, *_pi, *_pj, _nextp, *_apJ; \
190cf742fccSHong Zhang     PetscScalar *_aa, _valtmp, *_pa; \
191cf742fccSHong Zhang     _apJ = apj + api[i]; \
192cf742fccSHong Zhang     /* diagonal portion of A */ \
193cf742fccSHong Zhang     _ai  = ad->i; \
194cf742fccSHong Zhang     _anz = _ai[i + 1] - _ai[i]; \
195cf742fccSHong Zhang     _aj  = ad->j + _ai[i]; \
196cf742fccSHong Zhang     _aa  = ad->a + _ai[i]; \
197cf742fccSHong Zhang     for (_j = 0; _j < _anz; _j++) { \
198cf742fccSHong Zhang       _row = _aj[_j]; \
199cf742fccSHong Zhang       _pi  = p_loc->i; \
200cf742fccSHong Zhang       _pnz = _pi[_row + 1] - _pi[_row]; \
201cf742fccSHong Zhang       _pj  = p_loc->j + _pi[_row]; \
202cf742fccSHong Zhang       _pa  = p_loc->a + _pi[_row]; \
203cf742fccSHong Zhang       /* perform sparse axpy */ \
204cf742fccSHong Zhang       _valtmp = _aa[_j]; \
205cf742fccSHong Zhang       _nextp  = 0; \
206cf742fccSHong Zhang       for (_k = 0; _nextp < _pnz; _k++) { \
207cf742fccSHong Zhang         if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */ \
208cf742fccSHong Zhang           apa[_k] += _valtmp * _pa[_nextp++]; \
209cf742fccSHong Zhang         } \
210cf742fccSHong Zhang       } \
21146eb3cd7SSatish Balay       (void)PetscLogFlops(2.0 * _pnz); \
212cf742fccSHong Zhang     } \
213cf742fccSHong Zhang     /* off-diagonal portion of A */ \
21492ec70a1SHong Zhang     if (p_oth) { \
215cf742fccSHong Zhang       _ai  = ao->i; \
216cf742fccSHong Zhang       _anz = _ai[i + 1] - _ai[i]; \
217cf742fccSHong Zhang       _aj  = ao->j + _ai[i]; \
218cf742fccSHong Zhang       _aa  = ao->a + _ai[i]; \
219cf742fccSHong Zhang       for (_j = 0; _j < _anz; _j++) { \
220cf742fccSHong Zhang         _row = _aj[_j]; \
221cf742fccSHong Zhang         _pi  = p_oth->i; \
222cf742fccSHong Zhang         _pnz = _pi[_row + 1] - _pi[_row]; \
223cf742fccSHong Zhang         _pj  = p_oth->j + _pi[_row]; \
224cf742fccSHong Zhang         _pa  = p_oth->a + _pi[_row]; \
225cf742fccSHong Zhang         /* perform sparse axpy */ \
226cf742fccSHong Zhang         _valtmp = _aa[_j]; \
227cf742fccSHong Zhang         _nextp  = 0; \
228cf742fccSHong Zhang         for (_k = 0; _nextp < _pnz; _k++) { \
229cf742fccSHong Zhang           if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */ \
230cf742fccSHong Zhang             apa[_k] += _valtmp * _pa[_nextp++]; \
231cf742fccSHong Zhang           } \
232cf742fccSHong Zhang         } \
23346eb3cd7SSatish Balay         (void)PetscLogFlops(2.0 * _pnz); \
234cf742fccSHong Zhang       } \
23592ec70a1SHong Zhang     } \
236cf742fccSHong Zhang   }
237cf742fccSHong Zhang 
238e497d3c8SHong Zhang #define AProw_nonscalable(i, ad, ao, p_loc, p_oth, apa) \
239e497d3c8SHong Zhang   { \
240e497d3c8SHong Zhang     PetscInt     _anz, _pnz, _j, _k, *_ai, *_aj, _row, *_pi, *_pj; \
241e497d3c8SHong Zhang     PetscScalar *_aa, _valtmp, *_pa; \
242e497d3c8SHong Zhang     /* diagonal portion of A */ \
243e497d3c8SHong Zhang     _ai  = ad->i; \
244e497d3c8SHong Zhang     _anz = _ai[i + 1] - _ai[i]; \
245e497d3c8SHong Zhang     _aj  = ad->j + _ai[i]; \
246e497d3c8SHong Zhang     _aa  = ad->a + _ai[i]; \
247e497d3c8SHong Zhang     for (_j = 0; _j < _anz; _j++) { \
248e497d3c8SHong Zhang       _row = _aj[_j]; \
249e497d3c8SHong Zhang       _pi  = p_loc->i; \
250e497d3c8SHong Zhang       _pnz = _pi[_row + 1] - _pi[_row]; \
251e497d3c8SHong Zhang       _pj  = p_loc->j + _pi[_row]; \
252e497d3c8SHong Zhang       _pa  = p_loc->a + _pi[_row]; \
253e497d3c8SHong Zhang       /* perform dense axpy */ \
254e497d3c8SHong Zhang       _valtmp = _aa[_j]; \
255ad540459SPierre Jolivet       for (_k = 0; _k < _pnz; _k++) apa[_pj[_k]] += _valtmp * _pa[_k]; \
25646eb3cd7SSatish Balay       (void)PetscLogFlops(2.0 * _pnz); \
257e497d3c8SHong Zhang     } \
258e497d3c8SHong Zhang     /* off-diagonal portion of A */ \
25992ec70a1SHong Zhang     if (p_oth) { \
260e497d3c8SHong Zhang       _ai  = ao->i; \
261e497d3c8SHong Zhang       _anz = _ai[i + 1] - _ai[i]; \
262e497d3c8SHong Zhang       _aj  = ao->j + _ai[i]; \
263e497d3c8SHong Zhang       _aa  = ao->a + _ai[i]; \
264e497d3c8SHong Zhang       for (_j = 0; _j < _anz; _j++) { \
265e497d3c8SHong Zhang         _row = _aj[_j]; \
266e497d3c8SHong Zhang         _pi  = p_oth->i; \
267e497d3c8SHong Zhang         _pnz = _pi[_row + 1] - _pi[_row]; \
268e497d3c8SHong Zhang         _pj  = p_oth->j + _pi[_row]; \
269e497d3c8SHong Zhang         _pa  = p_oth->a + _pi[_row]; \
270e497d3c8SHong Zhang         /* perform dense axpy */ \
271e497d3c8SHong Zhang         _valtmp = _aa[_j]; \
272ad540459SPierre Jolivet         for (_k = 0; _k < _pnz; _k++) apa[_pj[_k]] += _valtmp * _pa[_k]; \
27346eb3cd7SSatish Balay         (void)PetscLogFlops(2.0 * _pnz); \
274e497d3c8SHong Zhang       } \
27592ec70a1SHong Zhang     } \
276e497d3c8SHong Zhang   }
277cf742fccSHong Zhang 
2783369ce9aSBarry Smith #endif
279