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