xref: /petsc/src/mat/impls/is/matis.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
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   PetscFunctionReturn(0);
57b4319ba4SBarry Smith }
58b4319ba4SBarry Smith 
59b4319ba4SBarry Smith #undef __FUNCT__
602e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
612e74eeadSLisandro Dalcin static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
622e74eeadSLisandro Dalcin {
63650997f4SStefano Zampini   Vec            temp_vec;
642e74eeadSLisandro Dalcin   PetscErrorCode ierr;
652e74eeadSLisandro Dalcin 
662e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
67650997f4SStefano Zampini   if (v3 != v2) {
68650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
69650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
70650997f4SStefano Zampini   } else {
71650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
72650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
73650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
74650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
75650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
76650997f4SStefano Zampini   }
772e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
782e74eeadSLisandro Dalcin }
792e74eeadSLisandro Dalcin 
802e74eeadSLisandro Dalcin #undef __FUNCT__
812e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
822e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
832e74eeadSLisandro Dalcin {
842e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
852e74eeadSLisandro Dalcin   PetscErrorCode ierr;
862e74eeadSLisandro Dalcin 
872e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  y = A' * x */
882e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
89ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
90ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
912e74eeadSLisandro Dalcin 
922e74eeadSLisandro Dalcin   /* multiply the local matrix */
932e74eeadSLisandro Dalcin   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
942e74eeadSLisandro Dalcin 
952e74eeadSLisandro Dalcin   /* scatter product back into global vector */
962e74eeadSLisandro Dalcin   ierr = VecSet(y,0);CHKERRQ(ierr);
97ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
98ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
992e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
1002e74eeadSLisandro Dalcin }
1012e74eeadSLisandro Dalcin 
1022e74eeadSLisandro Dalcin #undef __FUNCT__
1032e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
1042e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
1052e74eeadSLisandro Dalcin {
106650997f4SStefano Zampini   Vec            temp_vec;
1072e74eeadSLisandro Dalcin   PetscErrorCode ierr;
1082e74eeadSLisandro Dalcin 
1092e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
110650997f4SStefano Zampini   if (v3 != v2) {
111650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
112650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
113650997f4SStefano Zampini   } else {
114650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
115650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
116650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
117650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
118650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
119650997f4SStefano Zampini   }
1202e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
1212e74eeadSLisandro Dalcin }
1222e74eeadSLisandro Dalcin 
1232e74eeadSLisandro Dalcin #undef __FUNCT__
124b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
125dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
126b4319ba4SBarry Smith {
127b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
128dfbe8321SBarry Smith   PetscErrorCode ierr;
129b4319ba4SBarry Smith   PetscViewer    sviewer;
130b4319ba4SBarry Smith 
131b4319ba4SBarry Smith   PetscFunctionBegin;
132b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
133b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
134b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
135b4319ba4SBarry Smith   PetscFunctionReturn(0);
136b4319ba4SBarry Smith }
137b4319ba4SBarry Smith 
138b4319ba4SBarry Smith #undef __FUNCT__
139b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
140784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
141b4319ba4SBarry Smith {
142dfbe8321SBarry Smith   PetscErrorCode ierr;
1434e4c7dbeSStefano Zampini   PetscInt       n,bs;
144b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
145b4319ba4SBarry Smith   IS             from,to;
146b4319ba4SBarry Smith   Vec            global;
147b4319ba4SBarry Smith 
148b4319ba4SBarry Smith   PetscFunctionBegin;
149e7e72b3dSBarry Smith   if (is->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix");
150784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
151784ac674SJed Brown   if (rmapping != cmapping) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
152784ac674SJed Brown   ierr        = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
1536bf464f9SBarry Smith   ierr        = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
154784ac674SJed Brown   is->mapping = rmapping;
155b4319ba4SBarry Smith 
156b4319ba4SBarry Smith   /* Create the local matrix A */
157784ac674SJed Brown   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
1584e4c7dbeSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
159f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
160f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
1614e4c7dbeSStefano Zampini   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
162ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
163ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
164b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
165b4319ba4SBarry Smith 
166b4319ba4SBarry Smith   /* Create the local work vectors */
1674e4c7dbeSStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
1684e4c7dbeSStefano Zampini   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
1694e4c7dbeSStefano Zampini   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
170ff130e51SJed Brown   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
171ff130e51SJed Brown   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
1724e4c7dbeSStefano Zampini   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
173b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
174b4319ba4SBarry Smith 
175b4319ba4SBarry Smith   /* setup the global to local scatter */
176b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
177784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
1784e4c7dbeSStefano Zampini   ierr = MatGetVecs(A,&global,PETSC_NULL);CHKERRQ(ierr);
179b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
1806bf464f9SBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1816bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1826bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
183b4319ba4SBarry Smith   PetscFunctionReturn(0);
184b4319ba4SBarry Smith }
185b4319ba4SBarry Smith 
1862e74eeadSLisandro Dalcin #define ISG2LMapApply(mapping,n,in,out) 0; \
1872e74eeadSLisandro Dalcin   if (!(mapping)->globals) { \
1882e74eeadSLisandro Dalcin     PetscErrorCode _ierr = ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr); \
1892e74eeadSLisandro Dalcin   } \
1902e74eeadSLisandro Dalcin   { \
1912e74eeadSLisandro Dalcin     PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend; \
1922e74eeadSLisandro Dalcin     for (_i=0; _i<n; _i++) { \
1932e74eeadSLisandro Dalcin       if (in[_i] < 0)           out[_i] = in[_i]; \
1942e74eeadSLisandro Dalcin       else if (in[_i] < _start) out[_i] = -1; \
1952e74eeadSLisandro Dalcin       else if (in[_i] > _end)   out[_i] = -1; \
1962e74eeadSLisandro Dalcin       else                      out[_i] = _globals[in[_i] - _start]; \
1972e74eeadSLisandro Dalcin     } \
1982e74eeadSLisandro Dalcin   }
1992e74eeadSLisandro Dalcin 
2002e74eeadSLisandro Dalcin #undef __FUNCT__
2012e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
2022e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
2032e74eeadSLisandro Dalcin {
2042e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
2052e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
2062e74eeadSLisandro Dalcin   PetscErrorCode ierr;
2072e74eeadSLisandro Dalcin 
2082e74eeadSLisandro Dalcin   PetscFunctionBegin;
2092e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
210f23aa3ddSBarry Smith   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);
2112e74eeadSLisandro Dalcin #endif
2122e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
2132e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
2142e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
2152e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
2162e74eeadSLisandro Dalcin }
2172e74eeadSLisandro Dalcin 
2182e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
2192e74eeadSLisandro Dalcin #undef ISG2LMapApply
220b4319ba4SBarry Smith 
221b4319ba4SBarry Smith #undef __FUNCT__
222b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
22313f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
224b4319ba4SBarry Smith {
225dfbe8321SBarry Smith   PetscErrorCode ierr;
226b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
227b4319ba4SBarry Smith 
228b4319ba4SBarry Smith   PetscFunctionBegin;
229b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
230b4319ba4SBarry Smith   PetscFunctionReturn(0);
231b4319ba4SBarry Smith }
232b4319ba4SBarry Smith 
233b4319ba4SBarry Smith #undef __FUNCT__
2342e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
2352b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2362e74eeadSLisandro Dalcin {
2372e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
2382e74eeadSLisandro Dalcin   PetscInt       n_l = 0, *rows_l = PETSC_NULL;
2392e74eeadSLisandro Dalcin   PetscErrorCode ierr;
2402e74eeadSLisandro Dalcin 
2412e74eeadSLisandro Dalcin   PetscFunctionBegin;
24297b48c8fSBarry Smith   if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
2432e74eeadSLisandro Dalcin   if (n) {
2442e74eeadSLisandro Dalcin     ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr);
2452e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
2462e74eeadSLisandro Dalcin   }
2472b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
248c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
2492e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
2502e74eeadSLisandro Dalcin }
2512e74eeadSLisandro Dalcin 
2522e74eeadSLisandro Dalcin #undef __FUNCT__
253b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
2542b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
255b4319ba4SBarry Smith {
256b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
257dfbe8321SBarry Smith   PetscErrorCode ierr;
258f4df32b1SMatthew Knepley   PetscInt       i;
259b4319ba4SBarry Smith   PetscScalar    *array;
260b4319ba4SBarry Smith 
261b4319ba4SBarry Smith   PetscFunctionBegin;
2629c3e2b52SBarry Smith   if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support");
263b4319ba4SBarry Smith   {
264b4319ba4SBarry Smith     /*
265b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
266b4319ba4SBarry Smith        work properly in the interface nodes.
267b4319ba4SBarry Smith     */
268b4319ba4SBarry Smith     Vec         counter;
269b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
2704e4c7dbeSStefano Zampini     ierr = MatGetVecs(A,&counter,PETSC_NULL);CHKERRQ(ierr);
2712dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
2722dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
273ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
274ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
275ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
276ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2776bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
278b4319ba4SBarry Smith   }
279958c9bccSBarry Smith   if (!n) {
280b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
281b4319ba4SBarry Smith   } else {
282b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
283*2205254eSKarl Rupp 
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   PetscFunctionReturn(0);
294b4319ba4SBarry Smith }
295b4319ba4SBarry Smith 
296b4319ba4SBarry Smith #undef __FUNCT__
297b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
298dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
299b4319ba4SBarry Smith {
300b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
301dfbe8321SBarry Smith   PetscErrorCode ierr;
302b4319ba4SBarry Smith 
303b4319ba4SBarry Smith   PetscFunctionBegin;
304b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
305b4319ba4SBarry Smith   PetscFunctionReturn(0);
306b4319ba4SBarry Smith }
307b4319ba4SBarry Smith 
308b4319ba4SBarry Smith #undef __FUNCT__
309b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
310dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
311b4319ba4SBarry Smith {
312b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
313dfbe8321SBarry Smith   PetscErrorCode ierr;
314b4319ba4SBarry Smith 
315b4319ba4SBarry Smith   PetscFunctionBegin;
316b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
317b4319ba4SBarry Smith   PetscFunctionReturn(0);
318b4319ba4SBarry Smith }
319b4319ba4SBarry Smith 
320b4319ba4SBarry Smith EXTERN_C_BEGIN
321b4319ba4SBarry Smith #undef __FUNCT__
322b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
3237087cfbeSBarry Smith PetscErrorCode  MatISGetLocalMat_IS(Mat mat,Mat *local)
324b4319ba4SBarry Smith {
325b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
326b4319ba4SBarry Smith 
327b4319ba4SBarry Smith   PetscFunctionBegin;
328b4319ba4SBarry Smith   *local = is->A;
329b4319ba4SBarry Smith   PetscFunctionReturn(0);
330b4319ba4SBarry Smith }
331b4319ba4SBarry Smith EXTERN_C_END
332b4319ba4SBarry Smith 
333b4319ba4SBarry Smith #undef __FUNCT__
334b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
335b4319ba4SBarry Smith /*@
336b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
337b4319ba4SBarry Smith 
338b4319ba4SBarry Smith   Input Parameter:
339b4319ba4SBarry Smith .  mat - the matrix
340b4319ba4SBarry Smith 
341b4319ba4SBarry Smith   Output Parameter:
342b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
343b4319ba4SBarry Smith 
344b4319ba4SBarry Smith   Level: advanced
345b4319ba4SBarry Smith 
346b4319ba4SBarry Smith   Notes:
347b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
348b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
349b4319ba4SBarry Smith   of the MatSetValues() operation.
350b4319ba4SBarry Smith 
351b4319ba4SBarry Smith .seealso: MATIS
352b4319ba4SBarry Smith @*/
3537087cfbeSBarry Smith PetscErrorCode  MatISGetLocalMat(Mat mat,Mat *local)
354b4319ba4SBarry Smith {
3554ac538c5SBarry Smith   PetscErrorCode ierr;
356b4319ba4SBarry Smith 
357b4319ba4SBarry Smith   PetscFunctionBegin;
3580700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
359b4319ba4SBarry Smith   PetscValidPointer(local,2);
3604ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
361b4319ba4SBarry Smith   PetscFunctionReturn(0);
362b4319ba4SBarry Smith }
363b4319ba4SBarry Smith 
3643b03a366Sstefano_zampini EXTERN_C_BEGIN
3653b03a366Sstefano_zampini #undef __FUNCT__
3663b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
3673b03a366Sstefano_zampini PetscErrorCode  MatISSetLocalMat_IS(Mat mat,Mat local)
3683b03a366Sstefano_zampini {
3693b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
3703b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
3713b03a366Sstefano_zampini   PetscErrorCode ierr;
3723b03a366Sstefano_zampini 
3733b03a366Sstefano_zampini   PetscFunctionBegin;
3744e4c7dbeSStefano Zampini   if (is->A) {
3753b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
3763b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
377f23aa3ddSBarry Smith     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);
3784e4c7dbeSStefano Zampini   }
3793b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
3803b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
3813b03a366Sstefano_zampini   is->A = local;
3823b03a366Sstefano_zampini   PetscFunctionReturn(0);
3833b03a366Sstefano_zampini }
3843b03a366Sstefano_zampini EXTERN_C_END
3853b03a366Sstefano_zampini 
3863b03a366Sstefano_zampini #undef __FUNCT__
3873b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
3883b03a366Sstefano_zampini /*@
3893b03a366Sstefano_zampini     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
3903b03a366Sstefano_zampini 
3913b03a366Sstefano_zampini   Input Parameter:
3923b03a366Sstefano_zampini .  mat - the matrix
3933b03a366Sstefano_zampini .  local - the local matrix usually MATSEQAIJ
3943b03a366Sstefano_zampini 
3953b03a366Sstefano_zampini   Output Parameter:
3963b03a366Sstefano_zampini 
3973b03a366Sstefano_zampini   Level: advanced
3983b03a366Sstefano_zampini 
3993b03a366Sstefano_zampini   Notes:
4003b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
4013b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
4023b03a366Sstefano_zampini 
4033b03a366Sstefano_zampini .seealso: MATIS
4043b03a366Sstefano_zampini @*/
4053b03a366Sstefano_zampini PetscErrorCode  MatISSetLocalMat(Mat mat,Mat local)
4063b03a366Sstefano_zampini {
4073b03a366Sstefano_zampini   PetscErrorCode ierr;
4083b03a366Sstefano_zampini 
4093b03a366Sstefano_zampini   PetscFunctionBegin;
4103b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
4113b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
4123b03a366Sstefano_zampini   PetscFunctionReturn(0);
4133b03a366Sstefano_zampini }
4143b03a366Sstefano_zampini 
4156726f965SBarry Smith #undef __FUNCT__
4166726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
4176726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
4186726f965SBarry Smith {
4196726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
4206726f965SBarry Smith   PetscErrorCode ierr;
4216726f965SBarry Smith 
4226726f965SBarry Smith   PetscFunctionBegin;
4236726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
4246726f965SBarry Smith   PetscFunctionReturn(0);
4256726f965SBarry Smith }
4266726f965SBarry Smith 
4276726f965SBarry Smith #undef __FUNCT__
4282e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
4292e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
4302e74eeadSLisandro Dalcin {
4312e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
4322e74eeadSLisandro Dalcin   PetscErrorCode ierr;
4332e74eeadSLisandro Dalcin 
4342e74eeadSLisandro Dalcin   PetscFunctionBegin;
4352e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
4362e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4372e74eeadSLisandro Dalcin }
4382e74eeadSLisandro Dalcin 
4392e74eeadSLisandro Dalcin #undef __FUNCT__
4402e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
4412e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
4422e74eeadSLisandro Dalcin {
4432e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
4442e74eeadSLisandro Dalcin   PetscErrorCode ierr;
4452e74eeadSLisandro Dalcin 
4462e74eeadSLisandro Dalcin   PetscFunctionBegin;
4472e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
4482e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
4492e74eeadSLisandro Dalcin 
4502e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
4512e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
452ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
453ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4542e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4552e74eeadSLisandro Dalcin }
4562e74eeadSLisandro Dalcin 
4572e74eeadSLisandro Dalcin #undef __FUNCT__
4586726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
459ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
4606726f965SBarry Smith {
4616726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
4626726f965SBarry Smith   PetscErrorCode ierr;
4636726f965SBarry Smith 
4646726f965SBarry Smith   PetscFunctionBegin;
4654e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
4666726f965SBarry Smith   PetscFunctionReturn(0);
4676726f965SBarry Smith }
4686726f965SBarry Smith 
469284134d9SBarry Smith #undef __FUNCT__
470284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
471284134d9SBarry Smith /*@
472284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
473284134d9SBarry Smith        process but not across processes.
474284134d9SBarry Smith 
475284134d9SBarry Smith    Input Parameters:
476284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
4774e4c7dbeSStefano Zampini .     bs - local and global block size of the matrix
478284134d9SBarry Smith .     m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
479284134d9SBarry Smith -     map - mapping that defines the global number for each local number
480284134d9SBarry Smith 
481284134d9SBarry Smith    Output Parameter:
482284134d9SBarry Smith .    A - the resulting matrix
483284134d9SBarry Smith 
4848e6c10adSSatish Balay    Level: advanced
4858e6c10adSSatish Balay 
486284134d9SBarry Smith    Notes: See MATIS for more details
4878cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
4888cda6cd7SBarry Smith           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
4898cda6cd7SBarry Smith           plus the ghost points to global indices.
490284134d9SBarry Smith 
491284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
492284134d9SBarry Smith @*/
4934e4c7dbeSStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
494284134d9SBarry Smith {
495284134d9SBarry Smith   PetscErrorCode ierr;
496284134d9SBarry Smith 
497284134d9SBarry Smith   PetscFunctionBegin;
498284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
499d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
500284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
501284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
5023b03a366Sstefano_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 = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
57626283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
577b4319ba4SBarry Smith 
578b4319ba4SBarry Smith   b->A   = 0;
579b4319ba4SBarry Smith   b->ctx = 0;
580b4319ba4SBarry Smith   b->x   = 0;
581b4319ba4SBarry Smith   b->y   = 0;
582b4319ba4SBarry Smith   ierr   = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr);
5833b03a366Sstefano_zampini   ierr   = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISSetLocalMat_C","MatISSetLocalMat_IS",MatISSetLocalMat_IS);CHKERRQ(ierr);
58417667f90SBarry Smith   ierr   = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
585b4319ba4SBarry Smith   PetscFunctionReturn(0);
586b4319ba4SBarry Smith }
587b4319ba4SBarry Smith EXTERN_C_END
588