xref: /petsc/src/mat/impls/is/matis.c (revision 3b03a366a58d75015a64027c4058af9dfb8880e2)
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;
15013f74950SBarry Smith   PetscInt       n;
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);
165f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
166f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
1672e74eeadSLisandro Dalcin   ierr = MatSetOptionsPrefix(is->A,"is");CHKERRQ(ierr);
168b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
169b4319ba4SBarry Smith 
170b4319ba4SBarry Smith   /* Create the local work vectors */
171b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr);
172b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
173b4319ba4SBarry Smith 
174b4319ba4SBarry Smith   /* setup the global to local scatter */
175b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
176784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
177d0f46423SBarry Smith   ierr = VecCreateMPIWithArray(((PetscObject)A)->comm,A->cmap->n,A->cmap->N,PETSC_NULL,&global);CHKERRQ(ierr);
178b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
1796bf464f9SBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1806bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1816bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
182b4319ba4SBarry Smith   PetscFunctionReturn(0);
183b4319ba4SBarry Smith }
184b4319ba4SBarry Smith 
1852e74eeadSLisandro Dalcin #define ISG2LMapApply(mapping,n,in,out) 0;\
1862e74eeadSLisandro Dalcin   if (!(mapping)->globals) {\
1872e74eeadSLisandro Dalcin    PetscErrorCode _ierr = ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr);\
1882e74eeadSLisandro Dalcin   }\
1892e74eeadSLisandro Dalcin   {\
1902e74eeadSLisandro Dalcin     PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend;\
1912e74eeadSLisandro Dalcin     for (_i=0; _i<n; _i++) {\
1922e74eeadSLisandro Dalcin       if (in[_i] < 0)           out[_i] = in[_i];\
1932e74eeadSLisandro Dalcin       else if (in[_i] < _start) out[_i] = -1;\
1942e74eeadSLisandro Dalcin       else if (in[_i] > _end)   out[_i] = -1;\
1952e74eeadSLisandro Dalcin       else                      out[_i] = _globals[in[_i] - _start];\
1962e74eeadSLisandro Dalcin     }\
1972e74eeadSLisandro Dalcin   }
1982e74eeadSLisandro Dalcin 
1992e74eeadSLisandro Dalcin #undef __FUNCT__
2002e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
2012e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
2022e74eeadSLisandro Dalcin {
2032e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
2042e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
2052e74eeadSLisandro Dalcin   PetscErrorCode ierr;
2062e74eeadSLisandro Dalcin 
2072e74eeadSLisandro Dalcin   PetscFunctionBegin;
2082e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
2092e74eeadSLisandro Dalcin   if (m > 2048 || n > 2048) {
210e32f2f54SBarry Smith     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
2112e74eeadSLisandro Dalcin   }
2122e74eeadSLisandro Dalcin #endif
2132e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
2142e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
2152e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
2162e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
2172e74eeadSLisandro Dalcin }
2182e74eeadSLisandro Dalcin 
2192e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
2202e74eeadSLisandro Dalcin #undef ISG2LMapApply
221b4319ba4SBarry Smith 
222b4319ba4SBarry Smith #undef __FUNCT__
223b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
22413f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
225b4319ba4SBarry Smith {
226dfbe8321SBarry Smith   PetscErrorCode ierr;
227b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
228b4319ba4SBarry Smith 
229b4319ba4SBarry Smith   PetscFunctionBegin;
230b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
231b4319ba4SBarry Smith   PetscFunctionReturn(0);
232b4319ba4SBarry Smith }
233b4319ba4SBarry Smith 
234b4319ba4SBarry Smith #undef __FUNCT__
2352e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
2362b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2372e74eeadSLisandro Dalcin {
2382e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
2392e74eeadSLisandro Dalcin   PetscInt       n_l=0, *rows_l = PETSC_NULL;
2402e74eeadSLisandro Dalcin   PetscErrorCode ierr;
2412e74eeadSLisandro Dalcin 
2422e74eeadSLisandro Dalcin   PetscFunctionBegin;
24397b48c8fSBarry Smith   if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
2442e74eeadSLisandro Dalcin   if (n) {
2452e74eeadSLisandro Dalcin     ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr);
2462e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
2472e74eeadSLisandro Dalcin   }
2482b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
249c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
2502e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
2512e74eeadSLisandro Dalcin }
2522e74eeadSLisandro Dalcin 
2532e74eeadSLisandro Dalcin #undef __FUNCT__
254b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
2552b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
256b4319ba4SBarry Smith {
257b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
258dfbe8321SBarry Smith   PetscErrorCode ierr;
259f4df32b1SMatthew Knepley   PetscInt       i;
260b4319ba4SBarry Smith   PetscScalar    *array;
261b4319ba4SBarry Smith 
262b4319ba4SBarry Smith   PetscFunctionBegin;
2639c3e2b52SBarry Smith   if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
264b4319ba4SBarry Smith   {
265b4319ba4SBarry Smith     /*
266b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
267b4319ba4SBarry Smith        work properly in the interface nodes.
268b4319ba4SBarry Smith     */
269b4319ba4SBarry Smith     Vec         counter;
270b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
271d0f46423SBarry Smith     ierr = VecCreateMPI(((PetscObject)A)->comm,A->cmap->n,A->cmap->N,&counter);CHKERRQ(ierr);
2722dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
2732dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
274ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
275ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
276ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
277ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2786bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
279b4319ba4SBarry Smith   }
280958c9bccSBarry Smith   if (!n) {
281b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
282b4319ba4SBarry Smith   } else {
283b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
284b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
2852b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
286b4319ba4SBarry Smith     for (i=0; i<n; i++) {
287f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
288b4319ba4SBarry Smith     }
289b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
290b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
291b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
292b4319ba4SBarry Smith   }
293b4319ba4SBarry Smith 
294b4319ba4SBarry Smith   PetscFunctionReturn(0);
295b4319ba4SBarry Smith }
296b4319ba4SBarry Smith 
297b4319ba4SBarry Smith #undef __FUNCT__
298b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
299dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
300b4319ba4SBarry Smith {
301b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
302dfbe8321SBarry Smith   PetscErrorCode ierr;
303b4319ba4SBarry Smith 
304b4319ba4SBarry Smith   PetscFunctionBegin;
305b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
306b4319ba4SBarry Smith   PetscFunctionReturn(0);
307b4319ba4SBarry Smith }
308b4319ba4SBarry Smith 
309b4319ba4SBarry Smith #undef __FUNCT__
310b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
311dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
312b4319ba4SBarry Smith {
313b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
314dfbe8321SBarry Smith   PetscErrorCode ierr;
315b4319ba4SBarry Smith 
316b4319ba4SBarry Smith   PetscFunctionBegin;
317b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
318b4319ba4SBarry Smith   PetscFunctionReturn(0);
319b4319ba4SBarry Smith }
320b4319ba4SBarry Smith 
321b4319ba4SBarry Smith EXTERN_C_BEGIN
322b4319ba4SBarry Smith #undef __FUNCT__
323b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
3247087cfbeSBarry Smith PetscErrorCode  MatISGetLocalMat_IS(Mat mat,Mat *local)
325b4319ba4SBarry Smith {
326b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS *)mat->data;
327b4319ba4SBarry Smith 
328b4319ba4SBarry Smith   PetscFunctionBegin;
329b4319ba4SBarry Smith   *local = is->A;
330b4319ba4SBarry Smith   PetscFunctionReturn(0);
331b4319ba4SBarry Smith }
332b4319ba4SBarry Smith EXTERN_C_END
333b4319ba4SBarry Smith 
334b4319ba4SBarry Smith #undef __FUNCT__
335b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
336b4319ba4SBarry Smith /*@
337b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
338b4319ba4SBarry Smith 
339b4319ba4SBarry Smith   Input Parameter:
340b4319ba4SBarry Smith .  mat - the matrix
341b4319ba4SBarry Smith 
342b4319ba4SBarry Smith   Output Parameter:
343b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
344b4319ba4SBarry Smith 
345b4319ba4SBarry Smith   Level: advanced
346b4319ba4SBarry Smith 
347b4319ba4SBarry Smith   Notes:
348b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
349b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
350b4319ba4SBarry Smith   of the MatSetValues() operation.
351b4319ba4SBarry Smith 
352b4319ba4SBarry Smith .seealso: MATIS
353b4319ba4SBarry Smith @*/
3547087cfbeSBarry Smith PetscErrorCode  MatISGetLocalMat(Mat mat,Mat *local)
355b4319ba4SBarry Smith {
3564ac538c5SBarry Smith   PetscErrorCode ierr;
357b4319ba4SBarry Smith 
358b4319ba4SBarry Smith   PetscFunctionBegin;
3590700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
360b4319ba4SBarry Smith   PetscValidPointer(local,2);
3614ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat *),(mat,local));CHKERRQ(ierr);
362b4319ba4SBarry Smith   PetscFunctionReturn(0);
363b4319ba4SBarry Smith }
364b4319ba4SBarry Smith 
365*3b03a366Sstefano_zampini EXTERN_C_BEGIN
366*3b03a366Sstefano_zampini #undef __FUNCT__
367*3b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
368*3b03a366Sstefano_zampini PetscErrorCode  MatISSetLocalMat_IS(Mat mat,Mat local)
369*3b03a366Sstefano_zampini {
370*3b03a366Sstefano_zampini   Mat_IS *is = (Mat_IS *)mat->data;
371*3b03a366Sstefano_zampini   PetscInt nrows,ncols,orows,ocols;
372*3b03a366Sstefano_zampini   PetscErrorCode ierr;
373*3b03a366Sstefano_zampini 
374*3b03a366Sstefano_zampini   PetscFunctionBegin;
375*3b03a366Sstefano_zampini   ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
376*3b03a366Sstefano_zampini   ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
377*3b03a366Sstefano_zampini   if(orows != nrows || ocols != ncols ) {
378*3b03a366Sstefano_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);
379*3b03a366Sstefano_zampini   }
380*3b03a366Sstefano_zampini   ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
381*3b03a366Sstefano_zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
382*3b03a366Sstefano_zampini   is->A = local;
383*3b03a366Sstefano_zampini 
384*3b03a366Sstefano_zampini   PetscFunctionReturn(0);
385*3b03a366Sstefano_zampini }
386*3b03a366Sstefano_zampini EXTERN_C_END
387*3b03a366Sstefano_zampini 
388*3b03a366Sstefano_zampini #undef __FUNCT__
389*3b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
390*3b03a366Sstefano_zampini /*@
391*3b03a366Sstefano_zampini     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
392*3b03a366Sstefano_zampini 
393*3b03a366Sstefano_zampini   Input Parameter:
394*3b03a366Sstefano_zampini .  mat - the matrix
395*3b03a366Sstefano_zampini .  local - the local matrix usually MATSEQAIJ
396*3b03a366Sstefano_zampini 
397*3b03a366Sstefano_zampini   Output Parameter:
398*3b03a366Sstefano_zampini 
399*3b03a366Sstefano_zampini   Level: advanced
400*3b03a366Sstefano_zampini 
401*3b03a366Sstefano_zampini   Notes:
402*3b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
403*3b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
404*3b03a366Sstefano_zampini 
405*3b03a366Sstefano_zampini .seealso: MATIS
406*3b03a366Sstefano_zampini @*/
407*3b03a366Sstefano_zampini PetscErrorCode  MatISSetLocalMat(Mat mat,Mat local)
408*3b03a366Sstefano_zampini {
409*3b03a366Sstefano_zampini   PetscErrorCode ierr;
410*3b03a366Sstefano_zampini 
411*3b03a366Sstefano_zampini   PetscFunctionBegin;
412*3b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
413*3b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
414*3b03a366Sstefano_zampini   PetscFunctionReturn(0);
415*3b03a366Sstefano_zampini }
416*3b03a366Sstefano_zampini 
4176726f965SBarry Smith #undef __FUNCT__
4186726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
4196726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
4206726f965SBarry Smith {
4216726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
4226726f965SBarry Smith   PetscErrorCode ierr;
4236726f965SBarry Smith 
4246726f965SBarry Smith   PetscFunctionBegin;
4256726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
4266726f965SBarry Smith   PetscFunctionReturn(0);
4276726f965SBarry Smith }
4286726f965SBarry Smith 
4296726f965SBarry Smith #undef __FUNCT__
4302e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
4312e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
4322e74eeadSLisandro Dalcin {
4332e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
4342e74eeadSLisandro Dalcin   PetscErrorCode ierr;
4352e74eeadSLisandro Dalcin 
4362e74eeadSLisandro Dalcin   PetscFunctionBegin;
4372e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
4382e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4392e74eeadSLisandro Dalcin }
4402e74eeadSLisandro Dalcin 
4412e74eeadSLisandro Dalcin #undef __FUNCT__
4422e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
4432e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
4442e74eeadSLisandro Dalcin {
4452e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
4462e74eeadSLisandro Dalcin   PetscErrorCode ierr;
4472e74eeadSLisandro Dalcin 
4482e74eeadSLisandro Dalcin   PetscFunctionBegin;
4492e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
4502e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
4512e74eeadSLisandro Dalcin 
4522e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
4532e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
454ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
455ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4562e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4572e74eeadSLisandro Dalcin }
4582e74eeadSLisandro Dalcin 
4592e74eeadSLisandro Dalcin #undef __FUNCT__
4606726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
461ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool  flg)
4626726f965SBarry Smith {
4636726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
4646726f965SBarry Smith   PetscErrorCode ierr;
4656726f965SBarry Smith 
4666726f965SBarry Smith   PetscFunctionBegin;
4674e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
4686726f965SBarry Smith   PetscFunctionReturn(0);
4696726f965SBarry Smith }
4706726f965SBarry Smith 
471284134d9SBarry Smith #undef __FUNCT__
472284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
473284134d9SBarry Smith /*@
474284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
475284134d9SBarry Smith        process but not across processes.
476284134d9SBarry Smith 
477284134d9SBarry Smith    Input Parameters:
478284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
479284134d9SBarry Smith .     m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
480284134d9SBarry Smith -     map - mapping that defines the global number for each local number
481284134d9SBarry Smith 
482284134d9SBarry Smith    Output Parameter:
483284134d9SBarry Smith .    A - the resulting matrix
484284134d9SBarry Smith 
4858e6c10adSSatish Balay    Level: advanced
4868e6c10adSSatish Balay 
487284134d9SBarry Smith    Notes: See MATIS for more details
4888cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
4898cda6cd7SBarry Smith           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
4908cda6cd7SBarry Smith           plus the ghost points to global indices.
491284134d9SBarry Smith 
492284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
493284134d9SBarry Smith @*/
4947087cfbeSBarry Smith PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
495284134d9SBarry Smith {
496284134d9SBarry Smith   PetscErrorCode ierr;
497284134d9SBarry Smith 
498284134d9SBarry Smith   PetscFunctionBegin;
499284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
500284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
501284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
502*3b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
503784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
504284134d9SBarry Smith   PetscFunctionReturn(0);
505284134d9SBarry Smith }
506284134d9SBarry Smith 
507b4319ba4SBarry Smith /*MC
508b4319ba4SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
509b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
510b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
511b4319ba4SBarry Smith    product is handled "implicitly".
512b4319ba4SBarry Smith 
513b4319ba4SBarry Smith    Operations Provided:
5146726f965SBarry Smith +  MatMult()
5152e74eeadSLisandro Dalcin .  MatMultAdd()
5162e74eeadSLisandro Dalcin .  MatMultTranspose()
5172e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
5186726f965SBarry Smith .  MatZeroEntries()
5196726f965SBarry Smith .  MatSetOption()
5202e74eeadSLisandro Dalcin .  MatZeroRows()
5216726f965SBarry Smith .  MatZeroRowsLocal()
5222e74eeadSLisandro Dalcin .  MatSetValues()
5236726f965SBarry Smith .  MatSetValuesLocal()
5242e74eeadSLisandro Dalcin .  MatScale()
5252e74eeadSLisandro Dalcin .  MatGetDiagonal()
5266726f965SBarry Smith -  MatSetLocalToGlobalMapping()
527b4319ba4SBarry Smith 
528b4319ba4SBarry Smith    Options Database Keys:
529b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
530b4319ba4SBarry Smith 
531b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
532b4319ba4SBarry Smith 
533b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
534b4319ba4SBarry Smith 
535b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
536b4319ba4SBarry Smith           MatISGetLocalMat()
537b4319ba4SBarry Smith 
538b4319ba4SBarry Smith   Level: advanced
539b4319ba4SBarry Smith 
540b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
541b4319ba4SBarry Smith 
542b4319ba4SBarry Smith M*/
543b4319ba4SBarry Smith 
544b4319ba4SBarry Smith EXTERN_C_BEGIN
545b4319ba4SBarry Smith #undef __FUNCT__
546b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
5477087cfbeSBarry Smith PetscErrorCode  MatCreate_IS(Mat A)
548b4319ba4SBarry Smith {
549dfbe8321SBarry Smith   PetscErrorCode ierr;
550b4319ba4SBarry Smith   Mat_IS         *b;
551b4319ba4SBarry Smith 
552b4319ba4SBarry Smith   PetscFunctionBegin;
55338f2d2fdSLisandro Dalcin   ierr                = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr);
554b4319ba4SBarry Smith   A->data             = (void*)b;
555b4319ba4SBarry Smith   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
556b4319ba4SBarry Smith 
557b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
5582e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
5592e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
5602e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
561b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
562b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
5632e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
564b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
5652e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
566b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
567b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
568b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
569b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
5706726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
5712e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
5722e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
5736726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
574b4319ba4SBarry Smith 
57526283091SBarry Smith   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
57626283091SBarry Smith   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
57726283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
57826283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
579b4319ba4SBarry Smith 
580b4319ba4SBarry Smith   b->A          = 0;
581b4319ba4SBarry Smith   b->ctx        = 0;
582b4319ba4SBarry Smith   b->x          = 0;
583b4319ba4SBarry Smith   b->y          = 0;
584b4319ba4SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr);
585*3b03a366Sstefano_zampini   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISSetLocalMat_C","MatISSetLocalMat_IS",MatISSetLocalMat_IS);CHKERRQ(ierr);
58617667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
587b4319ba4SBarry Smith 
588b4319ba4SBarry Smith   PetscFunctionReturn(0);
589b4319ba4SBarry Smith }
590b4319ba4SBarry Smith EXTERN_C_END
591