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