xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision 394ed5eb51aeb2fd4d5be41f3bc40ab0c0d7d8f0)
13369ce9aSBarry Smith #if !defined(__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)
48f8487c73SHong Zhang   PetscTable 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 
70*394ed5ebSJunchao Zhang   /* MatSetValuesCOO() related stuff */
71*394ed5ebSJunchao Zhang   PetscCount   coo_n; /* Number of COOs passed to MatSetPreallocationCOO)() */
72*394ed5ebSJunchao Zhang   PetscSF      coo_sf; /* SF to send/recv remote values in MatSetValuesCOO() */
73*394ed5ebSJunchao Zhang   PetscCount   Atot1,Annz1,Btot1,Bnnz1; /* Number of local entries (tot1), number of uqique local entries(nnz1) belonging to diag (A) and offdiag (B) */
74*394ed5ebSJunchao Zhang   PetscCount   Atot2,Annz2,Btot2,Bnnz2; /* Number of received entries (tot2), number of uqique received entries(nnz2) belonging to diag (A) and offdiag (B) */
75*394ed5ebSJunchao Zhang   PetscCount   *Aimap1,*Ajmap1,*Aperm1; /* Lengths: [Annz1], [Annz1+1], [Atot1]. Local entries to diag */
76*394ed5ebSJunchao Zhang   PetscCount   *Bimap1,*Bjmap1,*Bperm1; /* Lengths: [Bnnz1], [Bnnz1+1], [Btot1]. Local entries to offdiag */
77*394ed5ebSJunchao Zhang   PetscCount   *Aimap2,*Ajmap2,*Aperm2; /* Lengths: [Annz2], [Annz2+1], [Atot2]. Remote entries to diag */
78*394ed5ebSJunchao Zhang   PetscCount   *Bimap2,*Bjmap2,*Bperm2; /* Lengths: [Bnnz2], [Bnnz2+1], [Btot2]. Remote entries to offdiag */
79*394ed5ebSJunchao Zhang   PetscCount   *Cperm1; /* [sendlen] Permutation to fill MPI send buffer. 'C' for communication */
80*394ed5ebSJunchao Zhang   PetscScalar  *sendbuf,*recvbuf; /* Buffers for remote values in MatSetValuesCOO() */
81*394ed5ebSJunchao Zhang   PetscInt     sendlen,recvlen; /* Lengths (in unit of PetscScalar) of send/recvbuf */
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 *[]);
98baa1ba78SHong Zhang PETSC_INTERN PetscErrorCode MatView_MPIAIJ(Mat,PetscViewer);
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);
10752f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ_Binary(Mat,PetscViewer);
108e3225e9dSHong Zhang PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
1094222ddf1SHong Zhang 
1104222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat);
1114e84afc0SStefano Zampini PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJBACKEND(Mat);
1124e84afc0SStefano Zampini PETSC_INTERN PetscErrorCode MatProductSymbolic_MPIAIJBACKEND(Mat);
1134222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat);
1144222ddf1SHong Zhang 
1154222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ(Mat);
1164222ddf1SHong Zhang 
1174222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_MPIAIJ_MPIAIJ(Mat);
1184222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductNumeric_RARt_MPIAIJ_MPIAIJ(Mat);
1194222ddf1SHong Zhang 
1204222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat);
1214222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat,Mat,PetscReal,Mat);
1224222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
1235a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
124f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
125f8487c73SHong Zhang 
1264222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat);
1275a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);
128f996eeb8SHong Zhang 
1294222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
1305a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
1310d3441aeSHong Zhang 
1324222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat,Mat,PetscReal,Mat);
1334222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat,Mat,PetscReal,Mat);
1344222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat,Mat,PetscReal,Mat);
135a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat,Mat,Mat);
1364a56b808SFande Kong PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat,Mat,Mat);
1374a56b808SFande Kong PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat,Mat,Mat);
138a9d7e0ccSandi selinger 
139a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE)
1404222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
141a4ffb656SHong Zhang #endif
1424e29119aSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPIDense(Mat,MatType,MatReuse,Mat*);
143d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
144d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
145d24d4204SJose E. Roman #endif
1460d3441aeSHong Zhang 
1475a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
1486718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(void*);
1496718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void*);
150c5af039cSHong Zhang 
1515a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
1525a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
153e9ede7d0Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
154904d1e70Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Mat,const PetscInt[],const PetscInt[]);
1555a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);
156187b3c17SHong Zhang 
1574222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat);
1584222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
1596da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
1602bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
161c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
1624222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
1635a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);
164c0d3702cSSatish Balay 
1654416b707SBarry Smith PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat);
16610e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1673a7fca6bSBarry Smith 
168570b7f6dSBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
1695a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*);
17097304618SKris Buschelman #endif
1715a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec);
1725a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*);
17397304618SKris Buschelman 
174001ddc4fSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
175001ddc4fSHong Zhang 
17611bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*);
1777087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);
1783a7fca6bSBarry Smith 
17953dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*);
18053dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat);
18153dd7562SDmitry Karpeev 
182*394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_MPIAIJ(Mat,PetscCount,const PetscInt[],const PetscInt[]);
183*394ed5ebSJunchao Zhang 
184cf742fccSHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
185cf742fccSHong Zhang #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \
186cf742fccSHong Zhang {\
187cf742fccSHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ;\
188cf742fccSHong Zhang   PetscScalar *_aa,_valtmp,*_pa;\
189cf742fccSHong Zhang   _apJ = apj + api[i];\
190cf742fccSHong Zhang   /* diagonal portion of A */\
191cf742fccSHong Zhang   _ai  = ad->i;\
192cf742fccSHong Zhang   _anz = _ai[i+1] - _ai[i];\
193cf742fccSHong Zhang   _aj  = ad->j + _ai[i];\
194cf742fccSHong Zhang   _aa  = ad->a + _ai[i];\
195cf742fccSHong Zhang   for (_j=0; _j<_anz; _j++) {\
196cf742fccSHong Zhang     _row = _aj[_j]; \
197cf742fccSHong Zhang     _pi  = p_loc->i;                             \
198cf742fccSHong Zhang     _pnz = _pi[_row+1] - _pi[_row];              \
199cf742fccSHong Zhang     _pj  = p_loc->j + _pi[_row];                 \
200cf742fccSHong Zhang     _pa  = p_loc->a + _pi[_row];                 \
201cf742fccSHong Zhang     /* perform sparse axpy */                    \
202cf742fccSHong Zhang     _valtmp = _aa[_j];                           \
203cf742fccSHong Zhang     _nextp  = 0; \
204cf742fccSHong Zhang     for (_k=0; _nextp<_pnz; _k++) {                    \
205cf742fccSHong Zhang       if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
206cf742fccSHong Zhang         apa[_k] += _valtmp*_pa[_nextp++];                             \
207cf742fccSHong Zhang       } \
208cf742fccSHong Zhang     }                                           \
20946eb3cd7SSatish Balay     (void)PetscLogFlops(2.0*_pnz);              \
210cf742fccSHong Zhang   }                                             \
211cf742fccSHong Zhang   /* off-diagonal portion of A */               \
21292ec70a1SHong Zhang   if (p_oth){ \
213cf742fccSHong Zhang     _ai  = ao->i;\
214cf742fccSHong Zhang     _anz = _ai[i+1] - _ai[i];                   \
215cf742fccSHong Zhang     _aj  = ao->j + _ai[i];                      \
216cf742fccSHong Zhang     _aa  = ao->a + _ai[i];                      \
217cf742fccSHong Zhang     for (_j=0; _j<_anz; _j++) {                 \
218cf742fccSHong Zhang       _row = _aj[_j];    \
219cf742fccSHong Zhang       _pi  = p_oth->i;                         \
220cf742fccSHong Zhang       _pnz = _pi[_row+1] - _pi[_row];          \
221cf742fccSHong Zhang       _pj  = p_oth->j + _pi[_row];             \
222cf742fccSHong Zhang       _pa  = p_oth->a + _pi[_row];             \
223cf742fccSHong Zhang       /* perform sparse axpy */                \
224cf742fccSHong Zhang       _valtmp = _aa[_j];                       \
225cf742fccSHong Zhang       _nextp  = 0; \
226cf742fccSHong Zhang       for (_k=0; _nextp<_pnz; _k++) {          \
227cf742fccSHong Zhang         if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
228cf742fccSHong Zhang           apa[_k] += _valtmp*_pa[_nextp++];    \
229cf742fccSHong Zhang         }                                      \
230cf742fccSHong Zhang       }                                        \
23146eb3cd7SSatish Balay       (void)PetscLogFlops(2.0*_pnz);           \
232cf742fccSHong Zhang     } \
23392ec70a1SHong Zhang   }\
234cf742fccSHong Zhang }
235cf742fccSHong Zhang 
236e497d3c8SHong Zhang #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
237e497d3c8SHong Zhang {\
238e497d3c8SHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;\
239e497d3c8SHong Zhang   PetscScalar *_aa,_valtmp,*_pa;                       \
240e497d3c8SHong Zhang   /* diagonal portion of A */\
241e497d3c8SHong Zhang   _ai  = ad->i;\
242e497d3c8SHong Zhang   _anz = _ai[i+1] - _ai[i];\
243e497d3c8SHong Zhang   _aj  = ad->j + _ai[i];\
244e497d3c8SHong Zhang   _aa  = ad->a + _ai[i];\
245e497d3c8SHong Zhang   for (_j=0; _j<_anz; _j++) {\
246e497d3c8SHong Zhang     _row = _aj[_j]; \
247e497d3c8SHong Zhang     _pi  = p_loc->i;                        \
248e497d3c8SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];         \
249e497d3c8SHong Zhang     _pj  = p_loc->j + _pi[_row];            \
250e497d3c8SHong Zhang     _pa  = p_loc->a + _pi[_row];            \
251e497d3c8SHong Zhang     /* perform dense axpy */                \
252e497d3c8SHong Zhang     _valtmp = _aa[_j];                      \
253e497d3c8SHong Zhang     for (_k=0; _k<_pnz; _k++) {             \
254e497d3c8SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];      \
255e497d3c8SHong Zhang     }                                       \
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];                    \
272e497d3c8SHong Zhang       for (_k=0; _k<_pnz; _k++) {           \
273e497d3c8SHong Zhang         apa[_pj[_k]] += _valtmp*_pa[_k];    \
274e497d3c8SHong Zhang       }                                     \
27546eb3cd7SSatish Balay       (void)PetscLogFlops(2.0*_pnz);        \
276e497d3c8SHong Zhang     }                                       \
27792ec70a1SHong Zhang   }\
278e497d3c8SHong Zhang }
279cf742fccSHong Zhang 
2803369ce9aSBarry Smith #endif
281