xref: /petsc/src/mat/impls/is/matis.c (revision f69a0ea3504314d029ee423e0de2950ece2dab86)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3b4319ba4SBarry Smith /*
4b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
5b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
6b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
7b4319ba4SBarry Smith    product is handled "implicitly".
8b4319ba4SBarry Smith 
9b4319ba4SBarry Smith      We provide:
10b4319ba4SBarry Smith          MatMult()
11b4319ba4SBarry Smith 
12b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
13b4319ba4SBarry Smith 
14b4319ba4SBarry Smith */
15b4319ba4SBarry Smith 
16b4319ba4SBarry Smith #include "src/mat/impls/is/matis.h"      /*I "petscmat.h" I*/
17b4319ba4SBarry Smith 
18b4319ba4SBarry Smith #undef __FUNCT__
19b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
20dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
21b4319ba4SBarry Smith {
22dfbe8321SBarry Smith   PetscErrorCode ierr;
23b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
24b4319ba4SBarry Smith 
25b4319ba4SBarry Smith   PetscFunctionBegin;
26b4319ba4SBarry Smith   if (b->A) {
27b4319ba4SBarry Smith     ierr = MatDestroy(b->A);CHKERRQ(ierr);
28b4319ba4SBarry Smith   }
29b4319ba4SBarry Smith   if (b->ctx) {
30b4319ba4SBarry Smith     ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr);
31b4319ba4SBarry Smith   }
32b4319ba4SBarry Smith   if (b->x) {
33b4319ba4SBarry Smith     ierr = VecDestroy(b->x);CHKERRQ(ierr);
34b4319ba4SBarry Smith   }
35b4319ba4SBarry Smith   if (b->y) {
36b4319ba4SBarry Smith     ierr = VecDestroy(b->y);CHKERRQ(ierr);
37b4319ba4SBarry Smith   }
38b4319ba4SBarry Smith   if (b->mapping) {
39b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr);
40b4319ba4SBarry Smith   }
41b4319ba4SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
42901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr);
43b4319ba4SBarry Smith   PetscFunctionReturn(0);
44b4319ba4SBarry Smith }
45b4319ba4SBarry Smith 
46b4319ba4SBarry Smith #undef __FUNCT__
47b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
48dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
49b4319ba4SBarry Smith {
50dfbe8321SBarry Smith   PetscErrorCode ierr;
51b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
52b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
53b4319ba4SBarry Smith 
54b4319ba4SBarry Smith   PetscFunctionBegin;
55b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
56b4319ba4SBarry Smith   ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
57b4319ba4SBarry Smith   ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
58b4319ba4SBarry Smith 
59b4319ba4SBarry Smith   /* multiply the local matrix */
60b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
61b4319ba4SBarry Smith 
62b4319ba4SBarry Smith   /* scatter product back into global memory */
632dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
64b4319ba4SBarry Smith   ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
65b4319ba4SBarry Smith   ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
66b4319ba4SBarry Smith 
67b4319ba4SBarry Smith   PetscFunctionReturn(0);
68b4319ba4SBarry Smith }
69b4319ba4SBarry Smith 
70b4319ba4SBarry Smith #undef __FUNCT__
71b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
72dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
73b4319ba4SBarry Smith {
74b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
75dfbe8321SBarry Smith   PetscErrorCode ierr;
76b4319ba4SBarry Smith   PetscViewer    sviewer;
77b4319ba4SBarry Smith 
78b4319ba4SBarry Smith   PetscFunctionBegin;
79b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
80b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
81b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
82b4319ba4SBarry Smith   PetscFunctionReturn(0);
83b4319ba4SBarry Smith }
84b4319ba4SBarry Smith 
85b4319ba4SBarry Smith #undef __FUNCT__
86b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
87dfbe8321SBarry Smith PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping)
88b4319ba4SBarry Smith {
89dfbe8321SBarry Smith   PetscErrorCode ierr;
9013f74950SBarry Smith   PetscInt       n;
91b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
92b4319ba4SBarry Smith   IS             from,to;
93b4319ba4SBarry Smith   Vec            global;
94b4319ba4SBarry Smith 
95b4319ba4SBarry Smith   PetscFunctionBegin;
96b4319ba4SBarry Smith   is->mapping = mapping;
97b4319ba4SBarry Smith   ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr);
98b4319ba4SBarry Smith 
99b4319ba4SBarry Smith   /* Create the local matrix A */
100b4319ba4SBarry Smith   ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr);
101*f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
102*f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
103b4319ba4SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");CHKERRQ(ierr);
104b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
105b4319ba4SBarry Smith 
106b4319ba4SBarry Smith   /* Create the local work vectors */
107b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr);
108b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
109b4319ba4SBarry Smith 
110b4319ba4SBarry Smith   /* setup the global to local scatter */
111b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
112b4319ba4SBarry Smith   ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr);
113b4319ba4SBarry Smith   ierr = VecCreateMPI(A->comm,A->n,A->N,&global);CHKERRQ(ierr);
114b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
115b4319ba4SBarry Smith   ierr = VecDestroy(global);CHKERRQ(ierr);
116b4319ba4SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
117b4319ba4SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
118b4319ba4SBarry Smith   PetscFunctionReturn(0);
119b4319ba4SBarry Smith }
120b4319ba4SBarry Smith 
121b4319ba4SBarry Smith 
122b4319ba4SBarry Smith #undef __FUNCT__
123b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
12413f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
125b4319ba4SBarry Smith {
126dfbe8321SBarry Smith   PetscErrorCode ierr;
127b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
128b4319ba4SBarry Smith 
129b4319ba4SBarry Smith   PetscFunctionBegin;
130b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
131b4319ba4SBarry Smith   PetscFunctionReturn(0);
132b4319ba4SBarry Smith }
133b4319ba4SBarry Smith 
134b4319ba4SBarry Smith #undef __FUNCT__
135b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
136dfbe8321SBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,IS isrows,const PetscScalar *diag)
137b4319ba4SBarry Smith {
138b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
139dfbe8321SBarry Smith   PetscErrorCode ierr;
14013f74950SBarry Smith   PetscInt       i,n,*rows;
141b4319ba4SBarry Smith   PetscScalar    *array;
142b4319ba4SBarry Smith 
143b4319ba4SBarry Smith   PetscFunctionBegin;
144b4319ba4SBarry Smith   {
145b4319ba4SBarry Smith     /*
146b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
147b4319ba4SBarry Smith        work properly in the interface nodes.
148b4319ba4SBarry Smith     */
149b4319ba4SBarry Smith     Vec         counter;
150b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
151b4319ba4SBarry Smith     ierr = VecCreateMPI(A->comm,A->n,A->N,&counter);CHKERRQ(ierr);
1522dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
1532dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
154b4319ba4SBarry Smith     ierr = VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
155b4319ba4SBarry Smith     ierr = VecScatterEnd  (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
156b4319ba4SBarry Smith     ierr = VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
157b4319ba4SBarry Smith     ierr = VecScatterEnd  (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
158b4319ba4SBarry Smith     ierr = VecDestroy(counter);CHKERRQ(ierr);
159b4319ba4SBarry Smith   }
160b4319ba4SBarry Smith   ierr = ISGetLocalSize(isrows,&n);CHKERRQ(ierr);
161958c9bccSBarry Smith   if (!n) {
162b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
163b4319ba4SBarry Smith   } else {
164b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
165b4319ba4SBarry Smith     ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
166b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
167b4319ba4SBarry Smith     ierr = MatZeroRows(is->A,isrows,diag);CHKERRQ(ierr);
168b4319ba4SBarry Smith     for (i=0; i<n; i++) {
169b4319ba4SBarry Smith       ierr = MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
170b4319ba4SBarry Smith     }
171b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
172b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
173b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
174b4319ba4SBarry Smith     ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
175b4319ba4SBarry Smith   }
176b4319ba4SBarry Smith 
177b4319ba4SBarry Smith   PetscFunctionReturn(0);
178b4319ba4SBarry Smith }
179b4319ba4SBarry Smith 
180b4319ba4SBarry Smith #undef __FUNCT__
181b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
182dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
183b4319ba4SBarry Smith {
184b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
185dfbe8321SBarry Smith   PetscErrorCode ierr;
186b4319ba4SBarry Smith 
187b4319ba4SBarry Smith   PetscFunctionBegin;
188b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
189b4319ba4SBarry Smith   PetscFunctionReturn(0);
190b4319ba4SBarry Smith }
191b4319ba4SBarry Smith 
192b4319ba4SBarry Smith #undef __FUNCT__
193b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
194dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
195b4319ba4SBarry Smith {
196b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
197dfbe8321SBarry Smith   PetscErrorCode ierr;
198b4319ba4SBarry Smith 
199b4319ba4SBarry Smith   PetscFunctionBegin;
200b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
201b4319ba4SBarry Smith   PetscFunctionReturn(0);
202b4319ba4SBarry Smith }
203b4319ba4SBarry Smith 
204b4319ba4SBarry Smith EXTERN_C_BEGIN
205b4319ba4SBarry Smith #undef __FUNCT__
206b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
207be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat_IS(Mat mat,Mat *local)
208b4319ba4SBarry Smith {
209b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS *)mat->data;
210b4319ba4SBarry Smith 
211b4319ba4SBarry Smith   PetscFunctionBegin;
212b4319ba4SBarry Smith   *local = is->A;
213b4319ba4SBarry Smith   PetscFunctionReturn(0);
214b4319ba4SBarry Smith }
215b4319ba4SBarry Smith EXTERN_C_END
216b4319ba4SBarry Smith 
217b4319ba4SBarry Smith #undef __FUNCT__
218b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
219b4319ba4SBarry Smith /*@
220b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
221b4319ba4SBarry Smith 
222b4319ba4SBarry Smith   Input Parameter:
223b4319ba4SBarry Smith .  mat - the matrix
224b4319ba4SBarry Smith 
225b4319ba4SBarry Smith   Output Parameter:
226b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
227b4319ba4SBarry Smith 
228b4319ba4SBarry Smith   Level: advanced
229b4319ba4SBarry Smith 
230b4319ba4SBarry Smith   Notes:
231b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
232b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
233b4319ba4SBarry Smith   of the MatSetValues() operation.
234b4319ba4SBarry Smith 
235b4319ba4SBarry Smith .seealso: MATIS
236b4319ba4SBarry Smith @*/
237be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat mat,Mat *local)
238b4319ba4SBarry Smith {
239dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,Mat *);
240b4319ba4SBarry Smith 
241b4319ba4SBarry Smith   PetscFunctionBegin;
242b4319ba4SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
243b4319ba4SBarry Smith   PetscValidPointer(local,2);
244b4319ba4SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr);
245b4319ba4SBarry Smith   if (f) {
246b4319ba4SBarry Smith     ierr = (*f)(mat,local);CHKERRQ(ierr);
247b4319ba4SBarry Smith   } else {
248b4319ba4SBarry Smith     local = 0;
249b4319ba4SBarry Smith   }
250b4319ba4SBarry Smith   PetscFunctionReturn(0);
251b4319ba4SBarry Smith }
252b4319ba4SBarry Smith 
2536726f965SBarry Smith #undef __FUNCT__
2546726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
2556726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
2566726f965SBarry Smith {
2576726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
2586726f965SBarry Smith   PetscErrorCode ierr;
2596726f965SBarry Smith 
2606726f965SBarry Smith   PetscFunctionBegin;
2616726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
2626726f965SBarry Smith   PetscFunctionReturn(0);
2636726f965SBarry Smith }
2646726f965SBarry Smith 
2656726f965SBarry Smith #undef __FUNCT__
2666726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
2676726f965SBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op)
2686726f965SBarry Smith {
2696726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
2706726f965SBarry Smith   PetscErrorCode ierr;
2716726f965SBarry Smith 
2726726f965SBarry Smith   PetscFunctionBegin;
2736726f965SBarry Smith   ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
2746726f965SBarry Smith   PetscFunctionReturn(0);
2756726f965SBarry Smith }
2766726f965SBarry Smith 
277b4319ba4SBarry Smith /*MC
278b4319ba4SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
279b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
280b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
281b4319ba4SBarry Smith    product is handled "implicitly".
282b4319ba4SBarry Smith 
283b4319ba4SBarry Smith    Operations Provided:
2846726f965SBarry Smith +  MatMult()
2856726f965SBarry Smith .  MatZeroEntries()
2866726f965SBarry Smith .  MatSetOption()
2876726f965SBarry Smith .  MatZeroRowsLocal()
2886726f965SBarry Smith .  MatSetValuesLocal()
2896726f965SBarry Smith -  MatSetLocalToGlobalMapping()
290b4319ba4SBarry Smith 
291b4319ba4SBarry Smith    Options Database Keys:
292b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
293b4319ba4SBarry Smith 
294b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
295b4319ba4SBarry Smith 
296b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
297b4319ba4SBarry Smith 
298b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
299b4319ba4SBarry Smith           MatISGetLocalMat()
300b4319ba4SBarry Smith 
301b4319ba4SBarry Smith   Level: advanced
302b4319ba4SBarry Smith 
303b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
304b4319ba4SBarry Smith 
305b4319ba4SBarry Smith M*/
306b4319ba4SBarry Smith 
307b4319ba4SBarry Smith EXTERN_C_BEGIN
308b4319ba4SBarry Smith #undef __FUNCT__
309b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
310be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_IS(Mat A)
311b4319ba4SBarry Smith {
312dfbe8321SBarry Smith   PetscErrorCode ierr;
313b4319ba4SBarry Smith   Mat_IS         *b;
314b4319ba4SBarry Smith 
315b4319ba4SBarry Smith   PetscFunctionBegin;
316b4319ba4SBarry Smith   ierr                = PetscNew(Mat_IS,&b);CHKERRQ(ierr);
317b4319ba4SBarry Smith   A->data             = (void*)b;
318b4319ba4SBarry Smith   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
319b4319ba4SBarry Smith   A->factor           = 0;
320b4319ba4SBarry Smith   A->mapping          = 0;
321b4319ba4SBarry Smith 
322b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
323b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
324b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
325b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
326b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
327b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
328b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
329b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
3306726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
3316726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
332b4319ba4SBarry Smith 
333b4319ba4SBarry Smith   ierr = PetscSplitOwnership(A->comm,&A->m,&A->M);CHKERRQ(ierr);
334b4319ba4SBarry Smith   ierr = PetscSplitOwnership(A->comm,&A->n,&A->N);CHKERRQ(ierr);
33513f74950SBarry Smith   ierr = MPI_Scan(&A->m,&b->rend,1,MPIU_INT,MPI_SUM,A->comm);CHKERRQ(ierr);
336b4319ba4SBarry Smith   b->rstart = b->rend - A->m;
337b4319ba4SBarry Smith 
338b4319ba4SBarry Smith   b->A          = 0;
339b4319ba4SBarry Smith   b->ctx        = 0;
340b4319ba4SBarry Smith   b->x          = 0;
341b4319ba4SBarry Smith   b->y          = 0;
342b4319ba4SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr);
343b4319ba4SBarry Smith 
344b4319ba4SBarry Smith   PetscFunctionReturn(0);
345b4319ba4SBarry Smith }
346b4319ba4SBarry Smith EXTERN_C_END
347b4319ba4SBarry Smith 
348b4319ba4SBarry Smith 
349b4319ba4SBarry Smith 
350b4319ba4SBarry Smith 
351b4319ba4SBarry Smith 
352b4319ba4SBarry Smith 
353b4319ba4SBarry Smith 
354b4319ba4SBarry Smith 
355b4319ba4SBarry Smith 
356b4319ba4SBarry Smith 
357b4319ba4SBarry Smith 
358b4319ba4SBarry Smith 
359b4319ba4SBarry Smith 
360