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