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 */ 27*a3bb6f32SFande Kong ISLocalToGlobalMapping ltog; /* mapping from local column indices to global column indices for A_loc */ 28ab1b013aSHong Zhang Mat Pt; /* used by MatTransposeMatMult(), Pt = P^T */ 297d0a89b7SHong Zhang PetscBool freestruct; /* flag for MatFreeIntermediateDataStructures() */ 309ce11a7cSHong Zhang Mat Rd,Ro,AP_loc,C_loc,C_oth; 31cf97cf3cSHong Zhang PetscInt algType; /* implementation algorithm */ 32f8487c73SHong Zhang 33f8487c73SHong Zhang Mat_Merge_SeqsToMPI *merge; 34f8487c73SHong Zhang PetscErrorCode (*destroy)(Mat); 354ae0847bSHong Zhang PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*); 36cf97cf3cSHong Zhang PetscErrorCode (*view)(Mat,PetscViewer); 373cdca5ebSHong Zhang } Mat_APMPI; 38f8487c73SHong Zhang 39f8487c73SHong Zhang typedef struct { 40f8487c73SHong Zhang Mat A,B; /* local submatrices: A (diag part), 41f8487c73SHong Zhang B (off-diag part) */ 42f8487c73SHong Zhang PetscMPIInt size; /* size of communicator */ 43f8487c73SHong Zhang PetscMPIInt rank; /* rank of proc in communicator */ 44f8487c73SHong Zhang 45f8487c73SHong Zhang /* The following variables are used for matrix assembly */ 46f8487c73SHong Zhang PetscBool donotstash; /* PETSC_TRUE if off processor entries dropped */ 47f8487c73SHong Zhang MPI_Request *send_waits; /* array of send requests */ 48f8487c73SHong Zhang MPI_Request *recv_waits; /* array of receive requests */ 49f8487c73SHong Zhang PetscInt nsends,nrecvs; /* numbers of sends and receives */ 50f8487c73SHong Zhang PetscScalar *svalues,*rvalues; /* sending and receiving data */ 51f8487c73SHong Zhang PetscInt rmax; /* maximum message length */ 52f8487c73SHong Zhang #if defined(PETSC_USE_CTABLE) 53f8487c73SHong Zhang PetscTable colmap; 54f8487c73SHong Zhang #else 55f8487c73SHong Zhang PetscInt *colmap; /* local col number of off-diag col */ 56f8487c73SHong Zhang #endif 57f8487c73SHong Zhang PetscInt *garray; /* global index of all off-processor columns */ 58f8487c73SHong Zhang 59f8487c73SHong Zhang /* The following variables are used for matrix-vector products */ 60f8487c73SHong Zhang Vec lvec; /* local vector */ 61f8487c73SHong Zhang Vec diag; 624b8d542aSHong Zhang VecScatter Mvctx,Mvctx_mpi1; /* scatter context for vector */ 63fa110396SHong Zhang PetscBool Mvctx_mpi1_flg; /* if true, additional Mvctx_mpi1 is requested for mat-mat ops, default false */ 64f8487c73SHong Zhang PetscBool roworiented; /* if true, row-oriented input, default true */ 65f8487c73SHong Zhang 66f8487c73SHong Zhang /* The following variables are for MatGetRow() */ 67f8487c73SHong Zhang PetscInt *rowindices; /* column indices for row */ 68f8487c73SHong Zhang PetscScalar *rowvalues; /* nonzero values in row */ 69f8487c73SHong Zhang PetscBool getrowactive; /* indicates MatGetRow(), not restored */ 70f8487c73SHong Zhang 71f8487c73SHong Zhang /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation and nonzero pattern */ 72f8487c73SHong Zhang PetscInt *ld; /* number of entries per row left of diagona block */ 73f8487c73SHong Zhang 74a1a4e84aSHong Zhang /* Used by MatMatMult() and MatPtAP() */ 753cdca5ebSHong Zhang Mat_APMPI *ap; 76ca45077fSPaul Mullowney 77f996eeb8SHong Zhang /* used by MatMatMatMult() */ 78f996eeb8SHong Zhang Mat_MatMatMatMult *matmatmatmult; 79f996eeb8SHong Zhang 80bbf3fe20SPaul Mullowney /* Used by MPICUSP and MPICUSPARSE classes */ 81bbf3fe20SPaul Mullowney void * spptr; 825cc03489SHong Zhang 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); 1105a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1115a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1120fc8cf34SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 113a9d7e0ccSandi selinger PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat,Mat,PetscReal,Mat*); 114b2405163SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 116f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 117f8487c73SHong Zhang 1185a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 1195a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*); 1205a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat); 121f996eeb8SHong Zhang 1225a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1235a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1245a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 1250d3441aeSHong Zhang 126a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat,Mat,PetscReal,Mat*); 127a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat,Mat,Mat); 1284624976aSHong Zhang PETSC_INTERN PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_AP(Mat); 1298a9b8717SHong Zhang PETSC_INTERN PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_BC(Mat); 130a9d7e0ccSandi selinger 131a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE) 132a4ffb656SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*); 133a4ffb656SHong Zhang #endif 1340d3441aeSHong Zhang 1355a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat); 1365a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 137c5af039cSHong Zhang 1388d45306eSHong Zhang PETSC_INTERN PetscErrorCode MatRARt_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1398d45306eSHong Zhang 1405a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 1415a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode); 142e9ede7d0Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 143904d1e70Sandi selinger PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Mat,const PetscInt[],const PetscInt[]); 1445a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat); 1455a576424SJed Brown PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*); 1465a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool); 147187b3c17SHong Zhang 1485a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1492bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 1506da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1516da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 1522bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 153c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat); 1548949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 1558949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 1568949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 1575a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*); 158c0d3702cSSatish Balay 1594416b707SBarry Smith PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat); 16010e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1613a7fca6bSBarry Smith 162570b7f6dSBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16) 1635a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*); 16497304618SKris Buschelman #endif 1655a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec); 1665a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*); 16797304618SKris Buschelman 168001ddc4fSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 169001ddc4fSHong Zhang 17011bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*); 1717087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec); 1723a7fca6bSBarry Smith 17353dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*); 17453dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat); 17553dd7562SDmitry Karpeev 176cf742fccSHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */ 177cf742fccSHong Zhang #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \ 178cf742fccSHong Zhang {\ 179cf742fccSHong Zhang PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ; \ 180cf742fccSHong Zhang PetscScalar *_aa,_valtmp,*_pa; \ 181cf742fccSHong Zhang _apJ = apj + api[i];\ 182cf742fccSHong Zhang /* diagonal portion of A */\ 183cf742fccSHong Zhang _ai = ad->i;\ 184cf742fccSHong Zhang _anz = _ai[i+1] - _ai[i];\ 185cf742fccSHong Zhang _aj = ad->j + _ai[i];\ 186cf742fccSHong Zhang _aa = ad->a + _ai[i];\ 187cf742fccSHong Zhang for (_j=0; _j<_anz; _j++) {\ 188cf742fccSHong Zhang _row = _aj[_j]; \ 189cf742fccSHong Zhang _pi = p_loc->i; \ 190cf742fccSHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 191cf742fccSHong Zhang _pj = p_loc->j + _pi[_row]; \ 192cf742fccSHong Zhang _pa = p_loc->a + _pi[_row]; \ 193cf742fccSHong Zhang /* perform sparse axpy */ \ 194cf742fccSHong Zhang _valtmp = _aa[_j]; \ 195cf742fccSHong Zhang _nextp = 0; \ 196cf742fccSHong Zhang for (_k=0; _nextp<_pnz; _k++) { \ 197cf742fccSHong Zhang if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */ \ 198cf742fccSHong Zhang apa[_k] += _valtmp*_pa[_nextp++]; \ 199cf742fccSHong Zhang } \ 200cf742fccSHong Zhang } \ 20146eb3cd7SSatish Balay (void)PetscLogFlops(2.0*_pnz); \ 202cf742fccSHong Zhang } \ 203cf742fccSHong Zhang /* off-diagonal portion of A */ \ 204cf742fccSHong Zhang _ai = ao->i;\ 205cf742fccSHong Zhang _anz = _ai[i+1] - _ai[i]; \ 206cf742fccSHong Zhang _aj = ao->j + _ai[i]; \ 207cf742fccSHong Zhang _aa = ao->a + _ai[i]; \ 208cf742fccSHong Zhang for (_j=0; _j<_anz; _j++) { \ 209cf742fccSHong Zhang _row = _aj[_j]; \ 210cf742fccSHong Zhang _pi = p_oth->i; \ 211cf742fccSHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 212cf742fccSHong Zhang _pj = p_oth->j + _pi[_row]; \ 213cf742fccSHong Zhang _pa = p_oth->a + _pi[_row]; \ 214cf742fccSHong Zhang /* perform sparse axpy */ \ 215cf742fccSHong Zhang _valtmp = _aa[_j]; \ 216cf742fccSHong Zhang _nextp = 0; \ 217cf742fccSHong Zhang for (_k=0; _nextp<_pnz; _k++) { \ 218cf742fccSHong Zhang if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\ 219cf742fccSHong Zhang apa[_k] += _valtmp*_pa[_nextp++]; \ 220cf742fccSHong Zhang } \ 221cf742fccSHong Zhang } \ 22246eb3cd7SSatish Balay (void)PetscLogFlops(2.0*_pnz); \ 223cf742fccSHong Zhang } \ 224cf742fccSHong Zhang } 225cf742fccSHong Zhang 226e497d3c8SHong Zhang #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \ 227e497d3c8SHong Zhang {\ 228e497d3c8SHong Zhang PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj; \ 229e497d3c8SHong Zhang PetscScalar *_aa,_valtmp,*_pa; \ 230e497d3c8SHong Zhang /* diagonal portion of A */\ 231e497d3c8SHong Zhang _ai = ad->i;\ 232e497d3c8SHong Zhang _anz = _ai[i+1] - _ai[i];\ 233e497d3c8SHong Zhang _aj = ad->j + _ai[i];\ 234e497d3c8SHong Zhang _aa = ad->a + _ai[i];\ 235e497d3c8SHong Zhang for (_j=0; _j<_anz; _j++) {\ 236e497d3c8SHong Zhang _row = _aj[_j]; \ 237e497d3c8SHong Zhang _pi = p_loc->i; \ 238e497d3c8SHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 239e497d3c8SHong Zhang _pj = p_loc->j + _pi[_row]; \ 240e497d3c8SHong Zhang _pa = p_loc->a + _pi[_row]; \ 241e497d3c8SHong Zhang /* perform dense axpy */ \ 242e497d3c8SHong Zhang _valtmp = _aa[_j]; \ 243e497d3c8SHong Zhang for (_k=0; _k<_pnz; _k++) { \ 244e497d3c8SHong Zhang apa[_pj[_k]] += _valtmp*_pa[_k]; \ 245e497d3c8SHong Zhang } \ 24646eb3cd7SSatish Balay (void)PetscLogFlops(2.0*_pnz); \ 247e497d3c8SHong Zhang } \ 248e497d3c8SHong Zhang /* off-diagonal portion of A */ \ 249e497d3c8SHong Zhang _ai = ao->i;\ 250e497d3c8SHong Zhang _anz = _ai[i+1] - _ai[i]; \ 251e497d3c8SHong Zhang _aj = ao->j + _ai[i]; \ 252e497d3c8SHong Zhang _aa = ao->a + _ai[i]; \ 253e497d3c8SHong Zhang for (_j=0; _j<_anz; _j++) { \ 254e497d3c8SHong Zhang _row = _aj[_j]; \ 255e497d3c8SHong Zhang _pi = p_oth->i; \ 256e497d3c8SHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 257e497d3c8SHong Zhang _pj = p_oth->j + _pi[_row]; \ 258e497d3c8SHong Zhang _pa = p_oth->a + _pi[_row]; \ 259e497d3c8SHong Zhang /* perform dense axpy */ \ 260e497d3c8SHong Zhang _valtmp = _aa[_j]; \ 261e497d3c8SHong Zhang for (_k=0; _k<_pnz; _k++) { \ 262e497d3c8SHong Zhang apa[_pj[_k]] += _valtmp*_pa[_k]; \ 263e497d3c8SHong Zhang } \ 26446eb3cd7SSatish Balay (void)PetscLogFlops(2.0*_pnz); \ 265e497d3c8SHong Zhang } \ 266e497d3c8SHong Zhang } 267cf742fccSHong Zhang 2683369ce9aSBarry Smith #endif 269