xref: /petsc/src/mat/tutorials/ex8.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Shows how to add a new MatOperation to AIJ MatType\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown #include <petscblaslapack.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown static PetscErrorCode MatScaleUserImpl_SeqAIJ(Mat inA,PetscScalar alpha)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(MatScale(inA,alpha));
11c4762a1bSJed Brown   PetscFunctionReturn(0);
12c4762a1bSJed Brown }
13c4762a1bSJed Brown 
14c4762a1bSJed Brown extern PetscErrorCode MatScaleUserImpl(Mat,PetscScalar);
15c4762a1bSJed Brown 
16c4762a1bSJed Brown static PetscErrorCode MatScaleUserImpl_MPIAIJ(Mat A,PetscScalar aa)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   Mat            AA,AB;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL));
229566063dSJacob Faibussowitsch   PetscCall(MatScaleUserImpl(AA,aa));
239566063dSJacob Faibussowitsch   PetscCall(MatScaleUserImpl(AB,aa));
24c4762a1bSJed Brown   PetscFunctionReturn(0);
25c4762a1bSJed Brown }
26c4762a1bSJed Brown 
27c4762a1bSJed Brown /* This routine registers MatScaleUserImpl_SeqAIJ() and
28c4762a1bSJed Brown    MatScaleUserImpl_MPIAIJ() as methods providing MatScaleUserImpl()
29c4762a1bSJed Brown    functionality for SeqAIJ and MPIAIJ matrix-types */
30c4762a1bSJed Brown PetscErrorCode RegisterMatScaleUserImpl(Mat mat)
31c4762a1bSJed Brown {
32c4762a1bSJed Brown   PetscMPIInt    size;
33c4762a1bSJed Brown 
34362febeeSStefano Zampini   PetscFunctionBegin;
359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
36c4762a1bSJed Brown   if (size == 1) { /* SeqAIJ Matrix */
379566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatScaleUserImpl_C",MatScaleUserImpl_SeqAIJ));
38c4762a1bSJed Brown   } else { /* MPIAIJ Matrix */
39c4762a1bSJed Brown     Mat AA,AB;
409566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJGetSeqAIJ(mat,&AA,&AB,NULL));
419566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatScaleUserImpl_C",MatScaleUserImpl_MPIAIJ));
429566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)AA,"MatScaleUserImpl_C",MatScaleUserImpl_SeqAIJ));
439566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)AB,"MatScaleUserImpl_C",MatScaleUserImpl_SeqAIJ));
44c4762a1bSJed Brown   }
45c4762a1bSJed Brown   PetscFunctionReturn(0);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48*2e956fe4SStefano Zampini /* This routine deregisters MatScaleUserImpl_SeqAIJ() and
49*2e956fe4SStefano Zampini    MatScaleUserImpl_MPIAIJ() as methods providing MatScaleUserImpl()
50*2e956fe4SStefano Zampini    functionality for SeqAIJ and MPIAIJ matrix-types */
51*2e956fe4SStefano Zampini PetscErrorCode DeRegisterMatScaleUserImpl(Mat mat)
52*2e956fe4SStefano Zampini {
53*2e956fe4SStefano Zampini   PetscMPIInt    size;
54*2e956fe4SStefano Zampini 
55*2e956fe4SStefano Zampini   PetscFunctionBegin;
56*2e956fe4SStefano Zampini   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
57*2e956fe4SStefano Zampini   if (size == 1) { /* SeqAIJ Matrix */
58*2e956fe4SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatScaleUserImpl_C",NULL));
59*2e956fe4SStefano Zampini   } else { /* MPIAIJ Matrix */
60*2e956fe4SStefano Zampini     Mat AA,AB;
61*2e956fe4SStefano Zampini     PetscCall(MatMPIAIJGetSeqAIJ(mat,&AA,&AB,NULL));
62*2e956fe4SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatScaleUserImpl_C",NULL));
63*2e956fe4SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)AA,"MatScaleUserImpl_C",NULL));
64*2e956fe4SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)AB,"MatScaleUserImpl_C",NULL));
65*2e956fe4SStefano Zampini   }
66*2e956fe4SStefano Zampini   PetscFunctionReturn(0);
67*2e956fe4SStefano Zampini }
68*2e956fe4SStefano Zampini 
69c4762a1bSJed Brown /* this routines queries the already registered MatScaleUserImp_XXX
70c4762a1bSJed Brown    implementations for the given matrix, and calls the correct
71c4762a1bSJed Brown    routine. i.e if MatType is SeqAIJ, MatScaleUserImpl_SeqAIJ() gets
72c4762a1bSJed Brown    called, and if MatType is MPIAIJ, MatScaleUserImpl_MPIAIJ() gets
73c4762a1bSJed Brown    called */
74c4762a1bSJed Brown PetscErrorCode MatScaleUserImpl(Mat mat,PetscScalar a)
75c4762a1bSJed Brown {
765f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(Mat,PetscScalar);
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   PetscFunctionBegin;
799566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)mat,"MatScaleUserImpl_C",&f));
809566063dSJacob Faibussowitsch   if (f) PetscCall((*f)(mat,a));
81c4762a1bSJed Brown   PetscFunctionReturn(0);
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown /* Main user code that uses MatScaleUserImpl() */
85c4762a1bSJed Brown 
86c4762a1bSJed Brown int main(int argc,char **args)
87c4762a1bSJed Brown {
88c4762a1bSJed Brown   Mat            mat;
89c4762a1bSJed Brown   PetscInt       i,j,m = 2,n,Ii,J;
90c4762a1bSJed Brown   PetscScalar    v,none = -1.0;
91c4762a1bSJed Brown   PetscMPIInt    rank,size;
92c4762a1bSJed Brown 
939566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
96c4762a1bSJed Brown   n    = 2*size;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* create the matrix */
999566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&mat));
1009566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
1019566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,MATAIJ));
1029566063dSJacob Faibussowitsch   PetscCall(MatSetUp(mat));
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* register user defined MatScaleUser() operation for both SeqAIJ
105c4762a1bSJed Brown      and MPIAIJ types */
1069566063dSJacob Faibussowitsch   PetscCall(RegisterMatScaleUserImpl(mat));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* assemble the matrix */
109c4762a1bSJed Brown   for (i=0; i<m; i++) {
110c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
111c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
1129566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(mat,1,&Ii,1,&J,&v,INSERT_VALUES));}
1139566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(mat,1,&Ii,1,&J,&v,INSERT_VALUES));}
1149566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(mat,1,&Ii,1,&J,&v,INSERT_VALUES));}
1159566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(mat,1,&Ii,1,&J,&v,INSERT_VALUES));}
1169566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(mat,1,&Ii,1,&Ii,&v,INSERT_VALUES));
117c4762a1bSJed Brown     }
118c4762a1bSJed Brown   }
1199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
1209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* check the matrix before and after scaling by -1.0 */
1239566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix _before_ MatScaleUserImpl() operation\n"));
1249566063dSJacob Faibussowitsch   PetscCall(MatView(mat,PETSC_VIEWER_STDOUT_WORLD));
1259566063dSJacob Faibussowitsch   PetscCall(MatScaleUserImpl(mat,none));
1269566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix _after_ MatScaleUserImpl() operation\n"));
1279566063dSJacob Faibussowitsch   PetscCall(MatView(mat,PETSC_VIEWER_STDOUT_WORLD));
128c4762a1bSJed Brown 
129*2e956fe4SStefano Zampini   /* deregister user defined MatScaleUser() operation for both SeqAIJ
130*2e956fe4SStefano Zampini      and MPIAIJ types */
131*2e956fe4SStefano Zampini   PetscCall(DeRegisterMatScaleUserImpl(mat));
1329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
1339566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
134b122ec5aSJacob Faibussowitsch   return 0;
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown /*TEST
138c4762a1bSJed Brown 
139c4762a1bSJed Brown    test:
140c4762a1bSJed Brown 
141c4762a1bSJed Brown    test:
142c4762a1bSJed Brown       suffix: 2
143c4762a1bSJed Brown       nsize: 2
144c4762a1bSJed Brown 
145c4762a1bSJed Brown TEST*/
146