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