xref: /petsc/src/mat/impls/scatter/mscatter.c (revision 38f2d2fdb3b6f522a3102c6eb796cebecf3224c0)
12a6744ebSBarry Smith #define PETSCMAT_DLL
22a6744ebSBarry Smith 
32a6744ebSBarry Smith /*
42a6744ebSBarry Smith    This provides a matrix that applies a VecScatter to a vector.
52a6744ebSBarry Smith */
62a6744ebSBarry Smith 
7b9147fbbSdalcinl #include "include/private/matimpl.h"        /*I "petscmat.h" I*/
82a6744ebSBarry Smith #include "private/vecimpl.h"
92a6744ebSBarry Smith 
102a6744ebSBarry Smith typedef struct {
112a6744ebSBarry Smith   VecScatter scatter;
122a6744ebSBarry Smith } Mat_Scatter;
132a6744ebSBarry Smith 
142a6744ebSBarry Smith #undef __FUNCT__
152a6744ebSBarry Smith #define __FUNCT__ "MatScatterGetVecScatter"
162a6744ebSBarry Smith /*@
172a6744ebSBarry Smith     MatScatterGetVecScatter - Returns the user-provided scatter set with MatScatterSetVecScatter()
182a6744ebSBarry Smith 
192a6744ebSBarry Smith     Not Collective, but not cannot use scatter if not used collectively on Mat
202a6744ebSBarry Smith 
212a6744ebSBarry Smith     Input Parameter:
222a6744ebSBarry Smith .   mat - the matrix, should have been created with MatCreateScatter() or have type MATSCATTER
232a6744ebSBarry Smith 
242a6744ebSBarry Smith     Output Parameter:
252a6744ebSBarry Smith .   scatter - the scatter context
262a6744ebSBarry Smith 
272a6744ebSBarry Smith     Level: intermediate
282a6744ebSBarry Smith 
292a6744ebSBarry Smith .keywords: matrix, scatter, get
302a6744ebSBarry Smith 
312a6744ebSBarry Smith .seealso: MatCreateScatter(), MatScatterSetVecScatter(), MATSCATTER
322a6744ebSBarry Smith @*/
332a6744ebSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat mat,VecScatter *scatter)
342a6744ebSBarry Smith {
352a6744ebSBarry Smith   Mat_Scatter    *mscatter;
362a6744ebSBarry Smith 
372a6744ebSBarry Smith   PetscFunctionBegin;
382a6744ebSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
392a6744ebSBarry Smith   PetscValidPointer(scatter,2);
402a6744ebSBarry Smith   mscatter = (Mat_Scatter*)mat->data;
412a6744ebSBarry Smith   *scatter = mscatter->scatter;
422a6744ebSBarry Smith   PetscFunctionReturn(0);
432a6744ebSBarry Smith }
442a6744ebSBarry Smith 
452a6744ebSBarry Smith #undef __FUNCT__
462a6744ebSBarry Smith #define __FUNCT__ "MatDestroy_Scatter"
472a6744ebSBarry Smith PetscErrorCode MatDestroy_Scatter(Mat mat)
482a6744ebSBarry Smith {
492a6744ebSBarry Smith   PetscErrorCode ierr;
502a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)mat->data;
512a6744ebSBarry Smith 
522a6744ebSBarry Smith   PetscFunctionBegin;
532a6744ebSBarry Smith   if (scatter->scatter) {ierr = VecScatterDestroy(scatter->scatter);CHKERRQ(ierr);}
542a6744ebSBarry Smith   ierr = PetscFree(scatter);CHKERRQ(ierr);
552a6744ebSBarry Smith   PetscFunctionReturn(0);
562a6744ebSBarry Smith }
572a6744ebSBarry Smith 
582a6744ebSBarry Smith #undef __FUNCT__
592a6744ebSBarry Smith #define __FUNCT__ "MatMult_Scatter"
602a6744ebSBarry Smith PetscErrorCode MatMult_Scatter(Mat A,Vec x,Vec y)
612a6744ebSBarry Smith {
622a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
632a6744ebSBarry Smith   PetscErrorCode ierr;
642a6744ebSBarry Smith 
652a6744ebSBarry Smith   PetscFunctionBegin;
662a6744ebSBarry Smith   if (!scatter->scatter) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
67ca9f406cSSatish Balay   ierr = VecScatterBegin(scatter->scatter,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
68ca9f406cSSatish Balay   ierr = VecScatterEnd(scatter->scatter,x,y,INSERT_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;
802a6744ebSBarry Smith   if (!scatter->scatter) SETERRQ(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;
952a6744ebSBarry Smith   if (!scatter->scatter) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
96ca9f406cSSatish Balay   ierr = VecScatterBegin(scatter->scatter,x,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
97ca9f406cSSatish Balay   ierr = VecScatterEnd(scatter->scatter,x,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
982a6744ebSBarry Smith   PetscFunctionReturn(0);
992a6744ebSBarry Smith }
1002a6744ebSBarry Smith 
1012a6744ebSBarry Smith #undef __FUNCT__
1022a6744ebSBarry Smith #define __FUNCT__ "MatMultTransposeAdd_Scatter"
1032a6744ebSBarry Smith PetscErrorCode MatMultTransposeAdd_Scatter(Mat A,Vec x,Vec y,Vec z)
1042a6744ebSBarry Smith {
1052a6744ebSBarry Smith   Mat_Scatter    *scatter = (Mat_Scatter*)A->data;
1062a6744ebSBarry Smith   PetscErrorCode ierr;
1072a6744ebSBarry Smith 
1082a6744ebSBarry Smith   PetscFunctionBegin;
1092a6744ebSBarry Smith   if (!scatter->scatter) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Need to first call MatScatterSetScatter()");
1102a6744ebSBarry Smith   if (z != y) {ierr = VecCopy(y,z);CHKERRQ(ierr);}
111ca9f406cSSatish Balay   ierr = VecScatterBegin(scatter->scatter,x,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
112ca9f406cSSatish Balay   ierr = VecScatterEnd(scatter->scatter,x,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1132a6744ebSBarry Smith   PetscFunctionReturn(0);
1142a6744ebSBarry Smith }
1152a6744ebSBarry Smith 
1162a6744ebSBarry Smith static struct _MatOps MatOps_Values = {0,
1172a6744ebSBarry Smith        0,
1182a6744ebSBarry Smith        0,
1192a6744ebSBarry Smith        MatMult_Scatter,
1202a6744ebSBarry Smith /* 4*/ MatMultAdd_Scatter,
1212a6744ebSBarry Smith        MatMultTranspose_Scatter,
1222a6744ebSBarry Smith        MatMultTransposeAdd_Scatter,
1232a6744ebSBarry Smith        0,
1242a6744ebSBarry Smith        0,
1252a6744ebSBarry Smith        0,
1262a6744ebSBarry Smith /*10*/ 0,
1272a6744ebSBarry Smith        0,
1282a6744ebSBarry Smith        0,
1292a6744ebSBarry Smith        0,
1302a6744ebSBarry Smith        0,
1312a6744ebSBarry Smith /*15*/ 0,
1322a6744ebSBarry Smith        0,
1332a6744ebSBarry Smith        0,
1342a6744ebSBarry Smith        0,
1352a6744ebSBarry Smith        0,
1362a6744ebSBarry Smith /*20*/ 0,
1372a6744ebSBarry Smith        0,
1382a6744ebSBarry Smith        0,
1392a6744ebSBarry Smith        0,
1402a6744ebSBarry Smith        0,
1412a6744ebSBarry Smith /*25*/ 0,
1422a6744ebSBarry Smith        0,
1432a6744ebSBarry Smith        0,
1442a6744ebSBarry Smith        0,
1452a6744ebSBarry Smith        0,
1462a6744ebSBarry Smith /*30*/ 0,
1472a6744ebSBarry Smith        0,
1482a6744ebSBarry Smith        0,
1492a6744ebSBarry Smith        0,
1502a6744ebSBarry Smith        0,
1512a6744ebSBarry Smith /*35*/ 0,
1522a6744ebSBarry Smith        0,
1532a6744ebSBarry Smith        0,
1542a6744ebSBarry Smith        0,
1552a6744ebSBarry Smith        0,
1562a6744ebSBarry Smith /*40*/ 0,
1572a6744ebSBarry Smith        0,
1582a6744ebSBarry Smith        0,
1592a6744ebSBarry Smith        0,
1602a6744ebSBarry Smith        0,
1612a6744ebSBarry Smith /*45*/ 0,
1622a6744ebSBarry Smith        0,
1632a6744ebSBarry Smith        0,
1642a6744ebSBarry Smith        0,
1652a6744ebSBarry Smith        0,
1662a6744ebSBarry Smith /*50*/ 0,
1672a6744ebSBarry Smith        0,
1682a6744ebSBarry Smith        0,
1692a6744ebSBarry Smith        0,
1702a6744ebSBarry Smith        0,
1712a6744ebSBarry Smith /*55*/ 0,
1722a6744ebSBarry Smith        0,
1732a6744ebSBarry Smith        0,
1742a6744ebSBarry Smith        0,
1752a6744ebSBarry Smith        0,
1762a6744ebSBarry Smith /*60*/ 0,
1772a6744ebSBarry Smith        MatDestroy_Scatter,
1782a6744ebSBarry Smith        0,
1792a6744ebSBarry Smith        0,
1802a6744ebSBarry Smith        0,
1812a6744ebSBarry Smith /*65*/ 0,
1822a6744ebSBarry Smith        0,
1832a6744ebSBarry Smith        0,
1842a6744ebSBarry Smith        0,
1852a6744ebSBarry Smith        0,
1862a6744ebSBarry Smith /*70*/ 0,
1872a6744ebSBarry Smith        0,
1882a6744ebSBarry Smith        0,
1892a6744ebSBarry Smith        0,
1902a6744ebSBarry Smith        0,
1912a6744ebSBarry Smith /*75*/ 0,
1922a6744ebSBarry Smith        0,
1932a6744ebSBarry Smith        0,
1942a6744ebSBarry Smith        0,
1952a6744ebSBarry Smith        0,
1962a6744ebSBarry Smith /*80*/ 0,
1972a6744ebSBarry Smith        0,
1982a6744ebSBarry Smith        0,
1992a6744ebSBarry Smith        0,
2002a6744ebSBarry Smith        0,
2012a6744ebSBarry Smith /*85*/ 0,
2022a6744ebSBarry Smith        0,
2032a6744ebSBarry Smith        0,
2042a6744ebSBarry Smith        0,
2052a6744ebSBarry Smith        0,
2062a6744ebSBarry Smith /*90*/ 0,
2072a6744ebSBarry Smith        0,
2082a6744ebSBarry Smith        0,
2092a6744ebSBarry Smith        0,
2102a6744ebSBarry Smith        0,
2112a6744ebSBarry Smith /*95*/ 0,
2122a6744ebSBarry Smith        0,
2132a6744ebSBarry Smith        0,
2142a6744ebSBarry Smith        0};
2152a6744ebSBarry Smith 
2162a6744ebSBarry Smith /*MC
2172a6744ebSBarry Smith    MATSCATTER - MATSCATTER = "scatter" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
2182a6744ebSBarry Smith 
2192a6744ebSBarry Smith   Level: advanced
2202a6744ebSBarry Smith 
2212a6744ebSBarry Smith .seealso: MatCreateScatter(), MatScatterSetVecScatter(), MatScatterGetVecScatter()
2222a6744ebSBarry Smith 
2232a6744ebSBarry Smith M*/
2242a6744ebSBarry Smith 
2252a6744ebSBarry Smith EXTERN_C_BEGIN
2262a6744ebSBarry Smith #undef __FUNCT__
2272a6744ebSBarry Smith #define __FUNCT__ "MatCreate_Scatter"
2282a6744ebSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Scatter(Mat A)
2292a6744ebSBarry Smith {
2302a6744ebSBarry Smith   Mat_Scatter    *b;
2312a6744ebSBarry Smith   PetscErrorCode ierr;
2322a6744ebSBarry Smith 
2332a6744ebSBarry Smith   PetscFunctionBegin;
2342a6744ebSBarry Smith   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
235*38f2d2fdSLisandro Dalcin   ierr = PetscNewLog(A,Mat_Scatter,&b);CHKERRQ(ierr);
2362a6744ebSBarry Smith 
2372a6744ebSBarry Smith   A->data = (void*)b;
2382a6744ebSBarry Smith 
2392813ca57SBarry Smith   ierr = PetscMapSetBlockSize(&A->rmap,1);CHKERRQ(ierr);
2402813ca57SBarry Smith   ierr = PetscMapSetBlockSize(&A->cmap,1);CHKERRQ(ierr);
2416148ca0dSBarry Smith   ierr = PetscMapSetUp(&A->rmap);CHKERRQ(ierr);
2426148ca0dSBarry Smith   ierr = PetscMapSetUp(&A->cmap);CHKERRQ(ierr);
2432a6744ebSBarry Smith 
2442a6744ebSBarry Smith   A->assembled     = PETSC_TRUE;
2452a6744ebSBarry Smith   A->preallocated  = PETSC_FALSE;
2462a6744ebSBarry Smith   PetscFunctionReturn(0);
2472a6744ebSBarry Smith }
2482a6744ebSBarry Smith EXTERN_C_END
2492a6744ebSBarry Smith 
2502a6744ebSBarry Smith #undef __FUNCT__
2512a6744ebSBarry Smith #define __FUNCT__ "MatCreateScatter"
2522a6744ebSBarry Smith /*@C
2532a6744ebSBarry Smith    MatCreateScatter - Creates a new matrix based on a VecScatter
2542a6744ebSBarry Smith 
2552a6744ebSBarry Smith   Collective on MPI_Comm
2562a6744ebSBarry Smith 
2572a6744ebSBarry Smith    Input Parameters:
2582a6744ebSBarry Smith +  comm - MPI communicator
2592a6744ebSBarry Smith -  scatter - a VecScatterContext
2602a6744ebSBarry Smith 
2612a6744ebSBarry Smith    Output Parameter:
2622a6744ebSBarry Smith .  A - the matrix
2632a6744ebSBarry Smith 
2642a6744ebSBarry Smith    Level: intermediate
2652a6744ebSBarry Smith 
2662a6744ebSBarry Smith    PETSc requires that matrices and vectors being used for certain
2672a6744ebSBarry Smith    operations are partitioned accordingly.  For example, when
2682a6744ebSBarry Smith    creating a scatter matrix, A, that supports parallel matrix-vector
2692a6744ebSBarry Smith    products using MatMult(A,x,y) the user should set the number
2702a6744ebSBarry Smith    of local matrix rows to be the number of local elements of the
2712a6744ebSBarry Smith    corresponding result vector, y. Note that this is information is
2722a6744ebSBarry Smith    required for use of the matrix interface routines, even though
2732a6744ebSBarry Smith    the scatter matrix may not actually be physically partitioned.
2742a6744ebSBarry Smith    For example,
2752a6744ebSBarry Smith 
2762a6744ebSBarry Smith .keywords: matrix, scatter, create
2772a6744ebSBarry Smith 
2782a6744ebSBarry Smith .seealso: MatScatterSetVecScatter(), MatScatterGetVecScatter(), MATSCATTER
2792a6744ebSBarry Smith @*/
2802a6744ebSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm comm,VecScatter scatter,Mat *A)
2812a6744ebSBarry Smith {
2822a6744ebSBarry Smith   PetscErrorCode ierr;
2832a6744ebSBarry Smith 
2842a6744ebSBarry Smith   PetscFunctionBegin;
2852a6744ebSBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2862a6744ebSBarry Smith   ierr = MatSetSizes(*A,scatter->to_n,scatter->from_n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
2872a6744ebSBarry Smith   ierr = MatSetType(*A,MATSCATTER);CHKERRQ(ierr);
2882a6744ebSBarry Smith   ierr = MatScatterSetVecScatter(*A,scatter);CHKERRQ(ierr);
2892a6744ebSBarry Smith   PetscFunctionReturn(0);
2902a6744ebSBarry Smith }
2912a6744ebSBarry Smith 
2922a6744ebSBarry Smith #undef __FUNCT__
2932a6744ebSBarry Smith #define __FUNCT__ "MatScatterSetVecScatter"
2942a6744ebSBarry Smith /*@
2952a6744ebSBarry Smith     MatScatterSetVecScatter - sets that scatter that the matrix is to apply as its linear operator
2962a6744ebSBarry Smith 
2972a6744ebSBarry Smith    Collective on Mat
2982a6744ebSBarry Smith 
2992a6744ebSBarry Smith     Input Parameters:
3002a6744ebSBarry Smith +   mat - the scatter matrix
3012a6744ebSBarry Smith -   scatter - the scatter context create with VecScatterCreate()
3022a6744ebSBarry Smith 
3032a6744ebSBarry Smith    Level: advanced
3042a6744ebSBarry Smith 
3052a6744ebSBarry Smith 
3062a6744ebSBarry Smith .seealso: MatCreateScatter(), MATSCATTER
3072a6744ebSBarry Smith @*/
3082a6744ebSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat mat,VecScatter scatter)
3092a6744ebSBarry Smith {
3102a6744ebSBarry Smith   Mat_Scatter    *mscatter = (Mat_Scatter*)mat->data;
3112a6744ebSBarry Smith   PetscErrorCode ierr;
3122a6744ebSBarry Smith 
3132a6744ebSBarry Smith   PetscFunctionBegin;
3142a6744ebSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
3152a6744ebSBarry Smith   PetscValidHeaderSpecific(scatter,VEC_SCATTER_COOKIE,2);
3164067f9b5SMatthew Knepley   PetscCheckSameComm((PetscObject)scatter,1,(PetscObject)mat,2);
3172a6744ebSBarry 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);
3182a6744ebSBarry 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);
3192a6744ebSBarry Smith 
320c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)scatter);CHKERRQ(ierr);
3212a6744ebSBarry Smith   if (mscatter->scatter) {ierr = VecScatterDestroy(mscatter->scatter);CHKERRQ(ierr);}
3222a6744ebSBarry Smith   mscatter->scatter = scatter;
3232a6744ebSBarry Smith   PetscFunctionReturn(0);
3242a6744ebSBarry Smith }
3252a6744ebSBarry Smith 
3262a6744ebSBarry Smith 
327