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