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