Actual source code: mcomposite.c
1: #define PETSCMAT_DLL
3: #include include/private/matimpl.h
5: typedef struct _Mat_CompositeLink *Mat_CompositeLink;
6: struct _Mat_CompositeLink {
7: Mat mat;
8: Mat_CompositeLink next;
9: };
10:
11: typedef struct {
12: Mat_CompositeLink head;
13: Vec work;
14: } Mat_Composite;
18: PetscErrorCode MatDestroy_Composite(Mat mat)
19: {
20: PetscErrorCode ierr;
21: Mat_Composite *shell = (Mat_Composite*)mat->data;
22: Mat_CompositeLink next = shell->head;
25: while (next) {
26: MatDestroy(next->mat);
27: next = next->next;
28: }
29: if (shell->work) {VecDestroy(shell->work);}
30: PetscFree(shell);
31: mat->data = 0;
32: return(0);
33: }
37: PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y)
38: {
39: Mat_Composite *shell = (Mat_Composite*)A->data;
40: Mat_CompositeLink next = shell->head;
41: PetscErrorCode ierr;
44: if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
45: MatMult(next->mat,x,y);
46: while ((next = next->next)) {
47: MatMultAdd(next->mat,x,y,y);
48: }
49: return(0);
50: }
54: PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y)
55: {
56: Mat_Composite *shell = (Mat_Composite*)A->data;
57: Mat_CompositeLink next = shell->head;
58: PetscErrorCode ierr;
61: if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
62: MatMultTranspose(next->mat,x,y);
63: while ((next = next->next)) {
64: MatMultTransposeAdd(next->mat,x,y,y);
65: }
66: return(0);
67: }
71: PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v)
72: {
73: Mat_Composite *shell = (Mat_Composite*)A->data;
74: Mat_CompositeLink next = shell->head;
75: PetscErrorCode ierr;
78: if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
79: MatGetDiagonal(next->mat,v);
80: if (next->next && !shell->work) {
81: VecDuplicate(v,&shell->work);
82: }
83: while ((next = next->next)) {
84: MatGetDiagonal(next->mat,shell->work);
85: VecAXPY(v,1.0,shell->work);
86: }
87: return(0);
88: }
92: PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t)
93: {
95: PetscTruth flg;
98: PetscOptionsHasName(Y->prefix,"-mat_composite_merge",&flg);
99: if (flg) {
100: MatCompositeMerge(Y);
101: }
102: return(0);
103: }
106: static struct _MatOps MatOps_Values = {0,
107: 0,
108: 0,
109: MatMult_Composite,
110: 0,
111: /* 5*/ MatMultTranspose_Composite,
112: 0,
113: 0,
114: 0,
115: 0,
116: /*10*/ 0,
117: 0,
118: 0,
119: 0,
120: 0,
121: /*15*/ 0,
122: 0,
123: MatGetDiagonal_Composite,
124: 0,
125: 0,
126: /*20*/ 0,
127: MatAssemblyEnd_Composite,
128: 0,
129: 0,
130: 0,
131: /*25*/ 0,
132: 0,
133: 0,
134: 0,
135: 0,
136: /*30*/ 0,
137: 0,
138: 0,
139: 0,
140: 0,
141: /*35*/ 0,
142: 0,
143: 0,
144: 0,
145: 0,
146: /*40*/ 0,
147: 0,
148: 0,
149: 0,
150: 0,
151: /*45*/ 0,
152: 0,
153: 0,
154: 0,
155: 0,
156: /*50*/ 0,
157: 0,
158: 0,
159: 0,
160: 0,
161: /*55*/ 0,
162: 0,
163: 0,
164: 0,
165: 0,
166: /*60*/ 0,
167: MatDestroy_Composite,
168: 0,
169: 0,
170: 0,
171: /*65*/ 0,
172: 0,
173: 0,
174: 0,
175: 0,
176: /*70*/ 0,
177: 0,
178: 0,
179: 0,
180: 0,
181: /*75*/ 0,
182: 0,
183: 0,
184: 0,
185: 0,
186: /*80*/ 0,
187: 0,
188: 0,
189: 0,
190: 0,
191: /*85*/ 0,
192: 0,
193: 0,
194: 0,
195: 0,
196: /*90*/ 0,
197: 0,
198: 0,
199: 0,
200: 0,
201: /*95*/ 0,
202: 0,
203: 0,
204: 0};
206: /*MC
207: MATCOMPOSITE - A matrix defined by the sum of one or more matrices (all matrices are of same size and parallel layout)
209: Level: advanced
211: .seealso: MatCreateComposite(), MatCompositeAddMat(), MatSetType(), MatCompositeMerge()
212: M*/
217: PetscErrorCode MatCreate_Composite(Mat A)
218: {
219: Mat_Composite *b;
223: PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));
225: PetscNew(Mat_Composite,&b);
226: A->data = (void*)b;
228: PetscMapSetBlockSize(&A->rmap,1);
229: PetscMapSetBlockSize(&A->cmap,1);
230: PetscMapSetUp(&A->rmap);
231: PetscMapSetUp(&A->cmap);
233: A->assembled = PETSC_TRUE;
234: A->preallocated = PETSC_FALSE;
235: PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);
236: return(0);
237: }
242: /*@C
243: MatCreateComposite - Creates a matrix as the sum of zero or more matrices
245: Collective on MPI_Comm
247: Input Parameters:
248: + comm - MPI communicator
249: . nmat - number of matrices to put in
250: - mats - the matrices
252: Output Parameter:
253: . A - the matrix
255: Level: advanced
257: Notes:
258: Alternative construction
259: $ MatCreate(comm,&mat);
260: $ MatSetSizes(mat,m,n,M,N);
261: $ MatSetType(mat,MATCOMPOSITE);
262: $ MatCompositeAddMat(mat,mats[0]);
263: $ ....
264: $ MatCompositeAddMat(mat,mats[nmat-1]);
265: $ MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
266: $ MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
268: .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCompositeMerge()
270: @*/
271: PetscErrorCode MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat)
272: {
274: PetscInt m,n,M,N,i;
277: if (nmat < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix");
280: MatGetLocalSize(mats[0],&m,&n);
281: MatGetSize(mats[0],&M,&N);
282: MatCreate(comm,mat);
283: MatSetSizes(*mat,m,n,M,N);
284: MatSetType(*mat,MATCOMPOSITE);
285: for (i=0; i<nmat; i++) {
286: MatCompositeAddMat(*mat,mats[i]);
287: }
288: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
289: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
290: return(0);
291: }
295: /*@
296: MatCompositeAddMat - add another matrix to a composite matrix
298: Collective on Mat
300: Input Parameters:
301: + mat - the composite matrix
302: - smat - the partial matrix
304: Level: advanced
306: .seealso: MatCreateComposite()
307: @*/
308: PetscErrorCode MatCompositeAddMat(Mat mat,Mat smat)
309: {
310: Mat_Composite *shell = (Mat_Composite*)mat->data;
311: PetscErrorCode ierr;
312: Mat_CompositeLink ilink,next;
317: PetscNew(struct _Mat_CompositeLink,&ilink);
318: ilink->next = 0;
319: PetscObjectReference((PetscObject)smat);
320: ilink->mat = smat;
322: next = shell->head;
323: if (!next) {
324: shell->head = ilink;
325: } else {
326: while (next->next) {
327: next = next->next;
328: }
329: next->next = ilink;
330: }
331: return(0);
332: }
336: /*@C
337: MatCompositeMerge - Given a composite matrix, replaces it with a "regular" matrix
338: by summing all the matrices inside the composite matrix.
340: Collective on MPI_Comm
342: Input Parameters:
343: . mat - the composite matrix
346: Options Database:
347: . -mat_composite_merge (you must call MatAssemblyBegin()/MatAssemblyEnd() to have this checked)
349: Level: advanced
351: Notes:
352: The MatType of the resulting matrix will be the same as the MatType of the FIRST
353: matrix in the composite matrix.
355: .seealso: MatDestroy(), MatMult(), MatCompositeAddMat(), MatCreateComposite(), MATCOMPOSITE
357: @*/
358: PetscErrorCode MatCompositeMerge(Mat mat)
359: {
360: Mat_Composite *shell = (Mat_Composite*)mat->data;
361: Mat_CompositeLink next = shell->head;
362: PetscErrorCode ierr;
363: Mat tmat;
366: if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
369: MatDuplicate(next->mat,MAT_COPY_VALUES,&tmat);
370: while ((next = next->next)) {
371: MatAXPY(tmat,1.0,next->mat,DIFFERENT_NONZERO_PATTERN);
372: }
373: MatDestroy_Composite(mat);
375: return(0);
376: }