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