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