xref: /petsc/src/mat/impls/is/matis.c (revision f347579b2c740aaf2706e3a53ee08f00957c2ec2)
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 
17 #undef __FUNCT__
18 #define __FUNCT__ "MatIsHermitian_IS"
19 PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
20 {
21   PetscErrorCode ierr;
22   Mat_IS         *matis = (Mat_IS*)A->data;
23   PetscBool      local_sym;
24 
25   PetscFunctionBegin;
26   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
27   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
28   PetscFunctionReturn(0);
29 }
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "MatIsSymmetric_IS"
33 PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
34 {
35   PetscErrorCode ierr;
36   Mat_IS         *matis = (Mat_IS*)A->data;
37   PetscBool      local_sym;
38 
39   PetscFunctionBegin;
40   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
41   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "MatDestroy_IS"
47 PetscErrorCode MatDestroy_IS(Mat A)
48 {
49   PetscErrorCode ierr;
50   Mat_IS         *b = (Mat_IS*)A->data;
51 
52   PetscFunctionBegin;
53   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
54   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
55   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
56   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
57   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
58   ierr = PetscFree(A->data);CHKERRQ(ierr);
59   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
60   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "MatMult_IS"
66 PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
67 {
68   PetscErrorCode ierr;
69   Mat_IS         *is  = (Mat_IS*)A->data;
70   PetscScalar    zero = 0.0;
71 
72   PetscFunctionBegin;
73   /*  scatter the global vector x into the local work vector */
74   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
75   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
76 
77   /* multiply the local matrix */
78   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
79 
80   /* scatter product back into global memory */
81   ierr = VecSet(y,zero);CHKERRQ(ierr);
82   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
83   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
84   PetscFunctionReturn(0);
85 }
86 
87 #undef __FUNCT__
88 #define __FUNCT__ "MatMultAdd_IS"
89 static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
90 {
91   Vec            temp_vec;
92   PetscErrorCode ierr;
93 
94   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
95   if (v3 != v2) {
96     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
97     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
98   } else {
99     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
100     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
101     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
102     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
103     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
104   }
105   PetscFunctionReturn(0);
106 }
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "MatMultTranspose_IS"
110 PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
111 {
112   Mat_IS         *is = (Mat_IS*)A->data;
113   PetscErrorCode ierr;
114 
115   PetscFunctionBegin; /*  y = A' * x */
116   /*  scatter the global vector x into the local work vector */
117   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119 
120   /* multiply the local matrix */
121   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
122 
123   /* scatter product back into global vector */
124   ierr = VecSet(y,0);CHKERRQ(ierr);
125   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
126   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
127   PetscFunctionReturn(0);
128 }
129 
130 #undef __FUNCT__
131 #define __FUNCT__ "MatMultTransposeAdd_IS"
132 PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
133 {
134   Vec            temp_vec;
135   PetscErrorCode ierr;
136 
137   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
138   if (v3 != v2) {
139     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
140     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
141   } else {
142     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
143     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
144     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
145     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
146     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
147   }
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "MatView_IS"
153 PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
154 {
155   Mat_IS         *a = (Mat_IS*)A->data;
156   PetscErrorCode ierr;
157   PetscViewer    sviewer;
158 
159   PetscFunctionBegin;
160   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
161   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
162   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
168 PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
169 {
170   PetscErrorCode ierr;
171   PetscInt       n,bs;
172   Mat_IS         *is = (Mat_IS*)A->data;
173   IS             from,to;
174   Vec            global;
175 
176   PetscFunctionBegin;
177   PetscCheckSameComm(A,1,rmapping,2);
178   if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
179   if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */
180     ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
181     ierr = VecDestroy(&is->x);CHKERRQ(ierr);
182     ierr = VecDestroy(&is->y);CHKERRQ(ierr);
183     ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr);
184     ierr = MatDestroy(&is->A);CHKERRQ(ierr);
185   }
186   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
187   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
188   is->mapping = rmapping;
189 /*
190   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
191   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
192 */
193 
194   /* Create the local matrix A */
195   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
196   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
197   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
198   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
199   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
200   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
201   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
202   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
203 
204   /* Create the local work vectors */
205   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
206   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
207   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
208   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
209   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
210   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
211   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
212 
213   /* setup the global to local scatter */
214   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
215   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
216   ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr);
217   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
218   ierr = VecDestroy(&global);CHKERRQ(ierr);
219   ierr = ISDestroy(&to);CHKERRQ(ierr);
220   ierr = ISDestroy(&from);CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 #undef __FUNCT__
225 #define __FUNCT__ "MatSetValues_IS"
226 PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
227 {
228   Mat_IS         *is = (Mat_IS*)mat->data;
229   PetscInt       rows_l[2048],cols_l[2048];
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233 #if defined(PETSC_USE_DEBUG)
234   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);
235 #endif
236   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
237   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
238   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
239   PetscFunctionReturn(0);
240 }
241 
242 #undef ISG2LMapSetUp
243 #undef ISG2LMapApply
244 
245 #undef __FUNCT__
246 #define __FUNCT__ "MatSetValuesLocal_IS"
247 PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
248 {
249   PetscErrorCode ierr;
250   Mat_IS         *is = (Mat_IS*)A->data;
251 
252   PetscFunctionBegin;
253   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
254   PetscFunctionReturn(0);
255 }
256 
257 #undef __FUNCT__
258 #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
259 PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
260 {
261   PetscErrorCode ierr;
262   Mat_IS         *is = (Mat_IS*)A->data;
263 
264   PetscFunctionBegin;
265   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "MatZeroRows_IS"
271 PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
272 {
273   Mat_IS         *is = (Mat_IS*)A->data;
274   PetscInt       n_l = 0, *rows_l = NULL;
275   PetscErrorCode ierr;
276 
277   PetscFunctionBegin;
278   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
279   if (n) {
280     ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr);
281     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
282   }
283   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
284   ierr = PetscFree(rows_l);CHKERRQ(ierr);
285   PetscFunctionReturn(0);
286 }
287 
288 #undef __FUNCT__
289 #define __FUNCT__ "MatZeroRowsLocal_IS"
290 PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
291 {
292   Mat_IS         *is = (Mat_IS*)A->data;
293   PetscErrorCode ierr;
294   PetscInt       i;
295   PetscScalar    *array;
296 
297   PetscFunctionBegin;
298   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
299   {
300     /*
301        Set up is->x as a "counting vector". This is in order to MatMult_IS
302        work properly in the interface nodes.
303     */
304     Vec         counter;
305     PetscScalar one=1.0, zero=0.0;
306     ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr);
307     ierr = VecSet(counter,zero);CHKERRQ(ierr);
308     ierr = VecSet(is->x,one);CHKERRQ(ierr);
309     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
310     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
311     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
312     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
313     ierr = VecDestroy(&counter);CHKERRQ(ierr);
314   }
315   if (!n) {
316     is->pure_neumann = PETSC_TRUE;
317   } else {
318     is->pure_neumann = PETSC_FALSE;
319 
320     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
321     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
322     for (i=0; i<n; i++) {
323       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
324     }
325     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
326     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
327     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
328   }
329   PetscFunctionReturn(0);
330 }
331 
332 #undef __FUNCT__
333 #define __FUNCT__ "MatAssemblyBegin_IS"
334 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
335 {
336   Mat_IS         *is = (Mat_IS*)A->data;
337   PetscErrorCode ierr;
338 
339   PetscFunctionBegin;
340   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
341   PetscFunctionReturn(0);
342 }
343 
344 #undef __FUNCT__
345 #define __FUNCT__ "MatAssemblyEnd_IS"
346 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
347 {
348   Mat_IS         *is = (Mat_IS*)A->data;
349   PetscErrorCode ierr;
350 
351   PetscFunctionBegin;
352   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
353   PetscFunctionReturn(0);
354 }
355 
356 #undef __FUNCT__
357 #define __FUNCT__ "MatISGetLocalMat_IS"
358 PetscErrorCode  MatISGetLocalMat_IS(Mat mat,Mat *local)
359 {
360   Mat_IS *is = (Mat_IS*)mat->data;
361 
362   PetscFunctionBegin;
363   *local = is->A;
364   PetscFunctionReturn(0);
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "MatISGetLocalMat"
369 /*@
370     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
371 
372   Input Parameter:
373 .  mat - the matrix
374 
375   Output Parameter:
376 .  local - the local matrix usually MATSEQAIJ
377 
378   Level: advanced
379 
380   Notes:
381     This can be called if you have precomputed the nonzero structure of the
382   matrix and want to provide it to the inner matrix object to improve the performance
383   of the MatSetValues() operation.
384 
385 .seealso: MATIS
386 @*/
387 PetscErrorCode  MatISGetLocalMat(Mat mat,Mat *local)
388 {
389   PetscErrorCode ierr;
390 
391   PetscFunctionBegin;
392   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
393   PetscValidPointer(local,2);
394   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNCT__
399 #define __FUNCT__ "MatISSetLocalMat_IS"
400 PetscErrorCode  MatISSetLocalMat_IS(Mat mat,Mat local)
401 {
402   Mat_IS         *is = (Mat_IS*)mat->data;
403   PetscInt       nrows,ncols,orows,ocols;
404   PetscErrorCode ierr;
405 
406   PetscFunctionBegin;
407   if (is->A) {
408     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
409     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
410     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);
411   }
412   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
413   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
414   is->A = local;
415   PetscFunctionReturn(0);
416 }
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "MatISSetLocalMat"
420 /*@
421     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
422 
423   Input Parameter:
424 .  mat - the matrix
425 .  local - the local matrix usually MATSEQAIJ
426 
427   Output Parameter:
428 
429   Level: advanced
430 
431   Notes:
432     This can be called if you have precomputed the local matrix and
433   want to provide it to the matrix object MATIS.
434 
435 .seealso: MATIS
436 @*/
437 PetscErrorCode  MatISSetLocalMat(Mat mat,Mat local)
438 {
439   PetscErrorCode ierr;
440 
441   PetscFunctionBegin;
442   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
443   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
444   PetscFunctionReturn(0);
445 }
446 
447 #undef __FUNCT__
448 #define __FUNCT__ "MatZeroEntries_IS"
449 PetscErrorCode MatZeroEntries_IS(Mat A)
450 {
451   Mat_IS         *a = (Mat_IS*)A->data;
452   PetscErrorCode ierr;
453 
454   PetscFunctionBegin;
455   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
456   PetscFunctionReturn(0);
457 }
458 
459 #undef __FUNCT__
460 #define __FUNCT__ "MatScale_IS"
461 PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
462 {
463   Mat_IS         *is = (Mat_IS*)A->data;
464   PetscErrorCode ierr;
465 
466   PetscFunctionBegin;
467   ierr = MatScale(is->A,a);CHKERRQ(ierr);
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "MatGetDiagonal_IS"
473 PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
474 {
475   Mat_IS         *is = (Mat_IS*)A->data;
476   PetscErrorCode ierr;
477 
478   PetscFunctionBegin;
479   /* get diagonal of the local matrix */
480   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
481 
482   /* scatter diagonal back into global vector */
483   ierr = VecSet(v,0);CHKERRQ(ierr);
484   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
485   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
486   PetscFunctionReturn(0);
487 }
488 
489 #undef __FUNCT__
490 #define __FUNCT__ "MatSetOption_IS"
491 PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
492 {
493   Mat_IS         *a = (Mat_IS*)A->data;
494   PetscErrorCode ierr;
495 
496   PetscFunctionBegin;
497   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 
501 #undef __FUNCT__
502 #define __FUNCT__ "MatCreateIS"
503 /*@
504     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
505        process but not across processes.
506 
507    Input Parameters:
508 +     comm - MPI communicator that will share the matrix
509 .     bs - local and global block size of the matrix
510 .     m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
511 -     map - mapping that defines the global number for each local number
512 
513    Output Parameter:
514 .    A - the resulting matrix
515 
516    Level: advanced
517 
518    Notes: See MATIS for more details
519           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
520           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
521           plus the ghost points to global indices.
522 
523 .seealso: MATIS, MatSetLocalToGlobalMapping()
524 @*/
525 PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
526 {
527   PetscErrorCode ierr;
528 
529   PetscFunctionBegin;
530   ierr = MatCreate(comm,A);CHKERRQ(ierr);
531   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
532   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
533   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
534   ierr = MatSetUp(*A);CHKERRQ(ierr);
535   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
536   PetscFunctionReturn(0);
537 }
538 
539 /*MC
540    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
541    This stores the matrices in globally unassembled form. Each processor
542    assembles only its local Neumann problem and the parallel matrix vector
543    product is handled "implicitly".
544 
545    Operations Provided:
546 +  MatMult()
547 .  MatMultAdd()
548 .  MatMultTranspose()
549 .  MatMultTransposeAdd()
550 .  MatZeroEntries()
551 .  MatSetOption()
552 .  MatZeroRows()
553 .  MatZeroRowsLocal()
554 .  MatSetValues()
555 .  MatSetValuesLocal()
556 .  MatScale()
557 .  MatGetDiagonal()
558 -  MatSetLocalToGlobalMapping()
559 
560    Options Database Keys:
561 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
562 
563    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
564 
565           You must call MatSetLocalToGlobalMapping() before using this matrix type.
566 
567           You can do matrix preallocation on the local matrix after you obtain it with
568           MatISGetLocalMat()
569 
570   Level: advanced
571 
572 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
573 
574 M*/
575 
576 #undef __FUNCT__
577 #define __FUNCT__ "MatCreate_IS"
578 PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
579 {
580   PetscErrorCode ierr;
581   Mat_IS         *b;
582 
583   PetscFunctionBegin;
584   ierr    = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr);
585   A->data = (void*)b;
586   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
587 
588   A->ops->mult                    = MatMult_IS;
589   A->ops->multadd                 = MatMultAdd_IS;
590   A->ops->multtranspose           = MatMultTranspose_IS;
591   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
592   A->ops->destroy                 = MatDestroy_IS;
593   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
594   A->ops->setvalues               = MatSetValues_IS;
595   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
596   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
597   A->ops->zerorows                = MatZeroRows_IS;
598   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
599   A->ops->assemblybegin           = MatAssemblyBegin_IS;
600   A->ops->assemblyend             = MatAssemblyEnd_IS;
601   A->ops->view                    = MatView_IS;
602   A->ops->zeroentries             = MatZeroEntries_IS;
603   A->ops->scale                   = MatScale_IS;
604   A->ops->getdiagonal             = MatGetDiagonal_IS;
605   A->ops->setoption               = MatSetOption_IS;
606   A->ops->ishermitian             = MatIsHermitian_IS;
607   A->ops->issymmetric             = MatIsSymmetric_IS;
608 
609   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
610   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
611 
612   b->A   = 0;
613   b->ctx = 0;
614   b->x   = 0;
615   b->y   = 0;
616   ierr   = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
617   ierr   = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
618   ierr   = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
619   PetscFunctionReturn(0);
620 }
621