xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision 8403a639082ce6decbea2cc545ff80e0678a8e23)
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>
682c86c8fSBarry Smith #include <petscctable.h>
78a729477SBarry Smith 
890431a8fSHong Zhang typedef struct { /* used by MatCreateMPIAIJSumSeqAIJ for reusing the merged matrix */
926283091SBarry Smith   PetscLayout rowmap;
10de0260b3SHong Zhang   PetscInt    **buf_ri,**buf_rj;
110febcb4bSHong Zhang   PetscMPIInt *len_s,*len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
120febcb4bSHong Zhang   PetscMPIInt nsend,nrecv;
13fc08c53fSHong Zhang   PetscInt    *bi,*bj;    /* i and j array of the local portion of mpi C (matrix product) - rename to ci, cj! */
14de0260b3SHong Zhang   PetscInt    *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
15dce485f0SBarry Smith   PetscErrorCode (*destroy)(Mat);
16dce485f0SBarry Smith   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
17b90dcfe3SHong Zhang } Mat_Merge_SeqsToMPI;
18b90dcfe3SHong Zhang 
19a1a4e84aSHong Zhang typedef struct { /* used by MatPtAP_MPIAIJ_MPIAIJ() and MatMatMult_MPIAIJ_MPIAIJ() */
20b7f45c76SHong Zhang   PetscInt    *startsj_s,*startsj_r;    /* used by MatGetBrowsOfAoCols_MPIAIJ */
21b7f45c76SHong Zhang   PetscScalar *bufa;                    /* used by MatGetBrowsOfAoCols_MPIAIJ */
22a1a4e84aSHong Zhang   Mat         P_loc,P_oth;     /* partial B_seq -- intend to replace B_seq */
23a1a4e84aSHong Zhang   PetscInt    *api,*apj;       /* symbolic i and j arrays of the local product A_loc*B_seq */
24445158ffSHong Zhang   PetscScalar *apv;
25d6ab1dc0SHong Zhang   PetscInt    rmax;            /* max num of nnz in a local row of the matrix product */
26b7f45c76SHong Zhang   MatReuse    reuse;           /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */
27598bc09dSHong Zhang   PetscScalar *apa;            /* tmp array for store a row of A*P used in MatMatMult() */
28c5af039cSHong Zhang   Mat         A_loc;           /* used by MatTransposeMatMult(), contains api and apj */
29ab1b013aSHong Zhang   Mat         Pt;              /* used by MatTransposeMatMult(), Pt = P^T */
30414894bdSHong Zhang   PetscBool   scalable;        /* flag determines scalable or non-scalable implementation */
319ce11a7cSHong Zhang   Mat         Rd,Ro,AP_loc,C_loc,C_oth;
32f8487c73SHong Zhang 
33f8487c73SHong Zhang   Mat_Merge_SeqsToMPI *merge;
34f8487c73SHong Zhang   PetscErrorCode (*destroy)(Mat);
354ae0847bSHong Zhang   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
36f8487c73SHong Zhang } Mat_PtAPMPI;
37f8487c73SHong Zhang 
38f8487c73SHong Zhang typedef struct {
39f8487c73SHong Zhang   Mat A,B;                             /* local submatrices: A (diag part),
40f8487c73SHong Zhang                                            B (off-diag part) */
41f8487c73SHong Zhang   PetscMPIInt size;                     /* size of communicator */
42f8487c73SHong Zhang   PetscMPIInt rank;                     /* rank of proc in communicator */
43f8487c73SHong Zhang 
44f8487c73SHong Zhang   /* The following variables are used for matrix assembly */
45f8487c73SHong Zhang   PetscBool   donotstash;               /* PETSC_TRUE if off processor entries dropped */
46f8487c73SHong Zhang   MPI_Request *send_waits;              /* array of send requests */
47f8487c73SHong Zhang   MPI_Request *recv_waits;              /* array of receive requests */
48f8487c73SHong Zhang   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
49f8487c73SHong Zhang   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
50f8487c73SHong Zhang   PetscInt    rmax;                     /* maximum message length */
51f8487c73SHong Zhang #if defined(PETSC_USE_CTABLE)
52f8487c73SHong Zhang   PetscTable colmap;
53f8487c73SHong Zhang #else
54f8487c73SHong Zhang   PetscInt *colmap;                     /* local col number of off-diag col */
55f8487c73SHong Zhang #endif
56f8487c73SHong Zhang   PetscInt *garray;                     /* global index of all off-processor columns */
57f8487c73SHong Zhang 
58f8487c73SHong Zhang   /* The following variables are used for matrix-vector products */
59f8487c73SHong Zhang   Vec        lvec;                 /* local vector */
60f8487c73SHong Zhang   Vec        diag;
61f8487c73SHong Zhang   VecScatter Mvctx;                /* scatter context for vector */
62f8487c73SHong Zhang   PetscBool  roworiented;          /* if true, row-oriented input, default true */
63f8487c73SHong Zhang 
64f8487c73SHong Zhang   /* The following variables are for MatGetRow() */
65f8487c73SHong Zhang   PetscInt    *rowindices;         /* column indices for row */
66f8487c73SHong Zhang   PetscScalar *rowvalues;          /* nonzero values in row */
67f8487c73SHong Zhang   PetscBool   getrowactive;        /* indicates MatGetRow(), not restored */
68f8487c73SHong Zhang 
69f8487c73SHong Zhang   /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation  and nonzero pattern */
70f8487c73SHong Zhang   PetscInt *ld;                    /* number of entries per row left of diagona block */
71f8487c73SHong Zhang 
72a1a4e84aSHong Zhang   /* Used by MatMatMult() and MatPtAP() */
73f8487c73SHong Zhang   Mat_PtAPMPI *ptap;
74ca45077fSPaul Mullowney 
75f996eeb8SHong Zhang   /* used by MatMatMatMult() */
76f996eeb8SHong Zhang   Mat_MatMatMatMult *matmatmatmult;
77f996eeb8SHong Zhang 
78bbf3fe20SPaul Mullowney   /* Used by MPICUSP and MPICUSPARSE classes */
79bbf3fe20SPaul Mullowney   void * spptr;
805cc03489SHong Zhang 
81f8487c73SHong Zhang } Mat_MPIAIJ;
82f8487c73SHong Zhang 
838cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
842327d61dSSatish Balay 
855a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetColoring_MPIAIJ(Mat,ISColoring);
865a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValuesAdifor_MPIAIJ(Mat,PetscInt,void*);
875a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
885a576424SJed Brown PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
895a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*);
905a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
9193dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
92f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
935a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
9453dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
955a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_All(Mat,MatGetSubMatrixOption,MatReuse,Mat *[]);
9653dd7562SDmitry Karpeev 
975494a064SHong Zhang 
985a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*);
995a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_Private (Mat,IS,IS,PetscInt,MatReuse,Mat*);
1003233f0ddSHong Zhang PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*);
101d6037b41SHong Zhang 
1025a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer);
103e3225e9dSHong Zhang PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1055a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1060fc8cf34SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
107b2405163SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
109f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
110f8487c73SHong Zhang 
1115a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
1125a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*);
1135a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);
114f996eeb8SHong Zhang 
1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1165a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1175a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
1180d3441aeSHong Zhang 
119*8403a639SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat,Mat,PetscReal,Mat*);
120*8403a639SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat,Mat,Mat);
1210d3441aeSHong Zhang 
1225a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat);
1235a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
124c5af039cSHong Zhang 
1255a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
1265a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
1275a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat);
1285a576424SJed Brown PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*);
1295a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);
130187b3c17SHong Zhang 
1315a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
1322bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
1336da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
1346da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
1352bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
136c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
1378949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
1388949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
1398949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
1405a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);
141c0d3702cSSatish Balay 
14210e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1433a7fca6bSBarry Smith 
144ce63c4c1SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
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 
158e497d3c8SHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth */
159e497d3c8SHong Zhang #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
160e497d3c8SHong Zhang {\
161e497d3c8SHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;      \
162e497d3c8SHong Zhang   PetscScalar *_aa,_valtmp,*_pa;                             \
163e497d3c8SHong Zhang   /* diagonal portion of A */\
164e497d3c8SHong Zhang   _ai  = ad->i;\
165e497d3c8SHong Zhang   _anz = _ai[i+1] - _ai[i];\
166e497d3c8SHong Zhang   _aj  = ad->j + _ai[i];\
167e497d3c8SHong Zhang   _aa  = ad->a + _ai[i];\
168e497d3c8SHong Zhang   for (_j=0; _j<_anz; _j++) {\
169e497d3c8SHong Zhang     _row = _aj[_j]; \
170e497d3c8SHong Zhang     _pi  = p_loc->i;                                 \
171e497d3c8SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];         \
172e497d3c8SHong Zhang     _pj  = p_loc->j + _pi[_row];                 \
173e497d3c8SHong Zhang     _pa  = p_loc->a + _pi[_row];                 \
174e497d3c8SHong Zhang     /* perform dense axpy */                    \
175e497d3c8SHong Zhang     _valtmp = _aa[_j];                           \
176e497d3c8SHong Zhang     for (_k=0; _k<_pnz; _k++) {                    \
177e497d3c8SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];               \
178e497d3c8SHong Zhang     }                                           \
179e497d3c8SHong Zhang     PetscLogFlops(2.0*_pnz);                    \
180e497d3c8SHong Zhang   }                                             \
181e497d3c8SHong Zhang   /* off-diagonal portion of A */               \
182e497d3c8SHong Zhang   _ai  = ao->i;\
183e497d3c8SHong Zhang   _anz = _ai[i+1] - _ai[i];                     \
184e497d3c8SHong Zhang   _aj  = ao->j + _ai[i];                         \
185e497d3c8SHong Zhang   _aa  = ao->a + _ai[i];                         \
186e497d3c8SHong Zhang   for (_j=0; _j<_anz; _j++) {                      \
187e497d3c8SHong Zhang     _row = _aj[_j];    \
188e497d3c8SHong Zhang     _pi  = p_oth->i;                         \
189e497d3c8SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];          \
190e497d3c8SHong Zhang     _pj  = p_oth->j + _pi[_row];                  \
191e497d3c8SHong Zhang     _pa  = p_oth->a + _pi[_row];                  \
192e497d3c8SHong Zhang     /* perform dense axpy */                     \
193e497d3c8SHong Zhang     _valtmp = _aa[_j];                             \
194e497d3c8SHong Zhang     for (_k=0; _k<_pnz; _k++) {                     \
195e497d3c8SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];                \
196e497d3c8SHong Zhang     }                                            \
197e497d3c8SHong Zhang     PetscLogFlops(2.0*_pnz);                     \
198e497d3c8SHong Zhang   } \
199e497d3c8SHong Zhang }
2003369ce9aSBarry Smith #endif
201