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