xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision 68d6bbd5ee09882b8e51c12fffe91e4d70783644)
1a4963045SJacob Faibussowitsch #pragma once
23369ce9aSBarry Smith 
3c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
48a729477SBarry Smith 
590431a8fSHong Zhang typedef struct { /* used by MatCreateMPIAIJSumSeqAIJ for reusing the merged matrix */
626283091SBarry Smith   PetscLayout  rowmap;
7de0260b3SHong Zhang   PetscInt   **buf_ri, **buf_rj;
80febcb4bSHong Zhang   PetscMPIInt *len_s, *len_r, *id_r; /* array of length of comm->size, store send/recv matrix values */
90febcb4bSHong Zhang   PetscMPIInt  nsend, nrecv;
10fc08c53fSHong Zhang   PetscInt    *bi, *bj;               /* i and j array of the local portion of mpi C (matrix product) - rename to ci, cj! */
11de0260b3SHong Zhang   PetscInt    *owners_co, *coi, *coj; /* i and j array of (p->B)^T*A*P - used in the communication */
12b90dcfe3SHong Zhang } Mat_Merge_SeqsToMPI;
13b90dcfe3SHong Zhang 
144222ddf1SHong Zhang typedef struct {                                /* used by MatPtAPXXX_MPIAIJ_MPIAIJ() and MatMatMultXXX_MPIAIJ_MPIAIJ() */
15b7f45c76SHong Zhang   PetscInt              *startsj_s, *startsj_r; /* used by MatGetBrowsOfAoCols_MPIAIJ */
16b7f45c76SHong Zhang   PetscScalar           *bufa;                  /* used by MatGetBrowsOfAoCols_MPIAIJ */
17a1a4e84aSHong Zhang   Mat                    P_loc, P_oth;          /* partial B_seq -- intend to replace B_seq */
18a1a4e84aSHong Zhang   PetscInt              *api, *apj;             /* symbolic i and j arrays of the local product A_loc*B_seq */
19445158ffSHong Zhang   PetscScalar           *apv;
20b7f45c76SHong Zhang   MatReuse               reuse; /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */
21598bc09dSHong Zhang   PetscScalar           *apa;   /* tmp array for store a row of A*P used in MatMatMult() */
22c5af039cSHong Zhang   Mat                    A_loc; /* used by MatTransposeMatMult(), contains api and apj */
23a3bb6f32SFande Kong   ISLocalToGlobalMapping ltog;  /* mapping from local column indices to global column indices for A_loc */
24ab1b013aSHong Zhang   Mat                    Pt;    /* used by MatTransposeMatMult(), Pt = P^T */
259ce11a7cSHong Zhang   Mat                    Rd, Ro, AP_loc, C_loc, C_oth;
26cf97cf3cSHong Zhang   PetscInt               algType; /* implementation algorithm */
274a56b808SFande Kong   PetscSF                sf;      /* use it to communicate remote part of C */
284a56b808SFande Kong   PetscInt              *c_othi, *c_rmti;
29f8487c73SHong Zhang 
30f8487c73SHong Zhang   Mat_Merge_SeqsToMPI *merge;
313cdca5ebSHong Zhang } Mat_APMPI;
32f8487c73SHong Zhang 
33f8487c73SHong Zhang #if defined(PETSC_USE_CTABLE)
34*68d6bbd5SJunchao Zhang   #define PETSCTABLE PetscHMapI
35f8487c73SHong Zhang #else
36*68d6bbd5SJunchao Zhang   #define PETSCTABLE PetscInt *
37f8487c73SHong Zhang #endif
38f8487c73SHong Zhang 
39*68d6bbd5SJunchao Zhang // Shared by MPIAIJ, MPIBAIJ, MPISBAIJ so that we can access common fields in the same way.
40*68d6bbd5SJunchao Zhang #define MPIAIJHEADER \
41*68d6bbd5SJunchao Zhang   Mat         A, B; /* local submatrices: A (diag part), B (off-diag part) */ \
42*68d6bbd5SJunchao Zhang   PetscMPIInt size; /* size of communicator */ \
43*68d6bbd5SJunchao Zhang   PetscMPIInt rank; /* rank of proc in communicator */ \
44*68d6bbd5SJunchao Zhang \
45*68d6bbd5SJunchao Zhang   /* The following variables are used for matrix assembly */ \
46*68d6bbd5SJunchao Zhang   PetscBool    donotstash;        /* if 1, off processor entries dropped */ \
47*68d6bbd5SJunchao Zhang   MPI_Request *send_waits;        /* array of send requests */ \
48*68d6bbd5SJunchao Zhang   MPI_Request *recv_waits;        /* array of receive requests */ \
49*68d6bbd5SJunchao Zhang   PetscInt     nsends, nrecvs;    /* numbers of sends and receives */ \
50*68d6bbd5SJunchao Zhang   MatScalar   *svalues, *rvalues; /* sending and receiving data */ \
51*68d6bbd5SJunchao Zhang   PetscInt     rmax;              /* maximum message length */ \
52*68d6bbd5SJunchao Zhang   PETSCTABLE   colmap;            /* local col number of off-diag col */ \
53*68d6bbd5SJunchao Zhang   PetscInt    *garray;            /* work array */ \
54*68d6bbd5SJunchao Zhang \
55*68d6bbd5SJunchao Zhang   /* The following variables are used for matrix-vector products */ \
56*68d6bbd5SJunchao Zhang   Vec        lvec;        /* local vector */ \
57*68d6bbd5SJunchao Zhang   VecScatter Mvctx;       /* scatter context for vector */ \
58*68d6bbd5SJunchao Zhang   PetscBool  roworiented; /* if true, row-oriented input, default true */ \
59*68d6bbd5SJunchao Zhang \
60*68d6bbd5SJunchao Zhang   /* The following variables are for MatGetRow() */ \
61*68d6bbd5SJunchao Zhang   PetscInt    *rowindices; /* column indices for row */ \
62*68d6bbd5SJunchao Zhang   PetscScalar *rowvalues;  /* nonzero values in row */ \
63*68d6bbd5SJunchao Zhang   PetscBool    getrowactive
64*68d6bbd5SJunchao Zhang 
65*68d6bbd5SJunchao Zhang typedef struct {
66*68d6bbd5SJunchao Zhang   MPIAIJHEADER;
67f8487c73SHong Zhang   Vec       diag;
68ce496241SStefano Zampini   PetscInt *ld; /* number of entries per row left of diagonal block */
69f8487c73SHong Zhang 
70ce496241SStefano Zampini   /* Used by device classes */
71bbf3fe20SPaul Mullowney   void *spptr;
725cc03489SHong Zhang 
732c4ab24aSJunchao Zhang   struct _MatOps cops;
742c4ab24aSJunchao Zhang } Mat_MPIAIJ;
752c4ab24aSJunchao Zhang 
762c4ab24aSJunchao Zhang typedef struct {
772c4ab24aSJunchao Zhang   PetscCount   n;                          /* Number of COOs passed to MatSetPreallocationCOO)() */
782c4ab24aSJunchao Zhang   PetscSF      sf;                         /* SF to send/recv remote values in MatSetValuesCOO() */
79158ec288SJunchao Zhang   PetscCount   Annz, Bnnz;                 /* Number of entries in diagonal A and off-diagonal B */
80158ec288SJunchao Zhang   PetscCount   Annz2, Bnnz2;               /* Number of unique remote entries belonging to A and B */
81158ec288SJunchao Zhang   PetscCount   Atot1, Atot2, Btot1, Btot2; /* Total local (tot1) and remote (tot2) entries (which might contain repeats) belonging to A and B */
82158ec288SJunchao Zhang   PetscCount  *Ajmap1, *Aperm1;            /* Lengths: [Annz+1], [Atot1]. Local entries to diag */
83158ec288SJunchao Zhang   PetscCount  *Bjmap1, *Bperm1;            /* Lengths: [Bnnz+1], [Btot1]. Local entries to offdiag */
84394ed5ebSJunchao Zhang   PetscCount  *Aimap2, *Ajmap2, *Aperm2;   /* Lengths: [Annz2], [Annz2+1], [Atot2]. Remote entries to diag */
85394ed5ebSJunchao Zhang   PetscCount  *Bimap2, *Bjmap2, *Bperm2;   /* Lengths: [Bnnz2], [Bnnz2+1], [Btot2]. Remote entries to offdiag */
86394ed5ebSJunchao Zhang   PetscCount  *Cperm1;                     /* [sendlen] Permutation to fill MPI send buffer. 'C' for communication */
87394ed5ebSJunchao Zhang   PetscScalar *sendbuf, *recvbuf;          /* Buffers for remote values in MatSetValuesCOO() */
88394ed5ebSJunchao Zhang   PetscInt     sendlen, recvlen;           /* Lengths (in unit of PetscScalar) of send/recvbuf */
892c4ab24aSJunchao Zhang } MatCOOStruct_MPIAIJ;
90f8487c73SHong Zhang 
918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
922327d61dSSatish Balay 
9393283c78SAlejandro Lamas Daviña PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat, MatAssemblyType);
9493283c78SAlejandro Lamas Daviña 
955a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
965a576424SJed Brown PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
975a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat, MatDuplicateOption, Mat *);
985a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat, PetscInt, IS[], PetscInt);
99b1b1104fSBarry Smith PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat, PetscInt, IS[], PetscInt);
10093dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat, ISColoring, MatFDColoring);
101f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat, ISColoring, MatFDColoring);
1027dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIAIJ(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
1037dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
1047dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat, MatCreateSubMatrixOption, MatReuse, Mat *[]);
105baa1ba78SHong Zhang PETSC_INTERN PetscErrorCode MatView_MPIAIJ(Mat, PetscViewer);
1065494a064SHong Zhang 
1077dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ(Mat, IS, IS, MatReuse, Mat *);
108618cbb4aSHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_nonscalable(Mat, IS, IS, PetscInt, MatReuse, Mat *);
1091358a193SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowDist(Mat, IS, IS, IS, MatReuse, Mat *);
1103b00a383SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowColDist(Mat, IS, IS, MatReuse, Mat *);
1113233f0ddSHong Zhang PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat, MPI_Comm, MatReuse, Mat *);
112d6037b41SHong Zhang 
1135a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat, PetscViewer);
11452f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ_Binary(Mat, PetscViewer);
115e3225e9dSHong Zhang PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
1164222ddf1SHong Zhang 
1174222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat);
1184e84afc0SStefano Zampini PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJBACKEND(Mat);
1194e84afc0SStefano Zampini PETSC_INTERN PetscErrorCode MatProductSymbolic_MPIAIJBACKEND(Mat);
1204222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat);
1214222ddf1SHong Zhang 
1224222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ(Mat);
1234222ddf1SHong Zhang 
1244222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_MPIAIJ_MPIAIJ(Mat);
1254222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductNumeric_RARt_MPIAIJ_MPIAIJ(Mat);
1264222ddf1SHong Zhang 
1274222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat, Mat, PetscReal, Mat);
1284222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat, Mat, PetscReal, Mat);
1294222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat, Mat, PetscReal, Mat);
1305a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat, Mat, Mat);
131f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat, Mat, Mat);
132f8487c73SHong Zhang 
1334222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat, Mat, Mat, PetscReal, Mat);
1345a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat, Mat, Mat, Mat);
135f996eeb8SHong Zhang 
1364222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat, Mat, PetscReal, Mat);
1375a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat, Mat, Mat);
1380d3441aeSHong Zhang 
1394222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat, Mat, PetscReal, Mat);
1404222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat, Mat, PetscReal, Mat);
1414222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat, Mat, PetscReal, Mat);
142a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat, Mat, Mat);
1434a56b808SFande Kong PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat, Mat, Mat);
1444a56b808SFande Kong PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat, Mat, Mat);
145a9d7e0ccSandi selinger 
146a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE)
1474222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
148a4ffb656SHong Zhang #endif
1494e29119aSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPIDense(Mat, MatType, MatReuse, Mat *);
150d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
151d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
152d24d4204SJose E. Roman #endif
1530d3441aeSHong Zhang 
1545a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
1556718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(void *);
1566718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void *);
157c5af039cSHong Zhang 
1585a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat, Mat, MatReuse, PetscInt **, PetscInt **, MatScalar **, Mat *);
1595a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
160e9ede7d0Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat(Mat, const PetscInt[], const PetscInt[], const PetscScalar[]);
161904d1e70Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Mat, const PetscInt[], const PetscInt[]);
1625a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat, MatOption, PetscBool);
163187b3c17SHong Zhang 
1644222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat, Mat, PetscReal, Mat);
1654222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat, Mat, PetscReal, Mat);
1666da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat, Mat, Mat);
1672bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat, Mat, Mat);
168c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat, Mat, Mat);
1694222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat, Mat, PetscReal, Mat);
1705a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat, Mat *);
171c0d3702cSSatish Balay 
172dbbe0bcdSBarry Smith PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(Mat, PetscOptionItems *);
17310e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]);
1743a7fca6bSBarry Smith 
175570b7f6dSBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
1765a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat, IS, IS, const MatFactorInfo *, Mat *);
17797304618SKris Buschelman #endif
1785a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat, Vec, Vec);
1795a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat, IS, IS, const MatFactorInfo *);
18097304618SKris Buschelman 
181001ddc4fSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *);
182001ddc4fSHong Zhang 
18311bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat, Mat *);
1847087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat, Vec);
1853a7fca6bSBarry Smith 
18653dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat, Mat *, Mat *);
18753dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat, IS, IS, IS, MatStructure, Mat, Mat);
18853dd7562SDmitry Karpeev 
189e8729f6fSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_MPIAIJ(Mat, PetscCount, PetscInt[], PetscInt[]);
190cbc6b225SStefano Zampini PETSC_INTERN PetscErrorCode MatResetPreallocationCOO_MPIAIJ(Mat);
191394ed5ebSJunchao Zhang 
192cf742fccSHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
193cf742fccSHong Zhang #define AProw_scalable(i, ad, ao, p_loc, p_oth, api, apj, apa) \
194a8f51744SPierre Jolivet   do { \
195cf742fccSHong Zhang     PetscInt     _anz, _pnz, _j, _k, *_ai, *_aj, _row, *_pi, *_pj, _nextp, *_apJ; \
196cf742fccSHong Zhang     PetscScalar *_aa, _valtmp, *_pa; \
197cf742fccSHong Zhang     _apJ = apj + api[i]; \
198cf742fccSHong Zhang     /* diagonal portion of A */ \
199cf742fccSHong Zhang     _ai  = ad->i; \
200cf742fccSHong Zhang     _anz = _ai[i + 1] - _ai[i]; \
201cf742fccSHong Zhang     _aj  = ad->j + _ai[i]; \
202cf742fccSHong Zhang     _aa  = ad->a + _ai[i]; \
203cf742fccSHong Zhang     for (_j = 0; _j < _anz; _j++) { \
204cf742fccSHong Zhang       _row = _aj[_j]; \
205cf742fccSHong Zhang       _pi  = p_loc->i; \
206cf742fccSHong Zhang       _pnz = _pi[_row + 1] - _pi[_row]; \
207cf742fccSHong Zhang       _pj  = p_loc->j + _pi[_row]; \
208cf742fccSHong Zhang       _pa  = p_loc->a + _pi[_row]; \
209cf742fccSHong Zhang       /* perform sparse axpy */ \
210cf742fccSHong Zhang       _valtmp = _aa[_j]; \
211cf742fccSHong Zhang       _nextp  = 0; \
212cf742fccSHong Zhang       for (_k = 0; _nextp < _pnz; _k++) { \
213cf742fccSHong Zhang         if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */ \
214cf742fccSHong Zhang           apa[_k] += _valtmp * _pa[_nextp++]; \
215cf742fccSHong Zhang         } \
216cf742fccSHong Zhang       } \
21746eb3cd7SSatish Balay       (void)PetscLogFlops(2.0 * _pnz); \
218cf742fccSHong Zhang     } \
219cf742fccSHong Zhang     /* off-diagonal portion of A */ \
22092ec70a1SHong Zhang     if (p_oth) { \
221cf742fccSHong Zhang       _ai  = ao->i; \
222cf742fccSHong Zhang       _anz = _ai[i + 1] - _ai[i]; \
223cf742fccSHong Zhang       _aj  = ao->j + _ai[i]; \
224cf742fccSHong Zhang       _aa  = ao->a + _ai[i]; \
225cf742fccSHong Zhang       for (_j = 0; _j < _anz; _j++) { \
226cf742fccSHong Zhang         _row = _aj[_j]; \
227cf742fccSHong Zhang         _pi  = p_oth->i; \
228cf742fccSHong Zhang         _pnz = _pi[_row + 1] - _pi[_row]; \
229cf742fccSHong Zhang         _pj  = p_oth->j + _pi[_row]; \
230cf742fccSHong Zhang         _pa  = p_oth->a + _pi[_row]; \
231cf742fccSHong Zhang         /* perform sparse axpy */ \
232cf742fccSHong Zhang         _valtmp = _aa[_j]; \
233cf742fccSHong Zhang         _nextp  = 0; \
234cf742fccSHong Zhang         for (_k = 0; _nextp < _pnz; _k++) { \
235cf742fccSHong Zhang           if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */ \
236cf742fccSHong Zhang             apa[_k] += _valtmp * _pa[_nextp++]; \
237cf742fccSHong Zhang           } \
238cf742fccSHong Zhang         } \
23946eb3cd7SSatish Balay         (void)PetscLogFlops(2.0 * _pnz); \
240cf742fccSHong Zhang       } \
24192ec70a1SHong Zhang     } \
242a8f51744SPierre Jolivet   } while (0)
243cf742fccSHong Zhang 
244e497d3c8SHong Zhang #define AProw_nonscalable(i, ad, ao, p_loc, p_oth, apa) \
245a8f51744SPierre Jolivet   do { \
246e497d3c8SHong Zhang     PetscInt     _anz, _pnz, _j, _k, *_ai, *_aj, _row, *_pi, *_pj; \
247e497d3c8SHong Zhang     PetscScalar *_aa, _valtmp, *_pa; \
248e497d3c8SHong Zhang     /* diagonal portion of A */ \
249e497d3c8SHong Zhang     _ai  = ad->i; \
250e497d3c8SHong Zhang     _anz = _ai[i + 1] - _ai[i]; \
251e497d3c8SHong Zhang     _aj  = ad->j + _ai[i]; \
252e497d3c8SHong Zhang     _aa  = ad->a + _ai[i]; \
253e497d3c8SHong Zhang     for (_j = 0; _j < _anz; _j++) { \
254e497d3c8SHong Zhang       _row = _aj[_j]; \
255e497d3c8SHong Zhang       _pi  = p_loc->i; \
256e497d3c8SHong Zhang       _pnz = _pi[_row + 1] - _pi[_row]; \
257e497d3c8SHong Zhang       _pj  = p_loc->j + _pi[_row]; \
258e497d3c8SHong Zhang       _pa  = p_loc->a + _pi[_row]; \
259e497d3c8SHong Zhang       /* perform dense axpy */ \
260e497d3c8SHong Zhang       _valtmp = _aa[_j]; \
261ad540459SPierre Jolivet       for (_k = 0; _k < _pnz; _k++) apa[_pj[_k]] += _valtmp * _pa[_k]; \
26246eb3cd7SSatish Balay       (void)PetscLogFlops(2.0 * _pnz); \
263e497d3c8SHong Zhang     } \
264e497d3c8SHong Zhang     /* off-diagonal portion of A */ \
26592ec70a1SHong Zhang     if (p_oth) { \
266e497d3c8SHong Zhang       _ai  = ao->i; \
267e497d3c8SHong Zhang       _anz = _ai[i + 1] - _ai[i]; \
268e497d3c8SHong Zhang       _aj  = ao->j + _ai[i]; \
269e497d3c8SHong Zhang       _aa  = ao->a + _ai[i]; \
270e497d3c8SHong Zhang       for (_j = 0; _j < _anz; _j++) { \
271e497d3c8SHong Zhang         _row = _aj[_j]; \
272e497d3c8SHong Zhang         _pi  = p_oth->i; \
273e497d3c8SHong Zhang         _pnz = _pi[_row + 1] - _pi[_row]; \
274e497d3c8SHong Zhang         _pj  = p_oth->j + _pi[_row]; \
275e497d3c8SHong Zhang         _pa  = p_oth->a + _pi[_row]; \
276e497d3c8SHong Zhang         /* perform dense axpy */ \
277e497d3c8SHong Zhang         _valtmp = _aa[_j]; \
278ad540459SPierre Jolivet         for (_k = 0; _k < _pnz; _k++) apa[_pj[_k]] += _valtmp * _pa[_k]; \
27946eb3cd7SSatish Balay         (void)PetscLogFlops(2.0 * _pnz); \
280e497d3c8SHong Zhang       } \
28192ec70a1SHong Zhang     } \
282a8f51744SPierre Jolivet   } while (0)
283