xref: /petsc/src/mat/impls/scatter/mscatter.c (revision f4259b3095b7e7c300c2b3c3729548b071684bd7)
12a6744ebSBarry Smith 
22a6744ebSBarry Smith /*
32a6744ebSBarry Smith    This provides a matrix that applies a VecScatter to a vector.
42a6744ebSBarry Smith */
52a6744ebSBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
7af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
82a6744ebSBarry Smith 
92a6744ebSBarry Smith typedef struct {
102a6744ebSBarry Smith   VecScatter scatter;
112a6744ebSBarry Smith } Mat_Scatter;
122a6744ebSBarry Smith 
132a6744ebSBarry Smith /*@
142a6744ebSBarry Smith     MatScatterGetVecScatter - Returns the user-provided scatter set with MatScatterSetVecScatter()
152a6744ebSBarry Smith 
162a6744ebSBarry Smith     Not Collective, but not cannot use scatter if not used collectively on Mat
172a6744ebSBarry Smith 
182a6744ebSBarry Smith     Input Parameter:
192a6744ebSBarry Smith .   mat - the matrix, should have been created with MatCreateScatter() or have type MATSCATTER
202a6744ebSBarry Smith 
212a6744ebSBarry Smith     Output Parameter:
222a6744ebSBarry Smith .   scatter - the scatter context
232a6744ebSBarry Smith 
242a6744ebSBarry Smith     Level: intermediate
252a6744ebSBarry Smith 
262a6744ebSBarry Smith .seealso: MatCreateScatter(), MatScatterSetVecScatter(), MATSCATTER
272a6744ebSBarry Smith @*/
287087cfbeSBarry Smith PetscErrorCode  MatScatterGetVecScatter(Mat mat,VecScatter *scatter)
292a6744ebSBarry Smith {
302a6744ebSBarry Smith   Mat_Scatter *mscatter;
312a6744ebSBarry Smith 
322a6744ebSBarry Smith   PetscFunctionBegin;
330700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
342a6744ebSBarry Smith   PetscValidPointer(scatter,2);
352a6744ebSBarry Smith   mscatter = (Mat_Scatter*)mat->data;
362a6744ebSBarry Smith   *scatter = mscatter->scatter;
372a6744ebSBarry Smith   PetscFunctionReturn(0);
382a6744ebSBarry Smith }
392a6744ebSBarry Smith 
402a6744ebSBarry Smith PetscErrorCode MatDestroy_Scatter(Mat mat)
412a6744ebSBarry Smith {
422a6744ebSBarry Smith   PetscErrorCode ierr;
432a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)mat->data;
442a6744ebSBarry Smith 
452a6744ebSBarry Smith   PetscFunctionBegin;
466bf464f9SBarry Smith   ierr = VecScatterDestroy(&scatter->scatter);CHKERRQ(ierr);
476bf464f9SBarry Smith   ierr = PetscFree(mat->data);CHKERRQ(ierr);
482a6744ebSBarry Smith   PetscFunctionReturn(0);
492a6744ebSBarry Smith }
502a6744ebSBarry Smith 
512a6744ebSBarry Smith PetscErrorCode MatMult_Scatter(Mat A,Vec x,Vec y)
522a6744ebSBarry Smith {
532a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
542a6744ebSBarry Smith   PetscErrorCode ierr;
552a6744ebSBarry Smith 
562a6744ebSBarry Smith   PetscFunctionBegin;
57ce94432eSBarry Smith   if (!scatter->scatter) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
5874d95942SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5974d95942SJed Brown   ierr = VecScatterBegin(scatter->scatter,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6074d95942SJed Brown   ierr = VecScatterEnd(scatter->scatter,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
612a6744ebSBarry Smith   PetscFunctionReturn(0);
622a6744ebSBarry Smith }
632a6744ebSBarry Smith 
642a6744ebSBarry Smith PetscErrorCode MatMultAdd_Scatter(Mat A,Vec x,Vec y,Vec z)
652a6744ebSBarry Smith {
662a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
672a6744ebSBarry Smith   PetscErrorCode ierr;
682a6744ebSBarry Smith 
692a6744ebSBarry Smith   PetscFunctionBegin;
70ce94432eSBarry Smith   if (!scatter->scatter) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
712a6744ebSBarry Smith   if (z != y) {ierr = VecCopy(y,z);CHKERRQ(ierr);}
72ca9f406cSSatish Balay   ierr = VecScatterBegin(scatter->scatter,x,z,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
73ca9f406cSSatish Balay   ierr = VecScatterEnd(scatter->scatter,x,z,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
742a6744ebSBarry Smith   PetscFunctionReturn(0);
752a6744ebSBarry Smith }
762a6744ebSBarry Smith 
772a6744ebSBarry Smith PetscErrorCode MatMultTranspose_Scatter(Mat A,Vec x,Vec y)
782a6744ebSBarry Smith {
792a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
802a6744ebSBarry Smith   PetscErrorCode ierr;
812a6744ebSBarry Smith 
822a6744ebSBarry Smith   PetscFunctionBegin;
83ce94432eSBarry Smith   if (!scatter->scatter) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
8474d95942SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
8574d95942SJed Brown   ierr = VecScatterBegin(scatter->scatter,x,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8674d95942SJed Brown   ierr = VecScatterEnd(scatter->scatter,x,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
872a6744ebSBarry Smith   PetscFunctionReturn(0);
882a6744ebSBarry Smith }
892a6744ebSBarry Smith 
902a6744ebSBarry Smith PetscErrorCode MatMultTransposeAdd_Scatter(Mat A,Vec x,Vec y,Vec z)
912a6744ebSBarry Smith {
922a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
932a6744ebSBarry Smith   PetscErrorCode ierr;
942a6744ebSBarry Smith 
952a6744ebSBarry Smith   PetscFunctionBegin;
96ce94432eSBarry Smith   if (!scatter->scatter) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
972a6744ebSBarry Smith   if (z != y) {ierr = VecCopy(y,z);CHKERRQ(ierr);}
98ca9f406cSSatish Balay   ierr = VecScatterBegin(scatter->scatter,x,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
99ca9f406cSSatish Balay   ierr = VecScatterEnd(scatter->scatter,x,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1002a6744ebSBarry Smith   PetscFunctionReturn(0);
1012a6744ebSBarry Smith }
1022a6744ebSBarry Smith 
103*f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL,
104*f4259b30SLisandro Dalcin                                        NULL,
105*f4259b30SLisandro Dalcin                                        NULL,
1062a6744ebSBarry Smith                                        MatMult_Scatter,
1072a6744ebSBarry Smith                                /*  4*/ MatMultAdd_Scatter,
1082a6744ebSBarry Smith                                        MatMultTranspose_Scatter,
1092a6744ebSBarry Smith                                        MatMultTransposeAdd_Scatter,
110*f4259b30SLisandro Dalcin                                        NULL,
111*f4259b30SLisandro Dalcin                                        NULL,
112*f4259b30SLisandro Dalcin                                        NULL,
113*f4259b30SLisandro Dalcin                                /* 10*/ NULL,
114*f4259b30SLisandro Dalcin                                        NULL,
115*f4259b30SLisandro Dalcin                                        NULL,
116*f4259b30SLisandro Dalcin                                        NULL,
117*f4259b30SLisandro Dalcin                                        NULL,
118*f4259b30SLisandro Dalcin                                /* 15*/ NULL,
119*f4259b30SLisandro Dalcin                                        NULL,
120*f4259b30SLisandro Dalcin                                        NULL,
121*f4259b30SLisandro Dalcin                                        NULL,
122*f4259b30SLisandro Dalcin                                        NULL,
123*f4259b30SLisandro Dalcin                                /* 20*/ NULL,
124*f4259b30SLisandro Dalcin                                        NULL,
125*f4259b30SLisandro Dalcin                                        NULL,
126*f4259b30SLisandro Dalcin                                        NULL,
127*f4259b30SLisandro Dalcin                                /* 24*/ NULL,
128*f4259b30SLisandro Dalcin                                        NULL,
129*f4259b30SLisandro Dalcin                                        NULL,
130*f4259b30SLisandro Dalcin                                        NULL,
131*f4259b30SLisandro Dalcin                                        NULL,
132*f4259b30SLisandro Dalcin                                /* 29*/ NULL,
133*f4259b30SLisandro Dalcin                                        NULL,
134*f4259b30SLisandro Dalcin                                        NULL,
135*f4259b30SLisandro Dalcin                                        NULL,
136*f4259b30SLisandro Dalcin                                        NULL,
137*f4259b30SLisandro Dalcin                                /* 34*/ NULL,
138*f4259b30SLisandro Dalcin                                        NULL,
139*f4259b30SLisandro Dalcin                                        NULL,
140*f4259b30SLisandro Dalcin                                        NULL,
141*f4259b30SLisandro Dalcin                                        NULL,
142*f4259b30SLisandro Dalcin                                /* 39*/ NULL,
143*f4259b30SLisandro Dalcin                                        NULL,
144*f4259b30SLisandro Dalcin                                        NULL,
145*f4259b30SLisandro Dalcin                                        NULL,
146*f4259b30SLisandro Dalcin                                        NULL,
147*f4259b30SLisandro Dalcin                                /* 44*/ NULL,
148*f4259b30SLisandro Dalcin                                        NULL,
1497d68702bSBarry Smith                                        MatShift_Basic,
150*f4259b30SLisandro Dalcin                                        NULL,
151*f4259b30SLisandro Dalcin                                        NULL,
152*f4259b30SLisandro Dalcin                                /* 49*/ NULL,
153*f4259b30SLisandro Dalcin                                        NULL,
154*f4259b30SLisandro Dalcin                                        NULL,
155*f4259b30SLisandro Dalcin                                        NULL,
156*f4259b30SLisandro Dalcin                                        NULL,
157*f4259b30SLisandro Dalcin                                /* 54*/ NULL,
158*f4259b30SLisandro Dalcin                                        NULL,
159*f4259b30SLisandro Dalcin                                        NULL,
160*f4259b30SLisandro Dalcin                                        NULL,
161*f4259b30SLisandro Dalcin                                        NULL,
162*f4259b30SLisandro Dalcin                                /* 59*/ NULL,
1632a6744ebSBarry Smith                                        MatDestroy_Scatter,
164*f4259b30SLisandro Dalcin                                        NULL,
165*f4259b30SLisandro Dalcin                                        NULL,
166*f4259b30SLisandro Dalcin                                        NULL,
167*f4259b30SLisandro Dalcin                                /* 64*/ NULL,
168*f4259b30SLisandro Dalcin                                        NULL,
169*f4259b30SLisandro Dalcin                                        NULL,
170*f4259b30SLisandro Dalcin                                        NULL,
171*f4259b30SLisandro Dalcin                                        NULL,
172*f4259b30SLisandro Dalcin                                /* 69*/ NULL,
173*f4259b30SLisandro Dalcin                                        NULL,
174*f4259b30SLisandro Dalcin                                        NULL,
175*f4259b30SLisandro Dalcin                                        NULL,
176*f4259b30SLisandro Dalcin                                        NULL,
177*f4259b30SLisandro Dalcin                                /* 74*/ NULL,
178*f4259b30SLisandro Dalcin                                        NULL,
179*f4259b30SLisandro Dalcin                                        NULL,
180*f4259b30SLisandro Dalcin                                        NULL,
181*f4259b30SLisandro Dalcin                                        NULL,
182*f4259b30SLisandro Dalcin                                /* 79*/ NULL,
183*f4259b30SLisandro Dalcin                                        NULL,
184*f4259b30SLisandro Dalcin                                        NULL,
185*f4259b30SLisandro Dalcin                                        NULL,
186*f4259b30SLisandro Dalcin                                        NULL,
187*f4259b30SLisandro Dalcin                                /* 84*/ NULL,
188*f4259b30SLisandro Dalcin                                        NULL,
189*f4259b30SLisandro Dalcin                                        NULL,
190*f4259b30SLisandro Dalcin                                        NULL,
191*f4259b30SLisandro Dalcin                                        NULL,
192*f4259b30SLisandro Dalcin                                /* 89*/ NULL,
193*f4259b30SLisandro Dalcin                                        NULL,
194*f4259b30SLisandro Dalcin                                        NULL,
195*f4259b30SLisandro Dalcin                                        NULL,
196*f4259b30SLisandro Dalcin                                        NULL,
197*f4259b30SLisandro Dalcin                                /* 94*/ NULL,
198*f4259b30SLisandro Dalcin                                        NULL,
199*f4259b30SLisandro Dalcin                                        NULL,
200*f4259b30SLisandro Dalcin                                        NULL,
201*f4259b30SLisandro Dalcin                                        NULL,
202*f4259b30SLisandro Dalcin                                 /*99*/ NULL,
203*f4259b30SLisandro Dalcin                                        NULL,
204*f4259b30SLisandro Dalcin                                        NULL,
205*f4259b30SLisandro Dalcin                                        NULL,
206*f4259b30SLisandro Dalcin                                        NULL,
207*f4259b30SLisandro Dalcin                                /*104*/ NULL,
208*f4259b30SLisandro Dalcin                                        NULL,
209*f4259b30SLisandro Dalcin                                        NULL,
210*f4259b30SLisandro Dalcin                                        NULL,
211*f4259b30SLisandro Dalcin                                        NULL,
212*f4259b30SLisandro Dalcin                                /*109*/ NULL,
213*f4259b30SLisandro Dalcin                                        NULL,
214*f4259b30SLisandro Dalcin                                        NULL,
215*f4259b30SLisandro Dalcin                                        NULL,
216*f4259b30SLisandro Dalcin                                        NULL,
217*f4259b30SLisandro Dalcin                                /*114*/ NULL,
218*f4259b30SLisandro Dalcin                                        NULL,
219*f4259b30SLisandro Dalcin                                        NULL,
220*f4259b30SLisandro Dalcin                                        NULL,
221*f4259b30SLisandro Dalcin                                        NULL,
222*f4259b30SLisandro Dalcin                                /*119*/ NULL,
223*f4259b30SLisandro Dalcin                                        NULL,
224*f4259b30SLisandro Dalcin                                        NULL,
225*f4259b30SLisandro Dalcin                                        NULL,
226*f4259b30SLisandro Dalcin                                        NULL,
227*f4259b30SLisandro Dalcin                                /*124*/ NULL,
228*f4259b30SLisandro Dalcin                                        NULL,
229*f4259b30SLisandro Dalcin                                        NULL,
230*f4259b30SLisandro Dalcin                                        NULL,
231*f4259b30SLisandro Dalcin                                        NULL,
232*f4259b30SLisandro Dalcin                                /*129*/ NULL,
233*f4259b30SLisandro Dalcin                                        NULL,
234*f4259b30SLisandro Dalcin                                        NULL,
235*f4259b30SLisandro Dalcin                                        NULL,
236*f4259b30SLisandro Dalcin                                        NULL,
237*f4259b30SLisandro Dalcin                                /*134*/ NULL,
238*f4259b30SLisandro Dalcin                                        NULL,
239*f4259b30SLisandro Dalcin                                        NULL,
240*f4259b30SLisandro Dalcin                                        NULL,
241*f4259b30SLisandro Dalcin                                        NULL,
242*f4259b30SLisandro Dalcin                                /*139*/ NULL,
243*f4259b30SLisandro Dalcin                                        NULL,
244*f4259b30SLisandro Dalcin                                        NULL
2453964eb88SJed Brown };
2462a6744ebSBarry Smith 
2472a6744ebSBarry Smith /*MC
2484a440b9bSBarry Smith    MATSCATTER - MATSCATTER = "scatter" - A matrix type that simply applies a VecScatterBegin/End()
2492a6744ebSBarry Smith 
2502a6744ebSBarry Smith   Level: advanced
2512a6744ebSBarry Smith 
2522a6744ebSBarry Smith .seealso: MatCreateScatter(), MatScatterSetVecScatter(), MatScatterGetVecScatter()
2532a6744ebSBarry Smith 
2542a6744ebSBarry Smith M*/
2552a6744ebSBarry Smith 
2568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Scatter(Mat A)
2572a6744ebSBarry Smith {
2582a6744ebSBarry Smith   Mat_Scatter    *b;
2592a6744ebSBarry Smith   PetscErrorCode ierr;
2602a6744ebSBarry Smith 
2612a6744ebSBarry Smith   PetscFunctionBegin;
2622a6744ebSBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
263b00a9115SJed Brown   ierr = PetscNewLog(A,&b);CHKERRQ(ierr);
2642a6744ebSBarry Smith 
2652a6744ebSBarry Smith   A->data = (void*)b;
2662a6744ebSBarry Smith 
26726283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
26826283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
2692a6744ebSBarry Smith 
2702a6744ebSBarry Smith   A->assembled    = PETSC_TRUE;
2712a6744ebSBarry Smith   A->preallocated = PETSC_FALSE;
27274d95942SJed Brown 
27374d95942SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSCATTER);CHKERRQ(ierr);
2742a6744ebSBarry Smith   PetscFunctionReturn(0);
2752a6744ebSBarry Smith }
2762a6744ebSBarry Smith 
2776eb45d04SBarry Smith #include <petsc/private/vecscatterimpl.h>
2782a6744ebSBarry Smith /*@C
2792a6744ebSBarry Smith    MatCreateScatter - Creates a new matrix based on a VecScatter
2802a6744ebSBarry Smith 
281d083f849SBarry Smith   Collective
2822a6744ebSBarry Smith 
2832a6744ebSBarry Smith    Input Parameters:
2842a6744ebSBarry Smith +  comm - MPI communicator
2852a6744ebSBarry Smith -  scatter - a VecScatterContext
2862a6744ebSBarry Smith 
2872a6744ebSBarry Smith    Output Parameter:
2882a6744ebSBarry Smith .  A - the matrix
2892a6744ebSBarry Smith 
2902a6744ebSBarry Smith    Level: intermediate
2912a6744ebSBarry Smith 
2922a6744ebSBarry Smith    PETSc requires that matrices and vectors being used for certain
2932a6744ebSBarry Smith    operations are partitioned accordingly.  For example, when
2942a6744ebSBarry Smith    creating a scatter matrix, A, that supports parallel matrix-vector
2952a6744ebSBarry Smith    products using MatMult(A,x,y) the user should set the number
2962a6744ebSBarry Smith    of local matrix rows to be the number of local elements of the
2972a6744ebSBarry Smith    corresponding result vector, y. Note that this is information is
2982a6744ebSBarry Smith    required for use of the matrix interface routines, even though
2992a6744ebSBarry Smith    the scatter matrix may not actually be physically partitioned.
3002a6744ebSBarry Smith 
3016eb45d04SBarry Smith   Developer Notes: This directly accesses information inside the VecScatter associated with the matrix-vector product
3026eb45d04SBarry Smith    for this matrix. This is not desirable..
3036eb45d04SBarry Smith 
3046eb45d04SBarry Smith 
3052a6744ebSBarry Smith .seealso: MatScatterSetVecScatter(), MatScatterGetVecScatter(), MATSCATTER
3062a6744ebSBarry Smith @*/
3077087cfbeSBarry Smith PetscErrorCode  MatCreateScatter(MPI_Comm comm,VecScatter scatter,Mat *A)
3082a6744ebSBarry Smith {
3092a6744ebSBarry Smith   PetscErrorCode ierr;
3102a6744ebSBarry Smith 
3112a6744ebSBarry Smith   PetscFunctionBegin;
3122a6744ebSBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
3132a6744ebSBarry Smith   ierr = MatSetSizes(*A,scatter->to_n,scatter->from_n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
3142a6744ebSBarry Smith   ierr = MatSetType(*A,MATSCATTER);CHKERRQ(ierr);
3152a6744ebSBarry Smith   ierr = MatScatterSetVecScatter(*A,scatter);CHKERRQ(ierr);
31637e1895aSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
3172a6744ebSBarry Smith   PetscFunctionReturn(0);
3182a6744ebSBarry Smith }
3192a6744ebSBarry Smith 
3202a6744ebSBarry Smith /*@
3212a6744ebSBarry Smith     MatScatterSetVecScatter - sets that scatter that the matrix is to apply as its linear operator
3222a6744ebSBarry Smith 
3232a6744ebSBarry Smith    Collective on Mat
3242a6744ebSBarry Smith 
3252a6744ebSBarry Smith     Input Parameters:
3262a6744ebSBarry Smith +   mat - the scatter matrix
3279448b7f1SJunchao Zhang -   scatter - the scatter context create with VecScatterCreate()
3282a6744ebSBarry Smith 
3292a6744ebSBarry Smith    Level: advanced
3302a6744ebSBarry Smith 
3312a6744ebSBarry Smith 
3322a6744ebSBarry Smith .seealso: MatCreateScatter(), MATSCATTER
3332a6744ebSBarry Smith @*/
3347087cfbeSBarry Smith PetscErrorCode  MatScatterSetVecScatter(Mat mat,VecScatter scatter)
3352a6744ebSBarry Smith {
3362a6744ebSBarry Smith   Mat_Scatter    *mscatter = (Mat_Scatter*)mat->data;
3372a6744ebSBarry Smith   PetscErrorCode ierr;
3382a6744ebSBarry Smith 
3392a6744ebSBarry Smith   PetscFunctionBegin;
3400700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
3410700a824SBarry Smith   PetscValidHeaderSpecific(scatter,VEC_SCATTER_CLASSID,2);
3424067f9b5SMatthew Knepley   PetscCheckSameComm((PetscObject)scatter,1,(PetscObject)mat,2);
343e32f2f54SBarry Smith   if (mat->rmap->n != scatter->to_n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number of local rows in matrix %D not equal local scatter size %D",mat->rmap->n,scatter->to_n);
344e32f2f54SBarry Smith   if (mat->cmap->n != scatter->from_n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number of local columns in matrix %D not equal local scatter size %D",mat->cmap->n,scatter->from_n);
3452a6744ebSBarry Smith 
346c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)scatter);CHKERRQ(ierr);
3476bf464f9SBarry Smith   ierr = VecScatterDestroy(&mscatter->scatter);CHKERRQ(ierr);
34826fbe8dcSKarl Rupp 
3492a6744ebSBarry Smith   mscatter->scatter = scatter;
3502a6744ebSBarry Smith   PetscFunctionReturn(0);
3512a6744ebSBarry Smith }
3522a6744ebSBarry Smith 
3532a6744ebSBarry Smith 
354