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