xref: /petsc/src/mat/impls/sell/mpi/mpisell.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
1d4002b98SHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>   /*I "petscmat.h" I*/
2d4002b98SHong Zhang #include <../src/mat/impls/sell/mpi/mpisell.h>   /*I "petscmat.h" I*/
3d4002b98SHong Zhang #include <petsc/private/vecimpl.h>
4d4002b98SHong Zhang #include <petsc/private/isimpl.h>
5d4002b98SHong Zhang #include <petscblaslapack.h>
6d4002b98SHong Zhang #include <petscsf.h>
7d4002b98SHong Zhang 
8d4002b98SHong Zhang /*MC
9d4002b98SHong Zhang    MATSELL - MATSELL = "sell" - A matrix type to be used for sparse matrices.
10d4002b98SHong Zhang 
11d4002b98SHong Zhang    This matrix type is identical to MATSEQSELL when constructed with a single process communicator,
12d4002b98SHong Zhang    and MATMPISELL otherwise.  As a result, for single process communicators,
13d4002b98SHong Zhang   MatSeqSELLSetPreallocation is supported, and similarly MatMPISELLSetPreallocation is supported
14d4002b98SHong Zhang   for communicators controlling multiple processes.  It is recommended that you call both of
15d4002b98SHong Zhang   the above preallocation routines for simplicity.
16d4002b98SHong Zhang 
17d4002b98SHong Zhang    Options Database Keys:
18d4002b98SHong Zhang . -mat_type sell - sets the matrix type to "sell" during a call to MatSetFromOptions()
19d4002b98SHong Zhang 
20d4002b98SHong Zhang   Level: beginner
21d4002b98SHong Zhang 
22d4002b98SHong Zhang .seealso: MatCreateSELL(), MatCreateSeqSELL(), MATSEQSELL, MATMPISELL
23d4002b98SHong Zhang M*/
24d4002b98SHong Zhang 
25d4002b98SHong Zhang PetscErrorCode MatDiagonalSet_MPISELL(Mat Y,Vec D,InsertMode is)
26d4002b98SHong Zhang {
27d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)Y->data;
28d4002b98SHong Zhang 
29d4002b98SHong Zhang   PetscFunctionBegin;
30d4002b98SHong Zhang   if (Y->assembled && Y->rmap->rstart == Y->cmap->rstart && Y->rmap->rend == Y->cmap->rend) {
319566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(sell->A,D,is));
32d4002b98SHong Zhang   } else {
339566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet_Default(Y,D,is));
34d4002b98SHong Zhang   }
35d4002b98SHong Zhang   PetscFunctionReturn(0);
36d4002b98SHong Zhang }
37d4002b98SHong Zhang 
38d4002b98SHong Zhang /*
39d4002b98SHong Zhang   Local utility routine that creates a mapping from the global column
40d4002b98SHong Zhang number to the local number in the off-diagonal part of the local
41d4002b98SHong Zhang storage of the matrix.  When PETSC_USE_CTABLE is used this is scalable at
42d4002b98SHong Zhang a slightly higher hash table cost; without it it is not scalable (each processor
43d4002b98SHong Zhang has an order N integer array but is fast to acess.
44d4002b98SHong Zhang */
45d4002b98SHong Zhang PetscErrorCode MatCreateColmap_MPISELL_Private(Mat mat)
46d4002b98SHong Zhang {
47d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
48d4002b98SHong Zhang   PetscInt       n=sell->B->cmap->n,i;
49d4002b98SHong Zhang 
50d4002b98SHong Zhang   PetscFunctionBegin;
5128b400f6SJacob Faibussowitsch   PetscCheck(sell->garray,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPISELL Matrix was assembled but is missing garray");
52d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
539566063dSJacob Faibussowitsch   PetscCall(PetscTableCreate(n,mat->cmap->N+1,&sell->colmap));
54d4002b98SHong Zhang   for (i=0; i<n; i++) {
559566063dSJacob Faibussowitsch     PetscCall(PetscTableAdd(sell->colmap,sell->garray[i]+1,i+1,INSERT_VALUES));
56d4002b98SHong Zhang   }
57d4002b98SHong Zhang #else
589566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(mat->cmap->N+1,&sell->colmap));
599566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)mat,(mat->cmap->N+1)*sizeof(PetscInt)));
60d4002b98SHong Zhang   for (i=0; i<n; i++) sell->colmap[sell->garray[i]] = i+1;
61d4002b98SHong Zhang #endif
62d4002b98SHong Zhang   PetscFunctionReturn(0);
63d4002b98SHong Zhang }
64d4002b98SHong Zhang 
65d4002b98SHong Zhang #define MatSetValues_SeqSELL_A_Private(row,col,value,addv,orow,ocol) \
66d4002b98SHong Zhang   { \
67d4002b98SHong Zhang     if (col <= lastcol1) low1 = 0; \
68d4002b98SHong Zhang     else                high1 = nrow1; \
69d4002b98SHong Zhang     lastcol1 = col; \
70d4002b98SHong Zhang     while (high1-low1 > 5) { \
71d4002b98SHong Zhang       t = (low1+high1)/2; \
72d4002b98SHong Zhang       if (*(cp1+8*t) > col) high1 = t; \
73d4002b98SHong Zhang       else                   low1 = t; \
74d4002b98SHong Zhang     } \
75d4002b98SHong Zhang     for (_i=low1; _i<high1; _i++) { \
76d4002b98SHong Zhang       if (*(cp1+8*_i) > col) break; \
77d4002b98SHong Zhang       if (*(cp1+8*_i) == col) { \
78d4002b98SHong Zhang         if (addv == ADD_VALUES) *(vp1+8*_i) += value;   \
79d4002b98SHong Zhang         else                     *(vp1+8*_i) = value; \
80d4002b98SHong Zhang         goto a_noinsert; \
81d4002b98SHong Zhang       } \
82d4002b98SHong Zhang     }  \
83d4002b98SHong Zhang     if (value == 0.0 && ignorezeroentries) {low1 = 0; high1 = nrow1;goto a_noinsert;} \
84d4002b98SHong Zhang     if (nonew == 1) {low1 = 0; high1 = nrow1; goto a_noinsert;} \
85*08401ef6SPierre Jolivet     PetscCheck(nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \
86d4002b98SHong Zhang     MatSeqXSELLReallocateSELL(A,am,1,nrow1,a->sliidx,row/8,row,col,a->colidx,a->val,cp1,vp1,nonew,MatScalar); \
87d4002b98SHong Zhang     /* shift up all the later entries in this row */ \
88d4002b98SHong Zhang     for (ii=nrow1-1; ii>=_i; ii--) { \
89d4002b98SHong Zhang       *(cp1+8*(ii+1)) = *(cp1+8*ii); \
90d4002b98SHong Zhang       *(vp1+8*(ii+1)) = *(vp1+8*ii); \
91d4002b98SHong Zhang     } \
92d4002b98SHong Zhang     *(cp1+8*_i) = col; \
93d4002b98SHong Zhang     *(vp1+8*_i) = value; \
94d4002b98SHong Zhang     a->nz++; nrow1++; A->nonzerostate++; \
95d4002b98SHong Zhang     a_noinsert: ; \
96d4002b98SHong Zhang     a->rlen[row] = nrow1; \
97d4002b98SHong Zhang   }
98d4002b98SHong Zhang 
99d4002b98SHong Zhang #define MatSetValues_SeqSELL_B_Private(row,col,value,addv,orow,ocol) \
100d4002b98SHong Zhang   { \
101d4002b98SHong Zhang     if (col <= lastcol2) low2 = 0; \
102d4002b98SHong Zhang     else                high2 = nrow2; \
103d4002b98SHong Zhang     lastcol2 = col; \
104d4002b98SHong Zhang     while (high2-low2 > 5) { \
105d4002b98SHong Zhang       t = (low2+high2)/2; \
106d4002b98SHong Zhang       if (*(cp2+8*t) > col) high2 = t; \
107d4002b98SHong Zhang       else low2  = t; \
108d4002b98SHong Zhang     } \
109d4002b98SHong Zhang     for (_i=low2; _i<high2; _i++) { \
110d4002b98SHong Zhang       if (*(cp2+8*_i) > col) break; \
111d4002b98SHong Zhang       if (*(cp2+8*_i) == col) { \
112d4002b98SHong Zhang         if (addv == ADD_VALUES) *(vp2+8*_i) += value; \
113d4002b98SHong Zhang         else                     *(vp2+8*_i) = value; \
114d4002b98SHong Zhang         goto b_noinsert; \
115d4002b98SHong Zhang       } \
116d4002b98SHong Zhang     } \
117d4002b98SHong Zhang     if (value == 0.0 && ignorezeroentries) {low2 = 0; high2 = nrow2; goto b_noinsert;} \
118d4002b98SHong Zhang     if (nonew == 1) {low2 = 0; high2 = nrow2; goto b_noinsert;} \
119*08401ef6SPierre Jolivet     PetscCheck(nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \
120d4002b98SHong Zhang     MatSeqXSELLReallocateSELL(B,bm,1,nrow2,b->sliidx,row/8,row,col,b->colidx,b->val,cp2,vp2,nonew,MatScalar); \
121d4002b98SHong Zhang     /* shift up all the later entries in this row */ \
122d4002b98SHong Zhang     for (ii=nrow2-1; ii>=_i; ii--) { \
123d4002b98SHong Zhang       *(cp2+8*(ii+1)) = *(cp2+8*ii); \
124d4002b98SHong Zhang       *(vp2+8*(ii+1)) = *(vp2+8*ii); \
125d4002b98SHong Zhang     } \
126d4002b98SHong Zhang     *(cp2+8*_i) = col; \
127d4002b98SHong Zhang     *(vp2+8*_i) = value; \
128d4002b98SHong Zhang     b->nz++; nrow2++; B->nonzerostate++; \
129d4002b98SHong Zhang     b_noinsert: ; \
130d4002b98SHong Zhang     b->rlen[row] = nrow2; \
131d4002b98SHong Zhang   }
132d4002b98SHong Zhang 
133d4002b98SHong Zhang PetscErrorCode MatSetValues_MPISELL(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
134d4002b98SHong Zhang {
135d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
136d4002b98SHong Zhang   PetscScalar    value;
137d4002b98SHong Zhang   PetscInt       i,j,rstart=mat->rmap->rstart,rend=mat->rmap->rend,shift1,shift2;
138d4002b98SHong Zhang   PetscInt       cstart=mat->cmap->rstart,cend=mat->cmap->rend,row,col;
139d4002b98SHong Zhang   PetscBool      roworiented=sell->roworiented;
140d4002b98SHong Zhang 
141d4002b98SHong Zhang   /* Some Variables required in the macro */
142d4002b98SHong Zhang   Mat            A=sell->A;
143d4002b98SHong Zhang   Mat_SeqSELL    *a=(Mat_SeqSELL*)A->data;
144d4002b98SHong Zhang   PetscBool      ignorezeroentries=a->ignorezeroentries,found;
145d4002b98SHong Zhang   Mat            B=sell->B;
146d4002b98SHong Zhang   Mat_SeqSELL    *b=(Mat_SeqSELL*)B->data;
147d4002b98SHong Zhang   PetscInt       *cp1,*cp2,ii,_i,nrow1,nrow2,low1,high1,low2,high2,t,lastcol1,lastcol2;
148d4002b98SHong Zhang   MatScalar      *vp1,*vp2;
149d4002b98SHong Zhang 
150d4002b98SHong Zhang   PetscFunctionBegin;
151d4002b98SHong Zhang   for (i=0; i<m; i++) {
152d4002b98SHong Zhang     if (im[i] < 0) continue;
1536bdcaf15SBarry Smith     PetscCheck(im[i] < mat->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,im[i],mat->rmap->N-1);
154d4002b98SHong Zhang     if (im[i] >= rstart && im[i] < rend) {
155d4002b98SHong Zhang       row      = im[i] - rstart;
156d4002b98SHong Zhang       lastcol1 = -1;
157d4002b98SHong Zhang       shift1   = a->sliidx[row>>3]+(row&0x07); /* starting index of the row */
158d4002b98SHong Zhang       cp1      = a->colidx+shift1;
159d4002b98SHong Zhang       vp1      = a->val+shift1;
160d4002b98SHong Zhang       nrow1    = a->rlen[row];
161d4002b98SHong Zhang       low1     = 0;
162d4002b98SHong Zhang       high1    = nrow1;
163d4002b98SHong Zhang       lastcol2 = -1;
164d4002b98SHong Zhang       shift2   = b->sliidx[row>>3]+(row&0x07); /* starting index of the row */
165d4002b98SHong Zhang       cp2      = b->colidx+shift2;
166d4002b98SHong Zhang       vp2      = b->val+shift2;
167d4002b98SHong Zhang       nrow2    = b->rlen[row];
168d4002b98SHong Zhang       low2     = 0;
169d4002b98SHong Zhang       high2    = nrow2;
170d4002b98SHong Zhang 
171d4002b98SHong Zhang       for (j=0; j<n; j++) {
172d4002b98SHong Zhang         if (roworiented) value = v[i*n+j];
173d4002b98SHong Zhang         else             value = v[i+j*m];
174d4002b98SHong Zhang         if (ignorezeroentries && value == 0.0 && (addv == ADD_VALUES)) continue;
175d4002b98SHong Zhang         if (in[j] >= cstart && in[j] < cend) {
176d4002b98SHong Zhang           col   = in[j] - cstart;
177d4002b98SHong Zhang           MatSetValue_SeqSELL_Private(A,row,col,value,addv,im[i],in[j],cp1,vp1,lastcol1,low1,high1); /* set one value */
178d4002b98SHong Zhang         } else if (in[j] < 0) continue;
1796bdcaf15SBarry Smith         else PetscCheck(in[j] < mat->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,in[j],mat->cmap->N-1);
180d4002b98SHong Zhang         else {
181d4002b98SHong Zhang           if (mat->was_assembled) {
182d4002b98SHong Zhang             if (!sell->colmap) {
1839566063dSJacob Faibussowitsch               PetscCall(MatCreateColmap_MPISELL_Private(mat));
184d4002b98SHong Zhang             }
185d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
1869566063dSJacob Faibussowitsch             PetscCall(PetscTableFind(sell->colmap,in[j]+1,&col));
187d4002b98SHong Zhang             col--;
188d4002b98SHong Zhang #else
189d4002b98SHong Zhang             col = sell->colmap[in[j]] - 1;
190d4002b98SHong Zhang #endif
191d4002b98SHong Zhang             if (col < 0 && !((Mat_SeqSELL*)(sell->B->data))->nonew) {
1929566063dSJacob Faibussowitsch               PetscCall(MatDisAssemble_MPISELL(mat));
193d4002b98SHong Zhang               col    = in[j];
194d4002b98SHong Zhang               /* Reinitialize the variables required by MatSetValues_SeqSELL_B_Private() */
195d4002b98SHong Zhang               B      = sell->B;
196d4002b98SHong Zhang               b      = (Mat_SeqSELL*)B->data;
197d4002b98SHong Zhang               shift2 = b->sliidx[row>>3]+(row&0x07); /* starting index of the row */
198d4002b98SHong Zhang               cp2    = b->colidx+shift2;
199d4002b98SHong Zhang               vp2    = b->val+shift2;
200d4002b98SHong Zhang               nrow2  = b->rlen[row];
201d4002b98SHong Zhang               low2   = 0;
202d4002b98SHong Zhang               high2  = nrow2;
203*08401ef6SPierre Jolivet             } else PetscCheck(col >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", im[i], in[j]);
204d4002b98SHong Zhang           } else col = in[j];
205d4002b98SHong Zhang           MatSetValue_SeqSELL_Private(B,row,col,value,addv,im[i],in[j],cp2,vp2,lastcol2,low2,high2); /* set one value */
206d4002b98SHong Zhang         }
207d4002b98SHong Zhang       }
208d4002b98SHong Zhang     } else {
20928b400f6SJacob Faibussowitsch       PetscCheck(!mat->nooffprocentries,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %" PetscInt_FMT " even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
210d4002b98SHong Zhang       if (!sell->donotstash) {
211d4002b98SHong Zhang         mat->assembled = PETSC_FALSE;
212d4002b98SHong Zhang         if (roworiented) {
2139566063dSJacob Faibussowitsch           PetscCall(MatStashValuesRow_Private(&mat->stash,im[i],n,in,v+i*n,(PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
214d4002b98SHong Zhang         } else {
2159566063dSJacob Faibussowitsch           PetscCall(MatStashValuesCol_Private(&mat->stash,im[i],n,in,v+i,m,(PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
216d4002b98SHong Zhang         }
217d4002b98SHong Zhang       }
218d4002b98SHong Zhang     }
219d4002b98SHong Zhang   }
220d4002b98SHong Zhang   PetscFunctionReturn(0);
221d4002b98SHong Zhang }
222d4002b98SHong Zhang 
223d4002b98SHong Zhang PetscErrorCode MatGetValues_MPISELL(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
224d4002b98SHong Zhang {
225d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
226d4002b98SHong Zhang   PetscInt       i,j,rstart=mat->rmap->rstart,rend=mat->rmap->rend;
227d4002b98SHong Zhang   PetscInt       cstart=mat->cmap->rstart,cend=mat->cmap->rend,row,col;
228d4002b98SHong Zhang 
229d4002b98SHong Zhang   PetscFunctionBegin;
230d4002b98SHong Zhang   for (i=0; i<m; i++) {
23154c59aa7SJacob Faibussowitsch     if (idxm[i] < 0) continue; /* negative row */
23254c59aa7SJacob Faibussowitsch     PetscCheck(idxm[i] < mat->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,idxm[i],mat->rmap->N-1);
233d4002b98SHong Zhang     if (idxm[i] >= rstart && idxm[i] < rend) {
234d4002b98SHong Zhang       row = idxm[i] - rstart;
235d4002b98SHong Zhang       for (j=0; j<n; j++) {
23654c59aa7SJacob Faibussowitsch         if (idxn[j] < 0) continue; /* negative column */
23754c59aa7SJacob Faibussowitsch         PetscCheck(idxn[j] < mat->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,idxn[j],mat->cmap->N-1);
238d4002b98SHong Zhang         if (idxn[j] >= cstart && idxn[j] < cend) {
239d4002b98SHong Zhang           col  = idxn[j] - cstart;
2409566063dSJacob Faibussowitsch           PetscCall(MatGetValues(sell->A,1,&row,1,&col,v+i*n+j));
241d4002b98SHong Zhang         } else {
242d4002b98SHong Zhang           if (!sell->colmap) {
2439566063dSJacob Faibussowitsch             PetscCall(MatCreateColmap_MPISELL_Private(mat));
244d4002b98SHong Zhang           }
245d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
2469566063dSJacob Faibussowitsch           PetscCall(PetscTableFind(sell->colmap,idxn[j]+1,&col));
247d4002b98SHong Zhang           col--;
248d4002b98SHong Zhang #else
249d4002b98SHong Zhang           col = sell->colmap[idxn[j]] - 1;
250d4002b98SHong Zhang #endif
251d4002b98SHong Zhang           if ((col < 0) || (sell->garray[col] != idxn[j])) *(v+i*n+j) = 0.0;
252d4002b98SHong Zhang           else {
2539566063dSJacob Faibussowitsch             PetscCall(MatGetValues(sell->B,1,&row,1,&col,v+i*n+j));
254d4002b98SHong Zhang           }
255d4002b98SHong Zhang         }
256d4002b98SHong Zhang       }
257d4002b98SHong Zhang     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
258d4002b98SHong Zhang   }
259d4002b98SHong Zhang   PetscFunctionReturn(0);
260d4002b98SHong Zhang }
261d4002b98SHong Zhang 
262d4002b98SHong Zhang extern PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat,Vec,Vec);
263d4002b98SHong Zhang 
264d4002b98SHong Zhang PetscErrorCode MatAssemblyBegin_MPISELL(Mat mat,MatAssemblyType mode)
265d4002b98SHong Zhang {
266d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
267d4002b98SHong Zhang   PetscInt       nstash,reallocs;
268d4002b98SHong Zhang 
269d4002b98SHong Zhang   PetscFunctionBegin;
270d4002b98SHong Zhang   if (sell->donotstash || mat->nooffprocentries) PetscFunctionReturn(0);
271d4002b98SHong Zhang 
2729566063dSJacob Faibussowitsch   PetscCall(MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range));
2739566063dSJacob Faibussowitsch   PetscCall(MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs));
2749566063dSJacob Faibussowitsch   PetscCall(PetscInfo(sell->A,"Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n",nstash,reallocs));
275d4002b98SHong Zhang   PetscFunctionReturn(0);
276d4002b98SHong Zhang }
277d4002b98SHong Zhang 
278d4002b98SHong Zhang PetscErrorCode MatAssemblyEnd_MPISELL(Mat mat,MatAssemblyType mode)
279d4002b98SHong Zhang {
280d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
281d4002b98SHong Zhang   PetscMPIInt    n;
282d4002b98SHong Zhang   PetscInt       i,flg;
283d4002b98SHong Zhang   PetscInt       *row,*col;
284d4002b98SHong Zhang   PetscScalar    *val;
285d4002b98SHong Zhang   PetscBool      other_disassembled;
286d4002b98SHong Zhang   /* do not use 'b = (Mat_SeqSELL*)sell->B->data' as B can be reset in disassembly */
287d4002b98SHong Zhang   PetscFunctionBegin;
288d4002b98SHong Zhang   if (!sell->donotstash && !mat->nooffprocentries) {
289d4002b98SHong Zhang     while (1) {
2909566063dSJacob Faibussowitsch       PetscCall(MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg));
291d4002b98SHong Zhang       if (!flg) break;
292d4002b98SHong Zhang 
293d4002b98SHong Zhang       for (i=0; i<n; i++) { /* assemble one by one */
2949566063dSJacob Faibussowitsch         PetscCall(MatSetValues_MPISELL(mat,1,row+i,1,col+i,val+i,mat->insertmode));
295d4002b98SHong Zhang       }
296d4002b98SHong Zhang     }
2979566063dSJacob Faibussowitsch     PetscCall(MatStashScatterEnd_Private(&mat->stash));
298d4002b98SHong Zhang   }
2999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sell->A,mode));
3009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sell->A,mode));
301d4002b98SHong Zhang 
302d4002b98SHong Zhang   /*
303d4002b98SHong Zhang      determine if any processor has disassembled, if so we must
304d4002b98SHong Zhang      also disassemble ourselfs, in order that we may reassemble.
305d4002b98SHong Zhang   */
306d4002b98SHong Zhang   /*
307d4002b98SHong Zhang      if nonzero structure of submatrix B cannot change then we know that
308d4002b98SHong Zhang      no processor disassembled thus we can skip this stuff
309d4002b98SHong Zhang   */
310d4002b98SHong Zhang   if (!((Mat_SeqSELL*)sell->B->data)->nonew) {
3111c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat)));
312*08401ef6SPierre Jolivet     PetscCheck(!mat->was_assembled || other_disassembled,PETSC_COMM_SELF,PETSC_ERR_SUP,"MatDisAssemble not implemented yet");
313d4002b98SHong Zhang   }
314d4002b98SHong Zhang   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
3159566063dSJacob Faibussowitsch     PetscCall(MatSetUpMultiply_MPISELL(mat));
316d4002b98SHong Zhang   }
317d4002b98SHong Zhang   /*
3189566063dSJacob Faibussowitsch   PetscCall(MatSetOption(sell->B,MAT_USE_INODES,PETSC_FALSE));
319d4002b98SHong Zhang   */
3209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sell->B,mode));
3219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sell->B,mode));
3229566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sell->rowvalues,sell->rowindices));
323f4259b30SLisandro Dalcin   sell->rowvalues = NULL;
3249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->diag));
325d4002b98SHong Zhang 
326d4002b98SHong Zhang   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
327d4002b98SHong Zhang   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqSELL*)(sell->A->data))->nonew) {
328d4002b98SHong Zhang     PetscObjectState state = sell->A->nonzerostate + sell->B->nonzerostate;
3291c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat)));
330d4002b98SHong Zhang   }
331d4002b98SHong Zhang   PetscFunctionReturn(0);
332d4002b98SHong Zhang }
333d4002b98SHong Zhang 
334d4002b98SHong Zhang PetscErrorCode MatZeroEntries_MPISELL(Mat A)
335d4002b98SHong Zhang {
336d4002b98SHong Zhang   Mat_MPISELL    *l=(Mat_MPISELL*)A->data;
337d4002b98SHong Zhang 
338d4002b98SHong Zhang   PetscFunctionBegin;
3399566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->A));
3409566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->B));
341d4002b98SHong Zhang   PetscFunctionReturn(0);
342d4002b98SHong Zhang }
343d4002b98SHong Zhang 
344d4002b98SHong Zhang PetscErrorCode MatMult_MPISELL(Mat A,Vec xx,Vec yy)
345d4002b98SHong Zhang {
346d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
347d4002b98SHong Zhang   PetscInt       nt;
348d4002b98SHong Zhang 
349d4002b98SHong Zhang   PetscFunctionBegin;
3509566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xx,&nt));
351*08401ef6SPierre Jolivet   PetscCheck(nt == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A (%" PetscInt_FMT ") and xx (%" PetscInt_FMT ")",A->cmap->n,nt);
3529566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD));
3539566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->mult)(a->A,xx,yy));
3549566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD));
3559566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B,a->lvec,yy,yy));
356d4002b98SHong Zhang   PetscFunctionReturn(0);
357d4002b98SHong Zhang }
358d4002b98SHong Zhang 
359d4002b98SHong Zhang PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat A,Vec bb,Vec xx)
360d4002b98SHong Zhang {
361d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
362d4002b98SHong Zhang 
363d4002b98SHong Zhang   PetscFunctionBegin;
3649566063dSJacob Faibussowitsch   PetscCall(MatMultDiagonalBlock(a->A,bb,xx));
365d4002b98SHong Zhang   PetscFunctionReturn(0);
366d4002b98SHong Zhang }
367d4002b98SHong Zhang 
368d4002b98SHong Zhang PetscErrorCode MatMultAdd_MPISELL(Mat A,Vec xx,Vec yy,Vec zz)
369d4002b98SHong Zhang {
370d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
371d4002b98SHong Zhang 
372d4002b98SHong Zhang   PetscFunctionBegin;
3739566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD));
3749566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multadd)(a->A,xx,yy,zz));
3759566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD));
3769566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B,a->lvec,zz,zz));
377d4002b98SHong Zhang   PetscFunctionReturn(0);
378d4002b98SHong Zhang }
379d4002b98SHong Zhang 
380d4002b98SHong Zhang PetscErrorCode MatMultTranspose_MPISELL(Mat A,Vec xx,Vec yy)
381d4002b98SHong Zhang {
382d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
383d4002b98SHong Zhang 
384d4002b98SHong Zhang   PetscFunctionBegin;
385d4002b98SHong Zhang   /* do nondiagonal part */
3869566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B,xx,a->lvec));
387d4002b98SHong Zhang   /* do local part */
3889566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multtranspose)(a->A,xx,yy));
389a29b4eb7SJunchao Zhang   /* add partial results together */
3909566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE));
3919566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE));
392d4002b98SHong Zhang   PetscFunctionReturn(0);
393d4002b98SHong Zhang }
394d4002b98SHong Zhang 
395d4002b98SHong Zhang PetscErrorCode MatIsTranspose_MPISELL(Mat Amat,Mat Bmat,PetscReal tol,PetscBool *f)
396d4002b98SHong Zhang {
397d4002b98SHong Zhang   MPI_Comm       comm;
398d4002b98SHong Zhang   Mat_MPISELL    *Asell=(Mat_MPISELL*)Amat->data,*Bsell;
399d4002b98SHong Zhang   Mat            Adia=Asell->A,Bdia,Aoff,Boff,*Aoffs,*Boffs;
400d4002b98SHong Zhang   IS             Me,Notme;
401d4002b98SHong Zhang   PetscInt       M,N,first,last,*notme,i;
402d4002b98SHong Zhang   PetscMPIInt    size;
403d4002b98SHong Zhang 
404d4002b98SHong Zhang   PetscFunctionBegin;
405d4002b98SHong Zhang   /* Easy test: symmetric diagonal block */
406d4002b98SHong Zhang   Bsell = (Mat_MPISELL*)Bmat->data; Bdia = Bsell->A;
4079566063dSJacob Faibussowitsch   PetscCall(MatIsTranspose(Adia,Bdia,tol,f));
408d4002b98SHong Zhang   if (!*f) PetscFunctionReturn(0);
4099566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
4109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
411d4002b98SHong Zhang   if (size == 1) PetscFunctionReturn(0);
412d4002b98SHong Zhang 
413d4002b98SHong Zhang   /* Hard test: off-diagonal block. This takes a MatCreateSubMatrix. */
4149566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Amat,&M,&N));
4159566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat,&first,&last));
4169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N-last+first,&notme));
417d4002b98SHong Zhang   for (i=0; i<first; i++) notme[i] = i;
418d4002b98SHong Zhang   for (i=last; i<M; i++) notme[i-last+first] = i;
4199566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(MPI_COMM_SELF,N-last+first,notme,PETSC_COPY_VALUES,&Notme));
4209566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(MPI_COMM_SELF,last-first,first,1,&Me));
4219566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Amat,1,&Me,&Notme,MAT_INITIAL_MATRIX,&Aoffs));
422d4002b98SHong Zhang   Aoff = Aoffs[0];
4239566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Bmat,1,&Notme,&Me,MAT_INITIAL_MATRIX,&Boffs));
424d4002b98SHong Zhang   Boff = Boffs[0];
4259566063dSJacob Faibussowitsch   PetscCall(MatIsTranspose(Aoff,Boff,tol,f));
4269566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1,&Aoffs));
4279566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1,&Boffs));
4289566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Me));
4299566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Notme));
4309566063dSJacob Faibussowitsch   PetscCall(PetscFree(notme));
431d4002b98SHong Zhang   PetscFunctionReturn(0);
432d4002b98SHong Zhang }
433d4002b98SHong Zhang 
434d4002b98SHong Zhang PetscErrorCode MatMultTransposeAdd_MPISELL(Mat A,Vec xx,Vec yy,Vec zz)
435d4002b98SHong Zhang {
436d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
437d4002b98SHong Zhang 
438d4002b98SHong Zhang   PetscFunctionBegin;
439d4002b98SHong Zhang   /* do nondiagonal part */
4409566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B,xx,a->lvec));
441d4002b98SHong Zhang   /* do local part */
4429566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multtransposeadd)(a->A,xx,yy,zz));
443e4a140f6SJunchao Zhang   /* add partial results together */
4449566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE));
4459566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE));
446d4002b98SHong Zhang   PetscFunctionReturn(0);
447d4002b98SHong Zhang }
448d4002b98SHong Zhang 
449d4002b98SHong Zhang /*
450d4002b98SHong Zhang   This only works correctly for square matrices where the subblock A->A is the
451d4002b98SHong Zhang    diagonal block
452d4002b98SHong Zhang */
453d4002b98SHong Zhang PetscErrorCode MatGetDiagonal_MPISELL(Mat A,Vec v)
454d4002b98SHong Zhang {
455d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
456d4002b98SHong Zhang 
457d4002b98SHong Zhang   PetscFunctionBegin;
458*08401ef6SPierre Jolivet   PetscCheck(A->rmap->N == A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block");
4592c71b3e2SJacob Faibussowitsch   PetscCheckFalse(A->rmap->rstart != A->cmap->rstart || A->rmap->rend != A->cmap->rend,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"row partition must equal col partition");
4609566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(a->A,v));
461d4002b98SHong Zhang   PetscFunctionReturn(0);
462d4002b98SHong Zhang }
463d4002b98SHong Zhang 
464d4002b98SHong Zhang PetscErrorCode MatScale_MPISELL(Mat A,PetscScalar aa)
465d4002b98SHong Zhang {
466d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
467d4002b98SHong Zhang 
468d4002b98SHong Zhang   PetscFunctionBegin;
4699566063dSJacob Faibussowitsch   PetscCall(MatScale(a->A,aa));
4709566063dSJacob Faibussowitsch   PetscCall(MatScale(a->B,aa));
471d4002b98SHong Zhang   PetscFunctionReturn(0);
472d4002b98SHong Zhang }
473d4002b98SHong Zhang 
474d4002b98SHong Zhang PetscErrorCode MatDestroy_MPISELL(Mat mat)
475d4002b98SHong Zhang {
476d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
477d4002b98SHong Zhang 
478d4002b98SHong Zhang   PetscFunctionBegin;
479d4002b98SHong Zhang #if defined(PETSC_USE_LOG)
480c0aa6a63SJacob Faibussowitsch   PetscLogObjectState((PetscObject)mat,"Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT,mat->rmap->N,mat->cmap->N);
481d4002b98SHong Zhang #endif
4829566063dSJacob Faibussowitsch   PetscCall(MatStashDestroy_Private(&mat->stash));
4839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->diag));
4849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sell->A));
4859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sell->B));
486d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
4879566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&sell->colmap));
488d4002b98SHong Zhang #else
4899566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->colmap));
490d4002b98SHong Zhang #endif
4919566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->garray));
4929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->lvec));
4939566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sell->Mvctx));
4949566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sell->rowvalues,sell->rowindices));
4959566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->ld));
4969566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
497d4002b98SHong Zhang 
4989566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)mat,NULL));
4999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL));
5009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL));
5019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatIsTranspose_C",NULL));
5029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPISELLSetPreallocation_C",NULL));
5039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisell_mpiaij_C",NULL));
5049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDiagonalScaleLocal_C",NULL));
505d4002b98SHong Zhang   PetscFunctionReturn(0);
506d4002b98SHong Zhang }
507d4002b98SHong Zhang 
508d4002b98SHong Zhang #include <petscdraw.h>
509d4002b98SHong Zhang PetscErrorCode MatView_MPISELL_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
510d4002b98SHong Zhang {
511d4002b98SHong Zhang   Mat_MPISELL       *sell=(Mat_MPISELL*)mat->data;
512d4002b98SHong Zhang   PetscErrorCode    ierr;
513d4002b98SHong Zhang   PetscMPIInt       rank=sell->rank,size=sell->size;
514d4002b98SHong Zhang   PetscBool         isdraw,iascii,isbinary;
515d4002b98SHong Zhang   PetscViewer       sviewer;
516d4002b98SHong Zhang   PetscViewerFormat format;
517d4002b98SHong Zhang 
518d4002b98SHong Zhang   PetscFunctionBegin;
5199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
522d4002b98SHong Zhang   if (iascii) {
5239566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
524d4002b98SHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
525d4002b98SHong Zhang       MatInfo   info;
5266335e310SSatish Balay       PetscInt *inodes;
527d4002b98SHong Zhang 
5289566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank));
5299566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(mat,MAT_LOCAL,&info));
5309566063dSJacob Faibussowitsch       PetscCall(MatInodeGetInodeSizes(sell->A,NULL,&inodes,NULL));
5319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
532d4002b98SHong Zhang       if (!inodes) {
533c0aa6a63SJacob Faibussowitsch         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %" PetscInt_FMT " nz %" PetscInt_FMT " nz alloced %" PetscInt_FMT " mem %" PetscInt_FMT ", not using I-node routines\n",
5349566063dSJacob Faibussowitsch                                                   rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);PetscCall(ierr);
535d4002b98SHong Zhang       } else {
536c0aa6a63SJacob Faibussowitsch         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %" PetscInt_FMT " nz %" PetscInt_FMT " nz alloced %" PetscInt_FMT " mem %" PetscInt_FMT ", using I-node routines\n",
5379566063dSJacob Faibussowitsch                                                   rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);PetscCall(ierr);
538d4002b98SHong Zhang       }
5399566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(sell->A,MAT_LOCAL,&info));
5409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %" PetscInt_FMT " \n",rank,(PetscInt)info.nz_used));
5419566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(sell->B,MAT_LOCAL,&info));
5429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %" PetscInt_FMT " \n",rank,(PetscInt)info.nz_used));
5439566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
5449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
5459566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n"));
5469566063dSJacob Faibussowitsch       PetscCall(VecScatterView(sell->Mvctx,viewer));
547d4002b98SHong Zhang       PetscFunctionReturn(0);
548d4002b98SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_INFO) {
549d4002b98SHong Zhang       PetscInt inodecount,inodelimit,*inodes;
5509566063dSJacob Faibussowitsch       PetscCall(MatInodeGetInodeSizes(sell->A,&inodecount,&inodes,&inodelimit));
551d4002b98SHong Zhang       if (inodes) {
5529566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"using I-node (on process 0) routines: found %" PetscInt_FMT " nodes, limit used is %" PetscInt_FMT "\n",inodecount,inodelimit));
553d4002b98SHong Zhang       } else {
5549566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"not using I-node (on process 0) routines\n"));
555d4002b98SHong Zhang       }
556d4002b98SHong Zhang       PetscFunctionReturn(0);
557d4002b98SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
558d4002b98SHong Zhang       PetscFunctionReturn(0);
559d4002b98SHong Zhang     }
560d4002b98SHong Zhang   } else if (isbinary) {
561d4002b98SHong Zhang     if (size == 1) {
5629566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)sell->A,((PetscObject)mat)->name));
5639566063dSJacob Faibussowitsch       PetscCall(MatView(sell->A,viewer));
564d4002b98SHong Zhang     } else {
5659566063dSJacob Faibussowitsch       /* PetscCall(MatView_MPISELL_Binary(mat,viewer)); */
566d4002b98SHong Zhang     }
567d4002b98SHong Zhang     PetscFunctionReturn(0);
568d4002b98SHong Zhang   } else if (isdraw) {
569d4002b98SHong Zhang     PetscDraw draw;
570d4002b98SHong Zhang     PetscBool isnull;
5719566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
5729566063dSJacob Faibussowitsch     PetscCall(PetscDrawIsNull(draw,&isnull));
573d4002b98SHong Zhang     if (isnull) PetscFunctionReturn(0);
574d4002b98SHong Zhang   }
575d4002b98SHong Zhang 
576d4002b98SHong Zhang   {
577d4002b98SHong Zhang     /* assemble the entire matrix onto first processor. */
578d4002b98SHong Zhang     Mat         A;
579d4002b98SHong Zhang     Mat_SeqSELL *Aloc;
580d4002b98SHong Zhang     PetscInt    M=mat->rmap->N,N=mat->cmap->N,*acolidx,row,col,i,j;
581d4002b98SHong Zhang     MatScalar   *aval;
582d4002b98SHong Zhang     PetscBool   isnonzero;
583d4002b98SHong Zhang 
5849566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat),&A));
585dd400576SPatrick Sanan     if (rank == 0) {
5869566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A,M,N,M,N));
587d4002b98SHong Zhang     } else {
5889566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A,0,0,M,N));
589d4002b98SHong Zhang     }
590d4002b98SHong Zhang     /* This is just a temporary matrix, so explicitly using MATMPISELL is probably best */
5919566063dSJacob Faibussowitsch     PetscCall(MatSetType(A,MATMPISELL));
5929566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(A,0,NULL,0,NULL));
5939566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE));
5949566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)A));
595d4002b98SHong Zhang 
596d4002b98SHong Zhang     /* copy over the A part */
597d4002b98SHong Zhang     Aloc = (Mat_SeqSELL*)sell->A->data;
598d4002b98SHong Zhang     acolidx = Aloc->colidx; aval = Aloc->val;
599d4002b98SHong Zhang     for (i=0; i<Aloc->totalslices; i++) { /* loop over slices */
600d4002b98SHong Zhang       for (j=Aloc->sliidx[i]; j<Aloc->sliidx[i+1]; j++) {
601d4002b98SHong Zhang         isnonzero = (PetscBool)((j-Aloc->sliidx[i])/8 < Aloc->rlen[(i<<3)+(j&0x07)]);
602d4002b98SHong Zhang         if (isnonzero) { /* check the mask bit */
603d4002b98SHong Zhang           row  = (i<<3)+(j&0x07) + mat->rmap->rstart; /* i<<3 is the starting row of this slice */
604d4002b98SHong Zhang           col  = *acolidx + mat->rmap->rstart;
6059566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A,1,&row,1,&col,aval,INSERT_VALUES));
606d4002b98SHong Zhang         }
607d4002b98SHong Zhang         aval++; acolidx++;
608d4002b98SHong Zhang       }
609d4002b98SHong Zhang     }
610d4002b98SHong Zhang 
611d4002b98SHong Zhang     /* copy over the B part */
612d4002b98SHong Zhang     Aloc = (Mat_SeqSELL*)sell->B->data;
613d4002b98SHong Zhang     acolidx = Aloc->colidx; aval = Aloc->val;
614d4002b98SHong Zhang     for (i=0; i<Aloc->totalslices; i++) {
615d4002b98SHong Zhang       for (j=Aloc->sliidx[i]; j<Aloc->sliidx[i+1]; j++) {
616d4002b98SHong Zhang         isnonzero = (PetscBool)((j-Aloc->sliidx[i])/8 < Aloc->rlen[(i<<3)+(j&0x07)]);
617d4002b98SHong Zhang         if (isnonzero) {
618d4002b98SHong Zhang           row  = (i<<3)+(j&0x07) + mat->rmap->rstart;
619d4002b98SHong Zhang           col  = sell->garray[*acolidx];
6209566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A,1,&row,1,&col,aval,INSERT_VALUES));
621d4002b98SHong Zhang         }
622d4002b98SHong Zhang         aval++; acolidx++;
623d4002b98SHong Zhang       }
624d4002b98SHong Zhang     }
625d4002b98SHong Zhang 
6269566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
6279566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
628d4002b98SHong Zhang     /*
629d4002b98SHong Zhang        Everyone has to call to draw the matrix since the graphics waits are
630d4002b98SHong Zhang        synchronized across all processors that share the PetscDraw object
631d4002b98SHong Zhang     */
6329566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
633dd400576SPatrick Sanan     if (rank == 0) {
6349566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)((Mat_MPISELL*)(A->data))->A,((PetscObject)mat)->name));
6359566063dSJacob Faibussowitsch       PetscCall(MatView_SeqSELL(((Mat_MPISELL*)(A->data))->A,sviewer));
636d4002b98SHong Zhang     }
6379566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
6389566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
6399566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
640d4002b98SHong Zhang   }
641d4002b98SHong Zhang   PetscFunctionReturn(0);
642d4002b98SHong Zhang }
643d4002b98SHong Zhang 
644d4002b98SHong Zhang PetscErrorCode MatView_MPISELL(Mat mat,PetscViewer viewer)
645d4002b98SHong Zhang {
646d4002b98SHong Zhang   PetscBool      iascii,isdraw,issocket,isbinary;
647d4002b98SHong Zhang 
648d4002b98SHong Zhang   PetscFunctionBegin;
6499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
6509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
6519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
6529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket));
653d4002b98SHong Zhang   if (iascii || isdraw || isbinary || issocket) {
6549566063dSJacob Faibussowitsch     PetscCall(MatView_MPISELL_ASCIIorDraworSocket(mat,viewer));
655d4002b98SHong Zhang   }
656d4002b98SHong Zhang   PetscFunctionReturn(0);
657d4002b98SHong Zhang }
658d4002b98SHong Zhang 
659d4002b98SHong Zhang PetscErrorCode MatGetGhosts_MPISELL(Mat mat,PetscInt *nghosts,const PetscInt *ghosts[])
660d4002b98SHong Zhang {
661d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
662d4002b98SHong Zhang 
663d4002b98SHong Zhang   PetscFunctionBegin;
6649566063dSJacob Faibussowitsch   PetscCall(MatGetSize(sell->B,NULL,nghosts));
665d4002b98SHong Zhang   if (ghosts) *ghosts = sell->garray;
666d4002b98SHong Zhang   PetscFunctionReturn(0);
667d4002b98SHong Zhang }
668d4002b98SHong Zhang 
669d4002b98SHong Zhang PetscErrorCode MatGetInfo_MPISELL(Mat matin,MatInfoType flag,MatInfo *info)
670d4002b98SHong Zhang {
671d4002b98SHong Zhang   Mat_MPISELL    *mat=(Mat_MPISELL*)matin->data;
672d4002b98SHong Zhang   Mat            A=mat->A,B=mat->B;
6733966268fSBarry Smith   PetscLogDouble isend[5],irecv[5];
674d4002b98SHong Zhang 
675d4002b98SHong Zhang   PetscFunctionBegin;
676d4002b98SHong Zhang   info->block_size = 1.0;
6779566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A,MAT_LOCAL,info));
678d4002b98SHong Zhang 
679d4002b98SHong Zhang   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
680d4002b98SHong Zhang   isend[3] = info->memory;  isend[4] = info->mallocs;
681d4002b98SHong Zhang 
6829566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(B,MAT_LOCAL,info));
683d4002b98SHong Zhang 
684d4002b98SHong Zhang   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
685d4002b98SHong Zhang   isend[3] += info->memory;  isend[4] += info->mallocs;
686d4002b98SHong Zhang   if (flag == MAT_LOCAL) {
687d4002b98SHong Zhang     info->nz_used      = isend[0];
688d4002b98SHong Zhang     info->nz_allocated = isend[1];
689d4002b98SHong Zhang     info->nz_unneeded  = isend[2];
690d4002b98SHong Zhang     info->memory       = isend[3];
691d4002b98SHong Zhang     info->mallocs      = isend[4];
692d4002b98SHong Zhang   } else if (flag == MAT_GLOBAL_MAX) {
6931c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_MAX,PetscObjectComm((PetscObject)matin)));
694d4002b98SHong Zhang 
695d4002b98SHong Zhang     info->nz_used      = irecv[0];
696d4002b98SHong Zhang     info->nz_allocated = irecv[1];
697d4002b98SHong Zhang     info->nz_unneeded  = irecv[2];
698d4002b98SHong Zhang     info->memory       = irecv[3];
699d4002b98SHong Zhang     info->mallocs      = irecv[4];
700d4002b98SHong Zhang   } else if (flag == MAT_GLOBAL_SUM) {
7011c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_SUM,PetscObjectComm((PetscObject)matin)));
702d4002b98SHong Zhang 
703d4002b98SHong Zhang     info->nz_used      = irecv[0];
704d4002b98SHong Zhang     info->nz_allocated = irecv[1];
705d4002b98SHong Zhang     info->nz_unneeded  = irecv[2];
706d4002b98SHong Zhang     info->memory       = irecv[3];
707d4002b98SHong Zhang     info->mallocs      = irecv[4];
708d4002b98SHong Zhang   }
709d4002b98SHong Zhang   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
710d4002b98SHong Zhang   info->fill_ratio_needed = 0;
711d4002b98SHong Zhang   info->factor_mallocs    = 0;
712d4002b98SHong Zhang   PetscFunctionReturn(0);
713d4002b98SHong Zhang }
714d4002b98SHong Zhang 
715d4002b98SHong Zhang PetscErrorCode MatSetOption_MPISELL(Mat A,MatOption op,PetscBool flg)
716d4002b98SHong Zhang {
717d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
718d4002b98SHong Zhang 
719d4002b98SHong Zhang   PetscFunctionBegin;
720d4002b98SHong Zhang   switch (op) {
721d4002b98SHong Zhang   case MAT_NEW_NONZERO_LOCATIONS:
722d4002b98SHong Zhang   case MAT_NEW_NONZERO_ALLOCATION_ERR:
723d4002b98SHong Zhang   case MAT_UNUSED_NONZERO_LOCATION_ERR:
724d4002b98SHong Zhang   case MAT_KEEP_NONZERO_PATTERN:
725d4002b98SHong Zhang   case MAT_NEW_NONZERO_LOCATION_ERR:
726d4002b98SHong Zhang   case MAT_USE_INODES:
727d4002b98SHong Zhang   case MAT_IGNORE_ZERO_ENTRIES:
728d4002b98SHong Zhang     MatCheckPreallocated(A,1);
7299566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
7309566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B,op,flg));
731d4002b98SHong Zhang     break;
732d4002b98SHong Zhang   case MAT_ROW_ORIENTED:
733d4002b98SHong Zhang     MatCheckPreallocated(A,1);
734d4002b98SHong Zhang     a->roworiented = flg;
735d4002b98SHong Zhang 
7369566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
7379566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B,op,flg));
738d4002b98SHong Zhang     break;
7398c78258cSHong Zhang   case MAT_FORCE_DIAGONAL_ENTRIES:
740071fcb05SBarry Smith   case MAT_SORTED_FULL:
7419566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Option %s ignored\n",MatOptions[op]));
742d4002b98SHong Zhang     break;
743d4002b98SHong Zhang   case MAT_IGNORE_OFF_PROC_ENTRIES:
744d4002b98SHong Zhang     a->donotstash = flg;
745d4002b98SHong Zhang     break;
746d4002b98SHong Zhang   case MAT_SPD:
747d4002b98SHong Zhang     A->spd_set = PETSC_TRUE;
748d4002b98SHong Zhang     A->spd     = flg;
749d4002b98SHong Zhang     if (flg) {
750d4002b98SHong Zhang       A->symmetric                  = PETSC_TRUE;
751d4002b98SHong Zhang       A->structurally_symmetric     = PETSC_TRUE;
752d4002b98SHong Zhang       A->symmetric_set              = PETSC_TRUE;
753d4002b98SHong Zhang       A->structurally_symmetric_set = PETSC_TRUE;
754d4002b98SHong Zhang     }
755d4002b98SHong Zhang     break;
756d4002b98SHong Zhang   case MAT_SYMMETRIC:
757d4002b98SHong Zhang     MatCheckPreallocated(A,1);
7589566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
759d4002b98SHong Zhang     break;
760d4002b98SHong Zhang   case MAT_STRUCTURALLY_SYMMETRIC:
761d4002b98SHong Zhang     MatCheckPreallocated(A,1);
7629566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
763d4002b98SHong Zhang     break;
764d4002b98SHong Zhang   case MAT_HERMITIAN:
765d4002b98SHong Zhang     MatCheckPreallocated(A,1);
7669566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
767d4002b98SHong Zhang     break;
768d4002b98SHong Zhang   case MAT_SYMMETRY_ETERNAL:
769d4002b98SHong Zhang     MatCheckPreallocated(A,1);
7709566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
771d4002b98SHong Zhang     break;
772d4002b98SHong Zhang   default:
77398921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
774d4002b98SHong Zhang   }
775d4002b98SHong Zhang   PetscFunctionReturn(0);
776d4002b98SHong Zhang }
777d4002b98SHong Zhang 
778d4002b98SHong Zhang PetscErrorCode MatDiagonalScale_MPISELL(Mat mat,Vec ll,Vec rr)
779d4002b98SHong Zhang {
780d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
781d4002b98SHong Zhang   Mat            a=sell->A,b=sell->B;
782d4002b98SHong Zhang   PetscInt       s1,s2,s3;
783d4002b98SHong Zhang 
784d4002b98SHong Zhang   PetscFunctionBegin;
7859566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mat,&s2,&s3));
786d4002b98SHong Zhang   if (rr) {
7879566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rr,&s1));
788*08401ef6SPierre Jolivet     PetscCheck(s1==s3,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"right vector non-conforming local size");
789d4002b98SHong Zhang     /* Overlap communication with computation. */
7909566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(sell->Mvctx,rr,sell->lvec,INSERT_VALUES,SCATTER_FORWARD));
791d4002b98SHong Zhang   }
792d4002b98SHong Zhang   if (ll) {
7939566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ll,&s1));
794*08401ef6SPierre Jolivet     PetscCheck(s1==s2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"left vector non-conforming local size");
7959566063dSJacob Faibussowitsch     PetscCall((*b->ops->diagonalscale)(b,ll,NULL));
796d4002b98SHong Zhang   }
797d4002b98SHong Zhang   /* scale  the diagonal block */
7989566063dSJacob Faibussowitsch   PetscCall((*a->ops->diagonalscale)(a,ll,rr));
799d4002b98SHong Zhang 
800d4002b98SHong Zhang   if (rr) {
801d4002b98SHong Zhang     /* Do a scatter end and then right scale the off-diagonal block */
8029566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(sell->Mvctx,rr,sell->lvec,INSERT_VALUES,SCATTER_FORWARD));
8039566063dSJacob Faibussowitsch     PetscCall((*b->ops->diagonalscale)(b,NULL,sell->lvec));
804d4002b98SHong Zhang   }
805d4002b98SHong Zhang   PetscFunctionReturn(0);
806d4002b98SHong Zhang }
807d4002b98SHong Zhang 
808d4002b98SHong Zhang PetscErrorCode MatSetUnfactored_MPISELL(Mat A)
809d4002b98SHong Zhang {
810d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
811d4002b98SHong Zhang 
812d4002b98SHong Zhang   PetscFunctionBegin;
8139566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(a->A));
814d4002b98SHong Zhang   PetscFunctionReturn(0);
815d4002b98SHong Zhang }
816d4002b98SHong Zhang 
817d4002b98SHong Zhang PetscErrorCode MatEqual_MPISELL(Mat A,Mat B,PetscBool  *flag)
818d4002b98SHong Zhang {
819d4002b98SHong Zhang   Mat_MPISELL    *matB=(Mat_MPISELL*)B->data,*matA=(Mat_MPISELL*)A->data;
820d4002b98SHong Zhang   Mat            a,b,c,d;
821d4002b98SHong Zhang   PetscBool      flg;
822d4002b98SHong Zhang 
823d4002b98SHong Zhang   PetscFunctionBegin;
824d4002b98SHong Zhang   a = matA->A; b = matA->B;
825d4002b98SHong Zhang   c = matB->A; d = matB->B;
826d4002b98SHong Zhang 
8279566063dSJacob Faibussowitsch   PetscCall(MatEqual(a,c,&flg));
828d4002b98SHong Zhang   if (flg) {
8299566063dSJacob Faibussowitsch     PetscCall(MatEqual(b,d,&flg));
830d4002b98SHong Zhang   }
8311c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A)));
832d4002b98SHong Zhang   PetscFunctionReturn(0);
833d4002b98SHong Zhang }
834d4002b98SHong Zhang 
835d4002b98SHong Zhang PetscErrorCode MatCopy_MPISELL(Mat A,Mat B,MatStructure str)
836d4002b98SHong Zhang {
837d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
838d4002b98SHong Zhang   Mat_MPISELL    *b=(Mat_MPISELL*)B->data;
839d4002b98SHong Zhang 
840d4002b98SHong Zhang   PetscFunctionBegin;
841d4002b98SHong Zhang   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
842d4002b98SHong Zhang   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
843d4002b98SHong Zhang     /* because of the column compression in the off-processor part of the matrix a->B,
844d4002b98SHong Zhang        the number of columns in a->B and b->B may be different, hence we cannot call
845d4002b98SHong Zhang        the MatCopy() directly on the two parts. If need be, we can provide a more
846d4002b98SHong Zhang        efficient copy than the MatCopy_Basic() by first uncompressing the a->B matrices
847d4002b98SHong Zhang        then copying the submatrices */
8489566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A,B,str));
849d4002b98SHong Zhang   } else {
8509566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->A,b->A,str));
8519566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->B,b->B,str));
852d4002b98SHong Zhang   }
853d4002b98SHong Zhang   PetscFunctionReturn(0);
854d4002b98SHong Zhang }
855d4002b98SHong Zhang 
856d4002b98SHong Zhang PetscErrorCode MatSetUp_MPISELL(Mat A)
857d4002b98SHong Zhang {
858d4002b98SHong Zhang   PetscFunctionBegin;
8599566063dSJacob Faibussowitsch   PetscCall(MatMPISELLSetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL));
860d4002b98SHong Zhang   PetscFunctionReturn(0);
861d4002b98SHong Zhang }
862d4002b98SHong Zhang 
863d4002b98SHong Zhang extern PetscErrorCode MatConjugate_SeqSELL(Mat);
864d4002b98SHong Zhang 
865d4002b98SHong Zhang PetscErrorCode MatConjugate_MPISELL(Mat mat)
866d4002b98SHong Zhang {
8675f80ce2aSJacob Faibussowitsch   PetscFunctionBegin;
8685f80ce2aSJacob Faibussowitsch   if (PetscDefined(USE_COMPLEX)) {
869d4002b98SHong Zhang     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
870d4002b98SHong Zhang 
8719566063dSJacob Faibussowitsch     PetscCall(MatConjugate_SeqSELL(sell->A));
8729566063dSJacob Faibussowitsch     PetscCall(MatConjugate_SeqSELL(sell->B));
8735f80ce2aSJacob Faibussowitsch   }
874d4002b98SHong Zhang   PetscFunctionReturn(0);
875d4002b98SHong Zhang }
876d4002b98SHong Zhang 
877d4002b98SHong Zhang PetscErrorCode MatRealPart_MPISELL(Mat A)
878d4002b98SHong Zhang {
879d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
880d4002b98SHong Zhang 
881d4002b98SHong Zhang   PetscFunctionBegin;
8829566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->A));
8839566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->B));
884d4002b98SHong Zhang   PetscFunctionReturn(0);
885d4002b98SHong Zhang }
886d4002b98SHong Zhang 
887d4002b98SHong Zhang PetscErrorCode MatImaginaryPart_MPISELL(Mat A)
888d4002b98SHong Zhang {
889d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
890d4002b98SHong Zhang 
891d4002b98SHong Zhang   PetscFunctionBegin;
8929566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->A));
8939566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->B));
894d4002b98SHong Zhang   PetscFunctionReturn(0);
895d4002b98SHong Zhang }
896d4002b98SHong Zhang 
897d4002b98SHong Zhang PetscErrorCode MatInvertBlockDiagonal_MPISELL(Mat A,const PetscScalar **values)
898d4002b98SHong Zhang {
899d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
900d4002b98SHong Zhang 
901d4002b98SHong Zhang   PetscFunctionBegin;
9029566063dSJacob Faibussowitsch   PetscCall(MatInvertBlockDiagonal(a->A,values));
903d4002b98SHong Zhang   A->factorerrortype = a->A->factorerrortype;
904d4002b98SHong Zhang   PetscFunctionReturn(0);
905d4002b98SHong Zhang }
906d4002b98SHong Zhang 
907d4002b98SHong Zhang static PetscErrorCode MatSetRandom_MPISELL(Mat x,PetscRandom rctx)
908d4002b98SHong Zhang {
909d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)x->data;
910d4002b98SHong Zhang 
911d4002b98SHong Zhang   PetscFunctionBegin;
9129566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(sell->A,rctx));
9139566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(sell->B,rctx));
9149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(x,MAT_FINAL_ASSEMBLY));
9159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(x,MAT_FINAL_ASSEMBLY));
916d4002b98SHong Zhang   PetscFunctionReturn(0);
917d4002b98SHong Zhang }
918d4002b98SHong Zhang 
919d4002b98SHong Zhang PetscErrorCode MatSetFromOptions_MPISELL(PetscOptionItems *PetscOptionsObject,Mat A)
920d4002b98SHong Zhang {
921d4002b98SHong Zhang   PetscFunctionBegin;
9229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"MPISELL options"));
9239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
924d4002b98SHong Zhang   PetscFunctionReturn(0);
925d4002b98SHong Zhang }
926d4002b98SHong Zhang 
927d4002b98SHong Zhang PetscErrorCode MatShift_MPISELL(Mat Y,PetscScalar a)
928d4002b98SHong Zhang {
929d4002b98SHong Zhang   Mat_MPISELL    *msell=(Mat_MPISELL*)Y->data;
930d4002b98SHong Zhang   Mat_SeqSELL    *sell=(Mat_SeqSELL*)msell->A->data;
931d4002b98SHong Zhang 
932d4002b98SHong Zhang   PetscFunctionBegin;
933d4002b98SHong Zhang   if (!Y->preallocated) {
9349566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(Y,1,NULL,0,NULL));
935d4002b98SHong Zhang   } else if (!sell->nz) {
936d4002b98SHong Zhang     PetscInt nonew = sell->nonew;
9379566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(msell->A,1,NULL));
938d4002b98SHong Zhang     sell->nonew = nonew;
939d4002b98SHong Zhang   }
9409566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y,a));
941d4002b98SHong Zhang   PetscFunctionReturn(0);
942d4002b98SHong Zhang }
943d4002b98SHong Zhang 
944d4002b98SHong Zhang PetscErrorCode MatMissingDiagonal_MPISELL(Mat A,PetscBool  *missing,PetscInt *d)
945d4002b98SHong Zhang {
946d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
947d4002b98SHong Zhang 
948d4002b98SHong Zhang   PetscFunctionBegin;
949*08401ef6SPierre Jolivet   PetscCheck(A->rmap->n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only works for square matrices");
9509566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(a->A,missing,d));
951d4002b98SHong Zhang   if (d) {
952d4002b98SHong Zhang     PetscInt rstart;
9539566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A,&rstart,NULL));
954d4002b98SHong Zhang     *d += rstart;
955d4002b98SHong Zhang 
956d4002b98SHong Zhang   }
957d4002b98SHong Zhang   PetscFunctionReturn(0);
958d4002b98SHong Zhang }
959d4002b98SHong Zhang 
960d4002b98SHong Zhang PetscErrorCode MatGetDiagonalBlock_MPISELL(Mat A,Mat *a)
961d4002b98SHong Zhang {
962d4002b98SHong Zhang   PetscFunctionBegin;
963d4002b98SHong Zhang   *a = ((Mat_MPISELL*)A->data)->A;
964d4002b98SHong Zhang   PetscFunctionReturn(0);
965d4002b98SHong Zhang }
966d4002b98SHong Zhang 
967d4002b98SHong Zhang /* -------------------------------------------------------------------*/
968d4002b98SHong Zhang static struct _MatOps MatOps_Values = {MatSetValues_MPISELL,
969f4259b30SLisandro Dalcin                                        NULL,
970f4259b30SLisandro Dalcin                                        NULL,
971d4002b98SHong Zhang                                        MatMult_MPISELL,
972d4002b98SHong Zhang                                 /* 4*/ MatMultAdd_MPISELL,
973d4002b98SHong Zhang                                        MatMultTranspose_MPISELL,
974d4002b98SHong Zhang                                        MatMultTransposeAdd_MPISELL,
975f4259b30SLisandro Dalcin                                        NULL,
976f4259b30SLisandro Dalcin                                        NULL,
977f4259b30SLisandro Dalcin                                        NULL,
978f4259b30SLisandro Dalcin                                 /*10*/ NULL,
979f4259b30SLisandro Dalcin                                        NULL,
980f4259b30SLisandro Dalcin                                        NULL,
981d4002b98SHong Zhang                                        MatSOR_MPISELL,
982f4259b30SLisandro Dalcin                                        NULL,
983d4002b98SHong Zhang                                 /*15*/ MatGetInfo_MPISELL,
984d4002b98SHong Zhang                                        MatEqual_MPISELL,
985d4002b98SHong Zhang                                        MatGetDiagonal_MPISELL,
986d4002b98SHong Zhang                                        MatDiagonalScale_MPISELL,
987f4259b30SLisandro Dalcin                                        NULL,
988d4002b98SHong Zhang                                 /*20*/ MatAssemblyBegin_MPISELL,
989d4002b98SHong Zhang                                        MatAssemblyEnd_MPISELL,
990d4002b98SHong Zhang                                        MatSetOption_MPISELL,
991d4002b98SHong Zhang                                        MatZeroEntries_MPISELL,
992f4259b30SLisandro Dalcin                                 /*24*/ NULL,
993f4259b30SLisandro Dalcin                                        NULL,
994f4259b30SLisandro Dalcin                                        NULL,
995f4259b30SLisandro Dalcin                                        NULL,
996f4259b30SLisandro Dalcin                                        NULL,
997d4002b98SHong Zhang                                 /*29*/ MatSetUp_MPISELL,
998f4259b30SLisandro Dalcin                                        NULL,
999f4259b30SLisandro Dalcin                                        NULL,
1000d4002b98SHong Zhang                                        MatGetDiagonalBlock_MPISELL,
1001f4259b30SLisandro Dalcin                                        NULL,
1002d4002b98SHong Zhang                                 /*34*/ MatDuplicate_MPISELL,
1003f4259b30SLisandro Dalcin                                        NULL,
1004f4259b30SLisandro Dalcin                                        NULL,
1005f4259b30SLisandro Dalcin                                        NULL,
1006f4259b30SLisandro Dalcin                                        NULL,
1007f4259b30SLisandro Dalcin                                 /*39*/ NULL,
1008f4259b30SLisandro Dalcin                                        NULL,
1009f4259b30SLisandro Dalcin                                        NULL,
1010d4002b98SHong Zhang                                        MatGetValues_MPISELL,
1011d4002b98SHong Zhang                                        MatCopy_MPISELL,
1012f4259b30SLisandro Dalcin                                 /*44*/ NULL,
1013d4002b98SHong Zhang                                        MatScale_MPISELL,
1014d4002b98SHong Zhang                                        MatShift_MPISELL,
1015d4002b98SHong Zhang                                        MatDiagonalSet_MPISELL,
1016f4259b30SLisandro Dalcin                                        NULL,
1017d4002b98SHong Zhang                                 /*49*/ MatSetRandom_MPISELL,
1018f4259b30SLisandro Dalcin                                        NULL,
1019f4259b30SLisandro Dalcin                                        NULL,
1020f4259b30SLisandro Dalcin                                        NULL,
1021f4259b30SLisandro Dalcin                                        NULL,
1022d4002b98SHong Zhang                                 /*54*/ MatFDColoringCreate_MPIXAIJ,
1023f4259b30SLisandro Dalcin                                        NULL,
1024d4002b98SHong Zhang                                        MatSetUnfactored_MPISELL,
1025f4259b30SLisandro Dalcin                                        NULL,
1026f4259b30SLisandro Dalcin                                        NULL,
1027f4259b30SLisandro Dalcin                                 /*59*/ NULL,
1028d4002b98SHong Zhang                                        MatDestroy_MPISELL,
1029d4002b98SHong Zhang                                        MatView_MPISELL,
1030f4259b30SLisandro Dalcin                                        NULL,
1031f4259b30SLisandro Dalcin                                        NULL,
1032f4259b30SLisandro Dalcin                                 /*64*/ NULL,
1033f4259b30SLisandro Dalcin                                        NULL,
1034f4259b30SLisandro Dalcin                                        NULL,
1035f4259b30SLisandro Dalcin                                        NULL,
1036f4259b30SLisandro Dalcin                                        NULL,
1037f4259b30SLisandro Dalcin                                 /*69*/ NULL,
1038f4259b30SLisandro Dalcin                                        NULL,
1039f4259b30SLisandro Dalcin                                        NULL,
1040f4259b30SLisandro Dalcin                                        NULL,
1041f4259b30SLisandro Dalcin                                        NULL,
1042f4259b30SLisandro Dalcin                                        NULL,
1043d4002b98SHong Zhang                                 /*75*/ MatFDColoringApply_AIJ, /* reuse AIJ function */
1044d4002b98SHong Zhang                                        MatSetFromOptions_MPISELL,
1045f4259b30SLisandro Dalcin                                        NULL,
1046f4259b30SLisandro Dalcin                                        NULL,
1047f4259b30SLisandro Dalcin                                        NULL,
1048f4259b30SLisandro Dalcin                                 /*80*/ NULL,
1049f4259b30SLisandro Dalcin                                        NULL,
1050f4259b30SLisandro Dalcin                                        NULL,
1051f4259b30SLisandro Dalcin                                 /*83*/ NULL,
1052f4259b30SLisandro Dalcin                                        NULL,
1053f4259b30SLisandro Dalcin                                        NULL,
1054f4259b30SLisandro Dalcin                                        NULL,
1055f4259b30SLisandro Dalcin                                        NULL,
1056f4259b30SLisandro Dalcin                                        NULL,
1057f4259b30SLisandro Dalcin                                 /*89*/ NULL,
1058f4259b30SLisandro Dalcin                                        NULL,
1059f4259b30SLisandro Dalcin                                        NULL,
1060f4259b30SLisandro Dalcin                                        NULL,
1061f4259b30SLisandro Dalcin                                        NULL,
1062f4259b30SLisandro Dalcin                                 /*94*/ NULL,
1063f4259b30SLisandro Dalcin                                        NULL,
1064f4259b30SLisandro Dalcin                                        NULL,
1065f4259b30SLisandro Dalcin                                        NULL,
1066f4259b30SLisandro Dalcin                                        NULL,
1067f4259b30SLisandro Dalcin                                 /*99*/ NULL,
1068f4259b30SLisandro Dalcin                                        NULL,
1069f4259b30SLisandro Dalcin                                        NULL,
1070d4002b98SHong Zhang                                        MatConjugate_MPISELL,
1071f4259b30SLisandro Dalcin                                        NULL,
1072f4259b30SLisandro Dalcin                                 /*104*/NULL,
1073d4002b98SHong Zhang                                        MatRealPart_MPISELL,
1074d4002b98SHong Zhang                                        MatImaginaryPart_MPISELL,
1075f4259b30SLisandro Dalcin                                        NULL,
1076f4259b30SLisandro Dalcin                                        NULL,
1077f4259b30SLisandro Dalcin                                 /*109*/NULL,
1078f4259b30SLisandro Dalcin                                        NULL,
1079f4259b30SLisandro Dalcin                                        NULL,
1080f4259b30SLisandro Dalcin                                        NULL,
1081d4002b98SHong Zhang                                        MatMissingDiagonal_MPISELL,
1082f4259b30SLisandro Dalcin                                 /*114*/NULL,
1083f4259b30SLisandro Dalcin                                        NULL,
1084d4002b98SHong Zhang                                        MatGetGhosts_MPISELL,
1085f4259b30SLisandro Dalcin                                        NULL,
1086f4259b30SLisandro Dalcin                                        NULL,
1087f4259b30SLisandro Dalcin                                 /*119*/NULL,
1088f4259b30SLisandro Dalcin                                        NULL,
1089f4259b30SLisandro Dalcin                                        NULL,
1090f4259b30SLisandro Dalcin                                        NULL,
1091f4259b30SLisandro Dalcin                                        NULL,
1092f4259b30SLisandro Dalcin                                 /*124*/NULL,
1093f4259b30SLisandro Dalcin                                        NULL,
1094d4002b98SHong Zhang                                        MatInvertBlockDiagonal_MPISELL,
1095f4259b30SLisandro Dalcin                                        NULL,
1096f4259b30SLisandro Dalcin                                        NULL,
1097f4259b30SLisandro Dalcin                                 /*129*/NULL,
1098f4259b30SLisandro Dalcin                                        NULL,
1099f4259b30SLisandro Dalcin                                        NULL,
1100f4259b30SLisandro Dalcin                                        NULL,
1101f4259b30SLisandro Dalcin                                        NULL,
1102f4259b30SLisandro Dalcin                                 /*134*/NULL,
1103f4259b30SLisandro Dalcin                                        NULL,
1104f4259b30SLisandro Dalcin                                        NULL,
1105f4259b30SLisandro Dalcin                                        NULL,
1106f4259b30SLisandro Dalcin                                        NULL,
1107f4259b30SLisandro Dalcin                                 /*139*/NULL,
1108f4259b30SLisandro Dalcin                                        NULL,
1109f4259b30SLisandro Dalcin                                        NULL,
1110d4002b98SHong Zhang                                        MatFDColoringSetUp_MPIXAIJ,
1111f4259b30SLisandro Dalcin                                        NULL,
1112d70f29a3SPierre Jolivet                                 /*144*/NULL,
1113d70f29a3SPierre Jolivet                                        NULL,
1114d70f29a3SPierre Jolivet                                        NULL,
1115d70f29a3SPierre Jolivet                                        NULL
1116d4002b98SHong Zhang };
1117d4002b98SHong Zhang 
1118d4002b98SHong Zhang /* ----------------------------------------------------------------------------------------*/
1119d4002b98SHong Zhang 
1120d4002b98SHong Zhang PetscErrorCode MatStoreValues_MPISELL(Mat mat)
1121d4002b98SHong Zhang {
1122d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
1123d4002b98SHong Zhang 
1124d4002b98SHong Zhang   PetscFunctionBegin;
11259566063dSJacob Faibussowitsch   PetscCall(MatStoreValues(sell->A));
11269566063dSJacob Faibussowitsch   PetscCall(MatStoreValues(sell->B));
1127d4002b98SHong Zhang   PetscFunctionReturn(0);
1128d4002b98SHong Zhang }
1129d4002b98SHong Zhang 
1130d4002b98SHong Zhang PetscErrorCode MatRetrieveValues_MPISELL(Mat mat)
1131d4002b98SHong Zhang {
1132d4002b98SHong Zhang   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
1133d4002b98SHong Zhang 
1134d4002b98SHong Zhang   PetscFunctionBegin;
11359566063dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(sell->A));
11369566063dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(sell->B));
1137d4002b98SHong Zhang   PetscFunctionReturn(0);
1138d4002b98SHong Zhang }
1139d4002b98SHong Zhang 
1140d4002b98SHong Zhang PetscErrorCode MatMPISELLSetPreallocation_MPISELL(Mat B,PetscInt d_rlenmax,const PetscInt d_rlen[],PetscInt o_rlenmax,const PetscInt o_rlen[])
1141d4002b98SHong Zhang {
1142d4002b98SHong Zhang   Mat_MPISELL    *b;
1143d4002b98SHong Zhang 
1144d4002b98SHong Zhang   PetscFunctionBegin;
11459566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
11469566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
1147d4002b98SHong Zhang   b = (Mat_MPISELL*)B->data;
1148d4002b98SHong Zhang 
1149d4002b98SHong Zhang   if (!B->preallocated) {
1150d4002b98SHong Zhang     /* Explicitly create 2 MATSEQSELL matrices. */
11519566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF,&b->A));
11529566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n));
11539566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(b->A,B,B));
11549566063dSJacob Faibussowitsch     PetscCall(MatSetType(b->A,MATSEQSELL));
11559566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)B,(PetscObject)b->A));
11569566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF,&b->B));
11579566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N));
11589566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(b->B,B,B));
11599566063dSJacob Faibussowitsch     PetscCall(MatSetType(b->B,MATSEQSELL));
11609566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)B,(PetscObject)b->B));
1161d4002b98SHong Zhang   }
1162d4002b98SHong Zhang 
11639566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(b->A,d_rlenmax,d_rlen));
11649566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(b->B,o_rlenmax,o_rlen));
1165d4002b98SHong Zhang   B->preallocated  = PETSC_TRUE;
1166d4002b98SHong Zhang   B->was_assembled = PETSC_FALSE;
1167d4002b98SHong Zhang   /*
1168d4002b98SHong Zhang     critical for MatAssemblyEnd to work.
1169d4002b98SHong Zhang     MatAssemblyBegin checks it to set up was_assembled
1170d4002b98SHong Zhang     and MatAssemblyEnd checks was_assembled to determine whether to build garray
1171d4002b98SHong Zhang   */
1172d4002b98SHong Zhang   B->assembled     = PETSC_FALSE;
1173d4002b98SHong Zhang   PetscFunctionReturn(0);
1174d4002b98SHong Zhang }
1175d4002b98SHong Zhang 
1176d4002b98SHong Zhang PetscErrorCode MatDuplicate_MPISELL(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
1177d4002b98SHong Zhang {
1178d4002b98SHong Zhang   Mat            mat;
1179d4002b98SHong Zhang   Mat_MPISELL    *a,*oldmat=(Mat_MPISELL*)matin->data;
1180d4002b98SHong Zhang 
1181d4002b98SHong Zhang   PetscFunctionBegin;
1182f4259b30SLisandro Dalcin   *newmat = NULL;
11839566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)matin),&mat));
11849566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N));
11859566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(mat,matin,matin));
11869566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,((PetscObject)matin)->type_name));
1187d4002b98SHong Zhang   a       = (Mat_MPISELL*)mat->data;
1188d4002b98SHong Zhang 
1189d4002b98SHong Zhang   mat->factortype   = matin->factortype;
1190d4002b98SHong Zhang   mat->assembled    = PETSC_TRUE;
1191d4002b98SHong Zhang   mat->insertmode   = NOT_SET_VALUES;
1192d4002b98SHong Zhang   mat->preallocated = PETSC_TRUE;
1193d4002b98SHong Zhang 
1194d4002b98SHong Zhang   a->size         = oldmat->size;
1195d4002b98SHong Zhang   a->rank         = oldmat->rank;
1196d4002b98SHong Zhang   a->donotstash   = oldmat->donotstash;
1197d4002b98SHong Zhang   a->roworiented  = oldmat->roworiented;
1198f4259b30SLisandro Dalcin   a->rowindices   = NULL;
1199f4259b30SLisandro Dalcin   a->rowvalues    = NULL;
1200d4002b98SHong Zhang   a->getrowactive = PETSC_FALSE;
1201d4002b98SHong Zhang 
12029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->rmap,&mat->rmap));
12039566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->cmap,&mat->cmap));
1204d4002b98SHong Zhang 
1205d4002b98SHong Zhang   if (oldmat->colmap) {
1206d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
12079566063dSJacob Faibussowitsch     PetscCall(PetscTableCreateCopy(oldmat->colmap,&a->colmap));
1208d4002b98SHong Zhang #else
12099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mat->cmap->N,&a->colmap));
12109566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)mat,(mat->cmap->N)*sizeof(PetscInt)));
12119566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->colmap,oldmat->colmap,mat->cmap->N));
1212d4002b98SHong Zhang #endif
1213f4259b30SLisandro Dalcin   } else a->colmap = NULL;
1214d4002b98SHong Zhang   if (oldmat->garray) {
1215d4002b98SHong Zhang     PetscInt len;
1216d4002b98SHong Zhang     len  = oldmat->B->cmap->n;
12179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len+1,&a->garray));
12189566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt)));
12199566063dSJacob Faibussowitsch     if (len) PetscCall(PetscArraycpy(a->garray,oldmat->garray,len));
1220f4259b30SLisandro Dalcin   } else a->garray = NULL;
1221d4002b98SHong Zhang 
12229566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(oldmat->lvec,&a->lvec));
12239566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec));
12249566063dSJacob Faibussowitsch   PetscCall(VecScatterCopy(oldmat->Mvctx,&a->Mvctx));
12259566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx));
12269566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->A,cpvalues,&a->A));
12279566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A));
12289566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->B,cpvalues,&a->B));
12299566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B));
12309566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist));
1231d4002b98SHong Zhang   *newmat = mat;
1232d4002b98SHong Zhang   PetscFunctionReturn(0);
1233d4002b98SHong Zhang }
1234d4002b98SHong Zhang 
1235d4002b98SHong Zhang /*@C
1236d4002b98SHong Zhang    MatMPISELLSetPreallocation - Preallocates memory for a sparse parallel matrix in sell format.
1237d4002b98SHong Zhang    For good matrix assembly performance the user should preallocate the matrix storage by
1238d4002b98SHong Zhang    setting the parameters d_nz (or d_nnz) and o_nz (or o_nnz).
1239d4002b98SHong Zhang 
1240d083f849SBarry Smith    Collective
1241d4002b98SHong Zhang 
1242d4002b98SHong Zhang    Input Parameters:
1243d4002b98SHong Zhang +  B - the matrix
1244d4002b98SHong Zhang .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1245d4002b98SHong Zhang            (same value is used for all local rows)
1246d4002b98SHong Zhang .  d_nnz - array containing the number of nonzeros in the various rows of the
1247d4002b98SHong Zhang            DIAGONAL portion of the local submatrix (possibly different for each row)
1248d4002b98SHong Zhang            or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
1249d4002b98SHong Zhang            The size of this array is equal to the number of local rows, i.e 'm'.
1250d4002b98SHong Zhang            For matrices that will be factored, you must leave room for (and set)
1251d4002b98SHong Zhang            the diagonal entry even if it is zero.
1252d4002b98SHong Zhang .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1253d4002b98SHong Zhang            submatrix (same value is used for all local rows).
1254d4002b98SHong Zhang -  o_nnz - array containing the number of nonzeros in the various rows of the
1255d4002b98SHong Zhang            OFF-DIAGONAL portion of the local submatrix (possibly different for
1256d4002b98SHong Zhang            each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
1257d4002b98SHong Zhang            structure. The size of this array is equal to the number
1258d4002b98SHong Zhang            of local rows, i.e 'm'.
1259d4002b98SHong Zhang 
1260d4002b98SHong Zhang    If the *_nnz parameter is given then the *_nz parameter is ignored
1261d4002b98SHong Zhang 
1262d4002b98SHong Zhang    The stored row and column indices begin with zero.
1263d4002b98SHong Zhang 
1264d4002b98SHong Zhang    The parallel matrix is partitioned such that the first m0 rows belong to
1265d4002b98SHong Zhang    process 0, the next m1 rows belong to process 1, the next m2 rows belong
1266d4002b98SHong Zhang    to process 2 etc.. where m0,m1,m2... are the input parameter 'm'.
1267d4002b98SHong Zhang 
1268d4002b98SHong Zhang    The DIAGONAL portion of the local submatrix of a processor can be defined
1269d4002b98SHong Zhang    as the submatrix which is obtained by extraction the part corresponding to
1270d4002b98SHong Zhang    the rows r1-r2 and columns c1-c2 of the global matrix, where r1 is the
1271d4002b98SHong Zhang    first row that belongs to the processor, r2 is the last row belonging to
1272d4002b98SHong Zhang    the this processor, and c1-c2 is range of indices of the local part of a
1273d4002b98SHong Zhang    vector suitable for applying the matrix to.  This is an mxn matrix.  In the
1274d4002b98SHong Zhang    common case of a square matrix, the row and column ranges are the same and
1275d4002b98SHong Zhang    the DIAGONAL part is also square. The remaining portion of the local
1276d4002b98SHong Zhang    submatrix (mxN) constitute the OFF-DIAGONAL portion.
1277d4002b98SHong Zhang 
1278d4002b98SHong Zhang    If o_nnz, d_nnz are specified, then o_nz, and d_nz are ignored.
1279d4002b98SHong Zhang 
1280d4002b98SHong Zhang    You can call MatGetInfo() to get information on how effective the preallocation was;
1281d4002b98SHong Zhang    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
1282d4002b98SHong Zhang    You can also run with the option -info and look for messages with the string
1283d4002b98SHong Zhang    malloc in them to see if additional memory allocation was needed.
1284d4002b98SHong Zhang 
1285d4002b98SHong Zhang    Example usage:
1286d4002b98SHong Zhang 
1287d4002b98SHong Zhang    Consider the following 8x8 matrix with 34 non-zero values, that is
1288d4002b98SHong Zhang    assembled across 3 processors. Lets assume that proc0 owns 3 rows,
1289d4002b98SHong Zhang    proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
1290d4002b98SHong Zhang    as follows:
1291d4002b98SHong Zhang 
1292d4002b98SHong Zhang .vb
1293d4002b98SHong Zhang             1  2  0  |  0  3  0  |  0  4
1294d4002b98SHong Zhang     Proc0   0  5  6  |  7  0  0  |  8  0
1295d4002b98SHong Zhang             9  0 10  | 11  0  0  | 12  0
1296d4002b98SHong Zhang     -------------------------------------
1297d4002b98SHong Zhang            13  0 14  | 15 16 17  |  0  0
1298d4002b98SHong Zhang     Proc1   0 18  0  | 19 20 21  |  0  0
1299d4002b98SHong Zhang             0  0  0  | 22 23  0  | 24  0
1300d4002b98SHong Zhang     -------------------------------------
1301d4002b98SHong Zhang     Proc2  25 26 27  |  0  0 28  | 29  0
1302d4002b98SHong Zhang            30  0  0  | 31 32 33  |  0 34
1303d4002b98SHong Zhang .ve
1304d4002b98SHong Zhang 
1305d4002b98SHong Zhang    This can be represented as a collection of submatrices as:
1306d4002b98SHong Zhang 
1307d4002b98SHong Zhang .vb
1308d4002b98SHong Zhang       A B C
1309d4002b98SHong Zhang       D E F
1310d4002b98SHong Zhang       G H I
1311d4002b98SHong Zhang .ve
1312d4002b98SHong Zhang 
1313d4002b98SHong Zhang    Where the submatrices A,B,C are owned by proc0, D,E,F are
1314d4002b98SHong Zhang    owned by proc1, G,H,I are owned by proc2.
1315d4002b98SHong Zhang 
1316d4002b98SHong Zhang    The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1317d4002b98SHong Zhang    The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1318d4002b98SHong Zhang    The 'M','N' parameters are 8,8, and have the same values on all procs.
1319d4002b98SHong Zhang 
1320d4002b98SHong Zhang    The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
1321d4002b98SHong Zhang    submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
1322d4002b98SHong Zhang    corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
1323d4002b98SHong Zhang    Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
1324d4002b98SHong Zhang    part as SeqSELL matrices. for eg: proc1 will store [E] as a SeqSELL
1325d4002b98SHong Zhang    matrix, ans [DF] as another SeqSELL matrix.
1326d4002b98SHong Zhang 
1327d4002b98SHong Zhang    When d_nz, o_nz parameters are specified, d_nz storage elements are
1328d4002b98SHong Zhang    allocated for every row of the local diagonal submatrix, and o_nz
1329d4002b98SHong Zhang    storage locations are allocated for every row of the OFF-DIAGONAL submat.
1330d4002b98SHong Zhang    One way to choose d_nz and o_nz is to use the max nonzerors per local
1331d4002b98SHong Zhang    rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
1332d4002b98SHong Zhang    In this case, the values of d_nz,o_nz are:
1333d4002b98SHong Zhang .vb
1334d4002b98SHong Zhang      proc0 : dnz = 2, o_nz = 2
1335d4002b98SHong Zhang      proc1 : dnz = 3, o_nz = 2
1336d4002b98SHong Zhang      proc2 : dnz = 1, o_nz = 4
1337d4002b98SHong Zhang .ve
1338d4002b98SHong Zhang    We are allocating m*(d_nz+o_nz) storage locations for every proc. This
1339d4002b98SHong Zhang    translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
1340d4002b98SHong Zhang    for proc3. i.e we are using 12+15+10=37 storage locations to store
1341d4002b98SHong Zhang    34 values.
1342d4002b98SHong Zhang 
1343d4002b98SHong Zhang    When d_nnz, o_nnz parameters are specified, the storage is specified
1344a5b23f4aSJose E. Roman    for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
1345d4002b98SHong Zhang    In the above case the values for d_nnz,o_nnz are:
1346d4002b98SHong Zhang .vb
1347d4002b98SHong Zhang      proc0: d_nnz = [2,2,2] and o_nnz = [2,2,2]
1348d4002b98SHong Zhang      proc1: d_nnz = [3,3,2] and o_nnz = [2,1,1]
1349d4002b98SHong Zhang      proc2: d_nnz = [1,1]   and o_nnz = [4,4]
1350d4002b98SHong Zhang .ve
1351d4002b98SHong Zhang    Here the space allocated is according to nz (or maximum values in the nnz
1352d4002b98SHong Zhang    if nnz is provided) for DIAGONAL and OFF-DIAGONAL submatrices, i.e (2+2+3+2)*3+(1+4)*2=37
1353d4002b98SHong Zhang 
1354d4002b98SHong Zhang    Level: intermediate
1355d4002b98SHong Zhang 
1356d4002b98SHong Zhang .seealso: MatCreate(), MatCreateSeqSELL(), MatSetValues(), MatCreatesell(),
1357d4002b98SHong Zhang           MATMPISELL, MatGetInfo(), PetscSplitOwnership()
1358d4002b98SHong Zhang @*/
1359d4002b98SHong Zhang PetscErrorCode MatMPISELLSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1360d4002b98SHong Zhang {
1361d4002b98SHong Zhang   PetscFunctionBegin;
1362d4002b98SHong Zhang   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
1363d4002b98SHong Zhang   PetscValidType(B,1);
1364cac4c232SBarry Smith   PetscTryMethod(B,"MatMPISELLSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));
1365d4002b98SHong Zhang   PetscFunctionReturn(0);
1366d4002b98SHong Zhang }
1367d4002b98SHong Zhang 
1368ed73aabaSBarry Smith /*MC
1369ed73aabaSBarry Smith    MATMPISELL - MATMPISELL = "mpisell" - A matrix type to be used for MPI sparse matrices,
1370ed73aabaSBarry Smith    based on the sliced Ellpack format
1371ed73aabaSBarry Smith 
1372ed73aabaSBarry Smith    Options Database Keys:
1373ed73aabaSBarry Smith . -mat_type sell - sets the matrix type to "seqsell" during a call to MatSetFromOptions()
1374ed73aabaSBarry Smith 
1375ed73aabaSBarry Smith    Level: beginner
1376ed73aabaSBarry Smith 
1377ed73aabaSBarry Smith .seealso: MatCreateSell(), MATSEQSELL, MATSELL, MATSEQAIJ, MATAIJ, MATMPIAIJ
1378ed73aabaSBarry Smith M*/
1379ed73aabaSBarry Smith 
1380d4002b98SHong Zhang /*@C
1381d4002b98SHong Zhang    MatCreateSELL - Creates a sparse parallel matrix in SELL format.
1382d4002b98SHong Zhang 
1383d083f849SBarry Smith    Collective
1384d4002b98SHong Zhang 
1385d4002b98SHong Zhang    Input Parameters:
1386d4002b98SHong Zhang +  comm - MPI communicator
1387d4002b98SHong Zhang .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1388d4002b98SHong Zhang            This value should be the same as the local size used in creating the
1389d4002b98SHong Zhang            y vector for the matrix-vector product y = Ax.
1390d4002b98SHong Zhang .  n - This value should be the same as the local size used in creating the
1391d4002b98SHong Zhang        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
1392d4002b98SHong Zhang        calculated if N is given) For square matrices n is almost always m.
1393d4002b98SHong Zhang .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
1394d4002b98SHong Zhang .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
1395d4002b98SHong Zhang .  d_rlenmax - max number of nonzeros per row in DIAGONAL portion of local submatrix
1396d4002b98SHong Zhang                (same value is used for all local rows)
1397d4002b98SHong Zhang .  d_rlen - array containing the number of nonzeros in the various rows of the
1398d4002b98SHong Zhang             DIAGONAL portion of the local submatrix (possibly different for each row)
1399d4002b98SHong Zhang             or NULL, if d_rlenmax is used to specify the nonzero structure.
1400d4002b98SHong Zhang             The size of this array is equal to the number of local rows, i.e 'm'.
1401d4002b98SHong Zhang .  o_rlenmax - max number of nonzeros per row in the OFF-DIAGONAL portion of local
1402d4002b98SHong Zhang                submatrix (same value is used for all local rows).
1403d4002b98SHong Zhang -  o_rlen - array containing the number of nonzeros in the various rows of the
1404d4002b98SHong Zhang             OFF-DIAGONAL portion of the local submatrix (possibly different for
1405d4002b98SHong Zhang             each row) or NULL, if o_rlenmax is used to specify the nonzero
1406d4002b98SHong Zhang             structure. The size of this array is equal to the number
1407d4002b98SHong Zhang             of local rows, i.e 'm'.
1408d4002b98SHong Zhang 
1409d4002b98SHong Zhang    Output Parameter:
1410d4002b98SHong Zhang .  A - the matrix
1411d4002b98SHong Zhang 
1412d4002b98SHong Zhang    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
1413f6f02116SRichard Tran Mills    MatXXXXSetPreallocation() paradigm instead of this routine directly.
1414d4002b98SHong Zhang    [MatXXXXSetPreallocation() is, for example, MatSeqSELLSetPreallocation]
1415d4002b98SHong Zhang 
1416d4002b98SHong Zhang    Notes:
1417d4002b98SHong Zhang    If the *_rlen parameter is given then the *_rlenmax parameter is ignored
1418d4002b98SHong Zhang 
1419d4002b98SHong Zhang    m,n,M,N parameters specify the size of the matrix, and its partitioning across
1420d4002b98SHong Zhang    processors, while d_rlenmax,d_rlen,o_rlenmax,o_rlen parameters specify the approximate
1421d4002b98SHong Zhang    storage requirements for this matrix.
1422d4002b98SHong Zhang 
1423d4002b98SHong Zhang    If PETSC_DECIDE or  PETSC_DETERMINE is used for a particular argument on one
1424d4002b98SHong Zhang    processor than it must be used on all processors that share the object for
1425d4002b98SHong Zhang    that argument.
1426d4002b98SHong Zhang 
1427d4002b98SHong Zhang    The user MUST specify either the local or global matrix dimensions
1428d4002b98SHong Zhang    (possibly both).
1429d4002b98SHong Zhang 
1430d4002b98SHong Zhang    The parallel matrix is partitioned across processors such that the
1431d4002b98SHong Zhang    first m0 rows belong to process 0, the next m1 rows belong to
1432d4002b98SHong Zhang    process 1, the next m2 rows belong to process 2 etc.. where
1433d4002b98SHong Zhang    m0,m1,m2,.. are the input parameter 'm'. i.e each processor stores
1434d4002b98SHong Zhang    values corresponding to [m x N] submatrix.
1435d4002b98SHong Zhang 
1436d4002b98SHong Zhang    The columns are logically partitioned with the n0 columns belonging
1437d4002b98SHong Zhang    to 0th partition, the next n1 columns belonging to the next
1438d4002b98SHong Zhang    partition etc.. where n0,n1,n2... are the input parameter 'n'.
1439d4002b98SHong Zhang 
1440d4002b98SHong Zhang    The DIAGONAL portion of the local submatrix on any given processor
1441d4002b98SHong Zhang    is the submatrix corresponding to the rows and columns m,n
1442d4002b98SHong Zhang    corresponding to the given processor. i.e diagonal matrix on
1443d4002b98SHong Zhang    process 0 is [m0 x n0], diagonal matrix on process 1 is [m1 x n1]
1444d4002b98SHong Zhang    etc. The remaining portion of the local submatrix [m x (N-n)]
1445d4002b98SHong Zhang    constitute the OFF-DIAGONAL portion. The example below better
1446d4002b98SHong Zhang    illustrates this concept.
1447d4002b98SHong Zhang 
1448d4002b98SHong Zhang    For a square global matrix we define each processor's diagonal portion
1449d4002b98SHong Zhang    to be its local rows and the corresponding columns (a square submatrix);
1450d4002b98SHong Zhang    each processor's off-diagonal portion encompasses the remainder of the
1451d4002b98SHong Zhang    local matrix (a rectangular submatrix).
1452d4002b98SHong Zhang 
1453d4002b98SHong Zhang    If o_rlen, d_rlen are specified, then o_rlenmax, and d_rlenmax are ignored.
1454d4002b98SHong Zhang 
1455d4002b98SHong Zhang    When calling this routine with a single process communicator, a matrix of
1456d4002b98SHong Zhang    type SEQSELL is returned.  If a matrix of type MATMPISELL is desired for this
1457d4002b98SHong Zhang    type of communicator, use the construction mechanism:
1458d4002b98SHong Zhang      MatCreate(...,&A); MatSetType(A,MATMPISELL); MatSetSizes(A, m,n,M,N); MatMPISELLSetPreallocation(A,...);
1459d4002b98SHong Zhang 
1460d4002b98SHong Zhang    Options Database Keys:
1461d4002b98SHong Zhang -  -mat_sell_oneindex - Internally use indexing starting at 1
1462d4002b98SHong Zhang         rather than 0.  Note that when calling MatSetValues(),
1463d4002b98SHong Zhang         the user still MUST index entries starting at 0!
1464d4002b98SHong Zhang 
1465d4002b98SHong Zhang    Example usage:
1466d4002b98SHong Zhang 
1467d4002b98SHong Zhang    Consider the following 8x8 matrix with 34 non-zero values, that is
1468d4002b98SHong Zhang    assembled across 3 processors. Lets assume that proc0 owns 3 rows,
1469d4002b98SHong Zhang    proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
1470d4002b98SHong Zhang    as follows:
1471d4002b98SHong Zhang 
1472d4002b98SHong Zhang .vb
1473d4002b98SHong Zhang             1  2  0  |  0  3  0  |  0  4
1474d4002b98SHong Zhang     Proc0   0  5  6  |  7  0  0  |  8  0
1475d4002b98SHong Zhang             9  0 10  | 11  0  0  | 12  0
1476d4002b98SHong Zhang     -------------------------------------
1477d4002b98SHong Zhang            13  0 14  | 15 16 17  |  0  0
1478d4002b98SHong Zhang     Proc1   0 18  0  | 19 20 21  |  0  0
1479d4002b98SHong Zhang             0  0  0  | 22 23  0  | 24  0
1480d4002b98SHong Zhang     -------------------------------------
1481d4002b98SHong Zhang     Proc2  25 26 27  |  0  0 28  | 29  0
1482d4002b98SHong Zhang            30  0  0  | 31 32 33  |  0 34
1483d4002b98SHong Zhang .ve
1484d4002b98SHong Zhang 
1485d4002b98SHong Zhang    This can be represented as a collection of submatrices as:
1486d4002b98SHong Zhang 
1487d4002b98SHong Zhang .vb
1488d4002b98SHong Zhang       A B C
1489d4002b98SHong Zhang       D E F
1490d4002b98SHong Zhang       G H I
1491d4002b98SHong Zhang .ve
1492d4002b98SHong Zhang 
1493d4002b98SHong Zhang    Where the submatrices A,B,C are owned by proc0, D,E,F are
1494d4002b98SHong Zhang    owned by proc1, G,H,I are owned by proc2.
1495d4002b98SHong Zhang 
1496d4002b98SHong Zhang    The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1497d4002b98SHong Zhang    The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1498d4002b98SHong Zhang    The 'M','N' parameters are 8,8, and have the same values on all procs.
1499d4002b98SHong Zhang 
1500d4002b98SHong Zhang    The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
1501d4002b98SHong Zhang    submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
1502d4002b98SHong Zhang    corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
1503d4002b98SHong Zhang    Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
1504d4002b98SHong Zhang    part as SeqSELL matrices. for eg: proc1 will store [E] as a SeqSELL
1505d4002b98SHong Zhang    matrix, ans [DF] as another SeqSELL matrix.
1506d4002b98SHong Zhang 
1507d4002b98SHong Zhang    When d_rlenmax, o_rlenmax parameters are specified, d_rlenmax storage elements are
1508d4002b98SHong Zhang    allocated for every row of the local diagonal submatrix, and o_rlenmax
1509d4002b98SHong Zhang    storage locations are allocated for every row of the OFF-DIAGONAL submat.
1510d4002b98SHong Zhang    One way to choose d_rlenmax and o_rlenmax is to use the max nonzerors per local
1511d4002b98SHong Zhang    rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
1512d4002b98SHong Zhang    In this case, the values of d_rlenmax,o_rlenmax are:
1513d4002b98SHong Zhang .vb
1514d4002b98SHong Zhang      proc0 : d_rlenmax = 2, o_rlenmax = 2
1515d4002b98SHong Zhang      proc1 : d_rlenmax = 3, o_rlenmax = 2
1516d4002b98SHong Zhang      proc2 : d_rlenmax = 1, o_rlenmax = 4
1517d4002b98SHong Zhang .ve
1518d4002b98SHong Zhang    We are allocating m*(d_rlenmax+o_rlenmax) storage locations for every proc. This
1519d4002b98SHong Zhang    translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
1520d4002b98SHong Zhang    for proc3. i.e we are using 12+15+10=37 storage locations to store
1521d4002b98SHong Zhang    34 values.
1522d4002b98SHong Zhang 
1523d4002b98SHong Zhang    When d_rlen, o_rlen parameters are specified, the storage is specified
1524a5b23f4aSJose E. Roman    for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
1525d4002b98SHong Zhang    In the above case the values for d_nnz,o_nnz are:
1526d4002b98SHong Zhang .vb
1527d4002b98SHong Zhang      proc0: d_nnz = [2,2,2] and o_nnz = [2,2,2]
1528d4002b98SHong Zhang      proc1: d_nnz = [3,3,2] and o_nnz = [2,1,1]
1529d4002b98SHong Zhang      proc2: d_nnz = [1,1]   and o_nnz = [4,4]
1530d4002b98SHong Zhang .ve
1531d4002b98SHong Zhang    Here the space allocated is still 37 though there are 34 nonzeros because
1532d4002b98SHong Zhang    the allocation is always done according to rlenmax.
1533d4002b98SHong Zhang 
1534d4002b98SHong Zhang    Level: intermediate
1535d4002b98SHong Zhang 
1536d4002b98SHong Zhang .seealso: MatCreate(), MatCreateSeqSELL(), MatSetValues(), MatMPISELLSetPreallocation(), MatMPISELLSetPreallocationSELL(),
1537d4002b98SHong Zhang           MATMPISELL, MatCreateMPISELLWithArrays()
1538d4002b98SHong Zhang @*/
1539d4002b98SHong Zhang PetscErrorCode MatCreateSELL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_rlenmax,const PetscInt d_rlen[],PetscInt o_rlenmax,const PetscInt o_rlen[],Mat *A)
1540d4002b98SHong Zhang {
1541d4002b98SHong Zhang   PetscMPIInt    size;
1542d4002b98SHong Zhang 
1543d4002b98SHong Zhang   PetscFunctionBegin;
15449566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,A));
15459566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A,m,n,M,N));
15469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
1547d4002b98SHong Zhang   if (size > 1) {
15489566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A,MATMPISELL));
15499566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(*A,d_rlenmax,d_rlen,o_rlenmax,o_rlen));
1550d4002b98SHong Zhang   } else {
15519566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A,MATSEQSELL));
15529566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(*A,d_rlenmax,d_rlen));
1553d4002b98SHong Zhang   }
1554d4002b98SHong Zhang   PetscFunctionReturn(0);
1555d4002b98SHong Zhang }
1556d4002b98SHong Zhang 
1557d4002b98SHong Zhang PetscErrorCode MatMPISELLGetSeqSELL(Mat A,Mat *Ad,Mat *Ao,const PetscInt *colmap[])
1558d4002b98SHong Zhang {
1559d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
1560d4002b98SHong Zhang   PetscBool      flg;
1561d4002b98SHong Zhang 
1562d4002b98SHong Zhang   PetscFunctionBegin;
15639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPISELL,&flg));
156428b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"This function requires a MATMPISELL matrix as input");
1565d4002b98SHong Zhang   if (Ad)     *Ad     = a->A;
1566d4002b98SHong Zhang   if (Ao)     *Ao     = a->B;
1567d4002b98SHong Zhang   if (colmap) *colmap = a->garray;
1568d4002b98SHong Zhang   PetscFunctionReturn(0);
1569d4002b98SHong Zhang }
1570d4002b98SHong Zhang 
1571d4002b98SHong Zhang /*@C
1572d4002b98SHong Zhang      MatMPISELLGetLocalMatCondensed - Creates a SeqSELL matrix from an MATMPISELL matrix by taking all its local rows and NON-ZERO columns
1573d4002b98SHong Zhang 
1574d4002b98SHong Zhang     Not Collective
1575d4002b98SHong Zhang 
1576d4002b98SHong Zhang    Input Parameters:
1577d4002b98SHong Zhang +    A - the matrix
1578d4002b98SHong Zhang .    scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
1579d4002b98SHong Zhang -    row, col - index sets of rows and columns to extract (or NULL)
1580d4002b98SHong Zhang 
1581d4002b98SHong Zhang    Output Parameter:
1582d4002b98SHong Zhang .    A_loc - the local sequential matrix generated
1583d4002b98SHong Zhang 
1584d4002b98SHong Zhang     Level: developer
1585d4002b98SHong Zhang 
1586d4002b98SHong Zhang .seealso: MatGetOwnershipRange(), MatMPISELLGetLocalMat()
1587d4002b98SHong Zhang 
1588d4002b98SHong Zhang @*/
1589d4002b98SHong Zhang PetscErrorCode MatMPISELLGetLocalMatCondensed(Mat A,MatReuse scall,IS *row,IS *col,Mat *A_loc)
1590d4002b98SHong Zhang {
1591d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
1592d4002b98SHong Zhang   PetscInt       i,start,end,ncols,nzA,nzB,*cmap,imark,*idx;
1593d4002b98SHong Zhang   IS             isrowa,iscola;
1594d4002b98SHong Zhang   Mat            *aloc;
1595d4002b98SHong Zhang   PetscBool      match;
1596d4002b98SHong Zhang 
1597d4002b98SHong Zhang   PetscFunctionBegin;
15989566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPISELL,&match));
159928b400f6SJacob Faibussowitsch   PetscCheck(match,PetscObjectComm((PetscObject)A), PETSC_ERR_SUP,"Requires MATMPISELL matrix as input");
16009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_Getlocalmatcondensed,A,0,0,0));
1601d4002b98SHong Zhang   if (!row) {
1602d4002b98SHong Zhang     start = A->rmap->rstart; end = A->rmap->rend;
16039566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&isrowa));
1604d4002b98SHong Zhang   } else {
1605d4002b98SHong Zhang     isrowa = *row;
1606d4002b98SHong Zhang   }
1607d4002b98SHong Zhang   if (!col) {
1608d4002b98SHong Zhang     start = A->cmap->rstart;
1609d4002b98SHong Zhang     cmap  = a->garray;
1610d4002b98SHong Zhang     nzA   = a->A->cmap->n;
1611d4002b98SHong Zhang     nzB   = a->B->cmap->n;
16129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nzA+nzB, &idx));
1613d4002b98SHong Zhang     ncols = 0;
1614d4002b98SHong Zhang     for (i=0; i<nzB; i++) {
1615d4002b98SHong Zhang       if (cmap[i] < start) idx[ncols++] = cmap[i];
1616d4002b98SHong Zhang       else break;
1617d4002b98SHong Zhang     }
1618d4002b98SHong Zhang     imark = i;
1619d4002b98SHong Zhang     for (i=0; i<nzA; i++) idx[ncols++] = start + i;
1620d4002b98SHong Zhang     for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i];
16219566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&iscola));
1622d4002b98SHong Zhang   } else {
1623d4002b98SHong Zhang     iscola = *col;
1624d4002b98SHong Zhang   }
1625d4002b98SHong Zhang   if (scall != MAT_INITIAL_MATRIX) {
16269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1,&aloc));
1627d4002b98SHong Zhang     aloc[0] = *A_loc;
1628d4002b98SHong Zhang   }
16299566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,1,&isrowa,&iscola,scall,&aloc));
1630d4002b98SHong Zhang   *A_loc = aloc[0];
16319566063dSJacob Faibussowitsch   PetscCall(PetscFree(aloc));
1632d4002b98SHong Zhang   if (!row) {
16339566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isrowa));
1634d4002b98SHong Zhang   }
1635d4002b98SHong Zhang   if (!col) {
16369566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iscola));
1637d4002b98SHong Zhang   }
16389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_Getlocalmatcondensed,A,0,0,0));
1639d4002b98SHong Zhang   PetscFunctionReturn(0);
1640d4002b98SHong Zhang }
1641d4002b98SHong Zhang 
1642d4002b98SHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>
1643d4002b98SHong Zhang 
1644d4002b98SHong Zhang PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1645d4002b98SHong Zhang {
1646d4002b98SHong Zhang   Mat_MPISELL    *a=(Mat_MPISELL*)A->data;
1647d4002b98SHong Zhang   Mat            B;
1648d4002b98SHong Zhang   Mat_MPIAIJ     *b;
1649d4002b98SHong Zhang 
1650d4002b98SHong Zhang   PetscFunctionBegin;
165128b400f6SJacob Faibussowitsch   PetscCheck(A->assembled,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Matrix must be assembled");
1652d4002b98SHong Zhang 
165394a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
165494a8b381SRichard Tran Mills     B = *newmat;
165594a8b381SRichard Tran Mills   } else {
16569566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
16579566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,MATMPIAIJ));
16589566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
16599566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B,A->rmap->bs,A->cmap->bs));
16609566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B,0,NULL));
16619566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(B,0,NULL,0,NULL));
166294a8b381SRichard Tran Mills   }
1663d4002b98SHong Zhang   b    = (Mat_MPIAIJ*) B->data;
166494a8b381SRichard Tran Mills 
166594a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
16669566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_REUSE_MATRIX, &b->A));
16679566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_REUSE_MATRIX, &b->B));
166894a8b381SRichard Tran Mills   } else {
16699566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->A));
16709566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->B));
16719566063dSJacob Faibussowitsch     PetscCall(MatDisAssemble_MPISELL(A));
16729566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->A));
16739566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->B));
16749566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
16759566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
16769566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
16779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
167894a8b381SRichard Tran Mills   }
1679d4002b98SHong Zhang 
1680d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
16819566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&B));
1682d4002b98SHong Zhang   } else {
1683d4002b98SHong Zhang     *newmat = B;
1684d4002b98SHong Zhang   }
1685d4002b98SHong Zhang   PetscFunctionReturn(0);
1686d4002b98SHong Zhang }
1687d4002b98SHong Zhang 
1688d4002b98SHong Zhang PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1689d4002b98SHong Zhang {
1690d4002b98SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
1691d4002b98SHong Zhang   Mat            B;
1692d4002b98SHong Zhang   Mat_MPISELL    *b;
1693d4002b98SHong Zhang 
1694d4002b98SHong Zhang   PetscFunctionBegin;
169528b400f6SJacob Faibussowitsch   PetscCheck(A->assembled,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Matrix must be assembled");
1696d4002b98SHong Zhang 
169794a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
169894a8b381SRichard Tran Mills     B = *newmat;
169994a8b381SRichard Tran Mills   } else {
17009566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
17019566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,MATMPISELL));
17029566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
17039566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B,A->rmap->bs,A->cmap->bs));
17049566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B,0,NULL));
17059566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(B,0,NULL,0,NULL));
170694a8b381SRichard Tran Mills   }
1707d4002b98SHong Zhang   b    = (Mat_MPISELL*) B->data;
170894a8b381SRichard Tran Mills 
170994a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
17109566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_REUSE_MATRIX, &b->A));
17119566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_REUSE_MATRIX, &b->B));
171294a8b381SRichard Tran Mills   } else {
17139566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->A));
17149566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->B));
17159566063dSJacob Faibussowitsch     PetscCall(MatDisAssemble_MPIAIJ(A));
17169566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_INITIAL_MATRIX, &b->A));
17179566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_INITIAL_MATRIX, &b->B));
17189566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
17199566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
17209566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
17219566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
172294a8b381SRichard Tran Mills   }
1723d4002b98SHong Zhang 
1724d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
17259566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&B));
1726d4002b98SHong Zhang   } else {
1727d4002b98SHong Zhang     *newmat = B;
1728d4002b98SHong Zhang   }
1729d4002b98SHong Zhang   PetscFunctionReturn(0);
1730d4002b98SHong Zhang }
1731d4002b98SHong Zhang 
1732d4002b98SHong Zhang PetscErrorCode MatSOR_MPISELL(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1733d4002b98SHong Zhang {
1734d4002b98SHong Zhang   Mat_MPISELL    *mat=(Mat_MPISELL*)matin->data;
1735f4259b30SLisandro Dalcin   Vec            bb1=NULL;
1736d4002b98SHong Zhang 
1737d4002b98SHong Zhang   PetscFunctionBegin;
1738d4002b98SHong Zhang   if (flag == SOR_APPLY_UPPER) {
17399566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx));
1740d4002b98SHong Zhang     PetscFunctionReturn(0);
1741d4002b98SHong Zhang   }
1742d4002b98SHong Zhang 
1743d4002b98SHong Zhang   if (its > 1 || ~flag & SOR_ZERO_INITIAL_GUESS || flag & SOR_EISENSTAT) {
17449566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(bb,&bb1));
1745d4002b98SHong Zhang   }
1746d4002b98SHong Zhang 
1747d4002b98SHong Zhang   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
1748d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17499566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx));
1750d4002b98SHong Zhang       its--;
1751d4002b98SHong Zhang     }
1752d4002b98SHong Zhang 
1753d4002b98SHong Zhang     while (its--) {
17549566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD));
17559566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD));
1756d4002b98SHong Zhang 
1757d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17589566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec,-1.0));
17599566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1));
1760d4002b98SHong Zhang 
1761d4002b98SHong Zhang       /* local sweep */
17629566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,1,xx));
1763d4002b98SHong Zhang     }
1764d4002b98SHong Zhang   } else if (flag & SOR_LOCAL_FORWARD_SWEEP) {
1765d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17669566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx));
1767d4002b98SHong Zhang       its--;
1768d4002b98SHong Zhang     }
1769d4002b98SHong Zhang     while (its--) {
17709566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD));
17719566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD));
1772d4002b98SHong Zhang 
1773d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17749566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec,-1.0));
17759566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1));
1776d4002b98SHong Zhang 
1777d4002b98SHong Zhang       /* local sweep */
17789566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A,bb1,omega,SOR_FORWARD_SWEEP,fshift,lits,1,xx));
1779d4002b98SHong Zhang     }
1780d4002b98SHong Zhang   } else if (flag & SOR_LOCAL_BACKWARD_SWEEP) {
1781d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17829566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx));
1783d4002b98SHong Zhang       its--;
1784d4002b98SHong Zhang     }
1785d4002b98SHong Zhang     while (its--) {
17869566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD));
17879566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD));
1788d4002b98SHong Zhang 
1789d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17909566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec,-1.0));
17919566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1));
1792d4002b98SHong Zhang 
1793d4002b98SHong Zhang       /* local sweep */
17949566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A,bb1,omega,SOR_BACKWARD_SWEEP,fshift,lits,1,xx));
1795d4002b98SHong Zhang     }
1796d4002b98SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)matin),PETSC_ERR_SUP,"Parallel SOR not supported");
1797d4002b98SHong Zhang 
17989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bb1));
1799d4002b98SHong Zhang 
1800d4002b98SHong Zhang   matin->factorerrortype = mat->A->factorerrortype;
1801d4002b98SHong Zhang   PetscFunctionReturn(0);
1802d4002b98SHong Zhang }
1803d4002b98SHong Zhang 
1804d4002b98SHong Zhang /*MC
1805d4002b98SHong Zhang    MATMPISELL - MATMPISELL = "MPISELL" - A matrix type to be used for parallel sparse matrices.
1806d4002b98SHong Zhang 
1807d4002b98SHong Zhang    Options Database Keys:
1808d4002b98SHong Zhang . -mat_type MPISELL - sets the matrix type to "MPISELL" during a call to MatSetFromOptions()
1809d4002b98SHong Zhang 
1810d4002b98SHong Zhang   Level: beginner
1811d4002b98SHong Zhang 
1812d4002b98SHong Zhang .seealso: MatCreateSELL()
1813d4002b98SHong Zhang M*/
1814d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat B)
1815d4002b98SHong Zhang {
1816d4002b98SHong Zhang   Mat_MPISELL    *b;
1817d4002b98SHong Zhang   PetscMPIInt    size;
1818d4002b98SHong Zhang 
1819d4002b98SHong Zhang   PetscFunctionBegin;
18209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B),&size));
18219566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&b));
1822d4002b98SHong Zhang   B->data       = (void*)b;
18239566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps)));
1824d4002b98SHong Zhang   B->assembled  = PETSC_FALSE;
1825d4002b98SHong Zhang   B->insertmode = NOT_SET_VALUES;
1826d4002b98SHong Zhang   b->size       = size;
18279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank));
1828d4002b98SHong Zhang   /* build cache for off array entries formed */
18299566063dSJacob Faibussowitsch   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash));
1830d4002b98SHong Zhang 
1831d4002b98SHong Zhang   b->donotstash  = PETSC_FALSE;
1832f4259b30SLisandro Dalcin   b->colmap      = NULL;
1833f4259b30SLisandro Dalcin   b->garray      = NULL;
1834d4002b98SHong Zhang   b->roworiented = PETSC_TRUE;
1835d4002b98SHong Zhang 
1836d4002b98SHong Zhang   /* stuff used for matrix vector multiply */
1837d4002b98SHong Zhang   b->lvec  = NULL;
1838d4002b98SHong Zhang   b->Mvctx = NULL;
1839d4002b98SHong Zhang 
1840d4002b98SHong Zhang   /* stuff for MatGetRow() */
1841f4259b30SLisandro Dalcin   b->rowindices   = NULL;
1842f4259b30SLisandro Dalcin   b->rowvalues    = NULL;
1843d4002b98SHong Zhang   b->getrowactive = PETSC_FALSE;
1844d4002b98SHong Zhang 
18459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISELL));
18469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISELL));
18479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatIsTranspose_C",MatIsTranspose_MPISELL));
18489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPISELLSetPreallocation_C",MatMPISELLSetPreallocation_MPISELL));
18499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisell_mpiaij_C",MatConvert_MPISELL_MPIAIJ));
18509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDiagonalScaleLocal_C",MatDiagonalScaleLocal_MPISELL));
18519566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B,MATMPISELL));
1852d4002b98SHong Zhang   PetscFunctionReturn(0);
1853d4002b98SHong Zhang }
1854