xref: /petsc/src/mat/impls/is/matis.c (revision 2e74eead718945e4cac82f5798e7ac58b88f4b2a)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3b4319ba4SBarry Smith /*
4b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
5b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
6b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
7b4319ba4SBarry Smith    product is handled "implicitly".
8b4319ba4SBarry Smith 
9b4319ba4SBarry Smith      We provide:
10b4319ba4SBarry Smith          MatMult()
11b4319ba4SBarry Smith 
12b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
13b4319ba4SBarry Smith 
14b4319ba4SBarry Smith */
15b4319ba4SBarry Smith 
16b4319ba4SBarry Smith #include "src/mat/impls/is/matis.h"      /*I "petscmat.h" I*/
17b4319ba4SBarry Smith 
18b4319ba4SBarry Smith #undef __FUNCT__
19b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
20dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
21b4319ba4SBarry Smith {
22dfbe8321SBarry Smith   PetscErrorCode ierr;
23b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
24b4319ba4SBarry Smith 
25b4319ba4SBarry Smith   PetscFunctionBegin;
26b4319ba4SBarry Smith   if (b->A) {
27b4319ba4SBarry Smith     ierr = MatDestroy(b->A);CHKERRQ(ierr);
28b4319ba4SBarry Smith   }
29b4319ba4SBarry Smith   if (b->ctx) {
30b4319ba4SBarry Smith     ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr);
31b4319ba4SBarry Smith   }
32b4319ba4SBarry Smith   if (b->x) {
33b4319ba4SBarry Smith     ierr = VecDestroy(b->x);CHKERRQ(ierr);
34b4319ba4SBarry Smith   }
35b4319ba4SBarry Smith   if (b->y) {
36b4319ba4SBarry Smith     ierr = VecDestroy(b->y);CHKERRQ(ierr);
37b4319ba4SBarry Smith   }
38b4319ba4SBarry Smith   if (b->mapping) {
39b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr);
40b4319ba4SBarry Smith   }
41b4319ba4SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
42dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
43901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr);
44b4319ba4SBarry Smith   PetscFunctionReturn(0);
45b4319ba4SBarry Smith }
46b4319ba4SBarry Smith 
47b4319ba4SBarry Smith #undef __FUNCT__
48b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
49dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
50b4319ba4SBarry Smith {
51dfbe8321SBarry Smith   PetscErrorCode ierr;
52b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
53b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
54b4319ba4SBarry Smith 
55b4319ba4SBarry Smith   PetscFunctionBegin;
56b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
57b4319ba4SBarry Smith   ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
58b4319ba4SBarry Smith   ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
59b4319ba4SBarry Smith 
60b4319ba4SBarry Smith   /* multiply the local matrix */
61b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
62b4319ba4SBarry Smith 
63b4319ba4SBarry Smith   /* scatter product back into global memory */
642dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
65b4319ba4SBarry Smith   ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
66b4319ba4SBarry Smith   ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
67b4319ba4SBarry Smith 
68b4319ba4SBarry Smith   PetscFunctionReturn(0);
69b4319ba4SBarry Smith }
70b4319ba4SBarry Smith 
71b4319ba4SBarry Smith #undef __FUNCT__
72*2e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
73*2e74eeadSLisandro Dalcin static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
74*2e74eeadSLisandro Dalcin {
75*2e74eeadSLisandro Dalcin   PetscErrorCode ierr;
76*2e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
77*2e74eeadSLisandro Dalcin 
78*2e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
79*2e74eeadSLisandro Dalcin   /*  scatter the global vector v1 into the local work vector */
80*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(v1,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
81*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd  (v1,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
82*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(v2,is->y,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
83*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd  (v2,is->y,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
84*2e74eeadSLisandro Dalcin 
85*2e74eeadSLisandro Dalcin   /* multiply the local matrix */
86*2e74eeadSLisandro Dalcin   ierr = MatMultAdd(is->A,is->x,is->y,is->y);CHKERRQ(ierr);
87*2e74eeadSLisandro Dalcin 
88*2e74eeadSLisandro Dalcin   /* scatter result back into global vector */
89*2e74eeadSLisandro Dalcin   ierr = VecSet(v3,0);CHKERRQ(ierr);
90*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(is->y,v3,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
91*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd  (is->y,v3,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
92*2e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
93*2e74eeadSLisandro Dalcin }
94*2e74eeadSLisandro Dalcin 
95*2e74eeadSLisandro Dalcin #undef __FUNCT__
96*2e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
97*2e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
98*2e74eeadSLisandro Dalcin {
99*2e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
100*2e74eeadSLisandro Dalcin   PetscErrorCode ierr;
101*2e74eeadSLisandro Dalcin 
102*2e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  y = A' * x */
103*2e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
104*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
105*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
106*2e74eeadSLisandro Dalcin 
107*2e74eeadSLisandro Dalcin   /* multiply the local matrix */
108*2e74eeadSLisandro Dalcin   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
109*2e74eeadSLisandro Dalcin 
110*2e74eeadSLisandro Dalcin   /* scatter product back into global vector */
111*2e74eeadSLisandro Dalcin   ierr = VecSet(y,0);CHKERRQ(ierr);
112*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
113*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
114*2e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
115*2e74eeadSLisandro Dalcin }
116*2e74eeadSLisandro Dalcin 
117*2e74eeadSLisandro Dalcin #undef __FUNCT__
118*2e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
119*2e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
120*2e74eeadSLisandro Dalcin {
121*2e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
122*2e74eeadSLisandro Dalcin   PetscErrorCode ierr;
123*2e74eeadSLisandro Dalcin 
124*2e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
125*2e74eeadSLisandro Dalcin   /*  scatter the global vectors v1 and v2  into the local work vectors */
126*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(v1,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
127*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd  (v1,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
128*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(v2,is->y,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
129*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd  (v2,is->y,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
130*2e74eeadSLisandro Dalcin 
131*2e74eeadSLisandro Dalcin   /* multiply the local matrix */
132*2e74eeadSLisandro Dalcin   ierr = MatMultTransposeAdd(is->A,is->x,is->y,is->y);CHKERRQ(ierr);
133*2e74eeadSLisandro Dalcin 
134*2e74eeadSLisandro Dalcin   /* scatter result back into global vector */
135*2e74eeadSLisandro Dalcin   ierr = VecSet(v3,0);CHKERRQ(ierr);
136*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(is->y,v3,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
137*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd  (is->y,v3,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
138*2e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
139*2e74eeadSLisandro Dalcin }
140*2e74eeadSLisandro Dalcin 
141*2e74eeadSLisandro Dalcin #undef __FUNCT__
142b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
143dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
144b4319ba4SBarry Smith {
145b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
146dfbe8321SBarry Smith   PetscErrorCode ierr;
147b4319ba4SBarry Smith   PetscViewer    sviewer;
148b4319ba4SBarry Smith 
149b4319ba4SBarry Smith   PetscFunctionBegin;
150b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
151b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
152b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
153b4319ba4SBarry Smith   PetscFunctionReturn(0);
154b4319ba4SBarry Smith }
155b4319ba4SBarry Smith 
156b4319ba4SBarry Smith #undef __FUNCT__
157b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
158dfbe8321SBarry Smith PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping)
159b4319ba4SBarry Smith {
160dfbe8321SBarry Smith   PetscErrorCode ierr;
16113f74950SBarry Smith   PetscInt       n;
162b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
163b4319ba4SBarry Smith   IS             from,to;
164b4319ba4SBarry Smith   Vec            global;
165b4319ba4SBarry Smith 
166b4319ba4SBarry Smith   PetscFunctionBegin;
167*2e74eeadSLisandro Dalcin   if (is->mapping) {
168*2e74eeadSLisandro Dalcin     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix");
169*2e74eeadSLisandro Dalcin   }
170*2e74eeadSLisandro Dalcin   PetscCheckSameComm(A,1,mapping,2);
171b4319ba4SBarry Smith   is->mapping = mapping;
172b4319ba4SBarry Smith   ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr);
173b4319ba4SBarry Smith 
174b4319ba4SBarry Smith   /* Create the local matrix A */
175b4319ba4SBarry Smith   ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr);
176f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
177f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
178*2e74eeadSLisandro Dalcin   ierr = MatSetOptionsPrefix(is->A,"is");CHKERRQ(ierr);
179b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
180b4319ba4SBarry Smith 
181b4319ba4SBarry Smith   /* Create the local work vectors */
182b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr);
183b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
184b4319ba4SBarry Smith 
185b4319ba4SBarry Smith   /* setup the global to local scatter */
186b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
187b4319ba4SBarry Smith   ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr);
188899cda47SBarry Smith   ierr = VecCreateMPI(A->comm,A->cmap.n,A->cmap.N,&global);CHKERRQ(ierr);
189b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
190b4319ba4SBarry Smith   ierr = VecDestroy(global);CHKERRQ(ierr);
191b4319ba4SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
192b4319ba4SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
193b4319ba4SBarry Smith   PetscFunctionReturn(0);
194b4319ba4SBarry Smith }
195b4319ba4SBarry Smith 
196*2e74eeadSLisandro Dalcin #define ISG2LMapApply(mapping,n,in,out) 0;\
197*2e74eeadSLisandro Dalcin   if (!(mapping)->globals) {\
198*2e74eeadSLisandro Dalcin    PetscErrorCode _ierr = ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr);\
199*2e74eeadSLisandro Dalcin   }\
200*2e74eeadSLisandro Dalcin   {\
201*2e74eeadSLisandro Dalcin     PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend;\
202*2e74eeadSLisandro Dalcin     for (_i=0; _i<n; _i++) {\
203*2e74eeadSLisandro Dalcin       if (in[_i] < 0)           out[_i] = in[_i];\
204*2e74eeadSLisandro Dalcin       else if (in[_i] < _start) out[_i] = -1;\
205*2e74eeadSLisandro Dalcin       else if (in[_i] > _end)   out[_i] = -1;\
206*2e74eeadSLisandro Dalcin       else                      out[_i] = _globals[in[_i] - _start];\
207*2e74eeadSLisandro Dalcin     }\
208*2e74eeadSLisandro Dalcin   }
209*2e74eeadSLisandro Dalcin 
210*2e74eeadSLisandro Dalcin #undef __FUNCT__
211*2e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
212*2e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
213*2e74eeadSLisandro Dalcin {
214*2e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
215*2e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
216*2e74eeadSLisandro Dalcin   PetscErrorCode ierr;
217*2e74eeadSLisandro Dalcin 
218*2e74eeadSLisandro Dalcin   PetscFunctionBegin;
219*2e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
220*2e74eeadSLisandro Dalcin   if (m > 2048 || n > 2048) {
221*2e74eeadSLisandro Dalcin     SETERRQ2(PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
222*2e74eeadSLisandro Dalcin   }
223*2e74eeadSLisandro Dalcin #endif
224*2e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
225*2e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
226*2e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
227*2e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
228*2e74eeadSLisandro Dalcin }
229*2e74eeadSLisandro Dalcin 
230*2e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
231*2e74eeadSLisandro Dalcin #undef ISG2LMapApply
232b4319ba4SBarry Smith 
233b4319ba4SBarry Smith #undef __FUNCT__
234b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
23513f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
236b4319ba4SBarry Smith {
237dfbe8321SBarry Smith   PetscErrorCode ierr;
238b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
239b4319ba4SBarry Smith 
240b4319ba4SBarry Smith   PetscFunctionBegin;
241b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
242b4319ba4SBarry Smith   PetscFunctionReturn(0);
243b4319ba4SBarry Smith }
244b4319ba4SBarry Smith 
245b4319ba4SBarry Smith #undef __FUNCT__
246*2e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
247*2e74eeadSLisandro Dalcin PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag)
248*2e74eeadSLisandro Dalcin {
249*2e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
250*2e74eeadSLisandro Dalcin   PetscInt       n_l=0, *rows_l = PETSC_NULL;
251*2e74eeadSLisandro Dalcin   PetscErrorCode ierr;
252*2e74eeadSLisandro Dalcin 
253*2e74eeadSLisandro Dalcin   PetscFunctionBegin;
254*2e74eeadSLisandro Dalcin   if (n) {
255*2e74eeadSLisandro Dalcin     ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr);
256*2e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
257*2e74eeadSLisandro Dalcin   }
258*2e74eeadSLisandro Dalcin   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag);CHKERRQ(ierr);
259*2e74eeadSLisandro Dalcin   if (rows_l) { ierr = PetscFree(rows_l);CHKERRQ(ierr); }
260*2e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
261*2e74eeadSLisandro Dalcin }
262*2e74eeadSLisandro Dalcin 
263*2e74eeadSLisandro Dalcin #undef __FUNCT__
264b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
265f4df32b1SMatthew Knepley PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag)
266b4319ba4SBarry Smith {
267b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
268dfbe8321SBarry Smith   PetscErrorCode ierr;
269f4df32b1SMatthew Knepley   PetscInt       i;
270b4319ba4SBarry Smith   PetscScalar    *array;
271b4319ba4SBarry Smith 
272b4319ba4SBarry Smith   PetscFunctionBegin;
273b4319ba4SBarry Smith   {
274b4319ba4SBarry Smith     /*
275b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
276b4319ba4SBarry Smith        work properly in the interface nodes.
277b4319ba4SBarry Smith     */
278b4319ba4SBarry Smith     Vec         counter;
279b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
280899cda47SBarry Smith     ierr = VecCreateMPI(A->comm,A->cmap.n,A->cmap.N,&counter);CHKERRQ(ierr);
2812dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
2822dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
283b4319ba4SBarry Smith     ierr = VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
284b4319ba4SBarry Smith     ierr = VecScatterEnd  (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
285b4319ba4SBarry Smith     ierr = VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
286b4319ba4SBarry Smith     ierr = VecScatterEnd  (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
287b4319ba4SBarry Smith     ierr = VecDestroy(counter);CHKERRQ(ierr);
288b4319ba4SBarry Smith   }
289958c9bccSBarry Smith   if (!n) {
290b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
291b4319ba4SBarry Smith   } else {
292b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
293b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
294f4df32b1SMatthew Knepley     ierr = MatZeroRows(is->A,n,rows,diag);CHKERRQ(ierr);
295b4319ba4SBarry Smith     for (i=0; i<n; i++) {
296f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
297b4319ba4SBarry Smith     }
298b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
299b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
300b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
301b4319ba4SBarry Smith   }
302b4319ba4SBarry Smith 
303b4319ba4SBarry Smith   PetscFunctionReturn(0);
304b4319ba4SBarry Smith }
305b4319ba4SBarry Smith 
306b4319ba4SBarry Smith #undef __FUNCT__
307b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
308dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
309b4319ba4SBarry Smith {
310b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
311dfbe8321SBarry Smith   PetscErrorCode ierr;
312b4319ba4SBarry Smith 
313b4319ba4SBarry Smith   PetscFunctionBegin;
314b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
315b4319ba4SBarry Smith   PetscFunctionReturn(0);
316b4319ba4SBarry Smith }
317b4319ba4SBarry Smith 
318b4319ba4SBarry Smith #undef __FUNCT__
319b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
320dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
321b4319ba4SBarry Smith {
322b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
323dfbe8321SBarry Smith   PetscErrorCode ierr;
324b4319ba4SBarry Smith 
325b4319ba4SBarry Smith   PetscFunctionBegin;
326b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
327b4319ba4SBarry Smith   PetscFunctionReturn(0);
328b4319ba4SBarry Smith }
329b4319ba4SBarry Smith 
330b4319ba4SBarry Smith EXTERN_C_BEGIN
331b4319ba4SBarry Smith #undef __FUNCT__
332b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
333be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat_IS(Mat mat,Mat *local)
334b4319ba4SBarry Smith {
335b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS *)mat->data;
336b4319ba4SBarry Smith 
337b4319ba4SBarry Smith   PetscFunctionBegin;
338b4319ba4SBarry Smith   *local = is->A;
339b4319ba4SBarry Smith   PetscFunctionReturn(0);
340b4319ba4SBarry Smith }
341b4319ba4SBarry Smith EXTERN_C_END
342b4319ba4SBarry Smith 
343b4319ba4SBarry Smith #undef __FUNCT__
344b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
345b4319ba4SBarry Smith /*@
346b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
347b4319ba4SBarry Smith 
348b4319ba4SBarry Smith   Input Parameter:
349b4319ba4SBarry Smith .  mat - the matrix
350b4319ba4SBarry Smith 
351b4319ba4SBarry Smith   Output Parameter:
352b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
353b4319ba4SBarry Smith 
354b4319ba4SBarry Smith   Level: advanced
355b4319ba4SBarry Smith 
356b4319ba4SBarry Smith   Notes:
357b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
358b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
359b4319ba4SBarry Smith   of the MatSetValues() operation.
360b4319ba4SBarry Smith 
361b4319ba4SBarry Smith .seealso: MATIS
362b4319ba4SBarry Smith @*/
363be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat mat,Mat *local)
364b4319ba4SBarry Smith {
365dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,Mat *);
366b4319ba4SBarry Smith 
367b4319ba4SBarry Smith   PetscFunctionBegin;
368b4319ba4SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
369b4319ba4SBarry Smith   PetscValidPointer(local,2);
370b4319ba4SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr);
371b4319ba4SBarry Smith   if (f) {
372b4319ba4SBarry Smith     ierr = (*f)(mat,local);CHKERRQ(ierr);
373b4319ba4SBarry Smith   } else {
374b4319ba4SBarry Smith     local = 0;
375b4319ba4SBarry Smith   }
376b4319ba4SBarry Smith   PetscFunctionReturn(0);
377b4319ba4SBarry Smith }
378b4319ba4SBarry Smith 
3796726f965SBarry Smith #undef __FUNCT__
3806726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
3816726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
3826726f965SBarry Smith {
3836726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
3846726f965SBarry Smith   PetscErrorCode ierr;
3856726f965SBarry Smith 
3866726f965SBarry Smith   PetscFunctionBegin;
3876726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
3886726f965SBarry Smith   PetscFunctionReturn(0);
3896726f965SBarry Smith }
3906726f965SBarry Smith 
3916726f965SBarry Smith #undef __FUNCT__
392*2e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
393*2e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
394*2e74eeadSLisandro Dalcin {
395*2e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
396*2e74eeadSLisandro Dalcin   PetscErrorCode ierr;
397*2e74eeadSLisandro Dalcin 
398*2e74eeadSLisandro Dalcin   PetscFunctionBegin;
399*2e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
400*2e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
401*2e74eeadSLisandro Dalcin }
402*2e74eeadSLisandro Dalcin 
403*2e74eeadSLisandro Dalcin #undef __FUNCT__
404*2e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
405*2e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
406*2e74eeadSLisandro Dalcin {
407*2e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
408*2e74eeadSLisandro Dalcin   PetscErrorCode ierr;
409*2e74eeadSLisandro Dalcin 
410*2e74eeadSLisandro Dalcin   PetscFunctionBegin;
411*2e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
412*2e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
413*2e74eeadSLisandro Dalcin 
414*2e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
415*2e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
416*2e74eeadSLisandro Dalcin   ierr = VecScatterBegin(is->x,v,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
417*2e74eeadSLisandro Dalcin   ierr = VecScatterEnd  (is->x,v,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
418*2e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
419*2e74eeadSLisandro Dalcin }
420*2e74eeadSLisandro Dalcin 
421*2e74eeadSLisandro Dalcin #undef __FUNCT__
4226726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
4236726f965SBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op)
4246726f965SBarry Smith {
4256726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
4266726f965SBarry Smith   PetscErrorCode ierr;
4276726f965SBarry Smith 
4286726f965SBarry Smith   PetscFunctionBegin;
4296726f965SBarry Smith   ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
4306726f965SBarry Smith   PetscFunctionReturn(0);
4316726f965SBarry Smith }
4326726f965SBarry Smith 
433284134d9SBarry Smith #undef __FUNCT__
434284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
435284134d9SBarry Smith /*@
436284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
437284134d9SBarry Smith        process but not across processes.
438284134d9SBarry Smith 
439284134d9SBarry Smith    Input Parameters:
440284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
441284134d9SBarry Smith .     m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
442284134d9SBarry Smith -     map - mapping that defines the global number for each local number
443284134d9SBarry Smith 
444284134d9SBarry Smith    Output Parameter:
445284134d9SBarry Smith .    A - the resulting matrix
446284134d9SBarry Smith 
4478e6c10adSSatish Balay    Level: advanced
4488e6c10adSSatish Balay 
449284134d9SBarry Smith    Notes: See MATIS for more details
450284134d9SBarry Smith           m and n are NOT related to the size of the map
451284134d9SBarry Smith 
452284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
453284134d9SBarry Smith @*/
454284134d9SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
455284134d9SBarry Smith {
456284134d9SBarry Smith   PetscErrorCode ierr;
457284134d9SBarry Smith 
458284134d9SBarry Smith   PetscFunctionBegin;
459284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
460284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
461284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
462284134d9SBarry Smith   ierr = MatSetLocalToGlobalMapping(*A,map);CHKERRQ(ierr);
463284134d9SBarry Smith   PetscFunctionReturn(0);
464284134d9SBarry Smith }
465284134d9SBarry Smith 
466b4319ba4SBarry Smith /*MC
467b4319ba4SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
468b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
469b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
470b4319ba4SBarry Smith    product is handled "implicitly".
471b4319ba4SBarry Smith 
472b4319ba4SBarry Smith    Operations Provided:
4736726f965SBarry Smith +  MatMult()
474*2e74eeadSLisandro Dalcin .  MatMultAdd()
475*2e74eeadSLisandro Dalcin .  MatMultTranspose()
476*2e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
4776726f965SBarry Smith .  MatZeroEntries()
4786726f965SBarry Smith .  MatSetOption()
479*2e74eeadSLisandro Dalcin .  MatZeroRows()
4806726f965SBarry Smith .  MatZeroRowsLocal()
481*2e74eeadSLisandro Dalcin .  MatSetValues()
4826726f965SBarry Smith .  MatSetValuesLocal()
483*2e74eeadSLisandro Dalcin .  MatScale()
484*2e74eeadSLisandro Dalcin .  MatGetDiagonal()
4856726f965SBarry Smith -  MatSetLocalToGlobalMapping()
486b4319ba4SBarry Smith 
487b4319ba4SBarry Smith    Options Database Keys:
488b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
489b4319ba4SBarry Smith 
490b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
491b4319ba4SBarry Smith 
492b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
493b4319ba4SBarry Smith 
494b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
495b4319ba4SBarry Smith           MatISGetLocalMat()
496b4319ba4SBarry Smith 
497b4319ba4SBarry Smith   Level: advanced
498b4319ba4SBarry Smith 
499b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
500b4319ba4SBarry Smith 
501b4319ba4SBarry Smith M*/
502b4319ba4SBarry Smith 
503b4319ba4SBarry Smith EXTERN_C_BEGIN
504b4319ba4SBarry Smith #undef __FUNCT__
505b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
506be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_IS(Mat A)
507b4319ba4SBarry Smith {
508dfbe8321SBarry Smith   PetscErrorCode ierr;
509b4319ba4SBarry Smith   Mat_IS         *b;
510b4319ba4SBarry Smith 
511b4319ba4SBarry Smith   PetscFunctionBegin;
512b4319ba4SBarry Smith   ierr                = PetscNew(Mat_IS,&b);CHKERRQ(ierr);
513b4319ba4SBarry Smith   A->data             = (void*)b;
514b4319ba4SBarry Smith   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
515b4319ba4SBarry Smith   A->factor           = 0;
516b4319ba4SBarry Smith   A->mapping          = 0;
517b4319ba4SBarry Smith 
518b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
519*2e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
520*2e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
521*2e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
522b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
523b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
524*2e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
525b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
526*2e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
527b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
528b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
529b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
530b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
5316726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
532*2e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
533*2e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
5346726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
535b4319ba4SBarry Smith 
536899cda47SBarry Smith   ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr);
537899cda47SBarry Smith   ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr);
538b4319ba4SBarry Smith 
539b4319ba4SBarry Smith   b->A          = 0;
540b4319ba4SBarry Smith   b->ctx        = 0;
541b4319ba4SBarry Smith   b->x          = 0;
542b4319ba4SBarry Smith   b->y          = 0;
543b4319ba4SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr);
54417667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
545b4319ba4SBarry Smith 
546b4319ba4SBarry Smith   PetscFunctionReturn(0);
547b4319ba4SBarry Smith }
548b4319ba4SBarry Smith EXTERN_C_END
549