xref: /petsc/src/mat/impls/is/matis.c (revision 285fb4e2b69b3de46a0633bd0adc6a7f684caa1e)
1 /*
2     Creates a matrix class for using the Neumann-Neumann type preconditioners.
3     This stores the matrices in globally unassembled form. Each processor
4     assembles only its local Neumann problem and the parallel matrix vector
5     product is handled "implicitly".
6 
7     Currently this allows for only one subdomain per processor.
8 */
9 
10 #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
11 #include <../src/mat/impls/aij/mpi/mpiaij.h>
12 #include <petsc/private/sfimpl.h>
13 
14 #define MATIS_MAX_ENTRIES_INSERTION 2048
15 static PetscErrorCode MatSetValuesLocal_IS(Mat,PetscInt,const PetscInt*,PetscInt,const PetscInt*,const PetscScalar*,InsertMode);
16 static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat,PetscInt,const PetscInt*,PetscInt,const PetscInt*,const PetscScalar*,InsertMode);
17 
18 static PetscErrorCode MatISContainerDestroyPtAP_Private(void *ptr)
19 {
20   MatISPtAP      ptap = (MatISPtAP)ptr;
21   PetscErrorCode ierr;
22 
23   PetscFunctionBegin;
24   ierr = MatDestroySubMatrices(ptap->ris1 ? 2 : 1,&ptap->lP);CHKERRQ(ierr);
25   ierr = ISDestroy(&ptap->cis0);CHKERRQ(ierr);
26   ierr = ISDestroy(&ptap->cis1);CHKERRQ(ierr);
27   ierr = ISDestroy(&ptap->ris0);CHKERRQ(ierr);
28   ierr = ISDestroy(&ptap->ris1);CHKERRQ(ierr);
29   ierr = PetscFree(ptap);CHKERRQ(ierr);
30   PetscFunctionReturn(0);
31 }
32 
33 static PetscErrorCode MatPtAPNumeric_IS_XAIJ(Mat A, Mat P, Mat C)
34 {
35   MatISPtAP      ptap;
36   Mat_IS         *matis = (Mat_IS*)A->data;
37   Mat            lA,lC;
38   MatReuse       reuse;
39   IS             ris[2],cis[2];
40   PetscContainer c;
41   PetscInt       n;
42   PetscErrorCode ierr;
43 
44   PetscFunctionBegin;
45   ierr = PetscObjectQuery((PetscObject)(C),"_MatIS_PtAP",(PetscObject*)&c);CHKERRQ(ierr);
46   if (!c) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing PtAP information");
47   ierr   = PetscContainerGetPointer(c,(void**)&ptap);CHKERRQ(ierr);
48   ris[0] = ptap->ris0;
49   ris[1] = ptap->ris1;
50   cis[0] = ptap->cis0;
51   cis[1] = ptap->cis1;
52   n      = ptap->ris1 ? 2 : 1;
53   reuse  = ptap->lP ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;
54   ierr   = MatCreateSubMatrices(P,n,ris,cis,reuse,&ptap->lP);CHKERRQ(ierr);
55 
56   ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr);
57   ierr = MatISGetLocalMat(C,&lC);CHKERRQ(ierr);
58   if (ptap->ris1) { /* unsymmetric A mapping */
59     Mat lPt;
60 
61     ierr = MatTranspose(ptap->lP[1],MAT_INITIAL_MATRIX,&lPt);CHKERRQ(ierr);
62     ierr = MatMatMatMult(lPt,lA,ptap->lP[0],reuse,ptap->fill,&lC);CHKERRQ(ierr);
63     if (matis->storel2l) {
64       ierr = PetscObjectCompose((PetscObject)(A),"_MatIS_PtAP_l2l",(PetscObject)lPt);CHKERRQ(ierr);
65     }
66     ierr = MatDestroy(&lPt);CHKERRQ(ierr);
67   } else {
68     ierr = MatPtAP(lA,ptap->lP[0],reuse,ptap->fill,&lC);CHKERRQ(ierr);
69     if (matis->storel2l) {
70      ierr = PetscObjectCompose((PetscObject)(C),"_MatIS_PtAP_l2l",(PetscObject)ptap->lP[0]);CHKERRQ(ierr);
71     }
72   }
73   if (reuse == MAT_INITIAL_MATRIX) {
74     ierr = MatISSetLocalMat(C,lC);CHKERRQ(ierr);
75     ierr = MatDestroy(&lC);CHKERRQ(ierr);
76   }
77   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79   PetscFunctionReturn(0);
80 }
81 
82 static PetscErrorCode MatGetNonzeroColumnsLocal_Private(Mat PT,IS *cis)
83 {
84   Mat            Po,Pd;
85   IS             zd,zo;
86   const PetscInt *garray;
87   PetscInt       *aux,i,bs;
88   PetscInt       dc,stc,oc,ctd,cto;
89   PetscBool      ismpiaij,ismpibaij,isseqaij,isseqbaij;
90   MPI_Comm       comm;
91   PetscErrorCode ierr;
92 
93   PetscFunctionBegin;
94   PetscValidHeaderSpecific(PT,MAT_CLASSID,1);
95   PetscValidPointer(cis,2);
96   ierr = PetscObjectGetComm((PetscObject)PT,&comm);CHKERRQ(ierr);
97   bs   = 1;
98   ierr = PetscObjectTypeCompare((PetscObject)PT,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
99   ierr = PetscObjectTypeCompare((PetscObject)PT,MATMPIBAIJ,&ismpibaij);CHKERRQ(ierr);
100   ierr = PetscObjectBaseTypeCompare((PetscObject)PT,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
101   ierr = PetscObjectTypeCompare((PetscObject)PT,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr);
102   if (isseqaij || isseqbaij) {
103     Pd = PT;
104     Po = NULL;
105     garray = NULL;
106   } else if (ismpiaij) {
107     ierr = MatMPIAIJGetSeqAIJ(PT,&Pd,&Po,&garray);CHKERRQ(ierr);
108   } else if (ismpibaij) {
109     ierr = MatMPIBAIJGetSeqBAIJ(PT,&Pd,&Po,&garray);CHKERRQ(ierr);
110     ierr = MatGetBlockSize(PT,&bs);CHKERRQ(ierr);
111   } else SETERRQ1(comm,PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)(PT))->type_name);
112 
113   /* identify any null columns in Pd or Po */
114   /* We use a tolerance comparison since it may happen that, with geometric multigrid,
115      some of the columns are not really zero, but very close to */
116   zo = zd = NULL;
117   if (Po) {
118     ierr = MatFindNonzeroRowsOrCols_Basic(Po,PETSC_TRUE,PETSC_SMALL,&zo);CHKERRQ(ierr);
119   }
120   ierr = MatFindNonzeroRowsOrCols_Basic(Pd,PETSC_TRUE,PETSC_SMALL,&zd);CHKERRQ(ierr);
121 
122   ierr = MatGetLocalSize(PT,NULL,&dc);CHKERRQ(ierr);
123   ierr = MatGetOwnershipRangeColumn(PT,&stc,NULL);CHKERRQ(ierr);
124   if (Po) { ierr = MatGetLocalSize(Po,NULL,&oc);CHKERRQ(ierr); }
125   else oc = 0;
126   ierr = PetscMalloc1((dc+oc)/bs,&aux);CHKERRQ(ierr);
127   if (zd) {
128     const PetscInt *idxs;
129     PetscInt       nz;
130 
131     /* this will throw an error if bs is not valid */
132     ierr = ISSetBlockSize(zd,bs);CHKERRQ(ierr);
133     ierr = ISGetLocalSize(zd,&nz);CHKERRQ(ierr);
134     ierr = ISGetIndices(zd,&idxs);CHKERRQ(ierr);
135     ctd  = nz/bs;
136     for (i=0; i<ctd; i++) aux[i] = (idxs[bs*i]+stc)/bs;
137     ierr = ISRestoreIndices(zd,&idxs);CHKERRQ(ierr);
138   } else {
139     ctd = dc/bs;
140     for (i=0; i<ctd; i++) aux[i] = i+stc/bs;
141   }
142   if (zo) {
143     const PetscInt *idxs;
144     PetscInt       nz;
145 
146     /* this will throw an error if bs is not valid */
147     ierr = ISSetBlockSize(zo,bs);CHKERRQ(ierr);
148     ierr = ISGetLocalSize(zo,&nz);CHKERRQ(ierr);
149     ierr = ISGetIndices(zo,&idxs);CHKERRQ(ierr);
150     cto  = nz/bs;
151     for (i=0; i<cto; i++) aux[i+ctd] = garray[idxs[bs*i]/bs];
152     ierr = ISRestoreIndices(zo,&idxs);CHKERRQ(ierr);
153   } else {
154     cto = oc/bs;
155     for (i=0; i<cto; i++) aux[i+ctd] = garray[i];
156   }
157   ierr = ISCreateBlock(comm,bs,ctd+cto,aux,PETSC_OWN_POINTER,cis);CHKERRQ(ierr);
158   ierr = ISDestroy(&zd);CHKERRQ(ierr);
159   ierr = ISDestroy(&zo);CHKERRQ(ierr);
160   PetscFunctionReturn(0);
161 }
162 
163 static PetscErrorCode MatPtAPSymbolic_IS_XAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
164 {
165   Mat                    PT;
166   MatISPtAP              ptap;
167   ISLocalToGlobalMapping Crl2g,Ccl2g,rl2g,cl2g;
168   PetscContainer         c;
169   const PetscInt         *garray;
170   PetscInt               ibs,N,dc;
171   MPI_Comm               comm;
172   PetscErrorCode         ierr;
173 
174   PetscFunctionBegin;
175   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
176   ierr = MatCreate(comm,C);CHKERRQ(ierr);
177   ierr = MatSetType(*C,MATIS);CHKERRQ(ierr);
178   ierr = MatGetSize(P,NULL,&N);CHKERRQ(ierr);
179   ierr = MatGetLocalSize(P,NULL,&dc);CHKERRQ(ierr);
180   ierr = MatSetSizes(*C,dc,dc,N,N);CHKERRQ(ierr);
181 /* Not sure about this
182   ierr = MatGetBlockSizes(P,NULL,&ibs);CHKERRQ(ierr);
183   ierr = MatSetBlockSize(*C,ibs);CHKERRQ(ierr);
184 */
185 
186   ierr = PetscNew(&ptap);CHKERRQ(ierr);
187   ierr = PetscContainerCreate(PETSC_COMM_SELF,&c);CHKERRQ(ierr);
188   ierr = PetscContainerSetPointer(c,ptap);CHKERRQ(ierr);
189   ierr = PetscContainerSetUserDestroy(c,MatISContainerDestroyPtAP_Private);CHKERRQ(ierr);
190   ierr = PetscObjectCompose((PetscObject)(*C),"_MatIS_PtAP",(PetscObject)c);CHKERRQ(ierr);
191   ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
192   ptap->fill = fill;
193 
194   ierr = MatGetLocalToGlobalMapping(A,&rl2g,&cl2g);CHKERRQ(ierr);
195 
196   ierr = ISLocalToGlobalMappingGetBlockSize(cl2g,&ibs);CHKERRQ(ierr);
197   ierr = ISLocalToGlobalMappingGetSize(cl2g,&N);CHKERRQ(ierr);
198   ierr = ISLocalToGlobalMappingGetBlockIndices(cl2g,&garray);CHKERRQ(ierr);
199   ierr = ISCreateBlock(comm,ibs,N/ibs,garray,PETSC_COPY_VALUES,&ptap->ris0);CHKERRQ(ierr);
200   ierr = ISLocalToGlobalMappingRestoreBlockIndices(cl2g,&garray);CHKERRQ(ierr);
201 
202   ierr = MatCreateSubMatrix(P,ptap->ris0,NULL,MAT_INITIAL_MATRIX,&PT);CHKERRQ(ierr);
203   ierr = MatGetNonzeroColumnsLocal_Private(PT,&ptap->cis0);CHKERRQ(ierr);
204   ierr = ISLocalToGlobalMappingCreateIS(ptap->cis0,&Ccl2g);CHKERRQ(ierr);
205   ierr = MatDestroy(&PT);CHKERRQ(ierr);
206 
207   Crl2g = NULL;
208   if (rl2g != cl2g) { /* unsymmetric A mapping */
209     PetscBool same,lsame = PETSC_FALSE;
210     PetscInt  N1,ibs1;
211 
212     ierr = ISLocalToGlobalMappingGetSize(rl2g,&N1);CHKERRQ(ierr);
213     ierr = ISLocalToGlobalMappingGetBlockSize(rl2g,&ibs1);CHKERRQ(ierr);
214     ierr = ISLocalToGlobalMappingGetBlockIndices(rl2g,&garray);CHKERRQ(ierr);
215     ierr = ISCreateBlock(comm,ibs,N/ibs,garray,PETSC_COPY_VALUES,&ptap->ris1);CHKERRQ(ierr);
216     ierr = ISLocalToGlobalMappingRestoreBlockIndices(rl2g,&garray);CHKERRQ(ierr);
217     if (ibs1 == ibs && N1 == N) { /* check if the l2gmaps are the same */
218       const PetscInt *i1,*i2;
219 
220       ierr = ISBlockGetIndices(ptap->ris0,&i1);CHKERRQ(ierr);
221       ierr = ISBlockGetIndices(ptap->ris1,&i2);CHKERRQ(ierr);
222       ierr = PetscMemcmp(i1,i2,N*sizeof(*i1),&lsame);CHKERRQ(ierr);
223     }
224     ierr = MPIU_Allreduce(&lsame,&same,1,MPIU_BOOL,MPI_LAND,comm);CHKERRQ(ierr);
225     if (same) {
226       ierr = ISDestroy(&ptap->ris1);CHKERRQ(ierr);
227     } else {
228       ierr = MatCreateSubMatrix(P,ptap->ris1,NULL,MAT_INITIAL_MATRIX,&PT);CHKERRQ(ierr);
229       ierr = MatGetNonzeroColumnsLocal_Private(PT,&ptap->cis1);CHKERRQ(ierr);
230       ierr = ISLocalToGlobalMappingCreateIS(ptap->cis1,&Crl2g);CHKERRQ(ierr);
231       ierr = MatDestroy(&PT);CHKERRQ(ierr);
232     }
233   }
234 /* Not sure about this
235   if (!Crl2g) {
236     ierr = MatGetBlockSize(*C,&ibs);CHKERRQ(ierr);
237     ierr = ISLocalToGlobalMappingSetBlockSize(Ccl2g,ibs);CHKERRQ(ierr);
238   }
239 */
240   ierr = MatSetLocalToGlobalMapping(*C,Crl2g ? Crl2g : Ccl2g,Ccl2g);CHKERRQ(ierr);
241   ierr = ISLocalToGlobalMappingDestroy(&Crl2g);CHKERRQ(ierr);
242   ierr = ISLocalToGlobalMappingDestroy(&Ccl2g);CHKERRQ(ierr);
243 
244   (*C)->ops->ptapnumeric = MatPtAPNumeric_IS_XAIJ;
245   PetscFunctionReturn(0);
246 }
247 
248 PETSC_INTERN PetscErrorCode MatPtAP_IS_XAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
249 {
250   PetscErrorCode ierr;
251 
252   PetscFunctionBegin;
253   if (scall == MAT_INITIAL_MATRIX) {
254     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
255     ierr = MatPtAPSymbolic_IS_XAIJ(A,P,fill,C);CHKERRQ(ierr);
256     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
257   }
258   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
259   ierr = ((*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
260   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
261   PetscFunctionReturn(0);
262 }
263 
264 static PetscErrorCode MatISContainerDestroyFields_Private(void *ptr)
265 {
266   MatISLocalFields lf = (MatISLocalFields)ptr;
267   PetscInt         i;
268   PetscErrorCode   ierr;
269 
270   PetscFunctionBegin;
271   for (i=0;i<lf->nr;i++) {
272     ierr = ISDestroy(&lf->rf[i]);CHKERRQ(ierr);
273   }
274   for (i=0;i<lf->nc;i++) {
275     ierr = ISDestroy(&lf->cf[i]);CHKERRQ(ierr);
276   }
277   ierr = PetscFree2(lf->rf,lf->cf);CHKERRQ(ierr);
278   ierr = PetscFree(lf);CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 
282 PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_IS(Mat A,MatType type,MatReuse reuse,Mat *newmat)
283 {
284   Mat_MPIAIJ             *aij  = (Mat_MPIAIJ*)A->data;
285   Mat_SeqAIJ             *diag = (Mat_SeqAIJ*)(aij->A->data);
286   Mat_SeqAIJ             *offd = (Mat_SeqAIJ*)(aij->B->data);
287   Mat                    lA;
288   ISLocalToGlobalMapping rl2g,cl2g;
289   IS                     is;
290   MPI_Comm               comm;
291   void                   *ptrs[2];
292   const char             *names[2] = {"_convert_csr_aux","_convert_csr_data"};
293   PetscScalar            *dd,*od,*aa,*data;
294   PetscInt               *di,*dj,*oi,*oj;
295   PetscInt               *aux,*ii,*jj;
296   PetscInt               lc,dr,dc,oc,str,stc,nnz,i,jd,jo,cum;
297   PetscErrorCode         ierr;
298 
299   PetscFunctionBegin;
300   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
301   if (!aij->garray) SETERRQ(comm,PETSC_ERR_SUP,"garray not present");
302 
303   /* access relevant information from MPIAIJ */
304   ierr = MatGetOwnershipRange(A,&str,NULL);CHKERRQ(ierr);
305   ierr = MatGetOwnershipRangeColumn(A,&stc,NULL);CHKERRQ(ierr);
306   ierr = MatGetLocalSize(A,&dr,&dc);CHKERRQ(ierr);
307   di   = diag->i;
308   dj   = diag->j;
309   dd   = diag->a;
310   oc   = aij->B->cmap->n;
311   oi   = offd->i;
312   oj   = offd->j;
313   od   = offd->a;
314   nnz  = diag->i[dr] + offd->i[dr];
315 
316   /* generate l2g maps for rows and cols */
317   ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
318   ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
319   ierr = ISDestroy(&is);CHKERRQ(ierr);
320   if (dr) {
321     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
322     for (i=0; i<dc; i++) aux[i]    = i+stc;
323     for (i=0; i<oc; i++) aux[i+dc] = aij->garray[i];
324     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
325     lc   = dc+oc;
326   } else {
327     ierr = ISCreateGeneral(comm,0,NULL,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
328     lc   = 0;
329   }
330   ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
331   ierr = ISDestroy(&is);CHKERRQ(ierr);
332 
333   /* create MATIS object */
334   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
335   ierr = MatSetSizes(*newmat,dr,dc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
336   ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr);
337   ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr);
338   ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
339   ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
340 
341   /* merge local matrices */
342   ierr = PetscMalloc1(nnz+dr+1,&aux);CHKERRQ(ierr);
343   ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
344   ii   = aux;
345   jj   = aux+dr+1;
346   aa   = data;
347   *ii  = *(di++) + *(oi++);
348   for (jd=0,jo=0,cum=0;*ii<nnz;cum++)
349   {
350      for (;jd<*di;jd++) { *jj++ = *dj++;      *aa++ = *dd++; }
351      for (;jo<*oi;jo++) { *jj++ = *oj++ + dc; *aa++ = *od++; }
352      *(++ii) = *(di++) + *(oi++);
353   }
354   for (;cum<dr;cum++) *(++ii) = nnz;
355   ii   = aux;
356   jj   = aux+dr+1;
357   aa   = data;
358   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,lc,ii,jj,aa,&lA);CHKERRQ(ierr);
359 
360   /* create containers to destroy the data */
361   ptrs[0] = aux;
362   ptrs[1] = data;
363   for (i=0; i<2; i++) {
364     PetscContainer c;
365 
366     ierr = PetscContainerCreate(PETSC_COMM_SELF,&c);CHKERRQ(ierr);
367     ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
368     ierr = PetscContainerSetUserDestroy(c,PetscContainerUserDestroyDefault);CHKERRQ(ierr);
369     ierr = PetscObjectCompose((PetscObject)lA,names[i],(PetscObject)c);CHKERRQ(ierr);
370     ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
371   }
372 
373   /* finalize matrix */
374   ierr = MatISSetLocalMat(*newmat,lA);CHKERRQ(ierr);
375   ierr = MatDestroy(&lA);CHKERRQ(ierr);
376   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
377   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
378   PetscFunctionReturn(0);
379 }
380 
381 /*@
382    MatISSetUpSF - Setup star forest objects used by MatIS.
383 
384    Collective on MPI_Comm
385 
386    Input Parameters:
387 +  A - the matrix
388 
389    Level: advanced
390 
391    Notes:
392     This function does not need to be called by the user.
393 
394 .keywords: matrix
395 
396 .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatISGetLocalMat()
397 @*/
398 PetscErrorCode  MatISSetUpSF(Mat A)
399 {
400   PetscErrorCode ierr;
401 
402   PetscFunctionBegin;
403   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
404   PetscValidType(A,1);
405   ierr = PetscTryMethod(A,"MatISSetUpSF_C",(Mat),(A));CHKERRQ(ierr);
406   PetscFunctionReturn(0);
407 }
408 
409 PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat A,MatType type,MatReuse reuse,Mat *newmat)
410 {
411   Mat                    **nest,*snest,**rnest,lA,B;
412   IS                     *iscol,*isrow,*islrow,*islcol;
413   ISLocalToGlobalMapping rl2g,cl2g;
414   MPI_Comm               comm;
415   PetscInt               *lr,*lc,*l2gidxs;
416   PetscInt               i,j,nr,nc,rbs,cbs;
417   PetscBool              convert,lreuse,*istrans;
418   PetscErrorCode         ierr;
419 
420   PetscFunctionBegin;
421   ierr   = MatNestGetSubMats(A,&nr,&nc,&nest);CHKERRQ(ierr);
422   lreuse = PETSC_FALSE;
423   rnest  = NULL;
424   if (reuse == MAT_REUSE_MATRIX) {
425     PetscBool ismatis,isnest;
426 
427     ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr);
428     if (!ismatis) SETERRQ1(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_USER,"Cannot reuse matrix of type %s",((PetscObject)(*newmat))->type);
429     ierr = MatISGetLocalMat(*newmat,&lA);CHKERRQ(ierr);
430     ierr = PetscObjectTypeCompare((PetscObject)lA,MATNEST,&isnest);CHKERRQ(ierr);
431     if (isnest) {
432       ierr   = MatNestGetSubMats(lA,&i,&j,&rnest);CHKERRQ(ierr);
433       lreuse = (PetscBool)(i == nr && j == nc);
434       if (!lreuse) rnest = NULL;
435     }
436   }
437   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
438   ierr = PetscCalloc2(nr,&lr,nc,&lc);CHKERRQ(ierr);
439   ierr = PetscCalloc6(nr,&isrow,nc,&iscol,
440                       nr,&islrow,nc,&islcol,
441                       nr*nc,&snest,nr*nc,&istrans);CHKERRQ(ierr);
442   ierr = MatNestGetISs(A,isrow,iscol);CHKERRQ(ierr);
443   for (i=0;i<nr;i++) {
444     for (j=0;j<nc;j++) {
445       PetscBool ismatis;
446       PetscInt  l1,l2,lb1,lb2,ij=i*nc+j;
447 
448       /* Null matrix pointers are allowed in MATNEST */
449       if (!nest[i][j]) continue;
450 
451       /* Nested matrices should be of type MATIS */
452       ierr = PetscObjectTypeCompare((PetscObject)nest[i][j],MATTRANSPOSEMAT,&istrans[ij]);CHKERRQ(ierr);
453       if (istrans[ij]) {
454         Mat T,lT;
455         ierr = MatTransposeGetMat(nest[i][j],&T);CHKERRQ(ierr);
456         ierr = PetscObjectTypeCompare((PetscObject)T,MATIS,&ismatis);CHKERRQ(ierr);
457         if (!ismatis) SETERRQ2(comm,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%D,%D) (transposed) is not of type MATIS",i,j);
458         ierr = MatISGetLocalMat(T,&lT);CHKERRQ(ierr);
459         ierr = MatCreateTranspose(lT,&snest[ij]);CHKERRQ(ierr);
460       } else {
461         ierr = PetscObjectTypeCompare((PetscObject)nest[i][j],MATIS,&ismatis);CHKERRQ(ierr);
462         if (!ismatis) SETERRQ2(comm,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%D,%D) is not of type MATIS",i,j);
463         ierr = MatISGetLocalMat(nest[i][j],&snest[ij]);CHKERRQ(ierr);
464       }
465 
466       /* Check compatibility of local sizes */
467       ierr = MatGetSize(snest[ij],&l1,&l2);CHKERRQ(ierr);
468       ierr = MatGetBlockSizes(snest[ij],&lb1,&lb2);CHKERRQ(ierr);
469       if (!l1 || !l2) continue;
470       if (lr[i] && l1 != lr[i]) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%D,%D) has invalid local size %D != %D",i,j,lr[i],l1);
471       if (lc[j] && l2 != lc[j]) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%D,%D) has invalid local size %D != %D",i,j,lc[j],l2);
472       lr[i] = l1;
473       lc[j] = l2;
474 
475       /* check compatibilty for local matrix reusage */
476       if (rnest && !rnest[i][j] != !snest[ij]) lreuse = PETSC_FALSE;
477     }
478   }
479 
480 #if defined (PETSC_USE_DEBUG)
481   /* Check compatibility of l2g maps for rows */
482   for (i=0;i<nr;i++) {
483     rl2g = NULL;
484     for (j=0;j<nc;j++) {
485       PetscInt n1,n2;
486 
487       if (!nest[i][j]) continue;
488       if (istrans[i*nc+j]) {
489         Mat T;
490 
491         ierr = MatTransposeGetMat(nest[i][j],&T);CHKERRQ(ierr);
492         ierr = MatGetLocalToGlobalMapping(T,NULL,&cl2g);CHKERRQ(ierr);
493       } else {
494         ierr = MatGetLocalToGlobalMapping(nest[i][j],&cl2g,NULL);CHKERRQ(ierr);
495       }
496       ierr = ISLocalToGlobalMappingGetSize(cl2g,&n1);CHKERRQ(ierr);
497       if (!n1) continue;
498       if (!rl2g) {
499         rl2g = cl2g;
500       } else {
501         const PetscInt *idxs1,*idxs2;
502         PetscBool      same;
503 
504         ierr = ISLocalToGlobalMappingGetSize(rl2g,&n2);CHKERRQ(ierr);
505         if (n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%D,%D) has invalid row l2gmap size %D != %D",i,j,n1,n2);
506         ierr = ISLocalToGlobalMappingGetIndices(cl2g,&idxs1);CHKERRQ(ierr);
507         ierr = ISLocalToGlobalMappingGetIndices(rl2g,&idxs2);CHKERRQ(ierr);
508         ierr = PetscMemcmp(idxs1,idxs2,n1*sizeof(PetscInt),&same);CHKERRQ(ierr);
509         ierr = ISLocalToGlobalMappingRestoreIndices(cl2g,&idxs1);CHKERRQ(ierr);
510         ierr = ISLocalToGlobalMappingRestoreIndices(rl2g,&idxs2);CHKERRQ(ierr);
511         if (!same) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%D,%D) has invalid row l2gmap",i,j);
512       }
513     }
514   }
515   /* Check compatibility of l2g maps for columns */
516   for (i=0;i<nc;i++) {
517     rl2g = NULL;
518     for (j=0;j<nr;j++) {
519       PetscInt n1,n2;
520 
521       if (!nest[j][i]) continue;
522       if (istrans[j*nc+i]) {
523         Mat T;
524 
525         ierr = MatTransposeGetMat(nest[j][i],&T);CHKERRQ(ierr);
526         ierr = MatGetLocalToGlobalMapping(T,&cl2g,NULL);CHKERRQ(ierr);
527       } else {
528         ierr = MatGetLocalToGlobalMapping(nest[j][i],NULL,&cl2g);CHKERRQ(ierr);
529       }
530       ierr = ISLocalToGlobalMappingGetSize(cl2g,&n1);CHKERRQ(ierr);
531       if (!n1) continue;
532       if (!rl2g) {
533         rl2g = cl2g;
534       } else {
535         const PetscInt *idxs1,*idxs2;
536         PetscBool      same;
537 
538         ierr = ISLocalToGlobalMappingGetSize(rl2g,&n2);CHKERRQ(ierr);
539         if (n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%D,%D) has invalid column l2gmap size %D != %D",j,i,n1,n2);
540         ierr = ISLocalToGlobalMappingGetIndices(cl2g,&idxs1);CHKERRQ(ierr);
541         ierr = ISLocalToGlobalMappingGetIndices(rl2g,&idxs2);CHKERRQ(ierr);
542         ierr = PetscMemcmp(idxs1,idxs2,n1*sizeof(PetscInt),&same);CHKERRQ(ierr);
543         ierr = ISLocalToGlobalMappingRestoreIndices(cl2g,&idxs1);CHKERRQ(ierr);
544         ierr = ISLocalToGlobalMappingRestoreIndices(rl2g,&idxs2);CHKERRQ(ierr);
545         if (!same) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot convert from MATNEST to MATIS! Matrix block (%D,%D) has invalid column l2gmap",j,i);
546       }
547     }
548   }
549 #endif
550 
551   B = NULL;
552   if (reuse != MAT_REUSE_MATRIX) {
553     PetscInt stl;
554 
555     /* Create l2g map for the rows of the new matrix and index sets for the local MATNEST */
556     for (i=0,stl=0;i<nr;i++) stl += lr[i];
557     ierr = PetscMalloc1(stl,&l2gidxs);CHKERRQ(ierr);
558     for (i=0,stl=0;i<nr;i++) {
559       Mat            usedmat;
560       Mat_IS         *matis;
561       const PetscInt *idxs;
562 
563       /* local IS for local NEST */
564       ierr  = ISCreateStride(PETSC_COMM_SELF,lr[i],stl,1,&islrow[i]);CHKERRQ(ierr);
565 
566       /* l2gmap */
567       j = 0;
568       usedmat = nest[i][j];
569       while (!usedmat && j < nc-1) usedmat = nest[i][++j];
570       if (!usedmat) SETERRQ(comm,PETSC_ERR_SUP,"Cannot find valid row mat");
571 
572       if (istrans[i*nc+j]) {
573         Mat T;
574         ierr    = MatTransposeGetMat(usedmat,&T);CHKERRQ(ierr);
575         usedmat = T;
576       }
577       ierr  = MatISSetUpSF(usedmat);CHKERRQ(ierr);
578       matis = (Mat_IS*)(usedmat->data);
579       ierr  = ISGetIndices(isrow[i],&idxs);CHKERRQ(ierr);
580       if (istrans[i*nc+j]) {
581         ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,idxs,l2gidxs+stl);CHKERRQ(ierr);
582         ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,idxs,l2gidxs+stl);CHKERRQ(ierr);
583       } else {
584         ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,idxs,l2gidxs+stl);CHKERRQ(ierr);
585         ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,idxs,l2gidxs+stl);CHKERRQ(ierr);
586       }
587       ierr = ISRestoreIndices(isrow[i],&idxs);CHKERRQ(ierr);
588       stl += lr[i];
589     }
590     ierr = ISLocalToGlobalMappingCreate(comm,1,stl,l2gidxs,PETSC_OWN_POINTER,&rl2g);CHKERRQ(ierr);
591 
592     /* Create l2g map for columns of the new matrix and index sets for the local MATNEST */
593     for (i=0,stl=0;i<nc;i++) stl += lc[i];
594     ierr = PetscMalloc1(stl,&l2gidxs);CHKERRQ(ierr);
595     for (i=0,stl=0;i<nc;i++) {
596       Mat            usedmat;
597       Mat_IS         *matis;
598       const PetscInt *idxs;
599 
600       /* local IS for local NEST */
601       ierr  = ISCreateStride(PETSC_COMM_SELF,lc[i],stl,1,&islcol[i]);CHKERRQ(ierr);
602 
603       /* l2gmap */
604       j = 0;
605       usedmat = nest[j][i];
606       while (!usedmat && j < nr-1) usedmat = nest[++j][i];
607       if (!usedmat) SETERRQ(comm,PETSC_ERR_SUP,"Cannot find valid column mat");
608       if (istrans[j*nc+i]) {
609         Mat T;
610         ierr    = MatTransposeGetMat(usedmat,&T);CHKERRQ(ierr);
611         usedmat = T;
612       }
613       ierr  = MatISSetUpSF(usedmat);CHKERRQ(ierr);
614       matis = (Mat_IS*)(usedmat->data);
615       ierr  = ISGetIndices(iscol[i],&idxs);CHKERRQ(ierr);
616       if (istrans[j*nc+i]) {
617         ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,idxs,l2gidxs+stl);CHKERRQ(ierr);
618         ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,idxs,l2gidxs+stl);CHKERRQ(ierr);
619       } else {
620         ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,idxs,l2gidxs+stl);CHKERRQ(ierr);
621         ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,idxs,l2gidxs+stl);CHKERRQ(ierr);
622       }
623       ierr = ISRestoreIndices(iscol[i],&idxs);CHKERRQ(ierr);
624       stl += lc[i];
625     }
626     ierr = ISLocalToGlobalMappingCreate(comm,1,stl,l2gidxs,PETSC_OWN_POINTER,&cl2g);CHKERRQ(ierr);
627 
628     /* Create MATIS */
629     ierr = MatCreate(comm,&B);CHKERRQ(ierr);
630     ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
631     ierr = MatGetBlockSizes(A,&rbs,&cbs);CHKERRQ(ierr);
632     ierr = MatSetBlockSizes(B,rbs,cbs);CHKERRQ(ierr);
633     ierr = MatSetType(B,MATIS);CHKERRQ(ierr);
634     ierr = MatSetLocalToGlobalMapping(B,rl2g,cl2g);CHKERRQ(ierr);
635     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
636     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
637     ierr = MatCreateNest(PETSC_COMM_SELF,nr,islrow,nc,islcol,snest,&lA);CHKERRQ(ierr);
638     for (i=0;i<nr*nc;i++) {
639       if (istrans[i]) {
640         ierr = MatDestroy(&snest[i]);CHKERRQ(ierr);
641       }
642     }
643     ierr = MatISSetLocalMat(B,lA);CHKERRQ(ierr);
644     ierr = MatDestroy(&lA);CHKERRQ(ierr);
645     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
646     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
647     if (reuse == MAT_INPLACE_MATRIX) {
648       ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
649     } else {
650       *newmat = B;
651     }
652   } else {
653     if (lreuse) {
654       ierr = MatISGetLocalMat(*newmat,&lA);CHKERRQ(ierr);
655       for (i=0;i<nr;i++) {
656         for (j=0;j<nc;j++) {
657           if (snest[i*nc+j]) {
658             ierr = MatNestSetSubMat(lA,i,j,snest[i*nc+j]);CHKERRQ(ierr);
659             if (istrans[i*nc+j]) {
660               ierr = MatDestroy(&snest[i*nc+j]);CHKERRQ(ierr);
661             }
662           }
663         }
664       }
665     } else {
666       PetscInt stl;
667       for (i=0,stl=0;i<nr;i++) {
668         ierr  = ISCreateStride(PETSC_COMM_SELF,lr[i],stl,1,&islrow[i]);CHKERRQ(ierr);
669         stl  += lr[i];
670       }
671       for (i=0,stl=0;i<nc;i++) {
672         ierr  = ISCreateStride(PETSC_COMM_SELF,lc[i],stl,1,&islcol[i]);CHKERRQ(ierr);
673         stl  += lc[i];
674       }
675       ierr = MatCreateNest(PETSC_COMM_SELF,nr,islrow,nc,islcol,snest,&lA);CHKERRQ(ierr);
676       for (i=0;i<nr*nc;i++) {
677         if (istrans[i]) {
678           ierr = MatDestroy(&snest[i]);CHKERRQ(ierr);
679         }
680       }
681       ierr = MatISSetLocalMat(*newmat,lA);CHKERRQ(ierr);
682       ierr = MatDestroy(&lA);CHKERRQ(ierr);
683     }
684     ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
685     ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
686   }
687 
688   /* Create local matrix in MATNEST format */
689   convert = PETSC_FALSE;
690   ierr = PetscOptionsGetBool(NULL,((PetscObject)A)->prefix,"-matis_convert_local_nest",&convert,NULL);CHKERRQ(ierr);
691   if (convert) {
692     Mat              M;
693     MatISLocalFields lf;
694     PetscContainer   c;
695 
696     ierr = MatISGetLocalMat(*newmat,&lA);CHKERRQ(ierr);
697     ierr = MatConvert(lA,MATAIJ,MAT_INITIAL_MATRIX,&M);CHKERRQ(ierr);
698     ierr = MatISSetLocalMat(*newmat,M);CHKERRQ(ierr);
699     ierr = MatDestroy(&M);CHKERRQ(ierr);
700 
701     /* attach local fields to the matrix */
702     ierr = PetscNew(&lf);CHKERRQ(ierr);
703     ierr = PetscCalloc2(nr,&lf->rf,nc,&lf->cf);CHKERRQ(ierr);
704     for (i=0;i<nr;i++) {
705       PetscInt n,st;
706 
707       ierr = ISGetLocalSize(islrow[i],&n);CHKERRQ(ierr);
708       ierr = ISStrideGetInfo(islrow[i],&st,NULL);CHKERRQ(ierr);
709       ierr = ISCreateStride(comm,n,st,1,&lf->rf[i]);CHKERRQ(ierr);
710     }
711     for (i=0;i<nc;i++) {
712       PetscInt n,st;
713 
714       ierr = ISGetLocalSize(islcol[i],&n);CHKERRQ(ierr);
715       ierr = ISStrideGetInfo(islcol[i],&st,NULL);CHKERRQ(ierr);
716       ierr = ISCreateStride(comm,n,st,1,&lf->cf[i]);CHKERRQ(ierr);
717     }
718     lf->nr = nr;
719     lf->nc = nc;
720     ierr = PetscContainerCreate(PetscObjectComm((PetscObject)(*newmat)),&c);CHKERRQ(ierr);
721     ierr = PetscContainerSetPointer(c,lf);CHKERRQ(ierr);
722     ierr = PetscContainerSetUserDestroy(c,MatISContainerDestroyFields_Private);CHKERRQ(ierr);
723     ierr = PetscObjectCompose((PetscObject)(*newmat),"_convert_nest_lfields",(PetscObject)c);CHKERRQ(ierr);
724     ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
725   }
726 
727   /* Free workspace */
728   for (i=0;i<nr;i++) {
729     ierr = ISDestroy(&islrow[i]);CHKERRQ(ierr);
730   }
731   for (i=0;i<nc;i++) {
732     ierr = ISDestroy(&islcol[i]);CHKERRQ(ierr);
733   }
734   ierr = PetscFree6(isrow,iscol,islrow,islcol,snest,istrans);CHKERRQ(ierr);
735   ierr = PetscFree2(lr,lc);CHKERRQ(ierr);
736   PetscFunctionReturn(0);
737 }
738 
739 static PetscErrorCode MatDiagonalScale_IS(Mat A, Vec l, Vec r)
740 {
741   Mat_IS            *matis = (Mat_IS*)A->data;
742   Vec               ll,rr;
743   const PetscScalar *Y,*X;
744   PetscScalar       *x,*y;
745   PetscErrorCode    ierr;
746 
747   PetscFunctionBegin;
748   ierr = MatISSetUpSF(A);CHKERRQ(ierr);
749   if (l) {
750     ll   = matis->y;
751     ierr = VecGetArrayRead(l,&Y);CHKERRQ(ierr);
752     ierr = VecGetArray(ll,&y);CHKERRQ(ierr);
753     ierr = PetscSFBcastBegin(matis->sf,MPIU_SCALAR,Y,y);CHKERRQ(ierr);
754   } else {
755     ll = NULL;
756   }
757   if (r) {
758     rr   = matis->x;
759     ierr = VecGetArrayRead(r,&X);CHKERRQ(ierr);
760     ierr = VecGetArray(rr,&x);CHKERRQ(ierr);
761     ierr = PetscSFBcastBegin(matis->csf,MPIU_SCALAR,X,x);CHKERRQ(ierr);
762   } else {
763     rr = NULL;
764   }
765   if (ll) {
766     ierr = PetscSFBcastEnd(matis->sf,MPIU_SCALAR,Y,y);CHKERRQ(ierr);
767     ierr = VecRestoreArrayRead(l,&Y);CHKERRQ(ierr);
768     ierr = VecRestoreArray(ll,&y);CHKERRQ(ierr);
769   }
770   if (rr) {
771     ierr = PetscSFBcastEnd(matis->csf,MPIU_SCALAR,X,x);CHKERRQ(ierr);
772     ierr = VecRestoreArrayRead(r,&X);CHKERRQ(ierr);
773     ierr = VecRestoreArray(rr,&x);CHKERRQ(ierr);
774   }
775   ierr = MatDiagonalScale(matis->A,ll,rr);CHKERRQ(ierr);
776   PetscFunctionReturn(0);
777 }
778 
779 static PetscErrorCode MatGetInfo_IS(Mat A,MatInfoType flag,MatInfo *ginfo)
780 {
781   Mat_IS         *matis = (Mat_IS*)A->data;
782   MatInfo        info;
783   PetscReal      isend[6],irecv[6];
784   PetscInt       bs;
785   PetscErrorCode ierr;
786 
787   PetscFunctionBegin;
788   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
789   if (matis->A->ops->getinfo) {
790     ierr     = MatGetInfo(matis->A,MAT_LOCAL,&info);CHKERRQ(ierr);
791     isend[0] = info.nz_used;
792     isend[1] = info.nz_allocated;
793     isend[2] = info.nz_unneeded;
794     isend[3] = info.memory;
795     isend[4] = info.mallocs;
796   } else {
797     isend[0] = 0.;
798     isend[1] = 0.;
799     isend[2] = 0.;
800     isend[3] = 0.;
801     isend[4] = 0.;
802   }
803   isend[5] = matis->A->num_ass;
804   if (flag == MAT_LOCAL) {
805     ginfo->nz_used      = isend[0];
806     ginfo->nz_allocated = isend[1];
807     ginfo->nz_unneeded  = isend[2];
808     ginfo->memory       = isend[3];
809     ginfo->mallocs      = isend[4];
810     ginfo->assemblies   = isend[5];
811   } else if (flag == MAT_GLOBAL_MAX) {
812     ierr = MPIU_Allreduce(isend,irecv,6,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
813 
814     ginfo->nz_used      = irecv[0];
815     ginfo->nz_allocated = irecv[1];
816     ginfo->nz_unneeded  = irecv[2];
817     ginfo->memory       = irecv[3];
818     ginfo->mallocs      = irecv[4];
819     ginfo->assemblies   = irecv[5];
820   } else if (flag == MAT_GLOBAL_SUM) {
821     ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
822 
823     ginfo->nz_used      = irecv[0];
824     ginfo->nz_allocated = irecv[1];
825     ginfo->nz_unneeded  = irecv[2];
826     ginfo->memory       = irecv[3];
827     ginfo->mallocs      = irecv[4];
828     ginfo->assemblies   = A->num_ass;
829   }
830   ginfo->block_size        = bs;
831   ginfo->fill_ratio_given  = 0;
832   ginfo->fill_ratio_needed = 0;
833   ginfo->factor_mallocs    = 0;
834   PetscFunctionReturn(0);
835 }
836 
837 PetscErrorCode MatTranspose_IS(Mat A,MatReuse reuse,Mat *B)
838 {
839   Mat                    C,lC,lA;
840   PetscErrorCode         ierr;
841 
842   PetscFunctionBegin;
843   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
844     ISLocalToGlobalMapping rl2g,cl2g;
845     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
846     ierr = MatSetSizes(C,A->cmap->n,A->rmap->n,A->cmap->N,A->rmap->N);CHKERRQ(ierr);
847     ierr = MatSetBlockSizes(C,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs));CHKERRQ(ierr);
848     ierr = MatSetType(C,MATIS);CHKERRQ(ierr);
849     ierr = MatGetLocalToGlobalMapping(A,&rl2g,&cl2g);CHKERRQ(ierr);
850     ierr = MatSetLocalToGlobalMapping(C,cl2g,rl2g);CHKERRQ(ierr);
851   } else {
852     C = *B;
853   }
854 
855   /* perform local transposition */
856   ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr);
857   ierr = MatTranspose(lA,MAT_INITIAL_MATRIX,&lC);CHKERRQ(ierr);
858   ierr = MatISSetLocalMat(C,lC);CHKERRQ(ierr);
859   ierr = MatDestroy(&lC);CHKERRQ(ierr);
860 
861   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
862     *B = C;
863   } else {
864     ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr);
865   }
866   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
867   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
868   PetscFunctionReturn(0);
869 }
870 
871 PetscErrorCode  MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode)
872 {
873   Mat_IS         *is = (Mat_IS*)A->data;
874   PetscErrorCode ierr;
875 
876   PetscFunctionBegin;
877   if (D) { /* MatShift_IS pass D = NULL */
878     ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
879     ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
880   }
881   ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr);
882   ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr);
883   PetscFunctionReturn(0);
884 }
885 
886 PetscErrorCode  MatShift_IS(Mat A,PetscScalar a)
887 {
888   Mat_IS         *is = (Mat_IS*)A->data;
889   PetscErrorCode ierr;
890 
891   PetscFunctionBegin;
892   ierr = VecSet(is->y,a);CHKERRQ(ierr);
893   ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr);
894   PetscFunctionReturn(0);
895 }
896 
897 static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
898 {
899   PetscErrorCode ierr;
900   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
901 
902   PetscFunctionBegin;
903 #if defined(PETSC_USE_DEBUG)
904   if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n);
905 #endif
906   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
907   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
908   ierr = MatSetValuesLocal_IS(A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
909   PetscFunctionReturn(0);
910 }
911 
912 static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
913 {
914   PetscErrorCode ierr;
915   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
916 
917   PetscFunctionBegin;
918 #if defined(PETSC_USE_DEBUG)
919   if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n);
920 #endif
921   ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
922   ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
923   ierr = MatSetValuesBlockedLocal_IS(A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
924   PetscFunctionReturn(0);
925 }
926 
927 static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs)
928 {
929   PetscInt      *owners = map->range;
930   PetscInt       n      = map->n;
931   PetscSF        sf;
932   PetscInt      *lidxs,*work = NULL;
933   PetscSFNode   *ridxs;
934   PetscMPIInt    rank;
935   PetscInt       r, p = 0, len = 0;
936   PetscErrorCode ierr;
937 
938   PetscFunctionBegin;
939   if (on) *on = 0;              /* squelch -Wmaybe-uninitialized */
940   /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */
941   ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr);
942   ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr);
943   for (r = 0; r < n; ++r) lidxs[r] = -1;
944   ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr);
945   for (r = 0; r < N; ++r) {
946     const PetscInt idx = idxs[r];
947     if (idx < 0 || map->N <= idx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D out of range [0,%D)",idx,map->N);
948     if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */
949       ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr);
950     }
951     ridxs[r].rank = p;
952     ridxs[r].index = idxs[r] - owners[p];
953   }
954   ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr);
955   ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr);
956   ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
957   ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
958   if (ogidxs) { /* communicate global idxs */
959     PetscInt cum = 0,start,*work2;
960 
961     ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
962     ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr);
963     for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++;
964     ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr);
965     start -= cum;
966     cum = 0;
967     for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++;
968     ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
969     ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
970     ierr = PetscFree(work2);CHKERRQ(ierr);
971   }
972   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
973   /* Compress and put in indices */
974   for (r = 0; r < n; ++r)
975     if (lidxs[r] >= 0) {
976       if (work) work[len] = work[r];
977       lidxs[len++] = r;
978     }
979   if (on) *on = len;
980   if (oidxs) *oidxs = lidxs;
981   if (ogidxs) *ogidxs = work;
982   PetscFunctionReturn(0);
983 }
984 
985 static PetscErrorCode MatCreateSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat)
986 {
987   Mat               locmat,newlocmat;
988   Mat_IS            *newmatis;
989 #if defined(PETSC_USE_DEBUG)
990   Vec               rtest,ltest;
991   const PetscScalar *array;
992 #endif
993   const PetscInt    *idxs;
994   PetscInt          i,m,n;
995   PetscErrorCode    ierr;
996 
997   PetscFunctionBegin;
998   if (scall == MAT_REUSE_MATRIX) {
999     PetscBool ismatis;
1000 
1001     ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr);
1002     if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type");
1003     newmatis = (Mat_IS*)(*newmat)->data;
1004     if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS");
1005     if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS");
1006   }
1007   /* irow and icol may not have duplicate entries */
1008 #if defined(PETSC_USE_DEBUG)
1009   ierr = MatCreateVecs(mat,&ltest,&rtest);CHKERRQ(ierr);
1010   ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr);
1011   ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
1012   for (i=0;i<n;i++) {
1013     ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
1014   }
1015   ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr);
1016   ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr);
1017   ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr);
1018   ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr);
1019   ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr);
1020   for (i=0;i<n;i++) if (array[i] != 0. && array[i] != 1.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Index %D counted %D times! Irow may not have duplicate entries",i+m,(PetscInt)PetscRealPart(array[i]));
1021   ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr);
1022   ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
1023   ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
1024   ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
1025   for (i=0;i<n;i++) {
1026     ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
1027   }
1028   ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr);
1029   ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr);
1030   ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr);
1031   ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr);
1032   ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr);
1033   for (i=0;i<n;i++) if (array[i] != 0. && array[i] != 1.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Index %D counted %D times! Icol may not have duplicate entries",i+m,(PetscInt)PetscRealPart(array[i]));
1034   ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr);
1035   ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
1036   ierr = VecDestroy(&rtest);CHKERRQ(ierr);
1037   ierr = VecDestroy(&ltest);CHKERRQ(ierr);
1038 #endif
1039   if (scall == MAT_INITIAL_MATRIX) {
1040     Mat_IS                 *matis = (Mat_IS*)mat->data;
1041     ISLocalToGlobalMapping rl2g;
1042     IS                     is;
1043     PetscInt               *lidxs,*lgidxs,*newgidxs;
1044     PetscInt               ll,newloc,irbs,icbs,arbs,acbs,rbs,cbs;
1045     MPI_Comm               comm;
1046 
1047     ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
1048     ierr = MatGetBlockSizes(mat,&arbs,&acbs);CHKERRQ(ierr);
1049     ierr = ISGetBlockSize(irow,&irbs);CHKERRQ(ierr);
1050     ierr = ISGetBlockSize(icol,&icbs);CHKERRQ(ierr);
1051     rbs  = arbs == irbs ? irbs : 1;
1052     cbs  = acbs == icbs ? icbs : 1;
1053     ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr);
1054     ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
1055     ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1056     ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr);
1057     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1058     ierr = MatSetBlockSizes(*newmat,rbs,cbs);CHKERRQ(ierr);
1059     /* communicate irow to their owners in the layout */
1060     ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
1061     ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
1062     ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
1063     ierr = MatISSetUpSF(mat);CHKERRQ(ierr);
1064     ierr = PetscMemzero(matis->sf_rootdata,matis->sf->nroots*sizeof(PetscInt));CHKERRQ(ierr);
1065     for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1;
1066     ierr = PetscFree(lidxs);CHKERRQ(ierr);
1067     ierr = PetscFree(lgidxs);CHKERRQ(ierr);
1068     ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
1069     ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
1070     for (i=0,newloc=0;i<matis->sf->nleaves;i++) if (matis->sf_leafdata[i]) newloc++;
1071     ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
1072     ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
1073     for (i=0,newloc=0;i<matis->sf->nleaves;i++)
1074       if (matis->sf_leafdata[i]) {
1075         lidxs[newloc] = i;
1076         newgidxs[newloc++] = matis->sf_leafdata[i]-1;
1077       }
1078     ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1079     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
1080     ierr = ISLocalToGlobalMappingSetBlockSize(rl2g,rbs);CHKERRQ(ierr);
1081     ierr = ISDestroy(&is);CHKERRQ(ierr);
1082     /* local is to extract local submatrix */
1083     newmatis = (Mat_IS*)(*newmat)->data;
1084     ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr);
1085     if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */
1086       PetscBool cong;
1087 
1088       ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr);
1089       if (cong) mat->congruentlayouts = 1;
1090       else      mat->congruentlayouts = 0;
1091     }
1092     if (mat->congruentlayouts && irow == icol && matis->csf == matis->sf) {
1093       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr);
1094       ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr);
1095       newmatis->getsub_cis = newmatis->getsub_ris;
1096     } else {
1097       ISLocalToGlobalMapping cl2g;
1098 
1099       /* communicate icol to their owners in the layout */
1100       ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
1101       ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
1102       ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
1103       ierr = PetscMemzero(matis->csf_rootdata,matis->csf->nroots*sizeof(PetscInt));CHKERRQ(ierr);
1104       for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1;
1105       ierr = PetscFree(lidxs);CHKERRQ(ierr);
1106       ierr = PetscFree(lgidxs);CHKERRQ(ierr);
1107       ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
1108       ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
1109       for (i=0,newloc=0;i<matis->csf->nleaves;i++) if (matis->csf_leafdata[i]) newloc++;
1110       ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
1111       ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
1112       for (i=0,newloc=0;i<matis->csf->nleaves;i++)
1113         if (matis->csf_leafdata[i]) {
1114           lidxs[newloc] = i;
1115           newgidxs[newloc++] = matis->csf_leafdata[i]-1;
1116         }
1117       ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1118       ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
1119       ierr = ISLocalToGlobalMappingSetBlockSize(cl2g,cbs);CHKERRQ(ierr);
1120       ierr = ISDestroy(&is);CHKERRQ(ierr);
1121       /* local is to extract local submatrix */
1122       ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr);
1123       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr);
1124       ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
1125     }
1126     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
1127   } else {
1128     ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr);
1129   }
1130   ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr);
1131   newmatis = (Mat_IS*)(*newmat)->data;
1132   ierr = MatCreateSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr);
1133   if (scall == MAT_INITIAL_MATRIX) {
1134     ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr);
1135     ierr = MatDestroy(&newlocmat);CHKERRQ(ierr);
1136   }
1137   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1138   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str)
1143 {
1144   Mat_IS         *a = (Mat_IS*)A->data,*b;
1145   PetscBool      ismatis;
1146   PetscErrorCode ierr;
1147 
1148   PetscFunctionBegin;
1149   ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr);
1150   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented");
1151   b = (Mat_IS*)B->data;
1152   ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr);
1153   ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr);
1154   PetscFunctionReturn(0);
1155 }
1156 
1157 static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool  *missing,PetscInt *d)
1158 {
1159   Vec               v;
1160   const PetscScalar *array;
1161   PetscInt          i,n;
1162   PetscErrorCode    ierr;
1163 
1164   PetscFunctionBegin;
1165   *missing = PETSC_FALSE;
1166   ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr);
1167   ierr = MatGetDiagonal(A,v);CHKERRQ(ierr);
1168   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1169   ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr);
1170   for (i=0;i<n;i++) if (array[i] == 0.) break;
1171   ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr);
1172   ierr = VecDestroy(&v);CHKERRQ(ierr);
1173   if (i != n) *missing = PETSC_TRUE;
1174   if (d) {
1175     *d = -1;
1176     if (*missing) {
1177       PetscInt rstart;
1178       ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
1179       *d = i+rstart;
1180     }
1181   }
1182   PetscFunctionReturn(0);
1183 }
1184 
1185 static PetscErrorCode MatISSetUpSF_IS(Mat B)
1186 {
1187   Mat_IS         *matis = (Mat_IS*)(B->data);
1188   const PetscInt *gidxs;
1189   PetscInt       nleaves;
1190   PetscErrorCode ierr;
1191 
1192   PetscFunctionBegin;
1193   if (matis->sf) PetscFunctionReturn(0);
1194   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
1195   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
1196   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
1197   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
1198   ierr = ISLocalToGlobalMappingGetSize(B->rmap->mapping,&nleaves);CHKERRQ(ierr);
1199   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
1200   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
1201   ierr = PetscMalloc2(matis->sf->nroots,&matis->sf_rootdata,matis->sf->nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
1202   if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */
1203     ierr = ISLocalToGlobalMappingGetSize(B->cmap->mapping,&nleaves);CHKERRQ(ierr);
1204     ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr);
1205     ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
1206     ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
1207     ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
1208     ierr = PetscMalloc2(matis->csf->nroots,&matis->csf_rootdata,matis->csf->nleaves,&matis->csf_leafdata);CHKERRQ(ierr);
1209   } else {
1210     matis->csf = matis->sf;
1211     matis->csf_leafdata = matis->sf_leafdata;
1212     matis->csf_rootdata = matis->sf_rootdata;
1213   }
1214   PetscFunctionReturn(0);
1215 }
1216 
1217 /*@
1218    MatISStoreL2L - Store local-to-local operators during the Galerkin process of MatPtAP.
1219 
1220    Collective on MPI_Comm
1221 
1222    Input Parameters:
1223 +  A - the matrix
1224 -  store - the boolean flag
1225 
1226    Level: advanced
1227 
1228    Notes:
1229 
1230 .keywords: matrix
1231 
1232 .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatPtAP()
1233 @*/
1234 PetscErrorCode  MatISStoreL2L(Mat A, PetscBool store)
1235 {
1236   PetscErrorCode ierr;
1237 
1238   PetscFunctionBegin;
1239   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1240   PetscValidType(A,1);
1241   PetscValidLogicalCollectiveBool(A,store,2);
1242   ierr = PetscTryMethod(A,"MatISStoreL2L_C",(Mat,PetscBool),(A,store));CHKERRQ(ierr);
1243   PetscFunctionReturn(0);
1244 }
1245 
1246 static PetscErrorCode MatISStoreL2L_IS(Mat A, PetscBool store)
1247 {
1248   Mat_IS         *matis = (Mat_IS*)(A->data);
1249   PetscErrorCode ierr;
1250 
1251   PetscFunctionBegin;
1252   matis->storel2l = store;
1253   if (!store) {
1254     ierr = PetscObjectCompose((PetscObject)(A),"_MatIS_PtAP_l2l",NULL);CHKERRQ(ierr);
1255   }
1256   PetscFunctionReturn(0);
1257 }
1258 
1259 /*@
1260    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
1261 
1262    Collective on MPI_Comm
1263 
1264    Input Parameters:
1265 +  B - the matrix
1266 .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1267            (same value is used for all local rows)
1268 .  d_nnz - array containing the number of nonzeros in the various rows of the
1269            DIAGONAL portion of the local submatrix (possibly different for each row)
1270            or NULL, if d_nz is used to specify the nonzero structure.
1271            The size of this array is equal to the number of local rows, i.e 'm'.
1272            For matrices that will be factored, you must leave room for (and set)
1273            the diagonal entry even if it is zero.
1274 .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1275            submatrix (same value is used for all local rows).
1276 -  o_nnz - array containing the number of nonzeros in the various rows of the
1277            OFF-DIAGONAL portion of the local submatrix (possibly different for
1278            each row) or NULL, if o_nz is used to specify the nonzero
1279            structure. The size of this array is equal to the number
1280            of local rows, i.e 'm'.
1281 
1282    If the *_nnz parameter is given then the *_nz parameter is ignored
1283 
1284    Level: intermediate
1285 
1286    Notes:
1287     This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
1288           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
1289           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
1290 
1291 .keywords: matrix
1292 
1293 .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS
1294 @*/
1295 PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1296 {
1297   PetscErrorCode ierr;
1298 
1299   PetscFunctionBegin;
1300   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
1301   PetscValidType(B,1);
1302   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 static PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1307 {
1308   Mat_IS         *matis = (Mat_IS*)(B->data);
1309   PetscInt       bs,i,nlocalcols;
1310   PetscErrorCode ierr;
1311 
1312   PetscFunctionBegin;
1313   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
1314   ierr = MatISSetUpSF(B);CHKERRQ(ierr);
1315 
1316   if (!d_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nz;
1317   else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
1318 
1319   if (!o_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nz;
1320   else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
1321 
1322   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
1323   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
1324   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
1325   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
1326 
1327   for (i=0;i<matis->sf->nleaves;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
1328   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
1329 #if defined(PETSC_HAVE_HYPRE)
1330   ierr = MatHYPRESetPreallocation(matis->A,0,matis->sf_leafdata,0,NULL);CHKERRQ(ierr);
1331 #endif
1332 
1333   for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
1334   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
1335 
1336   for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
1337   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
1338 
1339   /* for other matrix types */
1340   ierr = MatSetUp(matis->A);CHKERRQ(ierr);
1341   PetscFunctionReturn(0);
1342 }
1343 
1344 PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
1345 {
1346   Mat_IS          *matis = (Mat_IS*)(A->data);
1347   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
1348   const PetscInt  *global_indices_r,*global_indices_c;
1349   PetscInt        i,j,bs,rows,cols;
1350   PetscInt        lrows,lcols;
1351   PetscInt        local_rows,local_cols;
1352   PetscMPIInt     nsubdomains;
1353   PetscBool       isdense,issbaij;
1354   PetscErrorCode  ierr;
1355 
1356   PetscFunctionBegin;
1357   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
1358   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
1359   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1360   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
1361   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
1362   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
1363   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
1364   if (A->rmap->mapping != A->cmap->mapping) {
1365     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
1366   } else {
1367     global_indices_c = global_indices_r;
1368   }
1369 
1370   if (issbaij) {
1371     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
1372   }
1373   /*
1374      An SF reduce is needed to sum up properly on shared rows.
1375      Note that generally preallocation is not exact, since it overestimates nonzeros
1376   */
1377   ierr = MatISSetUpSF(A);CHKERRQ(ierr);
1378   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
1379   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
1380   /* All processes need to compute entire row ownership */
1381   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
1382   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1383   for (i=0;i<nsubdomains;i++) {
1384     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
1385       row_ownership[j] = i;
1386     }
1387   }
1388   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1389 
1390   /*
1391      my_dnz and my_onz contains exact contribution to preallocation from each local mat
1392      then, they will be summed up properly. This way, preallocation is always sufficient
1393   */
1394   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
1395   /* preallocation as a MATAIJ */
1396   if (isdense) { /* special case for dense local matrices */
1397     for (i=0;i<local_rows;i++) {
1398       PetscInt owner = row_ownership[global_indices_r[i]];
1399       for (j=0;j<local_cols;j++) {
1400         PetscInt index_col = global_indices_c[j];
1401         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
1402           my_dnz[i] += 1;
1403         } else { /* offdiag block */
1404           my_onz[i] += 1;
1405         }
1406       }
1407     }
1408   } else if (matis->A->ops->getrowij) {
1409     const PetscInt *ii,*jj,*jptr;
1410     PetscBool      done;
1411     ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
1412     if (!done) SETERRQ(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ");
1413     jptr = jj;
1414     for (i=0;i<local_rows;i++) {
1415       PetscInt index_row = global_indices_r[i];
1416       for (j=0;j<ii[i+1]-ii[i];j++,jptr++) {
1417         PetscInt owner = row_ownership[index_row];
1418         PetscInt index_col = global_indices_c[*jptr];
1419         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
1420           my_dnz[i] += 1;
1421         } else { /* offdiag block */
1422           my_onz[i] += 1;
1423         }
1424         /* same as before, interchanging rows and cols */
1425         if (issbaij && index_col != index_row) {
1426           owner = row_ownership[index_col];
1427           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
1428             my_dnz[*jptr] += 1;
1429           } else {
1430             my_onz[*jptr] += 1;
1431           }
1432         }
1433       }
1434     }
1435     ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
1436     if (!done) SETERRQ(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ");
1437   } else { /* loop over rows and use MatGetRow */
1438     for (i=0;i<local_rows;i++) {
1439       const PetscInt *cols;
1440       PetscInt       ncols,index_row = global_indices_r[i];
1441       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
1442       for (j=0;j<ncols;j++) {
1443         PetscInt owner = row_ownership[index_row];
1444         PetscInt index_col = global_indices_c[cols[j]];
1445         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
1446           my_dnz[i] += 1;
1447         } else { /* offdiag block */
1448           my_onz[i] += 1;
1449         }
1450         /* same as before, interchanging rows and cols */
1451         if (issbaij && index_col != index_row) {
1452           owner = row_ownership[index_col];
1453           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
1454             my_dnz[cols[j]] += 1;
1455           } else {
1456             my_onz[cols[j]] += 1;
1457           }
1458         }
1459       }
1460       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
1461     }
1462   }
1463   if (global_indices_c != global_indices_r) {
1464     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
1465   }
1466   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
1467   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
1468 
1469   /* Reduce my_dnz and my_onz */
1470   if (maxreduce) {
1471     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
1472     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
1473     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
1474     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
1475   } else {
1476     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
1477     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
1478     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
1479     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
1480   }
1481   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
1482 
1483   /* Resize preallocation if overestimated */
1484   for (i=0;i<lrows;i++) {
1485     dnz[i] = PetscMin(dnz[i],lcols);
1486     onz[i] = PetscMin(onz[i],cols-lcols);
1487   }
1488 
1489   /* Set preallocation */
1490   ierr = MatSeqAIJSetPreallocation(B,0,dnz);CHKERRQ(ierr);
1491   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
1492   for (i=0;i<lrows/bs;i++) {
1493     dnz[i] = dnz[i*bs]/bs;
1494     onz[i] = onz[i*bs]/bs;
1495   }
1496   ierr = MatSeqBAIJSetPreallocation(B,bs,0,dnz);CHKERRQ(ierr);
1497   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
1498   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
1499   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1500   if (issbaij) {
1501     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
1502   }
1503   PetscFunctionReturn(0);
1504 }
1505 
1506 static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
1507 {
1508   Mat_IS         *matis = (Mat_IS*)(mat->data);
1509   Mat            local_mat;
1510   /* info on mat */
1511   PetscInt       bs,rows,cols,lrows,lcols;
1512   PetscInt       local_rows,local_cols;
1513   PetscBool      isseqdense,isseqsbaij,isseqaij,isseqbaij;
1514 #if defined (PETSC_USE_DEBUG)
1515   PetscBool      lb[4],bb[4];
1516 #endif
1517   PetscMPIInt    nsubdomains;
1518   /* values insertion */
1519   PetscScalar    *array;
1520   /* work */
1521   PetscErrorCode ierr;
1522 
1523   PetscFunctionBegin;
1524   /* get info from mat */
1525   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
1526   if (nsubdomains == 1) {
1527     Mat            B;
1528     IS             rows,cols;
1529     IS             irows,icols;
1530     const PetscInt *ridxs,*cidxs;
1531 
1532     ierr = ISLocalToGlobalMappingGetIndices(mat->rmap->mapping,&ridxs);CHKERRQ(ierr);
1533     ierr = ISLocalToGlobalMappingGetIndices(mat->cmap->mapping,&cidxs);CHKERRQ(ierr);
1534     ierr = ISCreateGeneral(PETSC_COMM_SELF,mat->rmap->n,ridxs,PETSC_USE_POINTER,&rows);CHKERRQ(ierr);
1535     ierr = ISCreateGeneral(PETSC_COMM_SELF,mat->cmap->n,cidxs,PETSC_USE_POINTER,&cols);CHKERRQ(ierr);
1536     ierr = ISSetPermutation(rows);CHKERRQ(ierr);
1537     ierr = ISSetPermutation(cols);CHKERRQ(ierr);
1538     ierr = ISInvertPermutation(rows,mat->rmap->n,&irows);CHKERRQ(ierr);
1539     ierr = ISInvertPermutation(cols,mat->cmap->n,&icols);CHKERRQ(ierr);
1540     ierr = ISLocalToGlobalMappingRestoreIndices(mat->rmap->mapping,&ridxs);CHKERRQ(ierr);
1541     ierr = ISLocalToGlobalMappingRestoreIndices(mat->cmap->mapping,&cidxs);CHKERRQ(ierr);
1542     ierr = ISDestroy(&cols);CHKERRQ(ierr);
1543     ierr = ISDestroy(&rows);CHKERRQ(ierr);
1544     ierr = MatConvert(matis->A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
1545     ierr = MatCreateSubMatrix(B,irows,icols,reuse,M);CHKERRQ(ierr);
1546     ierr = MatDestroy(&B);CHKERRQ(ierr);
1547     ierr = ISDestroy(&icols);CHKERRQ(ierr);
1548     ierr = ISDestroy(&irows);CHKERRQ(ierr);
1549     PetscFunctionReturn(0);
1550   }
1551   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
1552   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1553   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
1554   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
1555   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isseqdense);CHKERRQ(ierr);
1556   ierr = PetscObjectBaseTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1557   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr);
1558   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
1559   if (!isseqdense && !isseqaij && !isseqbaij && !isseqsbaij) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)(matis->A))->type_name);
1560 #if defined (PETSC_USE_DEBUG)
1561   lb[0] = isseqdense;
1562   lb[1] = isseqaij;
1563   lb[2] = isseqbaij;
1564   lb[3] = isseqsbaij;
1565   ierr = MPIU_Allreduce(lb,bb,4,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
1566   if (!bb[0] && !bb[1] && !bb[2] && !bb[3]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local matrices must have the same type");
1567 #endif
1568 
1569   if (reuse == MAT_INITIAL_MATRIX) {
1570     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
1571     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
1572     if (!isseqsbaij) {
1573       ierr = MatSetType(*M,MATAIJ);CHKERRQ(ierr);
1574     } else {
1575       ierr = MatSetType(*M,MATSBAIJ);CHKERRQ(ierr);
1576     }
1577     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
1578     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
1579   } else {
1580     PetscInt mbs,mrows,mcols,mlrows,mlcols;
1581     /* some checks */
1582     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
1583     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
1584     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
1585     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
1586     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
1587     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
1588     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
1589     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
1590     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
1591   }
1592 
1593   if (isseqsbaij) {
1594     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
1595   } else {
1596     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
1597     local_mat = matis->A;
1598   }
1599 
1600   /* Set values */
1601   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
1602   if (isseqdense) { /* special case for dense local matrices */
1603     PetscInt i,*dummy;
1604 
1605     ierr = PetscMalloc1(PetscMax(local_rows,local_cols),&dummy);CHKERRQ(ierr);
1606     for (i=0;i<PetscMax(local_rows,local_cols);i++) dummy[i] = i;
1607     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
1608     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
1609     ierr = MatSetValuesLocal(*M,local_rows,dummy,local_cols,dummy,array,ADD_VALUES);CHKERRQ(ierr);
1610     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
1611     ierr = PetscFree(dummy);CHKERRQ(ierr);
1612   } else if (isseqaij) {
1613     PetscInt  i,nvtxs,*xadj,*adjncy;
1614     PetscBool done;
1615 
1616     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
1617     if (!done) SETERRQ(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ");
1618     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
1619     for (i=0;i<nvtxs;i++) {
1620       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
1621     }
1622     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
1623     if (!done) SETERRQ(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ");
1624     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
1625   } else { /* very basic values insertion for all other matrix types */
1626     PetscInt i;
1627 
1628     for (i=0;i<local_rows;i++) {
1629       PetscInt       j;
1630       const PetscInt *local_indices_cols;
1631 
1632       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
1633       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
1634       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
1635     }
1636   }
1637   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1638   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
1639   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1640   if (isseqdense) {
1641     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
1642   }
1643   PetscFunctionReturn(0);
1644 }
1645 
1646 /*@
1647     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
1648 
1649   Input Parameter:
1650 .  mat - the matrix (should be of type MATIS)
1651 .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
1652 
1653   Output Parameter:
1654 .  newmat - the matrix in AIJ format
1655 
1656   Level: developer
1657 
1658   Notes:
1659     mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
1660 
1661 .seealso: MATIS
1662 @*/
1663 PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
1664 {
1665   PetscErrorCode ierr;
1666 
1667   PetscFunctionBegin;
1668   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1669   PetscValidLogicalCollectiveEnum(mat,reuse,2);
1670   PetscValidPointer(newmat,3);
1671   if (reuse != MAT_INITIAL_MATRIX) {
1672     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
1673     PetscCheckSameComm(mat,1,*newmat,3);
1674     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
1675   }
1676   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
1677   PetscFunctionReturn(0);
1678 }
1679 
1680 PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
1681 {
1682   PetscErrorCode ierr;
1683   Mat_IS         *matis = (Mat_IS*)(mat->data);
1684   PetscInt       bs,m,n,M,N;
1685   Mat            B,localmat;
1686 
1687   PetscFunctionBegin;
1688   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1689   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
1690   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
1691   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
1692   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
1693   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
1694   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
1695   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1696   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1697   *newmat = B;
1698   PetscFunctionReturn(0);
1699 }
1700 
1701 static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
1702 {
1703   PetscErrorCode ierr;
1704   Mat_IS         *matis = (Mat_IS*)A->data;
1705   PetscBool      local_sym;
1706 
1707   PetscFunctionBegin;
1708   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
1709   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1710   PetscFunctionReturn(0);
1711 }
1712 
1713 static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
1714 {
1715   PetscErrorCode ierr;
1716   Mat_IS         *matis = (Mat_IS*)A->data;
1717   PetscBool      local_sym;
1718 
1719   PetscFunctionBegin;
1720   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
1721   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1722   PetscFunctionReturn(0);
1723 }
1724 
1725 static PetscErrorCode MatIsStructurallySymmetric_IS(Mat A,PetscBool  *flg)
1726 {
1727   PetscErrorCode ierr;
1728   Mat_IS         *matis = (Mat_IS*)A->data;
1729   PetscBool      local_sym;
1730 
1731   PetscFunctionBegin;
1732   if (A->rmap->mapping != A->cmap->mapping) {
1733     *flg = PETSC_FALSE;
1734     PetscFunctionReturn(0);
1735   }
1736   ierr = MatIsStructurallySymmetric(matis->A,&local_sym);CHKERRQ(ierr);
1737   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1738   PetscFunctionReturn(0);
1739 }
1740 
1741 static PetscErrorCode MatDestroy_IS(Mat A)
1742 {
1743   PetscErrorCode ierr;
1744   Mat_IS         *b = (Mat_IS*)A->data;
1745 
1746   PetscFunctionBegin;
1747   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
1748   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
1749   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
1750   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
1751   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
1752   ierr = VecDestroy(&b->counter);CHKERRQ(ierr);
1753   ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr);
1754   ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr);
1755   if (b->sf != b->csf) {
1756     ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr);
1757     ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr);
1758   }
1759   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
1760   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
1761   ierr = PetscFree(A->data);CHKERRQ(ierr);
1762   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
1763   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
1764   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
1765   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
1766   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
1767   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",NULL);CHKERRQ(ierr);
1768   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISStoreL2L_C",NULL);CHKERRQ(ierr);
1769   PetscFunctionReturn(0);
1770 }
1771 
1772 static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
1773 {
1774   PetscErrorCode ierr;
1775   Mat_IS         *is  = (Mat_IS*)A->data;
1776   PetscScalar    zero = 0.0;
1777 
1778   PetscFunctionBegin;
1779   /*  scatter the global vector x into the local work vector */
1780   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1781   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1782 
1783   /* multiply the local matrix */
1784   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
1785 
1786   /* scatter product back into global memory */
1787   ierr = VecSet(y,zero);CHKERRQ(ierr);
1788   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1789   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1790   PetscFunctionReturn(0);
1791 }
1792 
1793 static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
1794 {
1795   Vec            temp_vec;
1796   PetscErrorCode ierr;
1797 
1798   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
1799   if (v3 != v2) {
1800     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
1801     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
1802   } else {
1803     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
1804     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
1805     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
1806     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
1807     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
1808   }
1809   PetscFunctionReturn(0);
1810 }
1811 
1812 static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
1813 {
1814   Mat_IS         *is = (Mat_IS*)A->data;
1815   PetscErrorCode ierr;
1816 
1817   PetscFunctionBegin;
1818   /*  scatter the global vector x into the local work vector */
1819   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1820   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1821 
1822   /* multiply the local matrix */
1823   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
1824 
1825   /* scatter product back into global vector */
1826   ierr = VecSet(x,0);CHKERRQ(ierr);
1827   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1828   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1829   PetscFunctionReturn(0);
1830 }
1831 
1832 static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
1833 {
1834   Vec            temp_vec;
1835   PetscErrorCode ierr;
1836 
1837   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
1838   if (v3 != v2) {
1839     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
1840     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
1841   } else {
1842     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
1843     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
1844     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
1845     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
1846     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
1847   }
1848   PetscFunctionReturn(0);
1849 }
1850 
1851 static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
1852 {
1853   Mat_IS         *a = (Mat_IS*)A->data;
1854   PetscErrorCode ierr;
1855   PetscViewer    sviewer;
1856   PetscBool      isascii,view = PETSC_TRUE;
1857 
1858   PetscFunctionBegin;
1859   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1860   if (isascii)  {
1861     PetscViewerFormat format;
1862 
1863     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1864     if (format == PETSC_VIEWER_ASCII_INFO) view = PETSC_FALSE;
1865   }
1866   if (!view) PetscFunctionReturn(0);
1867   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
1868   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
1869   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
1870   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1871   PetscFunctionReturn(0);
1872 }
1873 
1874 static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
1875 {
1876   PetscErrorCode ierr;
1877   PetscInt       nr,rbs,nc,cbs;
1878   Mat_IS         *is = (Mat_IS*)A->data;
1879   IS             from,to;
1880   Vec            cglobal,rglobal;
1881 
1882   PetscFunctionBegin;
1883   PetscCheckSameComm(A,1,rmapping,2);
1884   PetscCheckSameComm(A,1,cmapping,3);
1885   /* Destroy any previous data */
1886   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
1887   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
1888   ierr = VecDestroy(&is->counter);CHKERRQ(ierr);
1889   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
1890   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
1891   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
1892   if (is->csf != is->sf) {
1893     ierr = PetscSFDestroy(&is->csf);CHKERRQ(ierr);
1894     ierr = PetscFree2(is->csf_rootdata,is->csf_leafdata);CHKERRQ(ierr);
1895   }
1896   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
1897   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
1898 
1899   /* Setup Layout and set local to global maps */
1900   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1901   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1902   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
1903   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
1904   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
1905   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
1906   /* check if the two mappings are actually the same for square matrices (DOLFIN passes 2 different objects) */
1907   if (rmapping != cmapping && A->rmap->N == A->cmap->N) {
1908     PetscBool same,gsame;
1909 
1910     same = PETSC_FALSE;
1911     if (nr == nc && cbs == rbs) {
1912       const PetscInt *idxs1,*idxs2;
1913 
1914       ierr = ISLocalToGlobalMappingGetBlockIndices(rmapping,&idxs1);CHKERRQ(ierr);
1915       ierr = ISLocalToGlobalMappingGetBlockIndices(cmapping,&idxs2);CHKERRQ(ierr);
1916       ierr = PetscMemcmp(idxs1,idxs2,(nr/rbs)*sizeof(PetscInt),&same);CHKERRQ(ierr);
1917       ierr = ISLocalToGlobalMappingRestoreBlockIndices(rmapping,&idxs1);CHKERRQ(ierr);
1918       ierr = ISLocalToGlobalMappingRestoreBlockIndices(cmapping,&idxs2);CHKERRQ(ierr);
1919     }
1920     ierr = MPIU_Allreduce(&same,&gsame,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1921     if (gsame) cmapping = rmapping;
1922   }
1923   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
1924   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
1925 
1926   /* Create the local matrix A */
1927   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
1928   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
1929   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
1930   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
1931   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
1932   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
1933   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
1934   ierr = PetscLayoutSetUp(is->A->rmap);CHKERRQ(ierr);
1935   ierr = PetscLayoutSetUp(is->A->cmap);CHKERRQ(ierr);
1936 
1937   if (!is->islocalref) { /* setup scatters and local vectors for MatMult */
1938     /* Create the local work vectors */
1939     ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
1940 
1941     /* setup the global to local scatters */
1942     ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
1943     ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
1944     ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
1945     ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
1946     if (rmapping != cmapping) {
1947       ierr = ISDestroy(&to);CHKERRQ(ierr);
1948       ierr = ISDestroy(&from);CHKERRQ(ierr);
1949       ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
1950       ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
1951       ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
1952     } else {
1953       ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
1954       is->cctx = is->rctx;
1955     }
1956 
1957     /* interface counter vector (local) */
1958     ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr);
1959     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
1960     ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1961     ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1962     ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1963     ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1964 
1965     /* free workspace */
1966     ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
1967     ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
1968     ierr = ISDestroy(&to);CHKERRQ(ierr);
1969     ierr = ISDestroy(&from);CHKERRQ(ierr);
1970   }
1971   ierr = MatSetUp(A);CHKERRQ(ierr);
1972   PetscFunctionReturn(0);
1973 }
1974 
1975 static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
1976 {
1977   Mat_IS         *is = (Mat_IS*)mat->data;
1978   PetscErrorCode ierr;
1979 #if defined(PETSC_USE_DEBUG)
1980   PetscInt       i,zm,zn;
1981 #endif
1982   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
1983 
1984   PetscFunctionBegin;
1985 #if defined(PETSC_USE_DEBUG)
1986   if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n);
1987   /* count negative indices */
1988   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
1989   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
1990 #endif
1991   ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
1992   ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
1993 #if defined(PETSC_USE_DEBUG)
1994   /* count negative indices (should be the same as before) */
1995   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
1996   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
1997   if (!is->A->rmap->mapping && zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS");
1998   if (!is->A->cmap->mapping && zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS");
1999 #endif
2000   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
2001   PetscFunctionReturn(0);
2002 }
2003 
2004 static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
2005 {
2006   Mat_IS         *is = (Mat_IS*)mat->data;
2007   PetscErrorCode ierr;
2008 #if defined(PETSC_USE_DEBUG)
2009   PetscInt       i,zm,zn;
2010 #endif
2011   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
2012 
2013   PetscFunctionBegin;
2014 #if defined(PETSC_USE_DEBUG)
2015   if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n);
2016   /* count negative indices */
2017   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
2018   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
2019 #endif
2020   ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
2021   ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
2022 #if defined(PETSC_USE_DEBUG)
2023   /* count negative indices (should be the same as before) */
2024   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
2025   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
2026   if (!is->A->rmap->mapping && zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS");
2027   if (!is->A->cmap->mapping && zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS");
2028 #endif
2029   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
2030   PetscFunctionReturn(0);
2031 }
2032 
2033 static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
2034 {
2035   PetscErrorCode ierr;
2036   Mat_IS         *is = (Mat_IS*)A->data;
2037 
2038   PetscFunctionBegin;
2039   if (is->A->rmap->mapping) {
2040     ierr = MatSetValuesLocal(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
2041   } else {
2042     ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
2043   }
2044   PetscFunctionReturn(0);
2045 }
2046 
2047 static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
2048 {
2049   PetscErrorCode ierr;
2050   Mat_IS         *is = (Mat_IS*)A->data;
2051 
2052   PetscFunctionBegin;
2053   if (is->A->rmap->mapping) {
2054 #if defined(PETSC_USE_DEBUG)
2055     PetscInt ibs,bs;
2056 
2057     ierr = ISLocalToGlobalMappingGetBlockSize(is->A->rmap->mapping,&ibs);CHKERRQ(ierr);
2058     ierr = MatGetBlockSize(is->A,&bs);CHKERRQ(ierr);
2059     if (ibs != bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Different block sizes! local mat %D, local l2g map %D",bs,ibs);
2060 #endif
2061     ierr = MatSetValuesBlockedLocal(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
2062   } else {
2063     ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
2064   }
2065   PetscFunctionReturn(0);
2066 }
2067 
2068 static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns)
2069 {
2070   Mat_IS         *is = (Mat_IS*)A->data;
2071   PetscErrorCode ierr;
2072 
2073   PetscFunctionBegin;
2074   if (!n) {
2075     is->pure_neumann = PETSC_TRUE;
2076   } else {
2077     PetscInt i;
2078     is->pure_neumann = PETSC_FALSE;
2079 
2080     if (columns) {
2081       ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
2082     } else {
2083       ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
2084     }
2085     if (diag != 0.) {
2086       const PetscScalar *array;
2087       ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr);
2088       for (i=0; i<n; i++) {
2089         ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
2090       }
2091       ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr);
2092     }
2093     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2094     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2095   }
2096   PetscFunctionReturn(0);
2097 }
2098 
2099 static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns)
2100 {
2101   Mat_IS         *matis = (Mat_IS*)A->data;
2102   PetscInt       nr,nl,len,i;
2103   PetscInt       *lrows;
2104   PetscErrorCode ierr;
2105 
2106   PetscFunctionBegin;
2107 #if defined(PETSC_USE_DEBUG)
2108   if (columns || diag != 0. || (x && b)) {
2109     PetscBool cong;
2110 
2111     ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr);
2112     cong = (PetscBool)(cong && matis->sf == matis->csf);
2113     if (!cong && columns) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Columns can be zeroed if and only if A->rmap and A->cmap are congruent and the l2g maps are the same for MATIS");
2114     if (!cong && diag != 0.) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Nonzero diagonal value supported if and only if A->rmap and A->cmap are congruent and the l2g maps are the same for MATIS");
2115     if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent, and the l2g maps be the same");
2116   }
2117 #endif
2118   /* get locally owned rows */
2119   ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr);
2120   /* fix right hand side if needed */
2121   if (x && b) {
2122     const PetscScalar *xx;
2123     PetscScalar       *bb;
2124 
2125     ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr);
2126     ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
2127     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
2128     ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr);
2129     ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
2130   }
2131   /* get rows associated to the local matrices */
2132   ierr = MatISSetUpSF(A);CHKERRQ(ierr);
2133   ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr);
2134   ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr);
2135   ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
2136   for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1;
2137   ierr = PetscFree(lrows);CHKERRQ(ierr);
2138   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
2139   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
2140   ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr);
2141   for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i;
2142   ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr);
2143   ierr = PetscFree(lrows);CHKERRQ(ierr);
2144   PetscFunctionReturn(0);
2145 }
2146 
2147 static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2148 {
2149   PetscErrorCode ierr;
2150 
2151   PetscFunctionBegin;
2152   ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr);
2153   PetscFunctionReturn(0);
2154 }
2155 
2156 static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2157 {
2158   PetscErrorCode ierr;
2159 
2160   PetscFunctionBegin;
2161   ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr);
2162   PetscFunctionReturn(0);
2163 }
2164 
2165 static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
2166 {
2167   Mat_IS         *is = (Mat_IS*)A->data;
2168   PetscErrorCode ierr;
2169 
2170   PetscFunctionBegin;
2171   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
2172   PetscFunctionReturn(0);
2173 }
2174 
2175 static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
2176 {
2177   Mat_IS         *is = (Mat_IS*)A->data;
2178   PetscErrorCode ierr;
2179 
2180   PetscFunctionBegin;
2181   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
2182   /* fix for local empty rows/cols */
2183   if (is->locempty && type == MAT_FINAL_ASSEMBLY) {
2184     Mat                    newlA;
2185     ISLocalToGlobalMapping l2g;
2186     IS                     tis;
2187     const PetscScalar      *v;
2188     PetscInt               i,n,cf,*loce,*locf;
2189     PetscBool              sym;
2190 
2191     ierr = MatGetRowMaxAbs(is->A,is->y,NULL);CHKERRQ(ierr);
2192     ierr = MatIsSymmetric(is->A,PETSC_SMALL,&sym);CHKERRQ(ierr);
2193     if (!sym) SETERRQ(PetscObjectComm((PetscObject)is->A),PETSC_ERR_SUP,"Not yet implemented for unsymmetric case");
2194     ierr = VecGetLocalSize(is->y,&n);CHKERRQ(ierr);
2195     ierr = PetscMalloc1(n,&loce);CHKERRQ(ierr);
2196     ierr = PetscMalloc1(n,&locf);CHKERRQ(ierr);
2197     ierr = VecGetArrayRead(is->y,&v);CHKERRQ(ierr);
2198     for (i=0,cf=0;i<n;i++) {
2199       if (v[i] == 0.0) {
2200         loce[i] = -1;
2201       } else {
2202         loce[i]    = cf;
2203         locf[cf++] = i;
2204       }
2205     }
2206     ierr = VecRestoreArrayRead(is->y,&v);CHKERRQ(ierr);
2207     /* extract valid submatrix */
2208     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is->A),cf,locf,PETSC_USE_POINTER,&tis);CHKERRQ(ierr);
2209     ierr = MatCreateSubMatrix(is->A,tis,tis,MAT_INITIAL_MATRIX,&newlA);CHKERRQ(ierr);
2210     ierr = ISDestroy(&tis);CHKERRQ(ierr);
2211     /* attach local l2g map for successive calls of MatSetValues */
2212     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is->A),1,n,loce,PETSC_OWN_POINTER,&l2g);CHKERRQ(ierr);
2213     ierr = MatSetLocalToGlobalMapping(newlA,l2g,l2g);CHKERRQ(ierr);
2214     ierr = ISLocalToGlobalMappingDestroy(&l2g);CHKERRQ(ierr);
2215     /* attach new global l2g map */
2216     ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,cf,locf,locf);CHKERRQ(ierr);
2217     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,cf,locf,PETSC_OWN_POINTER,&l2g);CHKERRQ(ierr);
2218     ierr = MatSetLocalToGlobalMapping(A,l2g,l2g);CHKERRQ(ierr);
2219     ierr = MatISSetLocalMat(A,newlA);CHKERRQ(ierr);
2220     ierr = MatDestroy(&newlA);CHKERRQ(ierr);
2221     ierr = ISLocalToGlobalMappingDestroy(&l2g);CHKERRQ(ierr);
2222   }
2223   is->locempty = PETSC_FALSE;
2224   PetscFunctionReturn(0);
2225 }
2226 
2227 static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
2228 {
2229   Mat_IS *is = (Mat_IS*)mat->data;
2230 
2231   PetscFunctionBegin;
2232   *local = is->A;
2233   PetscFunctionReturn(0);
2234 }
2235 
2236 static PetscErrorCode MatISRestoreLocalMat_IS(Mat mat,Mat *local)
2237 {
2238   PetscFunctionBegin;
2239   *local = NULL;
2240   PetscFunctionReturn(0);
2241 }
2242 
2243 /*@
2244     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
2245 
2246   Input Parameter:
2247 .  mat - the matrix
2248 
2249   Output Parameter:
2250 .  local - the local matrix
2251 
2252   Level: advanced
2253 
2254   Notes:
2255     This can be called if you have precomputed the nonzero structure of the
2256   matrix and want to provide it to the inner matrix object to improve the performance
2257   of the MatSetValues() operation.
2258 
2259   Call MatISRestoreLocalMat() when finished with the local matrix.
2260 
2261 .seealso: MATIS
2262 @*/
2263 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
2264 {
2265   PetscErrorCode ierr;
2266 
2267   PetscFunctionBegin;
2268   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2269   PetscValidPointer(local,2);
2270   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
2271   PetscFunctionReturn(0);
2272 }
2273 
2274 /*@
2275     MatISRestoreLocalMat - Restores the local matrix obtained with MatISGetLocalMat()
2276 
2277   Input Parameter:
2278 .  mat - the matrix
2279 
2280   Output Parameter:
2281 .  local - the local matrix
2282 
2283   Level: advanced
2284 
2285 .seealso: MATIS
2286 @*/
2287 PetscErrorCode MatISRestoreLocalMat(Mat mat,Mat *local)
2288 {
2289   PetscErrorCode ierr;
2290 
2291   PetscFunctionBegin;
2292   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2293   PetscValidPointer(local,2);
2294   ierr = PetscUseMethod(mat,"MatISRestoreLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
2295   PetscFunctionReturn(0);
2296 }
2297 
2298 static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
2299 {
2300   Mat_IS         *is = (Mat_IS*)mat->data;
2301   PetscInt       nrows,ncols,orows,ocols;
2302   PetscErrorCode ierr;
2303 
2304   PetscFunctionBegin;
2305   if (is->A) {
2306     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
2307     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
2308     if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %Dx%D (you passed a %Dx%D matrix)",orows,ocols,nrows,ncols);
2309   }
2310   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
2311   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
2312   is->A = local;
2313   PetscFunctionReturn(0);
2314 }
2315 
2316 /*@
2317     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
2318 
2319   Input Parameter:
2320 .  mat - the matrix
2321 .  local - the local matrix
2322 
2323   Output Parameter:
2324 
2325   Level: advanced
2326 
2327   Notes:
2328     This can be called if you have precomputed the local matrix and
2329   want to provide it to the matrix object MATIS.
2330 
2331 .seealso: MATIS
2332 @*/
2333 PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
2334 {
2335   PetscErrorCode ierr;
2336 
2337   PetscFunctionBegin;
2338   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
2339   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
2340   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
2341   PetscFunctionReturn(0);
2342 }
2343 
2344 static PetscErrorCode MatZeroEntries_IS(Mat A)
2345 {
2346   Mat_IS         *a = (Mat_IS*)A->data;
2347   PetscErrorCode ierr;
2348 
2349   PetscFunctionBegin;
2350   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
2351   PetscFunctionReturn(0);
2352 }
2353 
2354 static PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
2355 {
2356   Mat_IS         *is = (Mat_IS*)A->data;
2357   PetscErrorCode ierr;
2358 
2359   PetscFunctionBegin;
2360   ierr = MatScale(is->A,a);CHKERRQ(ierr);
2361   PetscFunctionReturn(0);
2362 }
2363 
2364 static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
2365 {
2366   Mat_IS         *is = (Mat_IS*)A->data;
2367   PetscErrorCode ierr;
2368 
2369   PetscFunctionBegin;
2370   /* get diagonal of the local matrix */
2371   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
2372 
2373   /* scatter diagonal back into global vector */
2374   ierr = VecSet(v,0);CHKERRQ(ierr);
2375   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2376   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2377   PetscFunctionReturn(0);
2378 }
2379 
2380 static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
2381 {
2382   Mat_IS         *a = (Mat_IS*)A->data;
2383   PetscErrorCode ierr;
2384 
2385   PetscFunctionBegin;
2386   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
2387   PetscFunctionReturn(0);
2388 }
2389 
2390 static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str)
2391 {
2392   Mat_IS         *y = (Mat_IS*)Y->data;
2393   Mat_IS         *x;
2394 #if defined(PETSC_USE_DEBUG)
2395   PetscBool      ismatis;
2396 #endif
2397   PetscErrorCode ierr;
2398 
2399   PetscFunctionBegin;
2400 #if defined(PETSC_USE_DEBUG)
2401   ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr);
2402   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS");
2403 #endif
2404   x = (Mat_IS*)X->data;
2405   ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr);
2406   PetscFunctionReturn(0);
2407 }
2408 
2409 static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat)
2410 {
2411   Mat                    lA;
2412   Mat_IS                 *matis;
2413   ISLocalToGlobalMapping rl2g,cl2g;
2414   IS                     is;
2415   const PetscInt         *rg,*rl;
2416   PetscInt               nrg;
2417   PetscInt               N,M,nrl,i,*idxs;
2418   PetscErrorCode         ierr;
2419 
2420   PetscFunctionBegin;
2421   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
2422   ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr);
2423   ierr = ISGetIndices(row,&rl);CHKERRQ(ierr);
2424   ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr);
2425 #if defined(PETSC_USE_DEBUG)
2426   for (i=0;i<nrl;i++) if (rl[i]>=nrg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row index %D -> %D greater than maximum possible %D",i,rl[i],nrg);
2427 #endif
2428   ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr);
2429   /* map from [0,nrl) to row */
2430   for (i=0;i<nrl;i++) idxs[i] = rl[i];
2431 #if defined(PETSC_USE_DEBUG)
2432   for (i=nrl;i<nrg;i++) idxs[i] = nrg;
2433 #else
2434   for (i=nrl;i<nrg;i++) idxs[i] = -1;
2435 #endif
2436   ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr);
2437   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
2438   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2439   ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
2440   ierr = ISDestroy(&is);CHKERRQ(ierr);
2441   /* compute new l2g map for columns */
2442   if (col != row || A->rmap->mapping != A->cmap->mapping) {
2443     const PetscInt *cg,*cl;
2444     PetscInt       ncg;
2445     PetscInt       ncl;
2446 
2447     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
2448     ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr);
2449     ierr = ISGetIndices(col,&cl);CHKERRQ(ierr);
2450     ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr);
2451 #if defined(PETSC_USE_DEBUG)
2452     for (i=0;i<ncl;i++) if (cl[i]>=ncg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local column index %D -> %D greater than maximum possible %D",i,cl[i],ncg);
2453 #endif
2454     ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr);
2455     /* map from [0,ncl) to col */
2456     for (i=0;i<ncl;i++) idxs[i] = cl[i];
2457 #if defined(PETSC_USE_DEBUG)
2458     for (i=ncl;i<ncg;i++) idxs[i] = ncg;
2459 #else
2460     for (i=ncl;i<ncg;i++) idxs[i] = -1;
2461 #endif
2462     ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr);
2463     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
2464     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2465     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
2466     ierr = ISDestroy(&is);CHKERRQ(ierr);
2467   } else {
2468     ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr);
2469     cl2g = rl2g;
2470   }
2471   /* create the MATIS submatrix */
2472   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
2473   ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr);
2474   ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
2475   ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr);
2476   matis = (Mat_IS*)((*submat)->data);
2477   matis->islocalref = PETSC_TRUE;
2478   ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr);
2479   ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr);
2480   ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr);
2481   ierr = MatSetUp(*submat);CHKERRQ(ierr);
2482   ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2483   ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2484   ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
2485   ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
2486   /* remove unsupported ops */
2487   ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
2488   (*submat)->ops->destroy               = MatDestroy_IS;
2489   (*submat)->ops->setvalueslocal        = MatSetValuesLocal_SubMat_IS;
2490   (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS;
2491   (*submat)->ops->assemblybegin         = MatAssemblyBegin_IS;
2492   (*submat)->ops->assemblyend           = MatAssemblyEnd_IS;
2493   PetscFunctionReturn(0);
2494 }
2495 
2496 static PetscErrorCode MatSetFromOptions_IS(PetscOptionItems *PetscOptionsObject, Mat A)
2497 {
2498   Mat_IS         *a = (Mat_IS*)A->data;
2499   PetscErrorCode ierr;
2500 
2501   PetscFunctionBegin;
2502   ierr = PetscOptionsHead(PetscOptionsObject,"MATIS options");CHKERRQ(ierr);
2503   ierr = PetscObjectOptionsBegin((PetscObject)A);
2504   ierr = PetscOptionsBool("-matis_fixempty","Fix local matrices in case of empty local rows/columns",NULL,a->locempty,&a->locempty,NULL);CHKERRQ(ierr);
2505   ierr = PetscOptionsBool("-matis_storel2l","Store local-to-local matrices generated from PtAP operations","MatISStoreL2L",a->storel2l,&a->storel2l,NULL);CHKERRQ(ierr);
2506   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2507   PetscFunctionReturn(0);
2508 }
2509 
2510 /*@
2511     MatCreateIS - Creates a "process" unassembled matrix, assembled on each
2512        process but not across processes.
2513 
2514    Input Parameters:
2515 +     comm    - MPI communicator that will share the matrix
2516 .     bs      - block size of the matrix
2517 .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
2518 .     rmap    - local to global map for rows
2519 -     cmap    - local to global map for cols
2520 
2521    Output Parameter:
2522 .    A - the resulting matrix
2523 
2524    Level: advanced
2525 
2526    Notes:
2527     See MATIS for more details.
2528           m and n are NOT related to the size of the map; they represent the size of the local parts of the vectors
2529           used in MatMult operations. The sizes of rmap and cmap define the size of the local matrices.
2530           If either rmap or cmap are NULL, then the matrix is assumed to be square.
2531 
2532 .seealso: MATIS, MatSetLocalToGlobalMapping()
2533 @*/
2534 PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
2535 {
2536   PetscErrorCode ierr;
2537 
2538   PetscFunctionBegin;
2539   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mappings");
2540   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2541   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
2542   if (bs > 0) {
2543     ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
2544   }
2545   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
2546   if (rmap && cmap) {
2547     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
2548   } else if (!rmap) {
2549     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
2550   } else {
2551     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
2552   }
2553   PetscFunctionReturn(0);
2554 }
2555 
2556 /*MC
2557    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP).
2558    This stores the matrices in globally unassembled form. Each processor
2559    assembles only its local Neumann problem and the parallel matrix vector
2560    product is handled "implicitly".
2561 
2562    Operations Provided:
2563 +  MatMult()
2564 .  MatMultAdd()
2565 .  MatMultTranspose()
2566 .  MatMultTransposeAdd()
2567 .  MatZeroEntries()
2568 .  MatSetOption()
2569 .  MatZeroRows()
2570 .  MatSetValues()
2571 .  MatSetValuesBlocked()
2572 .  MatSetValuesLocal()
2573 .  MatSetValuesBlockedLocal()
2574 .  MatScale()
2575 .  MatGetDiagonal()
2576 .  MatMissingDiagonal()
2577 .  MatDuplicate()
2578 .  MatCopy()
2579 .  MatAXPY()
2580 .  MatCreateSubMatrix()
2581 .  MatGetLocalSubMatrix()
2582 .  MatTranspose()
2583 .  MatPtAP() (with P of AIJ type)
2584 -  MatSetLocalToGlobalMapping()
2585 
2586    Options Database Keys:
2587 + -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
2588 . -matis_fixempty - Fixes local matrices in case of empty local rows/columns.
2589 - -matis_storel2l - stores the local-to-local operators generated by the Galerkin process of MatPtAP().
2590 
2591    Notes:
2592     Options prefix for the inner matrix are given by -is_mat_xxx
2593 
2594           You must call MatSetLocalToGlobalMapping() before using this matrix type.
2595 
2596           You can do matrix preallocation on the local matrix after you obtain it with
2597           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
2598 
2599   Level: advanced
2600 
2601 .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP
2602 
2603 M*/
2604 
2605 PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
2606 {
2607   PetscErrorCode ierr;
2608   Mat_IS         *b;
2609 
2610   PetscFunctionBegin;
2611   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
2612   A->data = (void*)b;
2613 
2614   /* matrix ops */
2615   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
2616   A->ops->mult                    = MatMult_IS;
2617   A->ops->multadd                 = MatMultAdd_IS;
2618   A->ops->multtranspose           = MatMultTranspose_IS;
2619   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
2620   A->ops->destroy                 = MatDestroy_IS;
2621   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
2622   A->ops->setvalues               = MatSetValues_IS;
2623   A->ops->setvaluesblocked        = MatSetValuesBlocked_IS;
2624   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
2625   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
2626   A->ops->zerorows                = MatZeroRows_IS;
2627   A->ops->zerorowscolumns         = MatZeroRowsColumns_IS;
2628   A->ops->assemblybegin           = MatAssemblyBegin_IS;
2629   A->ops->assemblyend             = MatAssemblyEnd_IS;
2630   A->ops->view                    = MatView_IS;
2631   A->ops->zeroentries             = MatZeroEntries_IS;
2632   A->ops->scale                   = MatScale_IS;
2633   A->ops->getdiagonal             = MatGetDiagonal_IS;
2634   A->ops->setoption               = MatSetOption_IS;
2635   A->ops->ishermitian             = MatIsHermitian_IS;
2636   A->ops->issymmetric             = MatIsSymmetric_IS;
2637   A->ops->isstructurallysymmetric = MatIsStructurallySymmetric_IS;
2638   A->ops->duplicate               = MatDuplicate_IS;
2639   A->ops->missingdiagonal         = MatMissingDiagonal_IS;
2640   A->ops->copy                    = MatCopy_IS;
2641   A->ops->getlocalsubmatrix       = MatGetLocalSubMatrix_IS;
2642   A->ops->createsubmatrix         = MatCreateSubMatrix_IS;
2643   A->ops->axpy                    = MatAXPY_IS;
2644   A->ops->diagonalset             = MatDiagonalSet_IS;
2645   A->ops->shift                   = MatShift_IS;
2646   A->ops->transpose               = MatTranspose_IS;
2647   A->ops->getinfo                 = MatGetInfo_IS;
2648   A->ops->diagonalscale           = MatDiagonalScale_IS;
2649   A->ops->setfromoptions          = MatSetFromOptions_IS;
2650 
2651   /* special MATIS functions */
2652   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
2653   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISRestoreLocalMat_C",MatISRestoreLocalMat_IS);CHKERRQ(ierr);
2654   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
2655   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
2656   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
2657   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",MatISSetUpSF_IS);CHKERRQ(ierr);
2658   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISStoreL2L_C",MatISStoreL2L_IS);CHKERRQ(ierr);
2659   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
2660   PetscFunctionReturn(0);
2661 }
2662