xref: /petsc/src/mat/impls/is/matis.c (revision 2e1947a5a08015068989a081e16c15ec463aefbc)
1 
2 /*
3     Creates a matrix class for using the Neumann-Neumann type preconditioners.
4    This stores the matrices in globally unassembled form. Each processor
5    assembles only its local Neumann problem and the parallel matrix vector
6    product is handled "implicitly".
7 
8      We provide:
9          MatMult()
10 
11     Currently this allows for only one subdomain per processor.
12 
13 */
14 
15 #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
16 #include <petscsf.h>
17 
18 #undef __FUNCT__
19 #define __FUNCT__ "MatISSetPreallocation"
20 PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
21 {
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
26   PetscValidType(B,1);
27   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
28   PetscFunctionReturn(0);
29 }
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "MatISSetPreallocation_IS"
33 PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
34 {
35   Mat_IS         *matis = (Mat_IS*)(B->data);
36   PetscSF        sf;
37   PetscInt       bs,i,nroots,*rootdata,nleaves,*leafdata,nlocalcols;
38   const PetscInt *gidxs;
39   PetscErrorCode ierr;
40 
41   PetscFunctionBegin;
42   if (!matis->A) {
43     SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
44   }
45   ierr = MatGetLocalSize(B,&nroots,NULL);CHKERRQ(ierr);
46   ierr = MatGetSize(matis->A,&nleaves,&nlocalcols);CHKERRQ(ierr);
47   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
48   ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr);
49   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&sf);CHKERRQ(ierr);
50   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
51   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr);
52   ierr = PetscSFSetGraphLayout(sf,B->rmap,nleaves,NULL,PETSC_COPY_VALUES,gidxs);CHKERRQ(ierr);
53   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr);
54   if (!d_nnz) {
55     for (i=0;i<nroots;i++) rootdata[i] += d_nz;
56   } else {
57     for (i=0;i<nroots;i++) rootdata[i] += d_nnz[i];
58   }
59   if (!o_nnz) {
60     for (i=0;i<nroots;i++) rootdata[i] += o_nz;
61   } else {
62     for (i=0;i<nroots;i++) rootdata[i] += o_nnz[i];
63   }
64   ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
65   ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
66   for (i=0;i<nleaves;i++) {
67     leafdata[i] = PetscMin(leafdata[i],nlocalcols);
68   }
69   ierr = MatSeqAIJSetPreallocation(matis->A,0,leafdata);CHKERRQ(ierr);
70   for (i=0;i<nleaves/bs;i++) {
71     leafdata[i] = leafdata[i*bs]/bs;
72   }
73   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr);
74   for (i=0;i<nleaves/bs;i++) {
75     leafdata[i] = leafdata[i]-i;
76   }
77   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr);
78   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
79   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
80   PetscFunctionReturn(0);
81 }
82 
83 #undef __FUNCT__
84 #define __FUNCT__ "MatISGetMPIXAIJ_IS"
85 PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
86 {
87   Mat_IS                 *matis = (Mat_IS*)(mat->data);
88   /* info on mat */
89   /* ISLocalToGlobalMapping rmapping,cmapping; */
90   PetscInt               bs,rows,cols;
91   PetscInt               lrows,lcols;
92   PetscInt               local_rows,local_cols;
93   PetscBool              isdense,issbaij,issbaij_red;
94   /* values insertion */
95   PetscScalar            *array;
96   PetscInt               *local_indices,*global_indices;
97   /* work */
98   PetscInt               i,j,index_row;
99   PetscErrorCode         ierr;
100 
101   PetscFunctionBegin;
102   /* MISSING CHECKS
103     - rectangular case not covered (it is not allowed by MATIS)
104   */
105   /* get info from mat */
106   /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */
107   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
108   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
109   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
110   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
111   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
112 
113   /* work */
114   ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr);
115   for (i=0;i<local_rows;i++) local_indices[i]=i;
116   /* map indices of local mat to global values */
117   ierr = PetscMalloc(PetscMax(local_cols,local_rows)*sizeof(*global_indices),&global_indices);CHKERRQ(ierr);
118   /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */
119   ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr);
120 
121   if (issbaij) {
122     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
123   }
124 
125   if (reuse == MAT_INITIAL_MATRIX) {
126     Mat         new_mat;
127     MatType     new_mat_type;
128     Vec         vec_dnz,vec_onz;
129     PetscScalar *my_dnz,*my_onz;
130     PetscInt    *dnz,*onz,*mat_ranges,*row_ownership;
131     PetscInt    index_col,owner;
132     PetscMPIInt nsubdomains;
133 
134     /* determining new matrix type */
135     ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
136     if (issbaij_red) {
137       new_mat_type = MATSBAIJ;
138     } else {
139       if (bs>1) {
140         new_mat_type = MATBAIJ;
141       } else {
142         new_mat_type = MATAIJ;
143       }
144     }
145 
146     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
147     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr);
148     ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr);
149     ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr);
150     ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr);
151     ierr = MatSetUp(new_mat);CHKERRQ(ierr);
152     ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr);
153 
154     /*
155       preallocation
156     */
157 
158     ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr);
159     /*
160        Some vectors are needed to sum up properly on shared interface dofs.
161        Preallocation macros cannot do the job.
162        Note that preallocation is not exact, since it overestimates nonzeros
163     */
164     ierr = MatGetVecs(new_mat,NULL,&vec_dnz);CHKERRQ(ierr);
165     /* ierr = VecSetLocalToGlobalMapping(vec_dnz,rmapping);CHKERRQ(ierr); */
166     ierr = VecSetLocalToGlobalMapping(vec_dnz,matis->mapping);CHKERRQ(ierr);
167     ierr = VecDuplicate(vec_dnz,&vec_onz);CHKERRQ(ierr);
168     /* All processes need to compute entire row ownership */
169     ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
170     ierr = MatGetOwnershipRanges(new_mat,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
171     for (i=0;i<nsubdomains;i++) {
172       for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
173         row_ownership[j]=i;
174       }
175     }
176 
177     /*
178        my_dnz and my_onz contains exact contribution to preallocation from each local mat
179        then, they will be summed up properly. This way, preallocation is always sufficient
180     */
181     ierr = PetscMalloc1(local_rows,&my_dnz);CHKERRQ(ierr);
182     ierr = PetscMalloc1(local_rows,&my_onz);CHKERRQ(ierr);
183     ierr = PetscMemzero(my_dnz,local_rows*sizeof(*my_dnz));CHKERRQ(ierr);
184     ierr = PetscMemzero(my_onz,local_rows*sizeof(*my_onz));CHKERRQ(ierr);
185     /* preallocation as a MATAIJ */
186     if (isdense) { /* special case for dense local matrices */
187       for (i=0;i<local_rows;i++) {
188         index_row = global_indices[i];
189         for (j=i;j<local_rows;j++) {
190           owner = row_ownership[index_row];
191           index_col = global_indices[j];
192           if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
193             my_dnz[i] += 1.0;
194           } else { /* offdiag block */
195             my_onz[i] += 1.0;
196           }
197           /* same as before, interchanging rows and cols */
198           if (i != j) {
199             owner = row_ownership[index_col];
200             if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
201               my_dnz[j] += 1.0;
202             } else {
203               my_onz[j] += 1.0;
204             }
205           }
206         }
207       }
208     } else {
209       for (i=0;i<local_rows;i++) {
210         PetscInt ncols;
211         const PetscInt *cols;
212         index_row = global_indices[i];
213         ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
214         for (j=0;j<ncols;j++) {
215           owner = row_ownership[index_row];
216           index_col = global_indices[cols[j]];
217           if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
218             my_dnz[i] += 1.0;
219           } else { /* offdiag block */
220             my_onz[i] += 1.0;
221           }
222           /* same as before, interchanging rows and cols */
223           if (issbaij) {
224             owner = row_ownership[index_col];
225             if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
226               my_dnz[j] += 1.0;
227             } else {
228               my_onz[j] += 1.0;
229             }
230           }
231         }
232         ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
233       }
234     }
235     ierr = VecSet(vec_dnz,0.0);CHKERRQ(ierr);
236     ierr = VecSet(vec_onz,0.0);CHKERRQ(ierr);
237     if (local_rows) { /* multilevel guard */
238       ierr = VecSetValuesLocal(vec_dnz,local_rows,local_indices,my_dnz,ADD_VALUES);CHKERRQ(ierr);
239       ierr = VecSetValuesLocal(vec_onz,local_rows,local_indices,my_onz,ADD_VALUES);CHKERRQ(ierr);
240     }
241     ierr = VecAssemblyBegin(vec_dnz);CHKERRQ(ierr);
242     ierr = VecAssemblyBegin(vec_onz);CHKERRQ(ierr);
243     ierr = VecAssemblyEnd(vec_dnz);CHKERRQ(ierr);
244     ierr = VecAssemblyEnd(vec_onz);CHKERRQ(ierr);
245     ierr = PetscFree(my_dnz);CHKERRQ(ierr);
246     ierr = PetscFree(my_onz);CHKERRQ(ierr);
247     ierr = PetscFree(row_ownership);CHKERRQ(ierr);
248 
249     /* set computed preallocation in dnz and onz */
250     ierr = VecGetArray(vec_dnz,&array);CHKERRQ(ierr);
251     for (i=0; i<lrows; i++) dnz[i] = (PetscInt)PetscRealPart(array[i]);
252     ierr = VecRestoreArray(vec_dnz,&array);CHKERRQ(ierr);
253     ierr = VecGetArray(vec_onz,&array);CHKERRQ(ierr);
254     for (i=0;i<lrows;i++) onz[i] = (PetscInt)PetscRealPart(array[i]);
255     ierr = VecRestoreArray(vec_onz,&array);CHKERRQ(ierr);
256     ierr = VecDestroy(&vec_dnz);CHKERRQ(ierr);
257     ierr = VecDestroy(&vec_onz);CHKERRQ(ierr);
258 
259     /* Resize preallocation if overestimated */
260     for (i=0;i<lrows;i++) {
261       dnz[i] = PetscMin(dnz[i],lcols);
262       onz[i] = PetscMin(onz[i],cols-lcols);
263     }
264     /* set preallocation */
265     ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr);
266     for (i=0;i<lrows/bs;i++) {
267       dnz[i] = dnz[i*bs]/bs;
268       onz[i] = onz[i*bs]/bs;
269     }
270     ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr);
271     for (i=0;i<lrows/bs;i++) {
272       dnz[i] = dnz[i]-i;
273     }
274     ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr);
275     ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
276     *M = new_mat;
277   } else {
278     PetscInt mbs,mrows,mcols;
279     /* some checks */
280     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
281     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
282     if (mrows != rows) {
283       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
284     }
285     if (mrows != rows) {
286       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
287     }
288     if (mbs != bs) {
289       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
290     }
291     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
292   }
293   /* set local to global mappings */
294   /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */
295   /* Set values */
296   if (isdense) { /* special case for dense local matrices */
297     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
298     ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr);
299     ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
300     ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr);
301     ierr = PetscFree(local_indices);CHKERRQ(ierr);
302     ierr = PetscFree(global_indices);CHKERRQ(ierr);
303   } else { /* very basic values insertion for all other matrix types */
304     ierr = PetscFree(local_indices);CHKERRQ(ierr);
305     for (i=0;i<local_rows;i++) {
306       ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
307       /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */
308       ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr);
309       ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr);
310       ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
311       ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
312     }
313     ierr = PetscFree(global_indices);CHKERRQ(ierr);
314   }
315   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
316   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
317   if (isdense) {
318     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
319   }
320   if (issbaij) {
321     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
322   }
323   PetscFunctionReturn(0);
324 }
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "MatISGetMPIXAIJ"
328 /*@
329     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
330 
331   Input Parameter:
332 .  mat - the matrix (should be of type MATIS)
333 .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
334 
335   Output Parameter:
336 .  newmat - the matrix in AIJ format
337 
338   Level: developer
339 
340   Notes:
341 
342 .seealso: MATIS
343 @*/
344 PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
345 {
346   PetscErrorCode ierr;
347 
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
350   PetscValidLogicalCollectiveEnum(mat,reuse,2);
351   PetscValidPointer(newmat,3);
352   if (reuse != MAT_INITIAL_MATRIX) {
353     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
354     PetscCheckSameComm(mat,1,*newmat,3);
355   }
356   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
357   PetscFunctionReturn(0);
358 }
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "MatDuplicate_IS"
362 PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
363 {
364   PetscErrorCode ierr;
365   Mat_IS         *matis = (Mat_IS*)(mat->data);
366   PetscInt       bs,m,n,M,N;
367   Mat            B,localmat;
368 
369   PetscFunctionBegin;
370   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
371   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
372   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
373   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr);
374   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
375   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
376   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
377   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
378   *newmat = B;
379   PetscFunctionReturn(0);
380 }
381 
382 #undef __FUNCT__
383 #define __FUNCT__ "MatIsHermitian_IS"
384 PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
385 {
386   PetscErrorCode ierr;
387   Mat_IS         *matis = (Mat_IS*)A->data;
388   PetscBool      local_sym;
389 
390   PetscFunctionBegin;
391   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
392   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
393   PetscFunctionReturn(0);
394 }
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "MatIsSymmetric_IS"
398 PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
399 {
400   PetscErrorCode ierr;
401   Mat_IS         *matis = (Mat_IS*)A->data;
402   PetscBool      local_sym;
403 
404   PetscFunctionBegin;
405   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
406   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
407   PetscFunctionReturn(0);
408 }
409 
410 #undef __FUNCT__
411 #define __FUNCT__ "MatDestroy_IS"
412 PetscErrorCode MatDestroy_IS(Mat A)
413 {
414   PetscErrorCode ierr;
415   Mat_IS         *b = (Mat_IS*)A->data;
416 
417   PetscFunctionBegin;
418   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
419   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
420   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
421   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
422   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
423   ierr = PetscFree(A->data);CHKERRQ(ierr);
424   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
425   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
426   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
427   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
428   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
429   PetscFunctionReturn(0);
430 }
431 
432 #undef __FUNCT__
433 #define __FUNCT__ "MatMult_IS"
434 PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
435 {
436   PetscErrorCode ierr;
437   Mat_IS         *is  = (Mat_IS*)A->data;
438   PetscScalar    zero = 0.0;
439 
440   PetscFunctionBegin;
441   /*  scatter the global vector x into the local work vector */
442   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
443   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
444 
445   /* multiply the local matrix */
446   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
447 
448   /* scatter product back into global memory */
449   ierr = VecSet(y,zero);CHKERRQ(ierr);
450   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
451   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
452   PetscFunctionReturn(0);
453 }
454 
455 #undef __FUNCT__
456 #define __FUNCT__ "MatMultAdd_IS"
457 PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
458 {
459   Vec            temp_vec;
460   PetscErrorCode ierr;
461 
462   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
463   if (v3 != v2) {
464     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
465     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
466   } else {
467     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
468     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
469     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
470     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
471     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
472   }
473   PetscFunctionReturn(0);
474 }
475 
476 #undef __FUNCT__
477 #define __FUNCT__ "MatMultTranspose_IS"
478 PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
479 {
480   Mat_IS         *is = (Mat_IS*)A->data;
481   PetscErrorCode ierr;
482 
483   PetscFunctionBegin; /*  y = A' * x */
484   /*  scatter the global vector x into the local work vector */
485   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
486   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
487 
488   /* multiply the local matrix */
489   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
490 
491   /* scatter product back into global vector */
492   ierr = VecSet(y,0);CHKERRQ(ierr);
493   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
494   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
495   PetscFunctionReturn(0);
496 }
497 
498 #undef __FUNCT__
499 #define __FUNCT__ "MatMultTransposeAdd_IS"
500 PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
501 {
502   Vec            temp_vec;
503   PetscErrorCode ierr;
504 
505   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
506   if (v3 != v2) {
507     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
508     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
509   } else {
510     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
511     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
512     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
513     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
514     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
515   }
516   PetscFunctionReturn(0);
517 }
518 
519 #undef __FUNCT__
520 #define __FUNCT__ "MatView_IS"
521 PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
522 {
523   Mat_IS         *a = (Mat_IS*)A->data;
524   PetscErrorCode ierr;
525   PetscViewer    sviewer;
526 
527   PetscFunctionBegin;
528   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
529   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
530   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
531   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
532   PetscFunctionReturn(0);
533 }
534 
535 #undef __FUNCT__
536 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
537 PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
538 {
539   PetscErrorCode ierr;
540   PetscInt       n,bs;
541   Mat_IS         *is = (Mat_IS*)A->data;
542   IS             from,to;
543   Vec            global;
544 
545   PetscFunctionBegin;
546   PetscCheckSameComm(A,1,rmapping,2);
547   if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
548   if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */
549     ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
550     ierr = VecDestroy(&is->x);CHKERRQ(ierr);
551     ierr = VecDestroy(&is->y);CHKERRQ(ierr);
552     ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr);
553     ierr = MatDestroy(&is->A);CHKERRQ(ierr);
554   }
555   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
556   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
557   is->mapping = rmapping;
558 /*
559   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
560   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
561 */
562 
563   /* Create the local matrix A */
564   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
565   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr);
566   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
567   if (bs > 1) {
568     ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr);
569   } else {
570     ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr);
571   }
572   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
573   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
574   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
575   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
576   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
577 
578   /* Create the local work vectors */
579   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
580   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
581   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
582   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
583   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
584   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
585   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
586 
587   /* setup the global to local scatter */
588   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
589   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
590   ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr);
591   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
592   ierr = VecDestroy(&global);CHKERRQ(ierr);
593   ierr = ISDestroy(&to);CHKERRQ(ierr);
594   ierr = ISDestroy(&from);CHKERRQ(ierr);
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "MatSetValues_IS"
600 PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
601 {
602   Mat_IS         *is = (Mat_IS*)mat->data;
603   PetscInt       rows_l[2048],cols_l[2048];
604   PetscErrorCode ierr;
605 
606   PetscFunctionBegin;
607 #if defined(PETSC_USE_DEBUG)
608   if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
609 #endif
610   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
611   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
612   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
613   PetscFunctionReturn(0);
614 }
615 
616 #undef ISG2LMapSetUp
617 #undef ISG2LMapApply
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "MatSetValuesLocal_IS"
621 PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
622 {
623   PetscErrorCode ierr;
624   Mat_IS         *is = (Mat_IS*)A->data;
625 
626   PetscFunctionBegin;
627   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
628   PetscFunctionReturn(0);
629 }
630 
631 #undef __FUNCT__
632 #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
633 PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
634 {
635   PetscErrorCode ierr;
636   Mat_IS         *is = (Mat_IS*)A->data;
637 
638   PetscFunctionBegin;
639   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
640   PetscFunctionReturn(0);
641 }
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "MatZeroRows_IS"
645 PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
646 {
647   Mat_IS         *is = (Mat_IS*)A->data;
648   PetscInt       n_l = 0, *rows_l = NULL;
649   PetscErrorCode ierr;
650 
651   PetscFunctionBegin;
652   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
653   if (n) {
654     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
655     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
656   }
657   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
658   ierr = PetscFree(rows_l);CHKERRQ(ierr);
659   PetscFunctionReturn(0);
660 }
661 
662 #undef __FUNCT__
663 #define __FUNCT__ "MatZeroRowsLocal_IS"
664 PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
665 {
666   Mat_IS         *is = (Mat_IS*)A->data;
667   PetscErrorCode ierr;
668   PetscInt       i;
669   PetscScalar    *array;
670 
671   PetscFunctionBegin;
672   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
673   {
674     /*
675        Set up is->x as a "counting vector". This is in order to MatMult_IS
676        work properly in the interface nodes.
677     */
678     Vec         counter;
679     PetscScalar one=1.0, zero=0.0;
680     ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr);
681     ierr = VecSet(counter,zero);CHKERRQ(ierr);
682     ierr = VecSet(is->x,one);CHKERRQ(ierr);
683     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
684     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
685     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
686     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
687     ierr = VecDestroy(&counter);CHKERRQ(ierr);
688   }
689   if (!n) {
690     is->pure_neumann = PETSC_TRUE;
691   } else {
692     is->pure_neumann = PETSC_FALSE;
693 
694     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
695     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
696     for (i=0; i<n; i++) {
697       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
698     }
699     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
700     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
701     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
702   }
703   PetscFunctionReturn(0);
704 }
705 
706 #undef __FUNCT__
707 #define __FUNCT__ "MatAssemblyBegin_IS"
708 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
709 {
710   Mat_IS         *is = (Mat_IS*)A->data;
711   PetscErrorCode ierr;
712 
713   PetscFunctionBegin;
714   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
715   PetscFunctionReturn(0);
716 }
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "MatAssemblyEnd_IS"
720 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
721 {
722   Mat_IS         *is = (Mat_IS*)A->data;
723   PetscErrorCode ierr;
724 
725   PetscFunctionBegin;
726   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
727   PetscFunctionReturn(0);
728 }
729 
730 #undef __FUNCT__
731 #define __FUNCT__ "MatISGetLocalMat_IS"
732 PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
733 {
734   Mat_IS *is = (Mat_IS*)mat->data;
735 
736   PetscFunctionBegin;
737   *local = is->A;
738   PetscFunctionReturn(0);
739 }
740 
741 #undef __FUNCT__
742 #define __FUNCT__ "MatISGetLocalMat"
743 /*@
744     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
745 
746   Input Parameter:
747 .  mat - the matrix
748 
749   Output Parameter:
750 .  local - the local matrix usually MATSEQAIJ
751 
752   Level: advanced
753 
754   Notes:
755     This can be called if you have precomputed the nonzero structure of the
756   matrix and want to provide it to the inner matrix object to improve the performance
757   of the MatSetValues() operation.
758 
759 .seealso: MATIS
760 @*/
761 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
762 {
763   PetscErrorCode ierr;
764 
765   PetscFunctionBegin;
766   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
767   PetscValidPointer(local,2);
768   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
769   PetscFunctionReturn(0);
770 }
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "MatISSetLocalMat_IS"
774 PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
775 {
776   Mat_IS         *is = (Mat_IS*)mat->data;
777   PetscInt       nrows,ncols,orows,ocols;
778   PetscErrorCode ierr;
779 
780   PetscFunctionBegin;
781   if (is->A) {
782     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
783     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
784     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)\n",orows,ocols,nrows,ncols);
785   }
786   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
787   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
788   is->A = local;
789   PetscFunctionReturn(0);
790 }
791 
792 #undef __FUNCT__
793 #define __FUNCT__ "MatISSetLocalMat"
794 /*@
795     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
796 
797   Input Parameter:
798 .  mat - the matrix
799 .  local - the local matrix usually MATSEQAIJ
800 
801   Output Parameter:
802 
803   Level: advanced
804 
805   Notes:
806     This can be called if you have precomputed the local matrix and
807   want to provide it to the matrix object MATIS.
808 
809 .seealso: MATIS
810 @*/
811 PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
812 {
813   PetscErrorCode ierr;
814 
815   PetscFunctionBegin;
816   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
817   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
818   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
819   PetscFunctionReturn(0);
820 }
821 
822 #undef __FUNCT__
823 #define __FUNCT__ "MatZeroEntries_IS"
824 PetscErrorCode MatZeroEntries_IS(Mat A)
825 {
826   Mat_IS         *a = (Mat_IS*)A->data;
827   PetscErrorCode ierr;
828 
829   PetscFunctionBegin;
830   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
831   PetscFunctionReturn(0);
832 }
833 
834 #undef __FUNCT__
835 #define __FUNCT__ "MatScale_IS"
836 PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
837 {
838   Mat_IS         *is = (Mat_IS*)A->data;
839   PetscErrorCode ierr;
840 
841   PetscFunctionBegin;
842   ierr = MatScale(is->A,a);CHKERRQ(ierr);
843   PetscFunctionReturn(0);
844 }
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "MatGetDiagonal_IS"
848 PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
849 {
850   Mat_IS         *is = (Mat_IS*)A->data;
851   PetscErrorCode ierr;
852 
853   PetscFunctionBegin;
854   /* get diagonal of the local matrix */
855   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
856 
857   /* scatter diagonal back into global vector */
858   ierr = VecSet(v,0);CHKERRQ(ierr);
859   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
860   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
861   PetscFunctionReturn(0);
862 }
863 
864 #undef __FUNCT__
865 #define __FUNCT__ "MatSetOption_IS"
866 PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
867 {
868   Mat_IS         *a = (Mat_IS*)A->data;
869   PetscErrorCode ierr;
870 
871   PetscFunctionBegin;
872   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
873   PetscFunctionReturn(0);
874 }
875 
876 #undef __FUNCT__
877 #define __FUNCT__ "MatCreateIS"
878 /*@
879     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
880        process but not across processes.
881 
882    Input Parameters:
883 +     comm - MPI communicator that will share the matrix
884 .     bs - local and global block size of the matrix
885 .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
886 -     map - mapping that defines the global number for each local number
887 
888    Output Parameter:
889 .    A - the resulting matrix
890 
891    Level: advanced
892 
893    Notes: See MATIS for more details
894           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
895           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
896           plus the ghost points to global indices.
897 
898 .seealso: MATIS, MatSetLocalToGlobalMapping()
899 @*/
900 PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
901 {
902   PetscErrorCode ierr;
903 
904   PetscFunctionBegin;
905   ierr = MatCreate(comm,A);CHKERRQ(ierr);
906   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
907   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
908   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
909   ierr = MatSetUp(*A);CHKERRQ(ierr);
910   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
911   PetscFunctionReturn(0);
912 }
913 
914 /*MC
915    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC.
916    This stores the matrices in globally unassembled form. Each processor
917    assembles only its local Neumann problem and the parallel matrix vector
918    product is handled "implicitly".
919 
920    Operations Provided:
921 +  MatMult()
922 .  MatMultAdd()
923 .  MatMultTranspose()
924 .  MatMultTransposeAdd()
925 .  MatZeroEntries()
926 .  MatSetOption()
927 .  MatZeroRows()
928 .  MatZeroRowsLocal()
929 .  MatSetValues()
930 .  MatSetValuesLocal()
931 .  MatScale()
932 .  MatGetDiagonal()
933 -  MatSetLocalToGlobalMapping()
934 
935    Options Database Keys:
936 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
937 
938    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
939 
940           You must call MatSetLocalToGlobalMapping() before using this matrix type.
941 
942           You can do matrix preallocation on the local matrix after you obtain it with
943           MatISGetLocalMat()
944 
945   Level: advanced
946 
947 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC
948 
949 M*/
950 
951 #undef __FUNCT__
952 #define __FUNCT__ "MatCreate_IS"
953 PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
954 {
955   PetscErrorCode ierr;
956   Mat_IS         *b;
957 
958   PetscFunctionBegin;
959   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
960   A->data = (void*)b;
961   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
962 
963   A->ops->mult                    = MatMult_IS;
964   A->ops->multadd                 = MatMultAdd_IS;
965   A->ops->multtranspose           = MatMultTranspose_IS;
966   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
967   A->ops->destroy                 = MatDestroy_IS;
968   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
969   A->ops->setvalues               = MatSetValues_IS;
970   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
971   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
972   A->ops->zerorows                = MatZeroRows_IS;
973   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
974   A->ops->assemblybegin           = MatAssemblyBegin_IS;
975   A->ops->assemblyend             = MatAssemblyEnd_IS;
976   A->ops->view                    = MatView_IS;
977   A->ops->zeroentries             = MatZeroEntries_IS;
978   A->ops->scale                   = MatScale_IS;
979   A->ops->getdiagonal             = MatGetDiagonal_IS;
980   A->ops->setoption               = MatSetOption_IS;
981   A->ops->ishermitian             = MatIsHermitian_IS;
982   A->ops->issymmetric             = MatIsSymmetric_IS;
983   A->ops->duplicate               = MatDuplicate_IS;
984 
985   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
986   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
987 
988   /* zeros pointer */
989   b->A       = 0;
990   b->ctx     = 0;
991   b->x       = 0;
992   b->y       = 0;
993   b->mapping = 0;
994 
995   /* special MATIS functions */
996   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
997   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
998   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
999   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
1000   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1001   PetscFunctionReturn(0);
1002 }
1003