xref: /petsc/src/mat/impls/is/matis.c (revision f0006bf2d9bd2c27ecf9006ac5be7569623c3ee0)
1be1d678aSKris Buschelman 
2b4319ba4SBarry Smith /*
3b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
4b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
5b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
6b4319ba4SBarry Smith    product is handled "implicitly".
7b4319ba4SBarry Smith 
8b4319ba4SBarry Smith      We provide:
9b4319ba4SBarry Smith          MatMult()
10b4319ba4SBarry Smith 
11b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
12b4319ba4SBarry Smith 
13b4319ba4SBarry Smith */
14b4319ba4SBarry Smith 
15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
16b4319ba4SBarry Smith 
17b4319ba4SBarry Smith #undef __FUNCT__
18b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
19dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
20b4319ba4SBarry Smith {
21dfbe8321SBarry Smith   PetscErrorCode ierr;
22b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
23b4319ba4SBarry Smith 
24b4319ba4SBarry Smith   PetscFunctionBegin;
256bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
266bf464f9SBarry Smith   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
276bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
286bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
296bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
30bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
31dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
32901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr);
33b4319ba4SBarry Smith   PetscFunctionReturn(0);
34b4319ba4SBarry Smith }
35b4319ba4SBarry Smith 
36b4319ba4SBarry Smith #undef __FUNCT__
37b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
38dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
39b4319ba4SBarry Smith {
40dfbe8321SBarry Smith   PetscErrorCode ierr;
41b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
42b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
43b4319ba4SBarry Smith 
44b4319ba4SBarry Smith   PetscFunctionBegin;
45b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
46ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
47ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
48b4319ba4SBarry Smith 
49b4319ba4SBarry Smith   /* multiply the local matrix */
50b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
51b4319ba4SBarry Smith 
52b4319ba4SBarry Smith   /* scatter product back into global memory */
532dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
54ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
55ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
56b4319ba4SBarry Smith 
57b4319ba4SBarry Smith   PetscFunctionReturn(0);
58b4319ba4SBarry Smith }
59b4319ba4SBarry Smith 
60b4319ba4SBarry Smith #undef __FUNCT__
612e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
622e74eeadSLisandro Dalcin static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
632e74eeadSLisandro Dalcin {
642e74eeadSLisandro Dalcin   PetscErrorCode ierr;
652e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
662e74eeadSLisandro Dalcin 
672e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
682e74eeadSLisandro Dalcin   /*  scatter the global vector v1 into the local work vector */
69ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
70ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
72ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
732e74eeadSLisandro Dalcin 
742e74eeadSLisandro Dalcin   /* multiply the local matrix */
752e74eeadSLisandro Dalcin   ierr = MatMultAdd(is->A,is->x,is->y,is->y);CHKERRQ(ierr);
762e74eeadSLisandro Dalcin 
772e74eeadSLisandro Dalcin   /* scatter result back into global vector */
782e74eeadSLisandro Dalcin   ierr = VecSet(v3,0);CHKERRQ(ierr);
79ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
80ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
812e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
822e74eeadSLisandro Dalcin }
832e74eeadSLisandro Dalcin 
842e74eeadSLisandro Dalcin #undef __FUNCT__
852e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
862e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
872e74eeadSLisandro Dalcin {
882e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
892e74eeadSLisandro Dalcin   PetscErrorCode ierr;
902e74eeadSLisandro Dalcin 
912e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  y = A' * x */
922e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
93ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
94ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
952e74eeadSLisandro Dalcin 
962e74eeadSLisandro Dalcin   /* multiply the local matrix */
972e74eeadSLisandro Dalcin   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
982e74eeadSLisandro Dalcin 
992e74eeadSLisandro Dalcin   /* scatter product back into global vector */
1002e74eeadSLisandro Dalcin   ierr = VecSet(y,0);CHKERRQ(ierr);
101ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
102ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1032e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
1042e74eeadSLisandro Dalcin }
1052e74eeadSLisandro Dalcin 
1062e74eeadSLisandro Dalcin #undef __FUNCT__
1072e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
1082e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
1092e74eeadSLisandro Dalcin {
1102e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
1112e74eeadSLisandro Dalcin   PetscErrorCode ierr;
1122e74eeadSLisandro Dalcin 
1132e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
1142e74eeadSLisandro Dalcin   /*  scatter the global vectors v1 and v2  into the local work vectors */
115ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
116ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
117ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1192e74eeadSLisandro Dalcin 
1202e74eeadSLisandro Dalcin   /* multiply the local matrix */
1212e74eeadSLisandro Dalcin   ierr = MatMultTransposeAdd(is->A,is->x,is->y,is->y);CHKERRQ(ierr);
1222e74eeadSLisandro Dalcin 
1232e74eeadSLisandro Dalcin   /* scatter result back into global vector */
1242e74eeadSLisandro Dalcin   ierr = VecSet(v3,0);CHKERRQ(ierr);
125ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
126ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1272e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
1282e74eeadSLisandro Dalcin }
1292e74eeadSLisandro Dalcin 
1302e74eeadSLisandro Dalcin #undef __FUNCT__
131b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
132dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
133b4319ba4SBarry Smith {
134b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
135dfbe8321SBarry Smith   PetscErrorCode ierr;
136b4319ba4SBarry Smith   PetscViewer    sviewer;
137b4319ba4SBarry Smith 
138b4319ba4SBarry Smith   PetscFunctionBegin;
139b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
140b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
141b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
142b4319ba4SBarry Smith   PetscFunctionReturn(0);
143b4319ba4SBarry Smith }
144b4319ba4SBarry Smith 
145b4319ba4SBarry Smith #undef __FUNCT__
146b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
147784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
148b4319ba4SBarry Smith {
149dfbe8321SBarry Smith   PetscErrorCode ierr;
1504e4c7dbeSStefano Zampini   PetscInt       n,bs;
151b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
152b4319ba4SBarry Smith   IS             from,to;
153b4319ba4SBarry Smith   Vec            global;
154b4319ba4SBarry Smith 
155b4319ba4SBarry Smith   PetscFunctionBegin;
156e7e72b3dSBarry Smith   if (is->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix");
157784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
158784ac674SJed Brown   if (rmapping != cmapping) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
159784ac674SJed Brown   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
1606bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
161784ac674SJed Brown   is->mapping = rmapping;
162b4319ba4SBarry Smith 
163b4319ba4SBarry Smith   /* Create the local matrix A */
164784ac674SJed Brown   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
1654e4c7dbeSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
166f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
167f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
1684e4c7dbeSStefano Zampini   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
169ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
170ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
171b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
172b4319ba4SBarry Smith 
173b4319ba4SBarry Smith   /* Create the local work vectors */
1744e4c7dbeSStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
1754e4c7dbeSStefano Zampini   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
1764e4c7dbeSStefano Zampini   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
177ff130e51SJed Brown   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
178ff130e51SJed Brown   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
1794e4c7dbeSStefano Zampini   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
180b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
181b4319ba4SBarry Smith 
182b4319ba4SBarry Smith   /* setup the global to local scatter */
183b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
184784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
1854e4c7dbeSStefano Zampini   ierr = MatGetVecs(A,&global,PETSC_NULL);CHKERRQ(ierr);
186b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
1876bf464f9SBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1886bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1896bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
190b4319ba4SBarry Smith   PetscFunctionReturn(0);
191b4319ba4SBarry Smith }
192b4319ba4SBarry Smith 
1932e74eeadSLisandro Dalcin #define ISG2LMapApply(mapping,n,in,out) 0;\
1942e74eeadSLisandro Dalcin   if (!(mapping)->globals) {\
1952e74eeadSLisandro Dalcin    PetscErrorCode _ierr = ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr);\
1962e74eeadSLisandro Dalcin   }\
1972e74eeadSLisandro Dalcin   {\
1982e74eeadSLisandro Dalcin     PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend;\
1992e74eeadSLisandro Dalcin     for (_i=0; _i<n; _i++) {\
2002e74eeadSLisandro Dalcin       if (in[_i] < 0)           out[_i] = in[_i];\
2012e74eeadSLisandro Dalcin       else if (in[_i] < _start) out[_i] = -1;\
2022e74eeadSLisandro Dalcin       else if (in[_i] > _end)   out[_i] = -1;\
2032e74eeadSLisandro Dalcin       else                      out[_i] = _globals[in[_i] - _start];\
2042e74eeadSLisandro Dalcin     }\
2052e74eeadSLisandro Dalcin   }
2062e74eeadSLisandro Dalcin 
2072e74eeadSLisandro Dalcin #undef __FUNCT__
2082e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
2092e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
2102e74eeadSLisandro Dalcin {
2112e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
2122e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
2132e74eeadSLisandro Dalcin   PetscErrorCode ierr;
2142e74eeadSLisandro Dalcin 
2152e74eeadSLisandro Dalcin   PetscFunctionBegin;
2162e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
2172e74eeadSLisandro Dalcin   if (m > 2048 || n > 2048) {
218e32f2f54SBarry Smith     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
2192e74eeadSLisandro Dalcin   }
2202e74eeadSLisandro Dalcin #endif
2212e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
2222e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
2232e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
2242e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
2252e74eeadSLisandro Dalcin }
2262e74eeadSLisandro Dalcin 
2272e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
2282e74eeadSLisandro Dalcin #undef ISG2LMapApply
229b4319ba4SBarry Smith 
230b4319ba4SBarry Smith #undef __FUNCT__
231b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
23213f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
233b4319ba4SBarry Smith {
234dfbe8321SBarry Smith   PetscErrorCode ierr;
235b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
236b4319ba4SBarry Smith 
237b4319ba4SBarry Smith   PetscFunctionBegin;
238b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
239b4319ba4SBarry Smith   PetscFunctionReturn(0);
240b4319ba4SBarry Smith }
241b4319ba4SBarry Smith 
242b4319ba4SBarry Smith #undef __FUNCT__
243*f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
244*f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
245*f0006bf2SLisandro Dalcin {
246*f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
247*f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
248*f0006bf2SLisandro Dalcin 
249*f0006bf2SLisandro Dalcin   PetscFunctionBegin;
250*f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
251*f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
252*f0006bf2SLisandro Dalcin }
253*f0006bf2SLisandro Dalcin 
254*f0006bf2SLisandro Dalcin #undef __FUNCT__
2552e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
2562b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2572e74eeadSLisandro Dalcin {
2582e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
2592e74eeadSLisandro Dalcin   PetscInt       n_l=0, *rows_l = PETSC_NULL;
2602e74eeadSLisandro Dalcin   PetscErrorCode ierr;
2612e74eeadSLisandro Dalcin 
2622e74eeadSLisandro Dalcin   PetscFunctionBegin;
26397b48c8fSBarry Smith   if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
2642e74eeadSLisandro Dalcin   if (n) {
2652e74eeadSLisandro Dalcin     ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr);
2662e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
2672e74eeadSLisandro Dalcin   }
2682b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
269c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
2702e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
2712e74eeadSLisandro Dalcin }
2722e74eeadSLisandro Dalcin 
2732e74eeadSLisandro Dalcin #undef __FUNCT__
274b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
2752b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
276b4319ba4SBarry Smith {
277b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
278dfbe8321SBarry Smith   PetscErrorCode ierr;
279f4df32b1SMatthew Knepley   PetscInt       i;
280b4319ba4SBarry Smith   PetscScalar    *array;
281b4319ba4SBarry Smith 
282b4319ba4SBarry Smith   PetscFunctionBegin;
2839c3e2b52SBarry Smith   if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
284b4319ba4SBarry Smith   {
285b4319ba4SBarry Smith     /*
286b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
287b4319ba4SBarry Smith        work properly in the interface nodes.
288b4319ba4SBarry Smith     */
289b4319ba4SBarry Smith     Vec         counter;
290b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
2914e4c7dbeSStefano Zampini     ierr = MatGetVecs(A,&counter,PETSC_NULL);CHKERRQ(ierr);
2922dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
2932dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
294ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
295ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
296ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
297ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2986bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
299b4319ba4SBarry Smith   }
300958c9bccSBarry Smith   if (!n) {
301b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
302b4319ba4SBarry Smith   } else {
303b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
304b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
3052b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
306b4319ba4SBarry Smith     for (i=0; i<n; i++) {
307f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
308b4319ba4SBarry Smith     }
309b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
310b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
311b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
312b4319ba4SBarry Smith   }
313b4319ba4SBarry Smith 
314b4319ba4SBarry Smith   PetscFunctionReturn(0);
315b4319ba4SBarry Smith }
316b4319ba4SBarry Smith 
317b4319ba4SBarry Smith #undef __FUNCT__
318b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
319dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
320b4319ba4SBarry Smith {
321b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
322dfbe8321SBarry Smith   PetscErrorCode ierr;
323b4319ba4SBarry Smith 
324b4319ba4SBarry Smith   PetscFunctionBegin;
325b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
326b4319ba4SBarry Smith   PetscFunctionReturn(0);
327b4319ba4SBarry Smith }
328b4319ba4SBarry Smith 
329b4319ba4SBarry Smith #undef __FUNCT__
330b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
331dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
332b4319ba4SBarry Smith {
333b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
334dfbe8321SBarry Smith   PetscErrorCode ierr;
335b4319ba4SBarry Smith 
336b4319ba4SBarry Smith   PetscFunctionBegin;
337b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
338b4319ba4SBarry Smith   PetscFunctionReturn(0);
339b4319ba4SBarry Smith }
340b4319ba4SBarry Smith 
341b4319ba4SBarry Smith EXTERN_C_BEGIN
342b4319ba4SBarry Smith #undef __FUNCT__
343b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
3447087cfbeSBarry Smith PetscErrorCode  MatISGetLocalMat_IS(Mat mat,Mat *local)
345b4319ba4SBarry Smith {
346b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS *)mat->data;
347b4319ba4SBarry Smith 
348b4319ba4SBarry Smith   PetscFunctionBegin;
349b4319ba4SBarry Smith   *local = is->A;
350b4319ba4SBarry Smith   PetscFunctionReturn(0);
351b4319ba4SBarry Smith }
352b4319ba4SBarry Smith EXTERN_C_END
353b4319ba4SBarry Smith 
354b4319ba4SBarry Smith #undef __FUNCT__
355b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
356b4319ba4SBarry Smith /*@
357b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
358b4319ba4SBarry Smith 
359b4319ba4SBarry Smith   Input Parameter:
360b4319ba4SBarry Smith .  mat - the matrix
361b4319ba4SBarry Smith 
362b4319ba4SBarry Smith   Output Parameter:
363b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
364b4319ba4SBarry Smith 
365b4319ba4SBarry Smith   Level: advanced
366b4319ba4SBarry Smith 
367b4319ba4SBarry Smith   Notes:
368b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
369b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
370b4319ba4SBarry Smith   of the MatSetValues() operation.
371b4319ba4SBarry Smith 
372b4319ba4SBarry Smith .seealso: MATIS
373b4319ba4SBarry Smith @*/
3747087cfbeSBarry Smith PetscErrorCode  MatISGetLocalMat(Mat mat,Mat *local)
375b4319ba4SBarry Smith {
3764ac538c5SBarry Smith   PetscErrorCode ierr;
377b4319ba4SBarry Smith 
378b4319ba4SBarry Smith   PetscFunctionBegin;
3790700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
380b4319ba4SBarry Smith   PetscValidPointer(local,2);
3814ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat *),(mat,local));CHKERRQ(ierr);
382b4319ba4SBarry Smith   PetscFunctionReturn(0);
383b4319ba4SBarry Smith }
384b4319ba4SBarry Smith 
3853b03a366Sstefano_zampini EXTERN_C_BEGIN
3863b03a366Sstefano_zampini #undef __FUNCT__
3873b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
3883b03a366Sstefano_zampini PetscErrorCode  MatISSetLocalMat_IS(Mat mat,Mat local)
3893b03a366Sstefano_zampini {
3903b03a366Sstefano_zampini   Mat_IS *is = (Mat_IS *)mat->data;
3913b03a366Sstefano_zampini   PetscInt nrows,ncols,orows,ocols;
3923b03a366Sstefano_zampini   PetscErrorCode ierr;
3933b03a366Sstefano_zampini 
3943b03a366Sstefano_zampini   PetscFunctionBegin;
3954e4c7dbeSStefano Zampini   if(is->A) {
3963b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
3973b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
3983b03a366Sstefano_zampini     if(orows != nrows || ocols != ncols ) {
3993b03a366Sstefano_zampini       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);
4003b03a366Sstefano_zampini     }
4014e4c7dbeSStefano Zampini   }
4023b03a366Sstefano_zampini   ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
4033b03a366Sstefano_zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
4043b03a366Sstefano_zampini   is->A = local;
4053b03a366Sstefano_zampini 
4063b03a366Sstefano_zampini   PetscFunctionReturn(0);
4073b03a366Sstefano_zampini }
4083b03a366Sstefano_zampini EXTERN_C_END
4093b03a366Sstefano_zampini 
4103b03a366Sstefano_zampini #undef __FUNCT__
4113b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
4123b03a366Sstefano_zampini /*@
4133b03a366Sstefano_zampini     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
4143b03a366Sstefano_zampini 
4153b03a366Sstefano_zampini   Input Parameter:
4163b03a366Sstefano_zampini .  mat - the matrix
4173b03a366Sstefano_zampini .  local - the local matrix usually MATSEQAIJ
4183b03a366Sstefano_zampini 
4193b03a366Sstefano_zampini   Output Parameter:
4203b03a366Sstefano_zampini 
4213b03a366Sstefano_zampini   Level: advanced
4223b03a366Sstefano_zampini 
4233b03a366Sstefano_zampini   Notes:
4243b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
4253b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
4263b03a366Sstefano_zampini 
4273b03a366Sstefano_zampini .seealso: MATIS
4283b03a366Sstefano_zampini @*/
4293b03a366Sstefano_zampini PetscErrorCode  MatISSetLocalMat(Mat mat,Mat local)
4303b03a366Sstefano_zampini {
4313b03a366Sstefano_zampini   PetscErrorCode ierr;
4323b03a366Sstefano_zampini 
4333b03a366Sstefano_zampini   PetscFunctionBegin;
4343b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
4353b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
4363b03a366Sstefano_zampini   PetscFunctionReturn(0);
4373b03a366Sstefano_zampini }
4383b03a366Sstefano_zampini 
4396726f965SBarry Smith #undef __FUNCT__
4406726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
4416726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
4426726f965SBarry Smith {
4436726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
4446726f965SBarry Smith   PetscErrorCode ierr;
4456726f965SBarry Smith 
4466726f965SBarry Smith   PetscFunctionBegin;
4476726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
4486726f965SBarry Smith   PetscFunctionReturn(0);
4496726f965SBarry Smith }
4506726f965SBarry Smith 
4516726f965SBarry Smith #undef __FUNCT__
4522e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
4532e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
4542e74eeadSLisandro Dalcin {
4552e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
4562e74eeadSLisandro Dalcin   PetscErrorCode ierr;
4572e74eeadSLisandro Dalcin 
4582e74eeadSLisandro Dalcin   PetscFunctionBegin;
4592e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
4602e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4612e74eeadSLisandro Dalcin }
4622e74eeadSLisandro Dalcin 
4632e74eeadSLisandro Dalcin #undef __FUNCT__
4642e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
4652e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
4662e74eeadSLisandro Dalcin {
4672e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
4682e74eeadSLisandro Dalcin   PetscErrorCode ierr;
4692e74eeadSLisandro Dalcin 
4702e74eeadSLisandro Dalcin   PetscFunctionBegin;
4712e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
4722e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
4732e74eeadSLisandro Dalcin 
4742e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
4752e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
476ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
477ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4782e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4792e74eeadSLisandro Dalcin }
4802e74eeadSLisandro Dalcin 
4812e74eeadSLisandro Dalcin #undef __FUNCT__
4826726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
483ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool  flg)
4846726f965SBarry Smith {
4856726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
4866726f965SBarry Smith   PetscErrorCode ierr;
4876726f965SBarry Smith 
4886726f965SBarry Smith   PetscFunctionBegin;
4894e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
4906726f965SBarry Smith   PetscFunctionReturn(0);
4916726f965SBarry Smith }
4926726f965SBarry Smith 
493284134d9SBarry Smith #undef __FUNCT__
494284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
495284134d9SBarry Smith /*@
496284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
497284134d9SBarry Smith        process but not across processes.
498284134d9SBarry Smith 
499284134d9SBarry Smith    Input Parameters:
500284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
5014e4c7dbeSStefano Zampini .     bs - local and global block size of the matrix
502284134d9SBarry Smith .     m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
503284134d9SBarry Smith -     map - mapping that defines the global number for each local number
504284134d9SBarry Smith 
505284134d9SBarry Smith    Output Parameter:
506284134d9SBarry Smith .    A - the resulting matrix
507284134d9SBarry Smith 
5088e6c10adSSatish Balay    Level: advanced
5098e6c10adSSatish Balay 
510284134d9SBarry Smith    Notes: See MATIS for more details
5118cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
5128cda6cd7SBarry Smith           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
5138cda6cd7SBarry Smith           plus the ghost points to global indices.
514284134d9SBarry Smith 
515284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
516284134d9SBarry Smith @*/
5174e4c7dbeSStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
518284134d9SBarry Smith {
519284134d9SBarry Smith   PetscErrorCode ierr;
520284134d9SBarry Smith 
521284134d9SBarry Smith   PetscFunctionBegin;
522284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
523d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
524284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
525284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
5263b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
527784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
528284134d9SBarry Smith   PetscFunctionReturn(0);
529284134d9SBarry Smith }
530284134d9SBarry Smith 
531b4319ba4SBarry Smith /*MC
532b4319ba4SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
533b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
534b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
535b4319ba4SBarry Smith    product is handled "implicitly".
536b4319ba4SBarry Smith 
537b4319ba4SBarry Smith    Operations Provided:
5386726f965SBarry Smith +  MatMult()
5392e74eeadSLisandro Dalcin .  MatMultAdd()
5402e74eeadSLisandro Dalcin .  MatMultTranspose()
5412e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
5426726f965SBarry Smith .  MatZeroEntries()
5436726f965SBarry Smith .  MatSetOption()
5442e74eeadSLisandro Dalcin .  MatZeroRows()
5456726f965SBarry Smith .  MatZeroRowsLocal()
5462e74eeadSLisandro Dalcin .  MatSetValues()
5476726f965SBarry Smith .  MatSetValuesLocal()
5482e74eeadSLisandro Dalcin .  MatScale()
5492e74eeadSLisandro Dalcin .  MatGetDiagonal()
5506726f965SBarry Smith -  MatSetLocalToGlobalMapping()
551b4319ba4SBarry Smith 
552b4319ba4SBarry Smith    Options Database Keys:
553b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
554b4319ba4SBarry Smith 
555b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
556b4319ba4SBarry Smith 
557b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
558b4319ba4SBarry Smith 
559b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
560b4319ba4SBarry Smith           MatISGetLocalMat()
561b4319ba4SBarry Smith 
562b4319ba4SBarry Smith   Level: advanced
563b4319ba4SBarry Smith 
564b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
565b4319ba4SBarry Smith 
566b4319ba4SBarry Smith M*/
567b4319ba4SBarry Smith 
568b4319ba4SBarry Smith EXTERN_C_BEGIN
569b4319ba4SBarry Smith #undef __FUNCT__
570b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
5717087cfbeSBarry Smith PetscErrorCode  MatCreate_IS(Mat A)
572b4319ba4SBarry Smith {
573dfbe8321SBarry Smith   PetscErrorCode ierr;
574b4319ba4SBarry Smith   Mat_IS         *b;
575b4319ba4SBarry Smith 
576b4319ba4SBarry Smith   PetscFunctionBegin;
57738f2d2fdSLisandro Dalcin   ierr                = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr);
578b4319ba4SBarry Smith   A->data             = (void*)b;
579b4319ba4SBarry Smith   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
580b4319ba4SBarry Smith 
581b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
5822e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
5832e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
5842e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
585b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
586b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
5872e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
588b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
589*f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
5902e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
591b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
592b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
593b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
594b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
5956726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
5962e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
5972e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
5986726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
599b4319ba4SBarry Smith 
60026283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
60126283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
602b4319ba4SBarry Smith 
603b4319ba4SBarry Smith   b->A          = 0;
604b4319ba4SBarry Smith   b->ctx        = 0;
605b4319ba4SBarry Smith   b->x          = 0;
606b4319ba4SBarry Smith   b->y          = 0;
607b4319ba4SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr);
6083b03a366Sstefano_zampini   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISSetLocalMat_C","MatISSetLocalMat_IS",MatISSetLocalMat_IS);CHKERRQ(ierr);
60917667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
610b4319ba4SBarry Smith 
611b4319ba4SBarry Smith   PetscFunctionReturn(0);
612b4319ba4SBarry Smith }
613b4319ba4SBarry Smith EXTERN_C_END
614