xref: /petsc/src/mat/impls/is/matis.c (revision 0fd2b57f38db073d333aeaf9bdf4aad09b7748ac)
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__ "MatISGetMPIXAIJ_IS"
19 PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
20 {
21   Mat_IS                 *matis = (Mat_IS*)(mat->data);
22   /* info on mat */
23   /* ISLocalToGlobalMapping rmapping,cmapping; */
24   PetscInt               bs,rows,cols;
25   PetscInt               lrows,lcols;
26   PetscInt               local_rows,local_cols;
27   PetscBool              isdense,issbaij,issbaij_red;
28   /* values insertion */
29   PetscScalar            *array;
30   PetscInt               *local_indices,*global_indices;
31   /* work */
32   PetscInt               i,j,index_row;
33   PetscErrorCode         ierr;
34 
35   PetscFunctionBegin;
36   /* MISSING CHECKS
37     - rectangular case not covered (it is not allowed by MATIS)
38   */
39   /* get info from mat */
40   /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */
41   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
42   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
43   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
44   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
45   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
46 
47   /* work */
48   ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr);
49   for (i=0;i<local_rows;i++) local_indices[i]=i;
50   /* map indices of local mat to global values */
51   ierr = PetscMalloc(PetscMax(local_cols,local_rows)*sizeof(*global_indices),&global_indices);CHKERRQ(ierr);
52   /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */
53   ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr);
54 
55   if (issbaij) {
56     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
57   }
58 
59   if (reuse == MAT_INITIAL_MATRIX) {
60     Mat         new_mat;
61     MatType     new_mat_type;
62     Vec         vec_dnz,vec_onz;
63     PetscScalar *my_dnz,*my_onz;
64     PetscInt    *dnz,*onz,*mat_ranges,*row_ownership;
65     PetscInt    index_col,owner;
66     PetscMPIInt nsubdomains;
67 
68     /* determining new matrix type */
69     ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
70     if (issbaij_red) {
71       new_mat_type = MATSBAIJ;
72     } else {
73       if (bs>1) {
74         new_mat_type = MATBAIJ;
75       } else {
76         new_mat_type = MATAIJ;
77       }
78     }
79 
80     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
81     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr);
82     ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr);
83     ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr);
84     ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr);
85     ierr = MatSetUp(new_mat);CHKERRQ(ierr);
86     ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr);
87 
88     /*
89       preallocation
90     */
91 
92     ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr);
93     /*
94        Some vectors are needed to sum up properly on shared interface dofs.
95        Preallocation macros cannot do the job.
96        Note that preallocation is not exact, since it overestimates nonzeros
97     */
98     ierr = MatGetVecs(new_mat,NULL,&vec_dnz);CHKERRQ(ierr);
99     /* ierr = VecSetLocalToGlobalMapping(vec_dnz,rmapping);CHKERRQ(ierr); */
100     ierr = VecSetLocalToGlobalMapping(vec_dnz,matis->mapping);CHKERRQ(ierr);
101     ierr = VecDuplicate(vec_dnz,&vec_onz);CHKERRQ(ierr);
102     /* All processes need to compute entire row ownership */
103     ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
104     ierr = MatGetOwnershipRanges(new_mat,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
105     for (i=0;i<nsubdomains;i++) {
106       for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
107         row_ownership[j]=i;
108       }
109     }
110 
111     /*
112        my_dnz and my_onz contains exact contribution to preallocation from each local mat
113        then, they will be summed up properly. This way, preallocation is always sufficient
114     */
115     ierr = PetscMalloc1(local_rows,&my_dnz);CHKERRQ(ierr);
116     ierr = PetscMalloc1(local_rows,&my_onz);CHKERRQ(ierr);
117     ierr = PetscMemzero(my_dnz,local_rows*sizeof(*my_dnz));CHKERRQ(ierr);
118     ierr = PetscMemzero(my_onz,local_rows*sizeof(*my_onz));CHKERRQ(ierr);
119     /* preallocation as a MATAIJ */
120     if (isdense) { /* special case for dense local matrices */
121       for (i=0;i<local_rows;i++) {
122         index_row = global_indices[i];
123         for (j=i;j<local_rows;j++) {
124           owner = row_ownership[index_row];
125           index_col = global_indices[j];
126           if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
127             my_dnz[i] += 1.0;
128           } else { /* offdiag block */
129             my_onz[i] += 1.0;
130           }
131           /* same as before, interchanging rows and cols */
132           if (i != j) {
133             owner = row_ownership[index_col];
134             if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
135               my_dnz[j] += 1.0;
136             } else {
137               my_onz[j] += 1.0;
138             }
139           }
140         }
141       }
142     } else {
143       for (i=0;i<local_rows;i++) {
144         PetscInt ncols;
145         const PetscInt *cols;
146         index_row = global_indices[i];
147         ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
148         for (j=0;j<ncols;j++) {
149           owner = row_ownership[index_row];
150           index_col = global_indices[cols[j]];
151           if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
152             my_dnz[i] += 1.0;
153           } else { /* offdiag block */
154             my_onz[i] += 1.0;
155           }
156           /* same as before, interchanging rows and cols */
157           if (issbaij) {
158             owner = row_ownership[index_col];
159             if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
160               my_dnz[j] += 1.0;
161             } else {
162               my_onz[j] += 1.0;
163             }
164           }
165         }
166         ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
167       }
168     }
169     ierr = VecSet(vec_dnz,0.0);CHKERRQ(ierr);
170     ierr = VecSet(vec_onz,0.0);CHKERRQ(ierr);
171     if (local_rows) { /* multilevel guard */
172       ierr = VecSetValuesLocal(vec_dnz,local_rows,local_indices,my_dnz,ADD_VALUES);CHKERRQ(ierr);
173       ierr = VecSetValuesLocal(vec_onz,local_rows,local_indices,my_onz,ADD_VALUES);CHKERRQ(ierr);
174     }
175     ierr = VecAssemblyBegin(vec_dnz);CHKERRQ(ierr);
176     ierr = VecAssemblyBegin(vec_onz);CHKERRQ(ierr);
177     ierr = VecAssemblyEnd(vec_dnz);CHKERRQ(ierr);
178     ierr = VecAssemblyEnd(vec_onz);CHKERRQ(ierr);
179     ierr = PetscFree(my_dnz);CHKERRQ(ierr);
180     ierr = PetscFree(my_onz);CHKERRQ(ierr);
181     ierr = PetscFree(row_ownership);CHKERRQ(ierr);
182 
183     /* set computed preallocation in dnz and onz */
184     ierr = VecGetArray(vec_dnz,&array);CHKERRQ(ierr);
185     for (i=0; i<lrows; i++) dnz[i] = (PetscInt)PetscRealPart(array[i]);
186     ierr = VecRestoreArray(vec_dnz,&array);CHKERRQ(ierr);
187     ierr = VecGetArray(vec_onz,&array);CHKERRQ(ierr);
188     for (i=0;i<lrows;i++) onz[i] = (PetscInt)PetscRealPart(array[i]);
189     ierr = VecRestoreArray(vec_onz,&array);CHKERRQ(ierr);
190     ierr = VecDestroy(&vec_dnz);CHKERRQ(ierr);
191     ierr = VecDestroy(&vec_onz);CHKERRQ(ierr);
192 
193     /* Resize preallocation if overestimated */
194     for (i=0;i<lrows;i++) {
195       dnz[i] = PetscMin(dnz[i],lcols);
196       onz[i] = PetscMin(onz[i],cols-lcols);
197     }
198     /* set preallocation */
199     ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr);
200     for (i=0;i<lrows/bs;i++) {
201       dnz[i] = dnz[i*bs]/bs;
202       onz[i] = onz[i*bs]/bs;
203     }
204     ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr);
205     for (i=0;i<lrows/bs;i++) {
206       dnz[i] = dnz[i]-i;
207     }
208     ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr);
209     ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
210     *M = new_mat;
211   } else {
212     PetscInt mbs,mrows,mcols;
213     /* some checks */
214     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
215     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
216     if (mrows != rows) {
217       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
218     }
219     if (mrows != rows) {
220       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
221     }
222     if (mbs != bs) {
223       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
224     }
225     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
226   }
227   /* set local to global mappings */
228   /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */
229   /* Set values */
230   if (isdense) { /* special case for dense local matrices */
231     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
232     ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr);
233     ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
234     ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr);
235     ierr = PetscFree(local_indices);CHKERRQ(ierr);
236     ierr = PetscFree(global_indices);CHKERRQ(ierr);
237   } else { /* very basic values insertion for all other matrix types */
238     ierr = PetscFree(local_indices);CHKERRQ(ierr);
239     for (i=0;i<local_rows;i++) {
240       ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
241       /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */
242       ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr);
243       ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr);
244       ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
245       ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
246     }
247     ierr = PetscFree(global_indices);CHKERRQ(ierr);
248   }
249   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
250   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251   if (isdense) {
252     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
253   }
254   if (issbaij) {
255     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
256   }
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "MatISGetMPIXAIJ"
262 /*@
263     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
264 
265   Input Parameter:
266 .  mat - the matrix (should be of type MATIS)
267 .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
268 
269   Output Parameter:
270 .  newmat - the matrix in AIJ format
271 
272   Level: developer
273 
274   Notes:
275 
276 .seealso: MATIS
277 @*/
278 PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
279 {
280   PetscErrorCode ierr;
281 
282   PetscFunctionBegin;
283   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
284   PetscValidLogicalCollectiveEnum(mat,reuse,2);
285   PetscValidPointer(newmat,3);
286   if (reuse != MAT_INITIAL_MATRIX) {
287     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
288     PetscCheckSameComm(mat,1,*newmat,3);
289   }
290   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "MatDuplicate_IS"
296 PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
297 {
298   PetscErrorCode ierr;
299   Mat_IS         *matis = (Mat_IS*)(mat->data);
300   PetscInt       bs,m,n,M,N;
301   Mat            B,localmat;
302 
303   PetscFunctionBegin;
304   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
305   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
306   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
307   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr);
308   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
309   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
310   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
311   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
312   *newmat = B;
313   PetscFunctionReturn(0);
314 }
315 
316 #undef __FUNCT__
317 #define __FUNCT__ "MatIsHermitian_IS"
318 PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
319 {
320   PetscErrorCode ierr;
321   Mat_IS         *matis = (Mat_IS*)A->data;
322   PetscBool      local_sym;
323 
324   PetscFunctionBegin;
325   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
326   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
327   PetscFunctionReturn(0);
328 }
329 
330 #undef __FUNCT__
331 #define __FUNCT__ "MatIsSymmetric_IS"
332 PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
333 {
334   PetscErrorCode ierr;
335   Mat_IS         *matis = (Mat_IS*)A->data;
336   PetscBool      local_sym;
337 
338   PetscFunctionBegin;
339   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
340   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
341   PetscFunctionReturn(0);
342 }
343 
344 #undef __FUNCT__
345 #define __FUNCT__ "MatDestroy_IS"
346 PetscErrorCode MatDestroy_IS(Mat A)
347 {
348   PetscErrorCode ierr;
349   Mat_IS         *b = (Mat_IS*)A->data;
350 
351   PetscFunctionBegin;
352   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
353   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
354   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
355   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
356   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
357   ierr = PetscFree(A->data);CHKERRQ(ierr);
358   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
359   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
360   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
361   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
362   PetscFunctionReturn(0);
363 }
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "MatMult_IS"
367 PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
368 {
369   PetscErrorCode ierr;
370   Mat_IS         *is  = (Mat_IS*)A->data;
371   PetscScalar    zero = 0.0;
372 
373   PetscFunctionBegin;
374   /*  scatter the global vector x into the local work vector */
375   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
376   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
377 
378   /* multiply the local matrix */
379   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
380 
381   /* scatter product back into global memory */
382   ierr = VecSet(y,zero);CHKERRQ(ierr);
383   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
384   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
385   PetscFunctionReturn(0);
386 }
387 
388 #undef __FUNCT__
389 #define __FUNCT__ "MatMultAdd_IS"
390 PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
391 {
392   Vec            temp_vec;
393   PetscErrorCode ierr;
394 
395   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
396   if (v3 != v2) {
397     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
398     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
399   } else {
400     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
401     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
402     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
403     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
404     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
405   }
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "MatMultTranspose_IS"
411 PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
412 {
413   Mat_IS         *is = (Mat_IS*)A->data;
414   PetscErrorCode ierr;
415 
416   PetscFunctionBegin; /*  y = A' * x */
417   /*  scatter the global vector x into the local work vector */
418   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
419   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
420 
421   /* multiply the local matrix */
422   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
423 
424   /* scatter product back into global vector */
425   ierr = VecSet(y,0);CHKERRQ(ierr);
426   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
427   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
428   PetscFunctionReturn(0);
429 }
430 
431 #undef __FUNCT__
432 #define __FUNCT__ "MatMultTransposeAdd_IS"
433 PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
434 {
435   Vec            temp_vec;
436   PetscErrorCode ierr;
437 
438   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
439   if (v3 != v2) {
440     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
441     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
442   } else {
443     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
444     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
445     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
446     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
447     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
448   }
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "MatView_IS"
454 PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
455 {
456   Mat_IS         *a = (Mat_IS*)A->data;
457   PetscErrorCode ierr;
458   PetscViewer    sviewer;
459 
460   PetscFunctionBegin;
461   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
462   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
463   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
464   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
465   PetscFunctionReturn(0);
466 }
467 
468 #undef __FUNCT__
469 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
470 PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
471 {
472   PetscErrorCode ierr;
473   PetscInt       n,bs;
474   Mat_IS         *is = (Mat_IS*)A->data;
475   IS             from,to;
476   Vec            global;
477 
478   PetscFunctionBegin;
479   PetscCheckSameComm(A,1,rmapping,2);
480   if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
481   if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */
482     ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
483     ierr = VecDestroy(&is->x);CHKERRQ(ierr);
484     ierr = VecDestroy(&is->y);CHKERRQ(ierr);
485     ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr);
486     ierr = MatDestroy(&is->A);CHKERRQ(ierr);
487   }
488   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
489   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
490   is->mapping = rmapping;
491 /*
492   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
493   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
494 */
495 
496   /* Create the local matrix A */
497   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
498   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
499   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
500   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
501   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
502   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
503   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
504   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
505 
506   /* Create the local work vectors */
507   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
508   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
509   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
510   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
511   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
512   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
513   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
514 
515   /* setup the global to local scatter */
516   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
517   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
518   ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr);
519   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
520   ierr = VecDestroy(&global);CHKERRQ(ierr);
521   ierr = ISDestroy(&to);CHKERRQ(ierr);
522   ierr = ISDestroy(&from);CHKERRQ(ierr);
523   PetscFunctionReturn(0);
524 }
525 
526 #undef __FUNCT__
527 #define __FUNCT__ "MatSetValues_IS"
528 PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
529 {
530   Mat_IS         *is = (Mat_IS*)mat->data;
531   PetscInt       rows_l[2048],cols_l[2048];
532   PetscErrorCode ierr;
533 
534   PetscFunctionBegin;
535 #if defined(PETSC_USE_DEBUG)
536   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);
537 #endif
538   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
539   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
540   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
541   PetscFunctionReturn(0);
542 }
543 
544 #undef ISG2LMapSetUp
545 #undef ISG2LMapApply
546 
547 #undef __FUNCT__
548 #define __FUNCT__ "MatSetValuesLocal_IS"
549 PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
550 {
551   PetscErrorCode ierr;
552   Mat_IS         *is = (Mat_IS*)A->data;
553 
554   PetscFunctionBegin;
555   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
556   PetscFunctionReturn(0);
557 }
558 
559 #undef __FUNCT__
560 #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
561 PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
562 {
563   PetscErrorCode ierr;
564   Mat_IS         *is = (Mat_IS*)A->data;
565 
566   PetscFunctionBegin;
567   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
568   PetscFunctionReturn(0);
569 }
570 
571 #undef __FUNCT__
572 #define __FUNCT__ "MatZeroRows_IS"
573 PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
574 {
575   Mat_IS         *is = (Mat_IS*)A->data;
576   PetscInt       n_l = 0, *rows_l = NULL;
577   PetscErrorCode ierr;
578 
579   PetscFunctionBegin;
580   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
581   if (n) {
582     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
583     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
584   }
585   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
586   ierr = PetscFree(rows_l);CHKERRQ(ierr);
587   PetscFunctionReturn(0);
588 }
589 
590 #undef __FUNCT__
591 #define __FUNCT__ "MatZeroRowsLocal_IS"
592 PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
593 {
594   Mat_IS         *is = (Mat_IS*)A->data;
595   PetscErrorCode ierr;
596   PetscInt       i;
597   PetscScalar    *array;
598 
599   PetscFunctionBegin;
600   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
601   {
602     /*
603        Set up is->x as a "counting vector". This is in order to MatMult_IS
604        work properly in the interface nodes.
605     */
606     Vec         counter;
607     PetscScalar one=1.0, zero=0.0;
608     ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr);
609     ierr = VecSet(counter,zero);CHKERRQ(ierr);
610     ierr = VecSet(is->x,one);CHKERRQ(ierr);
611     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
612     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
613     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
614     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
615     ierr = VecDestroy(&counter);CHKERRQ(ierr);
616   }
617   if (!n) {
618     is->pure_neumann = PETSC_TRUE;
619   } else {
620     is->pure_neumann = PETSC_FALSE;
621 
622     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
623     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
624     for (i=0; i<n; i++) {
625       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
626     }
627     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
628     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
629     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
630   }
631   PetscFunctionReturn(0);
632 }
633 
634 #undef __FUNCT__
635 #define __FUNCT__ "MatAssemblyBegin_IS"
636 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
637 {
638   Mat_IS         *is = (Mat_IS*)A->data;
639   PetscErrorCode ierr;
640 
641   PetscFunctionBegin;
642   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
643   PetscFunctionReturn(0);
644 }
645 
646 #undef __FUNCT__
647 #define __FUNCT__ "MatAssemblyEnd_IS"
648 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
649 {
650   Mat_IS         *is = (Mat_IS*)A->data;
651   PetscErrorCode ierr;
652 
653   PetscFunctionBegin;
654   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
655   PetscFunctionReturn(0);
656 }
657 
658 #undef __FUNCT__
659 #define __FUNCT__ "MatISGetLocalMat_IS"
660 PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
661 {
662   Mat_IS *is = (Mat_IS*)mat->data;
663 
664   PetscFunctionBegin;
665   *local = is->A;
666   PetscFunctionReturn(0);
667 }
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "MatISGetLocalMat"
671 /*@
672     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
673 
674   Input Parameter:
675 .  mat - the matrix
676 
677   Output Parameter:
678 .  local - the local matrix usually MATSEQAIJ
679 
680   Level: advanced
681 
682   Notes:
683     This can be called if you have precomputed the nonzero structure of the
684   matrix and want to provide it to the inner matrix object to improve the performance
685   of the MatSetValues() operation.
686 
687 .seealso: MATIS
688 @*/
689 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
690 {
691   PetscErrorCode ierr;
692 
693   PetscFunctionBegin;
694   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
695   PetscValidPointer(local,2);
696   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
697   PetscFunctionReturn(0);
698 }
699 
700 #undef __FUNCT__
701 #define __FUNCT__ "MatISSetLocalMat_IS"
702 PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
703 {
704   Mat_IS         *is = (Mat_IS*)mat->data;
705   PetscInt       nrows,ncols,orows,ocols;
706   PetscErrorCode ierr;
707 
708   PetscFunctionBegin;
709   if (is->A) {
710     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
711     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
712     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);
713   }
714   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
715   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
716   is->A = local;
717   PetscFunctionReturn(0);
718 }
719 
720 #undef __FUNCT__
721 #define __FUNCT__ "MatISSetLocalMat"
722 /*@
723     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
724 
725   Input Parameter:
726 .  mat - the matrix
727 .  local - the local matrix usually MATSEQAIJ
728 
729   Output Parameter:
730 
731   Level: advanced
732 
733   Notes:
734     This can be called if you have precomputed the local matrix and
735   want to provide it to the matrix object MATIS.
736 
737 .seealso: MATIS
738 @*/
739 PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
740 {
741   PetscErrorCode ierr;
742 
743   PetscFunctionBegin;
744   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
745   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
746   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
747   PetscFunctionReturn(0);
748 }
749 
750 #undef __FUNCT__
751 #define __FUNCT__ "MatZeroEntries_IS"
752 PetscErrorCode MatZeroEntries_IS(Mat A)
753 {
754   Mat_IS         *a = (Mat_IS*)A->data;
755   PetscErrorCode ierr;
756 
757   PetscFunctionBegin;
758   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
759   PetscFunctionReturn(0);
760 }
761 
762 #undef __FUNCT__
763 #define __FUNCT__ "MatScale_IS"
764 PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
765 {
766   Mat_IS         *is = (Mat_IS*)A->data;
767   PetscErrorCode ierr;
768 
769   PetscFunctionBegin;
770   ierr = MatScale(is->A,a);CHKERRQ(ierr);
771   PetscFunctionReturn(0);
772 }
773 
774 #undef __FUNCT__
775 #define __FUNCT__ "MatGetDiagonal_IS"
776 PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
777 {
778   Mat_IS         *is = (Mat_IS*)A->data;
779   PetscErrorCode ierr;
780 
781   PetscFunctionBegin;
782   /* get diagonal of the local matrix */
783   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
784 
785   /* scatter diagonal back into global vector */
786   ierr = VecSet(v,0);CHKERRQ(ierr);
787   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
788   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
789   PetscFunctionReturn(0);
790 }
791 
792 #undef __FUNCT__
793 #define __FUNCT__ "MatSetOption_IS"
794 PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
795 {
796   Mat_IS         *a = (Mat_IS*)A->data;
797   PetscErrorCode ierr;
798 
799   PetscFunctionBegin;
800   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
801   PetscFunctionReturn(0);
802 }
803 
804 #undef __FUNCT__
805 #define __FUNCT__ "MatCreateIS"
806 /*@
807     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
808        process but not across processes.
809 
810    Input Parameters:
811 +     comm - MPI communicator that will share the matrix
812 .     bs - local and global block size of the matrix
813 .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
814 -     map - mapping that defines the global number for each local number
815 
816    Output Parameter:
817 .    A - the resulting matrix
818 
819    Level: advanced
820 
821    Notes: See MATIS for more details
822           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
823           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
824           plus the ghost points to global indices.
825 
826 .seealso: MATIS, MatSetLocalToGlobalMapping()
827 @*/
828 PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
829 {
830   PetscErrorCode ierr;
831 
832   PetscFunctionBegin;
833   ierr = MatCreate(comm,A);CHKERRQ(ierr);
834   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
835   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
836   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
837   ierr = MatSetUp(*A);CHKERRQ(ierr);
838   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
839   PetscFunctionReturn(0);
840 }
841 
842 /*MC
843    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC.
844    This stores the matrices in globally unassembled form. Each processor
845    assembles only its local Neumann problem and the parallel matrix vector
846    product is handled "implicitly".
847 
848    Operations Provided:
849 +  MatMult()
850 .  MatMultAdd()
851 .  MatMultTranspose()
852 .  MatMultTransposeAdd()
853 .  MatZeroEntries()
854 .  MatSetOption()
855 .  MatZeroRows()
856 .  MatZeroRowsLocal()
857 .  MatSetValues()
858 .  MatSetValuesLocal()
859 .  MatScale()
860 .  MatGetDiagonal()
861 -  MatSetLocalToGlobalMapping()
862 
863    Options Database Keys:
864 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
865 
866    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
867 
868           You must call MatSetLocalToGlobalMapping() before using this matrix type.
869 
870           You can do matrix preallocation on the local matrix after you obtain it with
871           MatISGetLocalMat()
872 
873   Level: advanced
874 
875 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC
876 
877 M*/
878 
879 #undef __FUNCT__
880 #define __FUNCT__ "MatCreate_IS"
881 PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
882 {
883   PetscErrorCode ierr;
884   Mat_IS         *b;
885 
886   PetscFunctionBegin;
887   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
888   A->data = (void*)b;
889   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
890 
891   A->ops->mult                    = MatMult_IS;
892   A->ops->multadd                 = MatMultAdd_IS;
893   A->ops->multtranspose           = MatMultTranspose_IS;
894   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
895   A->ops->destroy                 = MatDestroy_IS;
896   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
897   A->ops->setvalues               = MatSetValues_IS;
898   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
899   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
900   A->ops->zerorows                = MatZeroRows_IS;
901   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
902   A->ops->assemblybegin           = MatAssemblyBegin_IS;
903   A->ops->assemblyend             = MatAssemblyEnd_IS;
904   A->ops->view                    = MatView_IS;
905   A->ops->zeroentries             = MatZeroEntries_IS;
906   A->ops->scale                   = MatScale_IS;
907   A->ops->getdiagonal             = MatGetDiagonal_IS;
908   A->ops->setoption               = MatSetOption_IS;
909   A->ops->ishermitian             = MatIsHermitian_IS;
910   A->ops->issymmetric             = MatIsSymmetric_IS;
911   A->ops->duplicate               = MatDuplicate_IS;
912 
913   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
914   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
915 
916   /* zeros pointer */
917   b->A       = 0;
918   b->ctx     = 0;
919   b->x       = 0;
920   b->y       = 0;
921   b->mapping = 0;
922 
923   /* special MATIS functions */
924   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
925   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
926   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
927   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
928   PetscFunctionReturn(0);
929 }
930