xref: /petsc/src/mat/impls/scatter/mscatter.c (revision 2a6744eb01855f5aa328eb8fdf4b0d01e72ad151)
1*2a6744ebSBarry Smith #define PETSCMAT_DLL
2*2a6744ebSBarry Smith 
3*2a6744ebSBarry Smith /*
4*2a6744ebSBarry Smith    This provides a matrix that applies a VecScatter to a vector.
5*2a6744ebSBarry Smith */
6*2a6744ebSBarry Smith 
7*2a6744ebSBarry Smith #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
8*2a6744ebSBarry Smith #include "private/vecimpl.h"
9*2a6744ebSBarry Smith 
10*2a6744ebSBarry Smith typedef struct {
11*2a6744ebSBarry Smith   VecScatter scatter;
12*2a6744ebSBarry Smith } Mat_Scatter;
13*2a6744ebSBarry Smith 
14*2a6744ebSBarry Smith #undef __FUNCT__
15*2a6744ebSBarry Smith #define __FUNCT__ "MatScatterGetVecScatter"
16*2a6744ebSBarry Smith /*@
17*2a6744ebSBarry Smith     MatScatterGetVecScatter - Returns the user-provided scatter set with MatScatterSetVecScatter()
18*2a6744ebSBarry Smith 
19*2a6744ebSBarry Smith     Not Collective, but not cannot use scatter if not used collectively on Mat
20*2a6744ebSBarry Smith 
21*2a6744ebSBarry Smith     Input Parameter:
22*2a6744ebSBarry Smith .   mat - the matrix, should have been created with MatCreateScatter() or have type MATSCATTER
23*2a6744ebSBarry Smith 
24*2a6744ebSBarry Smith     Output Parameter:
25*2a6744ebSBarry Smith .   scatter - the scatter context
26*2a6744ebSBarry Smith 
27*2a6744ebSBarry Smith     Level: intermediate
28*2a6744ebSBarry Smith 
29*2a6744ebSBarry Smith .keywords: matrix, scatter, get
30*2a6744ebSBarry Smith 
31*2a6744ebSBarry Smith .seealso: MatCreateScatter(), MatScatterSetVecScatter(), MATSCATTER
32*2a6744ebSBarry Smith @*/
33*2a6744ebSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat mat,VecScatter *scatter)
34*2a6744ebSBarry Smith {
35*2a6744ebSBarry Smith   Mat_Scatter    *mscatter;
36*2a6744ebSBarry Smith 
37*2a6744ebSBarry Smith   PetscFunctionBegin;
38*2a6744ebSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
39*2a6744ebSBarry Smith   PetscValidPointer(scatter,2);
40*2a6744ebSBarry Smith   mscatter = (Mat_Scatter*)mat->data;
41*2a6744ebSBarry Smith   *scatter = mscatter->scatter;
42*2a6744ebSBarry Smith   PetscFunctionReturn(0);
43*2a6744ebSBarry Smith }
44*2a6744ebSBarry Smith 
45*2a6744ebSBarry Smith #undef __FUNCT__
46*2a6744ebSBarry Smith #define __FUNCT__ "MatDestroy_Scatter"
47*2a6744ebSBarry Smith PetscErrorCode MatDestroy_Scatter(Mat mat)
48*2a6744ebSBarry Smith {
49*2a6744ebSBarry Smith   PetscErrorCode ierr;
50*2a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)mat->data;
51*2a6744ebSBarry Smith 
52*2a6744ebSBarry Smith   PetscFunctionBegin;
53*2a6744ebSBarry Smith   if (scatter->scatter) {ierr = VecScatterDestroy(scatter->scatter);CHKERRQ(ierr);}
54*2a6744ebSBarry Smith   ierr = PetscFree(scatter);CHKERRQ(ierr);
55*2a6744ebSBarry Smith   PetscFunctionReturn(0);
56*2a6744ebSBarry Smith }
57*2a6744ebSBarry Smith 
58*2a6744ebSBarry Smith #undef __FUNCT__
59*2a6744ebSBarry Smith #define __FUNCT__ "MatMult_Scatter"
60*2a6744ebSBarry Smith PetscErrorCode MatMult_Scatter(Mat A,Vec x,Vec y)
61*2a6744ebSBarry Smith {
62*2a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
63*2a6744ebSBarry Smith   PetscErrorCode ierr;
64*2a6744ebSBarry Smith 
65*2a6744ebSBarry Smith   PetscFunctionBegin;
66*2a6744ebSBarry Smith   if (!scatter->scatter) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
67*2a6744ebSBarry Smith   ierr = VecScatterBegin(x,y,INSERT_VALUES,SCATTER_FORWARD,scatter->scatter);CHKERRQ(ierr);
68*2a6744ebSBarry Smith   ierr = VecScatterEnd(x,y,INSERT_VALUES,SCATTER_FORWARD,scatter->scatter);CHKERRQ(ierr);
69*2a6744ebSBarry Smith   PetscFunctionReturn(0);
70*2a6744ebSBarry Smith }
71*2a6744ebSBarry Smith 
72*2a6744ebSBarry Smith #undef __FUNCT__
73*2a6744ebSBarry Smith #define __FUNCT__ "MatMultAdd_Scatter"
74*2a6744ebSBarry Smith PetscErrorCode MatMultAdd_Scatter(Mat A,Vec x,Vec y,Vec z)
75*2a6744ebSBarry Smith {
76*2a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
77*2a6744ebSBarry Smith   PetscErrorCode ierr;
78*2a6744ebSBarry Smith 
79*2a6744ebSBarry Smith   PetscFunctionBegin;
80*2a6744ebSBarry Smith   if (!scatter->scatter) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
81*2a6744ebSBarry Smith   if (z != y) {ierr = VecCopy(y,z);CHKERRQ(ierr);}
82*2a6744ebSBarry Smith   ierr = VecScatterBegin(x,z,ADD_VALUES,SCATTER_FORWARD,scatter->scatter);CHKERRQ(ierr);
83*2a6744ebSBarry Smith   ierr = VecScatterEnd(x,z,ADD_VALUES,SCATTER_FORWARD,scatter->scatter);CHKERRQ(ierr);
84*2a6744ebSBarry Smith   PetscFunctionReturn(0);
85*2a6744ebSBarry Smith }
86*2a6744ebSBarry Smith 
87*2a6744ebSBarry Smith #undef __FUNCT__
88*2a6744ebSBarry Smith #define __FUNCT__ "MatMultTranspose_Scatter"
89*2a6744ebSBarry Smith PetscErrorCode MatMultTranspose_Scatter(Mat A,Vec x,Vec y)
90*2a6744ebSBarry Smith {
91*2a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
92*2a6744ebSBarry Smith   PetscErrorCode ierr;
93*2a6744ebSBarry Smith 
94*2a6744ebSBarry Smith   PetscFunctionBegin;
95*2a6744ebSBarry Smith   if (!scatter->scatter) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
96*2a6744ebSBarry Smith   ierr = VecScatterBegin(x,y,INSERT_VALUES,SCATTER_REVERSE,scatter->scatter);CHKERRQ(ierr);
97*2a6744ebSBarry Smith   ierr = VecScatterEnd(x,y,INSERT_VALUES,SCATTER_REVERSE,scatter->scatter);CHKERRQ(ierr);
98*2a6744ebSBarry Smith   PetscFunctionReturn(0);
99*2a6744ebSBarry Smith }
100*2a6744ebSBarry Smith 
101*2a6744ebSBarry Smith #undef __FUNCT__
102*2a6744ebSBarry Smith #define __FUNCT__ "MatMultTransposeAdd_Scatter"
103*2a6744ebSBarry Smith PetscErrorCode MatMultTransposeAdd_Scatter(Mat A,Vec x,Vec y,Vec z)
104*2a6744ebSBarry Smith {
105*2a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
106*2a6744ebSBarry Smith   PetscErrorCode ierr;
107*2a6744ebSBarry Smith 
108*2a6744ebSBarry Smith   PetscFunctionBegin;
109*2a6744ebSBarry Smith   if (!scatter->scatter) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
110*2a6744ebSBarry Smith   if (z != y) {ierr = VecCopy(y,z);CHKERRQ(ierr);}
111*2a6744ebSBarry Smith   ierr = VecScatterBegin(x,z,ADD_VALUES,SCATTER_REVERSE,scatter->scatter);CHKERRQ(ierr);
112*2a6744ebSBarry Smith   ierr = VecScatterEnd(x,z,ADD_VALUES,SCATTER_REVERSE,scatter->scatter);CHKERRQ(ierr);
113*2a6744ebSBarry Smith   PetscFunctionReturn(0);
114*2a6744ebSBarry Smith }
115*2a6744ebSBarry Smith 
116*2a6744ebSBarry Smith static struct _MatOps MatOps_Values = {0,
117*2a6744ebSBarry Smith        0,
118*2a6744ebSBarry Smith        0,
119*2a6744ebSBarry Smith        MatMult_Scatter,
120*2a6744ebSBarry Smith /* 4*/ MatMultAdd_Scatter,
121*2a6744ebSBarry Smith        MatMultTranspose_Scatter,
122*2a6744ebSBarry Smith        MatMultTransposeAdd_Scatter,
123*2a6744ebSBarry Smith        0,
124*2a6744ebSBarry Smith        0,
125*2a6744ebSBarry Smith        0,
126*2a6744ebSBarry Smith /*10*/ 0,
127*2a6744ebSBarry Smith        0,
128*2a6744ebSBarry Smith        0,
129*2a6744ebSBarry Smith        0,
130*2a6744ebSBarry Smith        0,
131*2a6744ebSBarry Smith /*15*/ 0,
132*2a6744ebSBarry Smith        0,
133*2a6744ebSBarry Smith        0,
134*2a6744ebSBarry Smith        0,
135*2a6744ebSBarry Smith        0,
136*2a6744ebSBarry Smith /*20*/ 0,
137*2a6744ebSBarry Smith        0,
138*2a6744ebSBarry Smith        0,
139*2a6744ebSBarry Smith        0,
140*2a6744ebSBarry Smith        0,
141*2a6744ebSBarry Smith /*25*/ 0,
142*2a6744ebSBarry Smith        0,
143*2a6744ebSBarry Smith        0,
144*2a6744ebSBarry Smith        0,
145*2a6744ebSBarry Smith        0,
146*2a6744ebSBarry Smith /*30*/ 0,
147*2a6744ebSBarry Smith        0,
148*2a6744ebSBarry Smith        0,
149*2a6744ebSBarry Smith        0,
150*2a6744ebSBarry Smith        0,
151*2a6744ebSBarry Smith /*35*/ 0,
152*2a6744ebSBarry Smith        0,
153*2a6744ebSBarry Smith        0,
154*2a6744ebSBarry Smith        0,
155*2a6744ebSBarry Smith        0,
156*2a6744ebSBarry Smith /*40*/ 0,
157*2a6744ebSBarry Smith        0,
158*2a6744ebSBarry Smith        0,
159*2a6744ebSBarry Smith        0,
160*2a6744ebSBarry Smith        0,
161*2a6744ebSBarry Smith /*45*/ 0,
162*2a6744ebSBarry Smith        0,
163*2a6744ebSBarry Smith        0,
164*2a6744ebSBarry Smith        0,
165*2a6744ebSBarry Smith        0,
166*2a6744ebSBarry Smith /*50*/ 0,
167*2a6744ebSBarry Smith        0,
168*2a6744ebSBarry Smith        0,
169*2a6744ebSBarry Smith        0,
170*2a6744ebSBarry Smith        0,
171*2a6744ebSBarry Smith /*55*/ 0,
172*2a6744ebSBarry Smith        0,
173*2a6744ebSBarry Smith        0,
174*2a6744ebSBarry Smith        0,
175*2a6744ebSBarry Smith        0,
176*2a6744ebSBarry Smith /*60*/ 0,
177*2a6744ebSBarry Smith        MatDestroy_Scatter,
178*2a6744ebSBarry Smith        0,
179*2a6744ebSBarry Smith        0,
180*2a6744ebSBarry Smith        0,
181*2a6744ebSBarry Smith /*65*/ 0,
182*2a6744ebSBarry Smith        0,
183*2a6744ebSBarry Smith        0,
184*2a6744ebSBarry Smith        0,
185*2a6744ebSBarry Smith        0,
186*2a6744ebSBarry Smith /*70*/ 0,
187*2a6744ebSBarry Smith        0,
188*2a6744ebSBarry Smith        0,
189*2a6744ebSBarry Smith        0,
190*2a6744ebSBarry Smith        0,
191*2a6744ebSBarry Smith /*75*/ 0,
192*2a6744ebSBarry Smith        0,
193*2a6744ebSBarry Smith        0,
194*2a6744ebSBarry Smith        0,
195*2a6744ebSBarry Smith        0,
196*2a6744ebSBarry Smith /*80*/ 0,
197*2a6744ebSBarry Smith        0,
198*2a6744ebSBarry Smith        0,
199*2a6744ebSBarry Smith        0,
200*2a6744ebSBarry Smith        0,
201*2a6744ebSBarry Smith /*85*/ 0,
202*2a6744ebSBarry Smith        0,
203*2a6744ebSBarry Smith        0,
204*2a6744ebSBarry Smith        0,
205*2a6744ebSBarry Smith        0,
206*2a6744ebSBarry Smith /*90*/ 0,
207*2a6744ebSBarry Smith        0,
208*2a6744ebSBarry Smith        0,
209*2a6744ebSBarry Smith        0,
210*2a6744ebSBarry Smith        0,
211*2a6744ebSBarry Smith /*95*/ 0,
212*2a6744ebSBarry Smith        0,
213*2a6744ebSBarry Smith        0,
214*2a6744ebSBarry Smith        0};
215*2a6744ebSBarry Smith 
216*2a6744ebSBarry Smith /*MC
217*2a6744ebSBarry Smith    MATSCATTER - MATSCATTER = "scatter" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
218*2a6744ebSBarry Smith 
219*2a6744ebSBarry Smith   Level: advanced
220*2a6744ebSBarry Smith 
221*2a6744ebSBarry Smith .seealso: MatCreateScatter(), MatScatterSetVecScatter(), MatScatterGetVecScatter()
222*2a6744ebSBarry Smith 
223*2a6744ebSBarry Smith M*/
224*2a6744ebSBarry Smith 
225*2a6744ebSBarry Smith EXTERN_C_BEGIN
226*2a6744ebSBarry Smith #undef __FUNCT__
227*2a6744ebSBarry Smith #define __FUNCT__ "MatCreate_Scatter"
228*2a6744ebSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Scatter(Mat A)
229*2a6744ebSBarry Smith {
230*2a6744ebSBarry Smith   Mat_Scatter    *b;
231*2a6744ebSBarry Smith   PetscErrorCode ierr;
232*2a6744ebSBarry Smith 
233*2a6744ebSBarry Smith   PetscFunctionBegin;
234*2a6744ebSBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
235*2a6744ebSBarry Smith   ierr = PetscNew(Mat_Scatter,&b);CHKERRQ(ierr);
236*2a6744ebSBarry Smith 
237*2a6744ebSBarry Smith   A->data = (void*)b;
238*2a6744ebSBarry Smith 
239*2a6744ebSBarry Smith   ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr);
240*2a6744ebSBarry Smith   ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr);
241*2a6744ebSBarry Smith 
242*2a6744ebSBarry Smith   A->assembled     = PETSC_TRUE;
243*2a6744ebSBarry Smith   A->preallocated  = PETSC_FALSE;
244*2a6744ebSBarry Smith   PetscFunctionReturn(0);
245*2a6744ebSBarry Smith }
246*2a6744ebSBarry Smith EXTERN_C_END
247*2a6744ebSBarry Smith 
248*2a6744ebSBarry Smith #undef __FUNCT__
249*2a6744ebSBarry Smith #define __FUNCT__ "MatCreateScatter"
250*2a6744ebSBarry Smith /*@C
251*2a6744ebSBarry Smith    MatCreateScatter - Creates a new matrix based on a VecScatter
252*2a6744ebSBarry Smith 
253*2a6744ebSBarry Smith   Collective on MPI_Comm
254*2a6744ebSBarry Smith 
255*2a6744ebSBarry Smith    Input Parameters:
256*2a6744ebSBarry Smith +  comm - MPI communicator
257*2a6744ebSBarry Smith -  scatter - a VecScatterContext
258*2a6744ebSBarry Smith 
259*2a6744ebSBarry Smith    Output Parameter:
260*2a6744ebSBarry Smith .  A - the matrix
261*2a6744ebSBarry Smith 
262*2a6744ebSBarry Smith    Level: intermediate
263*2a6744ebSBarry Smith 
264*2a6744ebSBarry Smith    PETSc requires that matrices and vectors being used for certain
265*2a6744ebSBarry Smith    operations are partitioned accordingly.  For example, when
266*2a6744ebSBarry Smith    creating a scatter matrix, A, that supports parallel matrix-vector
267*2a6744ebSBarry Smith    products using MatMult(A,x,y) the user should set the number
268*2a6744ebSBarry Smith    of local matrix rows to be the number of local elements of the
269*2a6744ebSBarry Smith    corresponding result vector, y. Note that this is information is
270*2a6744ebSBarry Smith    required for use of the matrix interface routines, even though
271*2a6744ebSBarry Smith    the scatter matrix may not actually be physically partitioned.
272*2a6744ebSBarry Smith    For example,
273*2a6744ebSBarry Smith 
274*2a6744ebSBarry Smith .keywords: matrix, scatter, create
275*2a6744ebSBarry Smith 
276*2a6744ebSBarry Smith .seealso: MatScatterSetVecScatter(), MatScatterGetVecScatter(), MATSCATTER
277*2a6744ebSBarry Smith @*/
278*2a6744ebSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm comm,VecScatter scatter,Mat *A)
279*2a6744ebSBarry Smith {
280*2a6744ebSBarry Smith   PetscErrorCode ierr;
281*2a6744ebSBarry Smith 
282*2a6744ebSBarry Smith   PetscFunctionBegin;
283*2a6744ebSBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
284*2a6744ebSBarry Smith   ierr = MatSetSizes(*A,scatter->to_n,scatter->from_n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
285*2a6744ebSBarry Smith   ierr = MatSetType(*A,MATSCATTER);CHKERRQ(ierr);
286*2a6744ebSBarry Smith   ierr = MatScatterSetVecScatter(*A,scatter);CHKERRQ(ierr);
287*2a6744ebSBarry Smith   PetscFunctionReturn(0);
288*2a6744ebSBarry Smith }
289*2a6744ebSBarry Smith 
290*2a6744ebSBarry Smith #undef __FUNCT__
291*2a6744ebSBarry Smith #define __FUNCT__ "MatScatterSetVecScatter"
292*2a6744ebSBarry Smith /*@
293*2a6744ebSBarry Smith     MatScatterSetVecScatter - sets that scatter that the matrix is to apply as its linear operator
294*2a6744ebSBarry Smith 
295*2a6744ebSBarry Smith    Collective on Mat
296*2a6744ebSBarry Smith 
297*2a6744ebSBarry Smith     Input Parameters:
298*2a6744ebSBarry Smith +   mat - the scatter matrix
299*2a6744ebSBarry Smith -   scatter - the scatter context create with VecScatterCreate()
300*2a6744ebSBarry Smith 
301*2a6744ebSBarry Smith    Level: advanced
302*2a6744ebSBarry Smith 
303*2a6744ebSBarry Smith 
304*2a6744ebSBarry Smith .seealso: MatCreateScatter(), MATSCATTER
305*2a6744ebSBarry Smith @*/
306*2a6744ebSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat mat,VecScatter scatter)
307*2a6744ebSBarry Smith {
308*2a6744ebSBarry Smith   Mat_Scatter    *mscatter = (Mat_Scatter*)mat->data;
309*2a6744ebSBarry Smith   PetscErrorCode ierr;
310*2a6744ebSBarry Smith 
311*2a6744ebSBarry Smith   PetscFunctionBegin;
312*2a6744ebSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
313*2a6744ebSBarry Smith   PetscValidHeaderSpecific(scatter,VEC_SCATTER_COOKIE,2);
314*2a6744ebSBarry Smith   PetscCheckSameComm((PetscObject)scatter,(PetscObject)mat,1,2);
315*2a6744ebSBarry Smith   if (mat->rmap.n != scatter->to_n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Number of local rows in matrix %D not equal local scatter size %D",mat->rmap.n,scatter->to_n);
316*2a6744ebSBarry Smith   if (mat->cmap.n != scatter->from_n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Number of local columns in matrix %D not equal local scatter size %D",mat->cmap.n,scatter->from_n);
317*2a6744ebSBarry Smith 
318*2a6744ebSBarry Smith   if (mscatter->scatter) {ierr = VecScatterDestroy(mscatter->scatter);CHKERRQ(ierr);}
319*2a6744ebSBarry Smith   mscatter->scatter = scatter;
320*2a6744ebSBarry Smith   ierr = PetscObjectReference((PetscObject)scatter);CHKERRQ(ierr);
321*2a6744ebSBarry Smith   PetscFunctionReturn(0);
322*2a6744ebSBarry Smith }
323*2a6744ebSBarry Smith 
324*2a6744ebSBarry Smith 
325