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