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