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