xref: /petsc/src/mat/impls/composite/mcomposite.c (revision 793850ff306cb72643dd3e80f2be5928c3f41114)
1*793850ffSBarry Smith #define PETSCMAT_DLL
2*793850ffSBarry Smith 
3*793850ffSBarry Smith #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
4*793850ffSBarry Smith 
5*793850ffSBarry Smith typedef struct _Mat_CompositeLink *Mat_CompositeLink;
6*793850ffSBarry Smith struct _Mat_CompositeLink {
7*793850ffSBarry Smith   Mat               mat;
8*793850ffSBarry Smith   Mat_CompositeLink next;
9*793850ffSBarry Smith };
10*793850ffSBarry Smith 
11*793850ffSBarry Smith typedef struct {
12*793850ffSBarry Smith   Mat_CompositeLink head;
13*793850ffSBarry Smith   Vec               work;
14*793850ffSBarry Smith } Mat_Composite;
15*793850ffSBarry Smith 
16*793850ffSBarry Smith #undef __FUNCT__
17*793850ffSBarry Smith #define __FUNCT__ "MatDestroy_Composite"
18*793850ffSBarry Smith PetscErrorCode MatDestroy_Composite(Mat mat)
19*793850ffSBarry Smith {
20*793850ffSBarry Smith   PetscErrorCode   ierr;
21*793850ffSBarry Smith   Mat_Composite    *shell = (Mat_Composite*)mat->data;;
22*793850ffSBarry Smith   Mat_CompositeLink next = shell->head;
23*793850ffSBarry Smith 
24*793850ffSBarry Smith   PetscFunctionBegin;
25*793850ffSBarry Smith   while (next) {
26*793850ffSBarry Smith     ierr = MatDestroy(next->mat);CHKERRQ(ierr);
27*793850ffSBarry Smith     next = next->next;
28*793850ffSBarry Smith   }
29*793850ffSBarry Smith   if (shell->work) {ierr = VecDestroy(shell->work);CHKERRQ(ierr);}
30*793850ffSBarry Smith   ierr      = PetscFree(shell);CHKERRQ(ierr);
31*793850ffSBarry Smith   mat->data = 0;
32*793850ffSBarry Smith   PetscFunctionReturn(0);
33*793850ffSBarry Smith }
34*793850ffSBarry Smith 
35*793850ffSBarry Smith #undef __FUNCT__
36*793850ffSBarry Smith #define __FUNCT__ "MatMult_Composite"
37*793850ffSBarry Smith PetscErrorCode MatMult_Composite(Mat A,Vec x,Vec y)
38*793850ffSBarry Smith {
39*793850ffSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
40*793850ffSBarry Smith   Mat_CompositeLink next = shell->head;
41*793850ffSBarry Smith   PetscErrorCode    ierr;
42*793850ffSBarry Smith 
43*793850ffSBarry Smith   PetscFunctionBegin;
44*793850ffSBarry Smith   if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
45*793850ffSBarry Smith   ierr = MatMult(next->mat,x,y);CHKERRQ(ierr);
46*793850ffSBarry Smith   while ((next = next->next)) {
47*793850ffSBarry Smith     ierr = MatMultAdd(next->mat,x,y,y);CHKERRQ(ierr);
48*793850ffSBarry Smith   }
49*793850ffSBarry Smith   PetscFunctionReturn(0);
50*793850ffSBarry Smith }
51*793850ffSBarry Smith 
52*793850ffSBarry Smith #undef __FUNCT__
53*793850ffSBarry Smith #define __FUNCT__ "MatMultTranspose_Composite"
54*793850ffSBarry Smith PetscErrorCode MatMultTranspose_Composite(Mat A,Vec x,Vec y)
55*793850ffSBarry Smith {
56*793850ffSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
57*793850ffSBarry Smith   Mat_CompositeLink next = shell->head;
58*793850ffSBarry Smith   PetscErrorCode    ierr;
59*793850ffSBarry Smith 
60*793850ffSBarry Smith   PetscFunctionBegin;
61*793850ffSBarry Smith   if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
62*793850ffSBarry Smith   ierr = MatMultTranspose(next->mat,x,y);CHKERRQ(ierr);
63*793850ffSBarry Smith   while ((next = next->next)) {
64*793850ffSBarry Smith     ierr = MatMultTransposeAdd(next->mat,x,y,y);CHKERRQ(ierr);
65*793850ffSBarry Smith   }
66*793850ffSBarry Smith   PetscFunctionReturn(0);
67*793850ffSBarry Smith }
68*793850ffSBarry Smith 
69*793850ffSBarry Smith #undef __FUNCT__
70*793850ffSBarry Smith #define __FUNCT__ "MatGetDiagonal_Composite"
71*793850ffSBarry Smith PetscErrorCode MatGetDiagonal_Composite(Mat A,Vec v)
72*793850ffSBarry Smith {
73*793850ffSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)A->data;
74*793850ffSBarry Smith   Mat_CompositeLink next = shell->head;
75*793850ffSBarry Smith   PetscErrorCode    ierr;
76*793850ffSBarry Smith 
77*793850ffSBarry Smith   PetscFunctionBegin;
78*793850ffSBarry Smith   if (!next) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must provide at least one matrix with MatCompositeAddMat()");
79*793850ffSBarry Smith   ierr = MatGetDiagonal(next->mat,v);CHKERRQ(ierr);
80*793850ffSBarry Smith   if (next->next && !shell->work) {
81*793850ffSBarry Smith     ierr = VecDuplicate(v,&shell->work);CHKERRQ(ierr);
82*793850ffSBarry Smith   }
83*793850ffSBarry Smith   while ((next = next->next)) {
84*793850ffSBarry Smith     ierr = MatGetDiagonal(next->mat,shell->work);CHKERRQ(ierr);
85*793850ffSBarry Smith     ierr = VecAXPY(v,1.0,shell->work);CHKERRQ(ierr);
86*793850ffSBarry Smith   }
87*793850ffSBarry Smith   PetscFunctionReturn(0);
88*793850ffSBarry Smith }
89*793850ffSBarry Smith 
90*793850ffSBarry Smith #undef __FUNCT__
91*793850ffSBarry Smith #define __FUNCT__ "MatAssemblyEnd_Composite"
92*793850ffSBarry Smith PetscErrorCode MatAssemblyEnd_Composite(Mat Y,MatAssemblyType t)
93*793850ffSBarry Smith {
94*793850ffSBarry Smith   PetscFunctionBegin;
95*793850ffSBarry Smith   PetscFunctionReturn(0);
96*793850ffSBarry Smith }
97*793850ffSBarry Smith 
98*793850ffSBarry Smith 
99*793850ffSBarry Smith static struct _MatOps MatOps_Values = {0,
100*793850ffSBarry Smith        0,
101*793850ffSBarry Smith        0,
102*793850ffSBarry Smith        MatMult_Composite,
103*793850ffSBarry Smith        0,
104*793850ffSBarry Smith /* 5*/ MatMultTranspose_Composite,
105*793850ffSBarry Smith        0,
106*793850ffSBarry Smith        0,
107*793850ffSBarry Smith        0,
108*793850ffSBarry Smith        0,
109*793850ffSBarry Smith /*10*/ 0,
110*793850ffSBarry Smith        0,
111*793850ffSBarry Smith        0,
112*793850ffSBarry Smith        0,
113*793850ffSBarry Smith        0,
114*793850ffSBarry Smith /*15*/ 0,
115*793850ffSBarry Smith        0,
116*793850ffSBarry Smith        MatGetDiagonal_Composite,
117*793850ffSBarry Smith        0,
118*793850ffSBarry Smith        0,
119*793850ffSBarry Smith /*20*/ 0,
120*793850ffSBarry Smith        MatAssemblyEnd_Composite,
121*793850ffSBarry Smith        0,
122*793850ffSBarry Smith        0,
123*793850ffSBarry Smith        0,
124*793850ffSBarry Smith /*25*/ 0,
125*793850ffSBarry Smith        0,
126*793850ffSBarry Smith        0,
127*793850ffSBarry Smith        0,
128*793850ffSBarry Smith        0,
129*793850ffSBarry Smith /*30*/ 0,
130*793850ffSBarry Smith        0,
131*793850ffSBarry Smith        0,
132*793850ffSBarry Smith        0,
133*793850ffSBarry Smith        0,
134*793850ffSBarry Smith /*35*/ 0,
135*793850ffSBarry Smith        0,
136*793850ffSBarry Smith        0,
137*793850ffSBarry Smith        0,
138*793850ffSBarry Smith        0,
139*793850ffSBarry Smith /*40*/ 0,
140*793850ffSBarry Smith        0,
141*793850ffSBarry Smith        0,
142*793850ffSBarry Smith        0,
143*793850ffSBarry Smith        0,
144*793850ffSBarry Smith /*45*/ 0,
145*793850ffSBarry Smith        0,
146*793850ffSBarry Smith        0,
147*793850ffSBarry Smith        0,
148*793850ffSBarry Smith        0,
149*793850ffSBarry Smith /*50*/ 0,
150*793850ffSBarry Smith        0,
151*793850ffSBarry Smith        0,
152*793850ffSBarry Smith        0,
153*793850ffSBarry Smith        0,
154*793850ffSBarry Smith /*55*/ 0,
155*793850ffSBarry Smith        0,
156*793850ffSBarry Smith        0,
157*793850ffSBarry Smith        0,
158*793850ffSBarry Smith        0,
159*793850ffSBarry Smith /*60*/ 0,
160*793850ffSBarry Smith        MatDestroy_Composite,
161*793850ffSBarry Smith        0,
162*793850ffSBarry Smith        0,
163*793850ffSBarry Smith        0,
164*793850ffSBarry Smith /*65*/ 0,
165*793850ffSBarry Smith        0,
166*793850ffSBarry Smith        0,
167*793850ffSBarry Smith        0,
168*793850ffSBarry Smith        0,
169*793850ffSBarry Smith /*70*/ 0,
170*793850ffSBarry Smith        0,
171*793850ffSBarry Smith        0,
172*793850ffSBarry Smith        0,
173*793850ffSBarry Smith        0,
174*793850ffSBarry Smith /*75*/ 0,
175*793850ffSBarry Smith        0,
176*793850ffSBarry Smith        0,
177*793850ffSBarry Smith        0,
178*793850ffSBarry Smith        0,
179*793850ffSBarry Smith /*80*/ 0,
180*793850ffSBarry Smith        0,
181*793850ffSBarry Smith        0,
182*793850ffSBarry Smith        0,
183*793850ffSBarry Smith        0,
184*793850ffSBarry Smith /*85*/ 0,
185*793850ffSBarry Smith        0,
186*793850ffSBarry Smith        0,
187*793850ffSBarry Smith        0,
188*793850ffSBarry Smith        0,
189*793850ffSBarry Smith /*90*/ 0,
190*793850ffSBarry Smith        0,
191*793850ffSBarry Smith        0,
192*793850ffSBarry Smith        0,
193*793850ffSBarry Smith        0,
194*793850ffSBarry Smith /*95*/ 0,
195*793850ffSBarry Smith        0,
196*793850ffSBarry Smith        0,
197*793850ffSBarry Smith        0};
198*793850ffSBarry Smith 
199*793850ffSBarry Smith /*MC
200*793850ffSBarry Smith    MATCOMPOSITE - A matrix defined by the sum of one or more matrices (all matrices are of same size and parallel layout)
201*793850ffSBarry Smith 
202*793850ffSBarry Smith   Level: advanced
203*793850ffSBarry Smith 
204*793850ffSBarry Smith .seealso: MatCreateComposite
205*793850ffSBarry Smith M*/
206*793850ffSBarry Smith 
207*793850ffSBarry Smith EXTERN_C_BEGIN
208*793850ffSBarry Smith #undef __FUNCT__
209*793850ffSBarry Smith #define __FUNCT__ "MatCreate_Composite"
210*793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Composite(Mat A)
211*793850ffSBarry Smith {
212*793850ffSBarry Smith   Mat_Composite  *b;
213*793850ffSBarry Smith   PetscErrorCode ierr;
214*793850ffSBarry Smith 
215*793850ffSBarry Smith   PetscFunctionBegin;
216*793850ffSBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
217*793850ffSBarry Smith 
218*793850ffSBarry Smith   ierr = PetscNew(Mat_Composite,&b);CHKERRQ(ierr);
219*793850ffSBarry Smith   A->data = (void*)b;
220*793850ffSBarry Smith 
221*793850ffSBarry Smith   ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr);
222*793850ffSBarry Smith   ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr);
223*793850ffSBarry Smith 
224*793850ffSBarry Smith   A->assembled     = PETSC_TRUE;
225*793850ffSBarry Smith   A->preallocated  = PETSC_FALSE;
226*793850ffSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATCOMPOSITE);CHKERRQ(ierr);
227*793850ffSBarry Smith   PetscFunctionReturn(0);
228*793850ffSBarry Smith }
229*793850ffSBarry Smith EXTERN_C_END
230*793850ffSBarry Smith 
231*793850ffSBarry Smith #undef __FUNCT__
232*793850ffSBarry Smith #define __FUNCT__ "MatCreateComposite"
233*793850ffSBarry Smith /*@C
234*793850ffSBarry Smith    MatCreateComposite - Creates a matrix as the sum of zero or more matrices
235*793850ffSBarry Smith 
236*793850ffSBarry Smith   Collective on MPI_Comm
237*793850ffSBarry Smith 
238*793850ffSBarry Smith    Input Parameters:
239*793850ffSBarry Smith +  comm - MPI communicator
240*793850ffSBarry Smith .  nmat - number of matrices to put in
241*793850ffSBarry Smith -  mats - the matrices
242*793850ffSBarry Smith 
243*793850ffSBarry Smith    Output Parameter:
244*793850ffSBarry Smith .  A - the matrix
245*793850ffSBarry Smith 
246*793850ffSBarry Smith    Level: advanced
247*793850ffSBarry Smith 
248*793850ffSBarry Smith    Notes:
249*793850ffSBarry Smith      Alternative construction
250*793850ffSBarry Smith $       MatCreate(comm,&mat);
251*793850ffSBarry Smith $       MatSetSizes(mat,m,n,M,N);
252*793850ffSBarry Smith $       MatSetType(mat,MATCOMPOSITE);
253*793850ffSBarry Smith $       MatCompositeAddMat(mat,mats[0]);
254*793850ffSBarry Smith $       ....
255*793850ffSBarry Smith $       MatCompositeAddMat(mat,mats[nmat-1]);
256*793850ffSBarry Smith 
257*793850ffSBarry Smith .seealso: MatDestroy(), MatMult(), MatCompositeAddMat()
258*793850ffSBarry Smith 
259*793850ffSBarry Smith @*/
260*793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm comm,PetscInt nmat,const Mat *mats,Mat *mat)
261*793850ffSBarry Smith {
262*793850ffSBarry Smith   PetscErrorCode ierr;
263*793850ffSBarry Smith   PetscInt       m,n,M,N,i;
264*793850ffSBarry Smith 
265*793850ffSBarry Smith   PetscFunctionBegin;
266*793850ffSBarry Smith   if (nmat < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must pass in at least one matrix");
267*793850ffSBarry Smith   PetscValidPointer(mat,3);CHKERRQ(ierr);
268*793850ffSBarry Smith 
269*793850ffSBarry Smith   ierr = MatGetLocalSize(mats[0],&m,&n);CHKERRQ(ierr);
270*793850ffSBarry Smith   ierr = MatGetSize(mats[0],&M,&N);CHKERRQ(ierr);
271*793850ffSBarry Smith   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
272*793850ffSBarry Smith   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
273*793850ffSBarry Smith   ierr = MatSetType(*mat,MATCOMPOSITE);CHKERRQ(ierr);
274*793850ffSBarry Smith   for (i=0; i<nmat; i++) {
275*793850ffSBarry Smith     ierr = MatCompositeAddMat(*mat,mats[i]);CHKERRQ(ierr);
276*793850ffSBarry Smith   }
277*793850ffSBarry Smith   PetscFunctionReturn(0);
278*793850ffSBarry Smith }
279*793850ffSBarry Smith 
280*793850ffSBarry Smith #undef __FUNCT__
281*793850ffSBarry Smith #define __FUNCT__ "MatCompositeAddMat"
282*793850ffSBarry Smith /*@
283*793850ffSBarry Smith     MatCompositeAddMat - add another matrix to a composite matrix
284*793850ffSBarry Smith 
285*793850ffSBarry Smith    Collective on Mat
286*793850ffSBarry Smith 
287*793850ffSBarry Smith     Input Parameters:
288*793850ffSBarry Smith +   mat - the composite matrix
289*793850ffSBarry Smith -   smat - the partial matrix
290*793850ffSBarry Smith 
291*793850ffSBarry Smith    Level: advanced
292*793850ffSBarry Smith 
293*793850ffSBarry Smith .seealso: MatCreateComposite()
294*793850ffSBarry Smith @*/
295*793850ffSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat mat,Mat smat)
296*793850ffSBarry Smith {
297*793850ffSBarry Smith   Mat_Composite     *shell = (Mat_Composite*)mat->data;
298*793850ffSBarry Smith   PetscErrorCode    ierr;
299*793850ffSBarry Smith   Mat_CompositeLink ilink,next;
300*793850ffSBarry Smith 
301*793850ffSBarry Smith   PetscFunctionBegin;
302*793850ffSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
303*793850ffSBarry Smith   PetscValidHeaderSpecific(smat,MAT_COOKIE,2);
304*793850ffSBarry Smith   ierr       = PetscNew(struct _Mat_CompositeLink,&ilink);CHKERRQ(ierr);
305*793850ffSBarry Smith   ilink->next = 0;
306*793850ffSBarry Smith   ilink->mat  = smat;
307*793850ffSBarry Smith   ierr        = PetscObjectReference((PetscObject)smat);CHKERRQ(ierr);
308*793850ffSBarry Smith 
309*793850ffSBarry Smith   next = shell->head;
310*793850ffSBarry Smith   if (!next) {
311*793850ffSBarry Smith     shell->head  = ilink;
312*793850ffSBarry Smith   } else {
313*793850ffSBarry Smith     while (next->next) {
314*793850ffSBarry Smith       next = next->next;
315*793850ffSBarry Smith     }
316*793850ffSBarry Smith     next->next      = ilink;
317*793850ffSBarry Smith   }
318*793850ffSBarry Smith   PetscFunctionReturn(0);
319*793850ffSBarry Smith }
320*793850ffSBarry Smith 
321