xref: /petsc/src/mat/impls/scatter/mscatter.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
12a6744ebSBarry Smith 
22a6744ebSBarry Smith /*
32a6744ebSBarry Smith    This provides a matrix that applies a VecScatter to a vector.
42a6744ebSBarry Smith */
52a6744ebSBarry Smith 
6b45d2f2cSJed Brown #include <petsc-private/matimpl.h>        /*I "petscmat.h" I*/
7b45d2f2cSJed Brown #include <petsc-private/vecimpl.h>
82a6744ebSBarry Smith 
92a6744ebSBarry Smith typedef struct {
102a6744ebSBarry Smith   VecScatter scatter;
112a6744ebSBarry Smith } Mat_Scatter;
122a6744ebSBarry Smith 
132a6744ebSBarry Smith #undef __FUNCT__
142a6744ebSBarry Smith #define __FUNCT__ "MatScatterGetVecScatter"
152a6744ebSBarry Smith /*@
162a6744ebSBarry Smith     MatScatterGetVecScatter - Returns the user-provided scatter set with MatScatterSetVecScatter()
172a6744ebSBarry Smith 
182a6744ebSBarry Smith     Not Collective, but not cannot use scatter if not used collectively on Mat
192a6744ebSBarry Smith 
202a6744ebSBarry Smith     Input Parameter:
212a6744ebSBarry Smith .   mat - the matrix, should have been created with MatCreateScatter() or have type MATSCATTER
222a6744ebSBarry Smith 
232a6744ebSBarry Smith     Output Parameter:
242a6744ebSBarry Smith .   scatter - the scatter context
252a6744ebSBarry Smith 
262a6744ebSBarry Smith     Level: intermediate
272a6744ebSBarry Smith 
282a6744ebSBarry Smith .keywords: matrix, scatter, get
292a6744ebSBarry Smith 
302a6744ebSBarry Smith .seealso: MatCreateScatter(), MatScatterSetVecScatter(), MATSCATTER
312a6744ebSBarry Smith @*/
327087cfbeSBarry Smith PetscErrorCode  MatScatterGetVecScatter(Mat mat,VecScatter *scatter)
332a6744ebSBarry Smith {
342a6744ebSBarry Smith   Mat_Scatter *mscatter;
352a6744ebSBarry Smith 
362a6744ebSBarry Smith   PetscFunctionBegin;
370700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
382a6744ebSBarry Smith   PetscValidPointer(scatter,2);
392a6744ebSBarry Smith   mscatter = (Mat_Scatter*)mat->data;
402a6744ebSBarry Smith   *scatter = mscatter->scatter;
412a6744ebSBarry Smith   PetscFunctionReturn(0);
422a6744ebSBarry Smith }
432a6744ebSBarry Smith 
442a6744ebSBarry Smith #undef __FUNCT__
452a6744ebSBarry Smith #define __FUNCT__ "MatDestroy_Scatter"
462a6744ebSBarry Smith PetscErrorCode MatDestroy_Scatter(Mat mat)
472a6744ebSBarry Smith {
482a6744ebSBarry Smith   PetscErrorCode ierr;
492a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)mat->data;
502a6744ebSBarry Smith 
512a6744ebSBarry Smith   PetscFunctionBegin;
526bf464f9SBarry Smith   ierr = VecScatterDestroy(&scatter->scatter);CHKERRQ(ierr);
536bf464f9SBarry Smith   ierr = PetscFree(mat->data);CHKERRQ(ierr);
542a6744ebSBarry Smith   PetscFunctionReturn(0);
552a6744ebSBarry Smith }
562a6744ebSBarry Smith 
572a6744ebSBarry Smith #undef __FUNCT__
582a6744ebSBarry Smith #define __FUNCT__ "MatMult_Scatter"
592a6744ebSBarry Smith PetscErrorCode MatMult_Scatter(Mat A,Vec x,Vec y)
602a6744ebSBarry Smith {
612a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
622a6744ebSBarry Smith   PetscErrorCode ierr;
632a6744ebSBarry Smith 
642a6744ebSBarry Smith   PetscFunctionBegin;
65ce94432eSBarry Smith   if (!scatter->scatter) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
6674d95942SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
6774d95942SJed Brown   ierr = VecScatterBegin(scatter->scatter,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6874d95942SJed Brown   ierr = VecScatterEnd(scatter->scatter,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
692a6744ebSBarry Smith   PetscFunctionReturn(0);
702a6744ebSBarry Smith }
712a6744ebSBarry Smith 
722a6744ebSBarry Smith #undef __FUNCT__
732a6744ebSBarry Smith #define __FUNCT__ "MatMultAdd_Scatter"
742a6744ebSBarry Smith PetscErrorCode MatMultAdd_Scatter(Mat A,Vec x,Vec y,Vec z)
752a6744ebSBarry Smith {
762a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
772a6744ebSBarry Smith   PetscErrorCode ierr;
782a6744ebSBarry Smith 
792a6744ebSBarry Smith   PetscFunctionBegin;
80ce94432eSBarry Smith   if (!scatter->scatter) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
812a6744ebSBarry Smith   if (z != y) {ierr = VecCopy(y,z);CHKERRQ(ierr);}
82ca9f406cSSatish Balay   ierr = VecScatterBegin(scatter->scatter,x,z,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
83ca9f406cSSatish Balay   ierr = VecScatterEnd(scatter->scatter,x,z,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
842a6744ebSBarry Smith   PetscFunctionReturn(0);
852a6744ebSBarry Smith }
862a6744ebSBarry Smith 
872a6744ebSBarry Smith #undef __FUNCT__
882a6744ebSBarry Smith #define __FUNCT__ "MatMultTranspose_Scatter"
892a6744ebSBarry Smith PetscErrorCode MatMultTranspose_Scatter(Mat A,Vec x,Vec y)
902a6744ebSBarry Smith {
912a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
922a6744ebSBarry Smith   PetscErrorCode ierr;
932a6744ebSBarry Smith 
942a6744ebSBarry Smith   PetscFunctionBegin;
95ce94432eSBarry Smith   if (!scatter->scatter) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
9674d95942SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
9774d95942SJed Brown   ierr = VecScatterBegin(scatter->scatter,x,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9874d95942SJed Brown   ierr = VecScatterEnd(scatter->scatter,x,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
992a6744ebSBarry Smith   PetscFunctionReturn(0);
1002a6744ebSBarry Smith }
1012a6744ebSBarry Smith 
1022a6744ebSBarry Smith #undef __FUNCT__
1032a6744ebSBarry Smith #define __FUNCT__ "MatMultTransposeAdd_Scatter"
1042a6744ebSBarry Smith PetscErrorCode MatMultTransposeAdd_Scatter(Mat A,Vec x,Vec y,Vec z)
1052a6744ebSBarry Smith {
1062a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
1072a6744ebSBarry Smith   PetscErrorCode ierr;
1082a6744ebSBarry Smith 
1092a6744ebSBarry Smith   PetscFunctionBegin;
110ce94432eSBarry Smith   if (!scatter->scatter) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
1112a6744ebSBarry Smith   if (z != y) {ierr = VecCopy(y,z);CHKERRQ(ierr);}
112ca9f406cSSatish Balay   ierr = VecScatterBegin(scatter->scatter,x,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
113ca9f406cSSatish Balay   ierr = VecScatterEnd(scatter->scatter,x,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1142a6744ebSBarry Smith   PetscFunctionReturn(0);
1152a6744ebSBarry Smith }
1162a6744ebSBarry Smith 
1172a6744ebSBarry Smith static struct _MatOps MatOps_Values = {0,
1182a6744ebSBarry Smith         0,
1192a6744ebSBarry Smith         0,
1202a6744ebSBarry Smith         MatMult_Scatter,
1212a6744ebSBarry Smith /*  4*/ MatMultAdd_Scatter,
1222a6744ebSBarry Smith         MatMultTranspose_Scatter,
1232a6744ebSBarry Smith         MatMultTransposeAdd_Scatter,
1242a6744ebSBarry Smith         0,
1252a6744ebSBarry Smith         0,
1262a6744ebSBarry Smith         0,
1272a6744ebSBarry Smith /* 10*/ 0,
1282a6744ebSBarry Smith         0,
1292a6744ebSBarry Smith         0,
1302a6744ebSBarry Smith         0,
1312a6744ebSBarry Smith         0,
1322a6744ebSBarry Smith /* 15*/ 0,
1332a6744ebSBarry Smith         0,
1342a6744ebSBarry Smith         0,
1352a6744ebSBarry Smith         0,
1362a6744ebSBarry Smith         0,
1372a6744ebSBarry Smith /* 20*/ 0,
1382a6744ebSBarry Smith         0,
1392a6744ebSBarry Smith         0,
1402a6744ebSBarry Smith         0,
141d519adbfSMatthew Knepley /* 24*/ 0,
1422a6744ebSBarry Smith         0,
1432a6744ebSBarry Smith         0,
1442a6744ebSBarry Smith         0,
1452a6744ebSBarry Smith         0,
146d519adbfSMatthew Knepley /* 29*/ 0,
1472a6744ebSBarry Smith         0,
1482a6744ebSBarry Smith         0,
1492a6744ebSBarry Smith         0,
1502a6744ebSBarry Smith         0,
151d519adbfSMatthew Knepley /* 34*/ 0,
1522a6744ebSBarry Smith         0,
1532a6744ebSBarry Smith         0,
1542a6744ebSBarry Smith         0,
1552a6744ebSBarry Smith         0,
156d519adbfSMatthew Knepley /* 39*/ 0,
1572a6744ebSBarry Smith         0,
1582a6744ebSBarry Smith         0,
1592a6744ebSBarry Smith         0,
1602a6744ebSBarry Smith         0,
161d519adbfSMatthew Knepley /* 44*/ 0,
1622a6744ebSBarry Smith         0,
1632a6744ebSBarry Smith         0,
1642a6744ebSBarry Smith         0,
1652a6744ebSBarry Smith         0,
166d519adbfSMatthew Knepley /* 49*/ 0,
1672a6744ebSBarry Smith         0,
1682a6744ebSBarry Smith         0,
1692a6744ebSBarry Smith         0,
1702a6744ebSBarry Smith         0,
171d519adbfSMatthew Knepley /* 54*/ 0,
1722a6744ebSBarry Smith         0,
1732a6744ebSBarry Smith         0,
1742a6744ebSBarry Smith         0,
1752a6744ebSBarry Smith         0,
176d519adbfSMatthew Knepley /* 59*/ 0,
1772a6744ebSBarry Smith         MatDestroy_Scatter,
1782a6744ebSBarry Smith         0,
1792a6744ebSBarry Smith         0,
1802a6744ebSBarry Smith         0,
181d519adbfSMatthew Knepley /* 64*/ 0,
1822a6744ebSBarry Smith         0,
1832a6744ebSBarry Smith         0,
1842a6744ebSBarry Smith         0,
1852a6744ebSBarry Smith         0,
186d519adbfSMatthew Knepley /* 69*/ 0,
1872a6744ebSBarry Smith         0,
1882a6744ebSBarry Smith         0,
1892a6744ebSBarry Smith         0,
1902a6744ebSBarry Smith         0,
191d519adbfSMatthew Knepley /* 74*/ 0,
1922a6744ebSBarry Smith         0,
1932a6744ebSBarry Smith         0,
1942a6744ebSBarry Smith         0,
1952a6744ebSBarry Smith         0,
196d519adbfSMatthew Knepley /* 79*/ 0,
1972a6744ebSBarry Smith         0,
1982a6744ebSBarry Smith         0,
1992a6744ebSBarry Smith         0,
2002a6744ebSBarry Smith         0,
201d519adbfSMatthew Knepley /* 84*/ 0,
2022a6744ebSBarry Smith         0,
2032a6744ebSBarry Smith         0,
2042a6744ebSBarry Smith         0,
2052a6744ebSBarry Smith         0,
206d519adbfSMatthew Knepley /* 89*/ 0,
2072a6744ebSBarry Smith         0,
2082a6744ebSBarry Smith         0,
2092a6744ebSBarry Smith         0,
2102a6744ebSBarry Smith         0,
211d519adbfSMatthew Knepley /* 94*/ 0,
2122a6744ebSBarry Smith         0,
2132a6744ebSBarry Smith         0,
2142a6744ebSBarry Smith         0};
2152a6744ebSBarry Smith 
2162a6744ebSBarry Smith /*MC
2174a440b9bSBarry Smith    MATSCATTER - MATSCATTER = "scatter" - A matrix type that simply applies a VecScatterBegin/End()
2182a6744ebSBarry Smith 
2192a6744ebSBarry Smith   Level: advanced
2202a6744ebSBarry Smith 
2212a6744ebSBarry Smith .seealso: MatCreateScatter(), MatScatterSetVecScatter(), MatScatterGetVecScatter()
2222a6744ebSBarry Smith 
2232a6744ebSBarry Smith M*/
2242a6744ebSBarry Smith 
2252a6744ebSBarry Smith #undef __FUNCT__
2262a6744ebSBarry Smith #define __FUNCT__ "MatCreate_Scatter"
227*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Scatter(Mat A)
2282a6744ebSBarry Smith {
2292a6744ebSBarry Smith   Mat_Scatter    *b;
2302a6744ebSBarry Smith   PetscErrorCode ierr;
2312a6744ebSBarry Smith 
2322a6744ebSBarry Smith   PetscFunctionBegin;
2332a6744ebSBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
23438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(A,Mat_Scatter,&b);CHKERRQ(ierr);
2352a6744ebSBarry Smith 
2362a6744ebSBarry Smith   A->data = (void*)b;
2372a6744ebSBarry Smith 
23826283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
23926283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
2402a6744ebSBarry Smith 
2412a6744ebSBarry Smith   A->assembled    = PETSC_TRUE;
2422a6744ebSBarry Smith   A->preallocated = PETSC_FALSE;
24374d95942SJed Brown 
24474d95942SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSCATTER);CHKERRQ(ierr);
2452a6744ebSBarry Smith   PetscFunctionReturn(0);
2462a6744ebSBarry Smith }
2472a6744ebSBarry Smith 
2482a6744ebSBarry Smith #undef __FUNCT__
2492a6744ebSBarry Smith #define __FUNCT__ "MatCreateScatter"
2502a6744ebSBarry Smith /*@C
2512a6744ebSBarry Smith    MatCreateScatter - Creates a new matrix based on a VecScatter
2522a6744ebSBarry Smith 
2532a6744ebSBarry Smith   Collective on MPI_Comm
2542a6744ebSBarry Smith 
2552a6744ebSBarry Smith    Input Parameters:
2562a6744ebSBarry Smith +  comm - MPI communicator
2572a6744ebSBarry Smith -  scatter - a VecScatterContext
2582a6744ebSBarry Smith 
2592a6744ebSBarry Smith    Output Parameter:
2602a6744ebSBarry Smith .  A - the matrix
2612a6744ebSBarry Smith 
2622a6744ebSBarry Smith    Level: intermediate
2632a6744ebSBarry Smith 
2642a6744ebSBarry Smith    PETSc requires that matrices and vectors being used for certain
2652a6744ebSBarry Smith    operations are partitioned accordingly.  For example, when
2662a6744ebSBarry Smith    creating a scatter matrix, A, that supports parallel matrix-vector
2672a6744ebSBarry Smith    products using MatMult(A,x,y) the user should set the number
2682a6744ebSBarry Smith    of local matrix rows to be the number of local elements of the
2692a6744ebSBarry Smith    corresponding result vector, y. Note that this is information is
2702a6744ebSBarry Smith    required for use of the matrix interface routines, even though
2712a6744ebSBarry Smith    the scatter matrix may not actually be physically partitioned.
2722a6744ebSBarry Smith 
2732a6744ebSBarry Smith .keywords: matrix, scatter, create
2742a6744ebSBarry Smith 
2752a6744ebSBarry Smith .seealso: MatScatterSetVecScatter(), MatScatterGetVecScatter(), MATSCATTER
2762a6744ebSBarry Smith @*/
2777087cfbeSBarry Smith PetscErrorCode  MatCreateScatter(MPI_Comm comm,VecScatter scatter,Mat *A)
2782a6744ebSBarry Smith {
2792a6744ebSBarry Smith   PetscErrorCode ierr;
2802a6744ebSBarry Smith 
2812a6744ebSBarry Smith   PetscFunctionBegin;
2822a6744ebSBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2832a6744ebSBarry Smith   ierr = MatSetSizes(*A,scatter->to_n,scatter->from_n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
2842a6744ebSBarry Smith   ierr = MatSetType(*A,MATSCATTER);CHKERRQ(ierr);
2852a6744ebSBarry Smith   ierr = MatScatterSetVecScatter(*A,scatter);CHKERRQ(ierr);
28637e1895aSJed Brown   ierr = MatSetUp(*A);CHKERRQ(ierr);
2872a6744ebSBarry Smith   PetscFunctionReturn(0);
2882a6744ebSBarry Smith }
2892a6744ebSBarry Smith 
2902a6744ebSBarry Smith #undef __FUNCT__
2912a6744ebSBarry Smith #define __FUNCT__ "MatScatterSetVecScatter"
2922a6744ebSBarry Smith /*@
2932a6744ebSBarry Smith     MatScatterSetVecScatter - sets that scatter that the matrix is to apply as its linear operator
2942a6744ebSBarry Smith 
2952a6744ebSBarry Smith    Collective on Mat
2962a6744ebSBarry Smith 
2972a6744ebSBarry Smith     Input Parameters:
2982a6744ebSBarry Smith +   mat - the scatter matrix
2992a6744ebSBarry Smith -   scatter - the scatter context create with VecScatterCreate()
3002a6744ebSBarry Smith 
3012a6744ebSBarry Smith    Level: advanced
3022a6744ebSBarry Smith 
3032a6744ebSBarry Smith 
3042a6744ebSBarry Smith .seealso: MatCreateScatter(), MATSCATTER
3052a6744ebSBarry Smith @*/
3067087cfbeSBarry Smith PetscErrorCode  MatScatterSetVecScatter(Mat mat,VecScatter scatter)
3072a6744ebSBarry Smith {
3082a6744ebSBarry Smith   Mat_Scatter    *mscatter = (Mat_Scatter*)mat->data;
3092a6744ebSBarry Smith   PetscErrorCode ierr;
3102a6744ebSBarry Smith 
3112a6744ebSBarry Smith   PetscFunctionBegin;
3120700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
3130700a824SBarry Smith   PetscValidHeaderSpecific(scatter,VEC_SCATTER_CLASSID,2);
3144067f9b5SMatthew Knepley   PetscCheckSameComm((PetscObject)scatter,1,(PetscObject)mat,2);
315e32f2f54SBarry 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);
316e32f2f54SBarry 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);
3172a6744ebSBarry Smith 
318c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)scatter);CHKERRQ(ierr);
3196bf464f9SBarry Smith   ierr = VecScatterDestroy(&mscatter->scatter);CHKERRQ(ierr);
32026fbe8dcSKarl Rupp 
3212a6744ebSBarry Smith   mscatter->scatter = scatter;
3222a6744ebSBarry Smith   PetscFunctionReturn(0);
3232a6744ebSBarry Smith }
3242a6744ebSBarry Smith 
3252a6744ebSBarry Smith 
326