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; 30f8487c73SHong Zhang 31f8487c73SHong Zhang Mat_Merge_SeqsToMPI *merge; 32f8487c73SHong Zhang PetscErrorCode (*destroy)(Mat); 334ae0847bSHong Zhang PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*); 34f8487c73SHong Zhang } Mat_PtAPMPI; 35f8487c73SHong Zhang 36f8487c73SHong Zhang typedef struct { 37f8487c73SHong Zhang Mat A,B; /* local submatrices: A (diag part), 38f8487c73SHong Zhang B (off-diag part) */ 39f8487c73SHong Zhang PetscMPIInt size; /* size of communicator */ 40f8487c73SHong Zhang PetscMPIInt rank; /* rank of proc in communicator */ 41f8487c73SHong Zhang 42f8487c73SHong Zhang /* The following variables are used for matrix assembly */ 43f8487c73SHong Zhang PetscBool donotstash; /* PETSC_TRUE if off processor entries dropped */ 44f8487c73SHong Zhang MPI_Request *send_waits; /* array of send requests */ 45f8487c73SHong Zhang MPI_Request *recv_waits; /* array of receive requests */ 46f8487c73SHong Zhang PetscInt nsends,nrecvs; /* numbers of sends and receives */ 47f8487c73SHong Zhang PetscScalar *svalues,*rvalues; /* sending and receiving data */ 48f8487c73SHong Zhang PetscInt rmax; /* maximum message length */ 49f8487c73SHong Zhang #if defined(PETSC_USE_CTABLE) 50f8487c73SHong Zhang PetscTable colmap; 51f8487c73SHong Zhang #else 52f8487c73SHong Zhang PetscInt *colmap; /* local col number of off-diag col */ 53f8487c73SHong Zhang #endif 54f8487c73SHong Zhang PetscInt *garray; /* global index of all off-processor columns */ 55f8487c73SHong Zhang 56f8487c73SHong Zhang /* The following variables are used for matrix-vector products */ 57f8487c73SHong Zhang Vec lvec; /* local vector */ 58f8487c73SHong Zhang Vec diag; 59f8487c73SHong Zhang VecScatter Mvctx; /* scatter context for vector */ 60f8487c73SHong Zhang PetscBool roworiented; /* if true, row-oriented input, default true */ 61f8487c73SHong Zhang 62f8487c73SHong Zhang /* The following variables are for MatGetRow() */ 63f8487c73SHong Zhang PetscInt *rowindices; /* column indices for row */ 64f8487c73SHong Zhang PetscScalar *rowvalues; /* nonzero values in row */ 65f8487c73SHong Zhang PetscBool getrowactive; /* indicates MatGetRow(), not restored */ 66f8487c73SHong Zhang 67f8487c73SHong Zhang /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation and nonzero pattern */ 68f8487c73SHong Zhang PetscInt *ld; /* number of entries per row left of diagona block */ 69f8487c73SHong Zhang 70a1a4e84aSHong Zhang /* Used by MatMatMult() and MatPtAP() */ 71f8487c73SHong Zhang Mat_PtAPMPI *ptap; 72ca45077fSPaul Mullowney 73f996eeb8SHong Zhang /* used by MatMatMatMult() */ 74f996eeb8SHong Zhang Mat_MatMatMatMult *matmatmatmult; 75f996eeb8SHong Zhang 76bbf3fe20SPaul Mullowney /* Used by MPICUSP and MPICUSPARSE classes */ 77bbf3fe20SPaul Mullowney void * spptr; 785cc03489SHong Zhang 79f8487c73SHong Zhang } Mat_MPIAIJ; 80f8487c73SHong Zhang 818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat); 822327d61dSSatish Balay 835a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat); 845a576424SJed Brown PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat); 855a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*); 865a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt); 87b1b1104fSBarry Smith PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat,PetscInt,IS [],PetscInt); 8893dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring); 89f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring); 907dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 917dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 927dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat,MatCreateSubMatrixOption,MatReuse,Mat *[]); 9353dd7562SDmitry Karpeev 945494a064SHong Zhang 957dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*); 96618cbb4aSHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_nonscalable(Mat,IS,IS,PetscInt,MatReuse,Mat*); 97*3b00a383SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowDist(Mat,IS,IS,MatReuse,Mat*); 98*3b00a383SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowColDist(Mat,IS,IS,MatReuse,Mat*); 993233f0ddSHong Zhang PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*); 100d6037b41SHong Zhang 1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer); 102e3225e9dSHong Zhang PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat); 1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1050fc8cf34SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 106b2405163SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 108f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 109f8487c73SHong Zhang 1105a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 1115a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*); 1125a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat); 113f996eeb8SHong Zhang 1145a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1165a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 1170d3441aeSHong Zhang 1188403a639SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat,Mat,PetscReal,Mat*); 1198403a639SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat,Mat,Mat); 1200d3441aeSHong Zhang 1215a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat); 1225a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 123c5af039cSHong Zhang 1245a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 1255a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode); 1265a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat); 1275a576424SJed Brown PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*); 1285a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool); 129187b3c17SHong Zhang 1305a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1312bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 1326da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1336da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 1342bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 135c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat); 1368949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 1378949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 1388949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 1395a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*); 140c0d3702cSSatish Balay 1414416b707SBarry Smith PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat); 14210e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1433a7fca6bSBarry Smith 144570b7f6dSBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16) 1455a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*); 14697304618SKris Buschelman #endif 1475a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec); 1485a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*); 14997304618SKris Buschelman 150001ddc4fSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 151001ddc4fSHong Zhang 15211bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*); 1537087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec); 1543a7fca6bSBarry Smith 15553dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*); 15653dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat); 15753dd7562SDmitry Karpeev 158cf742fccSHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */ 159cf742fccSHong Zhang #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \ 160cf742fccSHong Zhang {\ 161cf742fccSHong Zhang PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ; \ 162cf742fccSHong Zhang PetscScalar *_aa,_valtmp,*_pa; \ 163cf742fccSHong Zhang _apJ = apj + api[i];\ 164cf742fccSHong Zhang /* diagonal portion of A */\ 165cf742fccSHong Zhang _ai = ad->i;\ 166cf742fccSHong Zhang _anz = _ai[i+1] - _ai[i];\ 167cf742fccSHong Zhang _aj = ad->j + _ai[i];\ 168cf742fccSHong Zhang _aa = ad->a + _ai[i];\ 169cf742fccSHong Zhang for (_j=0; _j<_anz; _j++) {\ 170cf742fccSHong Zhang _row = _aj[_j]; \ 171cf742fccSHong Zhang _pi = p_loc->i; \ 172cf742fccSHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 173cf742fccSHong Zhang _pj = p_loc->j + _pi[_row]; \ 174cf742fccSHong Zhang _pa = p_loc->a + _pi[_row]; \ 175cf742fccSHong Zhang /* perform sparse axpy */ \ 176cf742fccSHong Zhang _valtmp = _aa[_j]; \ 177cf742fccSHong Zhang _nextp = 0; \ 178cf742fccSHong Zhang for (_k=0; _nextp<_pnz; _k++) { \ 179cf742fccSHong Zhang if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */ \ 180cf742fccSHong Zhang apa[_k] += _valtmp*_pa[_nextp++]; \ 181cf742fccSHong Zhang } \ 182cf742fccSHong Zhang } \ 183cf742fccSHong Zhang PetscLogFlops(2.0*_pnz); \ 184cf742fccSHong Zhang } \ 185cf742fccSHong Zhang /* off-diagonal portion of A */ \ 186cf742fccSHong Zhang _ai = ao->i;\ 187cf742fccSHong Zhang _anz = _ai[i+1] - _ai[i]; \ 188cf742fccSHong Zhang _aj = ao->j + _ai[i]; \ 189cf742fccSHong Zhang _aa = ao->a + _ai[i]; \ 190cf742fccSHong Zhang for (_j=0; _j<_anz; _j++) { \ 191cf742fccSHong Zhang _row = _aj[_j]; \ 192cf742fccSHong Zhang _pi = p_oth->i; \ 193cf742fccSHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 194cf742fccSHong Zhang _pj = p_oth->j + _pi[_row]; \ 195cf742fccSHong Zhang _pa = p_oth->a + _pi[_row]; \ 196cf742fccSHong Zhang /* perform sparse axpy */ \ 197cf742fccSHong Zhang _valtmp = _aa[_j]; \ 198cf742fccSHong Zhang _nextp = 0; \ 199cf742fccSHong Zhang for (_k=0; _nextp<_pnz; _k++) { \ 200cf742fccSHong Zhang if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\ 201cf742fccSHong Zhang apa[_k] += _valtmp*_pa[_nextp++]; \ 202cf742fccSHong Zhang } \ 203cf742fccSHong Zhang } \ 204cf742fccSHong Zhang PetscLogFlops(2.0*_pnz); \ 205cf742fccSHong Zhang } \ 206cf742fccSHong Zhang } 207cf742fccSHong Zhang 208e497d3c8SHong Zhang #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \ 209e497d3c8SHong Zhang {\ 210e497d3c8SHong Zhang PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj; \ 211e497d3c8SHong Zhang PetscScalar *_aa,_valtmp,*_pa; \ 212e497d3c8SHong Zhang /* diagonal portion of A */\ 213e497d3c8SHong Zhang _ai = ad->i;\ 214e497d3c8SHong Zhang _anz = _ai[i+1] - _ai[i];\ 215e497d3c8SHong Zhang _aj = ad->j + _ai[i];\ 216e497d3c8SHong Zhang _aa = ad->a + _ai[i];\ 217e497d3c8SHong Zhang for (_j=0; _j<_anz; _j++) {\ 218e497d3c8SHong Zhang _row = _aj[_j]; \ 219e497d3c8SHong Zhang _pi = p_loc->i; \ 220e497d3c8SHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 221e497d3c8SHong Zhang _pj = p_loc->j + _pi[_row]; \ 222e497d3c8SHong Zhang _pa = p_loc->a + _pi[_row]; \ 223e497d3c8SHong Zhang /* perform dense axpy */ \ 224e497d3c8SHong Zhang _valtmp = _aa[_j]; \ 225e497d3c8SHong Zhang for (_k=0; _k<_pnz; _k++) { \ 226e497d3c8SHong Zhang apa[_pj[_k]] += _valtmp*_pa[_k]; \ 227e497d3c8SHong Zhang } \ 228e497d3c8SHong Zhang PetscLogFlops(2.0*_pnz); \ 229e497d3c8SHong Zhang } \ 230e497d3c8SHong Zhang /* off-diagonal portion of A */ \ 231e497d3c8SHong Zhang _ai = ao->i;\ 232e497d3c8SHong Zhang _anz = _ai[i+1] - _ai[i]; \ 233e497d3c8SHong Zhang _aj = ao->j + _ai[i]; \ 234e497d3c8SHong Zhang _aa = ao->a + _ai[i]; \ 235e497d3c8SHong Zhang for (_j=0; _j<_anz; _j++) { \ 236e497d3c8SHong Zhang _row = _aj[_j]; \ 237e497d3c8SHong Zhang _pi = p_oth->i; \ 238e497d3c8SHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 239e497d3c8SHong Zhang _pj = p_oth->j + _pi[_row]; \ 240e497d3c8SHong Zhang _pa = p_oth->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 } \ 246e497d3c8SHong Zhang PetscLogFlops(2.0*_pnz); \ 247e497d3c8SHong Zhang } \ 248e497d3c8SHong Zhang } 249cf742fccSHong Zhang 2503369ce9aSBarry Smith #endif 251