xref: /petsc/src/mat/impls/composite/mcomposite.c (revision 8cbdbec6f8317ddf7886f91eb9c6bd083b543c50)
1 
2 #include <petsc-private/matimpl.h>        /*I "petscmat.h" I*/
3 
4 typedef struct _Mat_CompositeLink *Mat_CompositeLink;
5 struct _Mat_CompositeLink {
6   Mat               mat;
7   Vec               work;
8   Mat_CompositeLink next,prev;
9 };
10 
11 typedef struct {
12   MatCompositeType  type;
13   Mat_CompositeLink head,tail;
14   Vec               work;
15   PetscScalar       scale;        /* scale factor supplied with MatScale() */
16   Vec               left,right;   /* left and right diagonal scaling provided with MatDiagonalScale() */
17   Vec               leftwork,rightwork;
18 } Mat_Composite;
19 
20 #undef __FUNCT__
21 #define __FUNCT__ "MatDestroy_Composite"
22 PetscErrorCode MatDestroy_Composite(Mat mat)
23 {
24   PetscErrorCode   ierr;
25   Mat_Composite    *shell = (Mat_Composite*)mat->data;
26   Mat_CompositeLink next = shell->head,oldnext;
27 
28   PetscFunctionBegin;
29   while (next) {
30     ierr = MatDestroy(&next->mat);CHKERRQ(ierr);
31     if (next->work && (!next->next || next->work != next->next->work)) {
32       ierr = VecDestroy(&next->work);CHKERRQ(ierr);
33     }
34     oldnext = next;
35     next     = next->next;
36     ierr     = PetscFree(oldnext);CHKERRQ(ierr);
37   }
38   ierr = VecDestroy(&shell->work);CHKERRQ(ierr);
39   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
40   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
41   ierr = VecDestroy(&shell->leftwork);CHKERRQ(ierr);
42   ierr = VecDestroy(&shell->rightwork);CHKERRQ(ierr);
43   ierr = PetscFree(mat->data);CHKERRQ(ierr);
44   PetscFunctionReturn(0);
45 }
46 
47 #undef __FUNCT__
48 #define __FUNCT__ "MatMult_Composite_Multiplicative"
49 PetscErrorCode MatMult_Composite_Multiplicative(Mat A,Vec x,Vec y)
50 {
51   Mat_Composite     *shell = (Mat_Composite*)A->data;
52   Mat_CompositeLink next = shell->head;
53   PetscErrorCode    ierr;
54   Vec               in,out;
55 
56   PetscFunctionBegin;
57   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
58   in = x;
59   if (shell->right) {
60     if (!shell->rightwork) {
61       ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr);
62     }
63     ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr);
64     in   = shell->rightwork;
65   }
66   while (next->next) {
67     if (!next->work) { /* should reuse previous work if the same size */
68       ierr = MatGetVecs(next->mat,PETSC_NULL,&next->work);CHKERRQ(ierr);
69     }
70     out = next->work;
71     ierr = MatMult(next->mat,in,out);CHKERRQ(ierr);
72     in   = out;
73     next = next->next;
74   }
75   ierr = MatMult(next->mat,in,y);CHKERRQ(ierr);
76   if (shell->left) {
77     ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr);
78   }
79   ierr = VecScale(y,shell->scale);CHKERRQ(ierr);
80   PetscFunctionReturn(0);
81 }
82 
83 #undef __FUNCT__
84 #define __FUNCT__ "MatMultTranspose_Composite_Multiplicative"
85 PetscErrorCode MatMultTranspose_Composite_Multiplicative(Mat A,Vec x,Vec y)
86 {
87   Mat_Composite     *shell = (Mat_Composite*)A->data;
88   Mat_CompositeLink tail = shell->tail;
89   PetscErrorCode    ierr;
90   Vec               in,out;
91 
92   PetscFunctionBegin;
93   if (!tail) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
94   in = x;
95   if (shell->left) {
96     if (!shell->leftwork) {
97       ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr);
98     }
99     ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr);
100     in   = shell->leftwork;
101   }
102   while (tail->prev) {
103     if (!tail->prev->work) { /* should reuse previous work if the same size */
104       ierr = MatGetVecs(tail->mat,PETSC_NULL,&tail->prev->work);CHKERRQ(ierr);
105     }
106     out = tail->prev->work;
107     ierr = MatMultTranspose(tail->mat,in,out);CHKERRQ(ierr);
108     in   = out;
109     tail = tail->prev;
110   }
111   ierr = MatMultTranspose(tail->mat,in,y);CHKERRQ(ierr);
112   if (shell->right) {
113     ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr);
114   }
115   ierr = VecScale(y,shell->scale);CHKERRQ(ierr);
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "MatMult_Composite"
121 PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y)
122 {
123   Mat_Composite     *shell = (Mat_Composite*)A->data;
124   Mat_CompositeLink next = shell->head;
125   PetscErrorCode    ierr;
126   Vec               in;
127 
128   PetscFunctionBegin;
129   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
130   in = x;
131   if (shell->right) {
132     if (!shell->rightwork) {
133       ierr = VecDuplicate(shell->right,&shell->rightwork);CHKERRQ(ierr);
134     }
135     ierr = VecPointwiseMult(shell->rightwork,shell->right,in);CHKERRQ(ierr);
136     in   = shell->rightwork;
137   }
138   ierr = MatMult(next->mat,in,y);CHKERRQ(ierr);
139   while ((next = next->next)) {
140     ierr = MatMultAdd(next->mat,in,y,y);CHKERRQ(ierr);
141   }
142   if (shell->left) {
143     ierr = VecPointwiseMult(y,shell->left,y);CHKERRQ(ierr);
144   }
145   ierr = VecScale(y,shell->scale);CHKERRQ(ierr);
146   PetscFunctionReturn(0);
147 }
148 
149 #undef __FUNCT__
150 #define __FUNCT__ "MatMultTranspose_Composite"
151 PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y)
152 {
153   Mat_Composite     *shell = (Mat_Composite*)A->data;
154   Mat_CompositeLink next = shell->head;
155   PetscErrorCode    ierr;
156   Vec               in;
157 
158   PetscFunctionBegin;
159   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
160   in = x;
161   if (shell->left) {
162     if (!shell->leftwork) {
163       ierr = VecDuplicate(shell->left,&shell->leftwork);CHKERRQ(ierr);
164     }
165     ierr = VecPointwiseMult(shell->leftwork,shell->left,in);CHKERRQ(ierr);
166     in   = shell->leftwork;
167   }
168   ierr = MatMultTranspose(next->mat,in,y);CHKERRQ(ierr);
169   while ((next = next->next)) {
170     ierr = MatMultTransposeAdd(next->mat,in,y,y);CHKERRQ(ierr);
171   }
172   if (shell->right) {
173     ierr = VecPointwiseMult(y,shell->right,y);CHKERRQ(ierr);
174   }
175   ierr = VecScale(y,shell->scale);CHKERRQ(ierr);
176   PetscFunctionReturn(0);
177 }
178 
179 #undef __FUNCT__
180 #define __FUNCT__ "MatGetDiagonal_Composite"
181 PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v)
182 {
183   Mat_Composite     *shell = (Mat_Composite*)A->data;
184   Mat_CompositeLink next = shell->head;
185   PetscErrorCode    ierr;
186 
187   PetscFunctionBegin;
188   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
189   if (shell->right || shell->left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get diagonal if left or right scaling");
190 
191   ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr);
192   if (next->next && !shell->work) {
193     ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr);
194   }
195   while ((next = next->next)) {
196     ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr);
197     ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr);
198   }
199   ierr = VecScale(v,shell->scale);CHKERRQ(ierr);
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "MatAssemblyEnd_Composite"
205 PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t)
206 {
207   PetscErrorCode ierr;
208   PetscBool      flg = PETSC_FALSE;
209 
210   PetscFunctionBegin;
211   ierr = PetscOptionsGetBool(((PetscObject)Y)->prefix,"-mat_composite_merge",&flg,PETSC_NULL);CHKERRQ(ierr);
212   if (flg) {
213     ierr = MatCompositeMerge(Y);CHKERRQ(ierr);
214   }
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "MatScale_Composite"
220 PetscErrorCode MatScale_Composite(Mat inA,PetscScalar alpha)
221 {
222   Mat_Composite  *a = (Mat_Composite*)inA->data;
223 
224   PetscFunctionBegin;
225   a->scale *= alpha;
226   PetscFunctionReturn(0);
227 }
228 
229 #undef __FUNCT__
230 #define __FUNCT__ "MatDiagonalScale_Composite"
231 PetscErrorCode MatDiagonalScale_Composite(Mat inA,Vec left,Vec right)
232 {
233   Mat_Composite  *a = (Mat_Composite*)inA->data;
234   PetscErrorCode ierr;
235 
236   PetscFunctionBegin;
237   if (left) {
238     if (!a->left) {
239       ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr);
240       ierr = VecCopy(left,a->left);CHKERRQ(ierr);
241     } else {
242       ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr);
243     }
244   }
245   if (right) {
246     if (!a->right) {
247       ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr);
248       ierr = VecCopy(right,a->right);CHKERRQ(ierr);
249     } else {
250       ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr);
251     }
252   }
253   PetscFunctionReturn(0);
254 }
255 
256 static struct _MatOps MatOps_Values = {0,
257        0,
258        0,
259        MatMult_Composite,
260        0,
261 /* 5*/ MatMultTranspose_Composite,
262        0,
263        0,
264        0,
265        0,
266 /*10*/ 0,
267        0,
268        0,
269        0,
270        0,
271 /*15*/ 0,
272        0,
273        MatGetDiagonal_Composite,
274        MatDiagonalScale_Composite,
275        0,
276 /*20*/ 0,
277        MatAssemblyEnd_Composite,
278        0,
279        0,
280 /*24*/ 0,
281        0,
282        0,
283        0,
284        0,
285 /*29*/ 0,
286        0,
287        0,
288        0,
289        0,
290 /*34*/ 0,
291        0,
292        0,
293        0,
294        0,
295 /*39*/ 0,
296        0,
297        0,
298        0,
299        0,
300 /*44*/ 0,
301        MatScale_Composite,
302        0,
303        0,
304        0,
305 /*49*/ 0,
306        0,
307        0,
308        0,
309        0,
310 /*54*/ 0,
311        0,
312        0,
313        0,
314        0,
315 /*59*/ 0,
316        MatDestroy_Composite,
317        0,
318        0,
319        0,
320 /*64*/ 0,
321        0,
322        0,
323        0,
324        0,
325 /*69*/ 0,
326        0,
327        0,
328        0,
329        0,
330 /*74*/ 0,
331        0,
332        0,
333        0,
334        0,
335 /*79*/ 0,
336        0,
337        0,
338        0,
339        0,
340 /*84*/ 0,
341        0,
342        0,
343        0,
344        0,
345 /*89*/ 0,
346        0,
347        0,
348        0,
349        0,
350 /*94*/ 0,
351        0,
352        0,
353        0};
354 
355 /*MC
356    MATCOMPOSITE - A matrix defined by the sum (or product) of one or more matrices (all matrices are of same size and parallel layout).
357 
358    Notes: to use the product of the matrices call MatCompositeSetType(mat,MAT_COMPOSITE_MULTIPLICATIVE);
359 
360   Level: advanced
361 
362 .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType
363 M*/
364 
365 EXTERN_C_BEGIN
366 #undef __FUNCT__
367 #define __FUNCT__ "MatCreate_Composite"
368 PetscErrorCode  MatCreate_Composite(Mat A)
369 {
370   Mat_Composite  *b;
371   PetscErrorCode ierr;
372 
373   PetscFunctionBegin;
374   ierr = PetscNewLog(A,Mat_Composite,&b);CHKERRQ(ierr);
375   A->data = (void*)b;
376   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
377 
378   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
379   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
380 
381   A->assembled     = PETSC_TRUE;
382   A->preallocated  = PETSC_TRUE;
383   b->type          = MAT_COMPOSITE_ADDITIVE;
384   b->scale         = 1.0;
385   ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr);
386   PetscFunctionReturn(0);
387 }
388 EXTERN_C_END
389 
390 #undef __FUNCT__
391 #define __FUNCT__ "MatCreateComposite"
392 /*@C
393    MatCreateComposite - Creates a matrix as the sum of zero or more matrices
394 
395   Collective on MPI_Comm
396 
397    Input Parameters:
398 +  comm - MPI communicator
399 .  nmat - number of matrices to put in
400 -  mats - the matrices
401 
402    Output Parameter:
403 .  A - the matrix
404 
405    Level: advanced
406 
407    Notes:
408      Alternative construction
409 $       MatCreate(comm,&mat);
410 $       MatSetSizes(mat,m,n,M,N);
411 $       MatSetType(mat,MATCOMPOSITE);
412 $       MatCompositeAddMat(mat,mats[0]);
413 $       ....
414 $       MatCompositeAddMat(mat,mats[nmat-1]);
415 $       MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
416 $       MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
417 
418      For the multiplicative form the product is mat[nmat-1]*mat[nmat-2]*....*mat[0]
419 
420 .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeMerge(), MatCompositeSetType(), MatCompositeType
421 
422 @*/
423 PetscErrorCode  MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat)
424 {
425   PetscErrorCode ierr;
426   PetscInt       m,n,M,N,i;
427 
428   PetscFunctionBegin;
429   if (nmat < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix");
430   PetscValidPointer(mat,3);
431 
432   ierr = MatGetLocalSize(mats[0],&m,&n);CHKERRQ(ierr);
433   ierr = MatGetSize(mats[0],&M,&N);CHKERRQ(ierr);
434   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
435   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
436   ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr);
437   for (i=0; i<nmat; i++) {
438     ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr);
439   }
440   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
441   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
442   PetscFunctionReturn(0);
443 }
444 
445 #undef __FUNCT__
446 #define __FUNCT__ "MatCompositeAddMat"
447 /*@
448     MatCompositeAddMat - add another matrix to a composite matrix
449 
450    Collective on Mat
451 
452     Input Parameters:
453 +   mat - the composite matrix
454 -   smat - the partial matrix
455 
456    Level: advanced
457 
458 .seealso: MatCreateComposite()
459 @*/
460 PetscErrorCode  MatCompositeAddMat(Mat mat,Mat smat)
461 {
462   Mat_Composite     *shell;
463   PetscErrorCode    ierr;
464   Mat_CompositeLink ilink,next;
465 
466   PetscFunctionBegin;
467   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
468   PetscValidHeaderSpecific(smat,MAT_CLASSID,2);
469   ierr        = PetscNewLog(mat,struct _Mat_CompositeLink,&ilink);CHKERRQ(ierr);
470   ilink->next = 0;
471   ierr        = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr);
472   ilink->mat  = smat;
473 
474   shell = (Mat_Composite*)mat->data;
475   next = shell->head;
476   if (!next) {
477     shell->head  = ilink;
478   } else {
479     while (next->next) {
480       next = next->next;
481     }
482     next->next      = ilink;
483     ilink->prev     = next;
484   }
485   shell->tail = ilink;
486   PetscFunctionReturn(0);
487 }
488 
489 #undef __FUNCT__
490 #define __FUNCT__ "MatCompositeSetType"
491 /*@C
492    MatCompositeSetType - Indicates if the matrix is defined as the sum of a set of matrices or the product
493 
494   Collective on MPI_Comm
495 
496    Input Parameters:
497 .  mat - the composite matrix
498 
499 
500    Level: advanced
501 
502    Notes:
503       The MatType of the resulting matrix will be the same as the MatType of the FIRST
504     matrix in the composite matrix.
505 
506 .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE
507 
508 @*/
509 PetscErrorCode  MatCompositeSetType(Mat mat,MatCompositeType type)
510 {
511   Mat_Composite  *b = (Mat_Composite*)mat->data;
512   PetscBool      flg;
513   PetscErrorCode ierr;
514 
515   PetscFunctionBegin;
516   ierr = PetscObjectTypeCompare((PetscObject)mat,MATCOMPOSITE,&flg);CHKERRQ(ierr);
517   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Can only use with composite matrix");
518   if (type == MAT_COMPOSITE_MULTIPLICATIVE) {
519     mat->ops->getdiagonal   = 0;
520     mat->ops->mult          = MatMult_Composite_Multiplicative;
521     mat->ops->multtranspose = MatMultTranspose_Composite_Multiplicative;
522     b->type                 = MAT_COMPOSITE_MULTIPLICATIVE;
523   } else {
524     mat->ops->getdiagonal   = MatGetDiagonal_Composite;
525     mat->ops->mult          = MatMult_Composite;
526     mat->ops->multtranspose = MatMultTranspose_Composite;
527     b->type                 = MAT_COMPOSITE_ADDITIVE;
528   }
529   PetscFunctionReturn(0);
530 }
531 
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "MatCompositeMerge"
535 /*@C
536    MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix
537      by summing all the matrices inside the composite matrix.
538 
539   Collective on MPI_Comm
540 
541    Input Parameters:
542 .  mat - the composite matrix
543 
544 
545    Options Database:
546 .  -mat_composite_merge  (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked)
547 
548    Level: advanced
549 
550    Notes:
551       The MatType of the resulting matrix will be the same as the MatType of the FIRST
552     matrix in the composite matrix.
553 
554 .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE
555 
556 @*/
557 PetscErrorCode  MatCompositeMerge(Mat mat)
558 {
559   Mat_Composite     *shell = (Mat_Composite*)mat->data;
560   Mat_CompositeLink next = shell->head, prev = shell->tail;
561   PetscErrorCode    ierr;
562   Mat               tmat,newmat;
563 
564   PetscFunctionBegin;
565   if (!next) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
566 
567   PetscFunctionBegin;
568   if (shell->type == MAT_COMPOSITE_ADDITIVE) {
569     ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr);
570     while ((next = next->next)) {
571       ierr = MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
572     }
573   } else {
574     ierr = MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);CHKERRQ(ierr);
575     while ((prev = prev->prev)) {
576       ierr = MatMatMult(tmat,prev->mat,MAT_INITIAL_MATRIX,PETSC_DECIDE,&newmat);CHKERRQ(ierr);
577       ierr = MatDestroy(&tmat);CHKERRQ(ierr);
578       tmat = newmat;
579     }
580   }
581   ierr = MatHeaderReplace(mat,tmat);CHKERRQ(ierr);
582   ierr = MatDiagonalScale(mat,shell->left,shell->right);CHKERRQ(ierr);
583   ierr = MatScale(mat,shell->scale);CHKERRQ(ierr);
584   PetscFunctionReturn(0);
585 }
586