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