xref: /petsc/src/mat/tutorials/ex8.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
7*9371c9d4SSatish Balay static PetscErrorCode MatScaleUserImpl_SeqAIJ(Mat inA, PetscScalar alpha) {
8c4762a1bSJed Brown   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(MatScale(inA, alpha));
10c4762a1bSJed Brown   PetscFunctionReturn(0);
11c4762a1bSJed Brown }
12c4762a1bSJed Brown 
13c4762a1bSJed Brown extern PetscErrorCode MatScaleUserImpl(Mat, PetscScalar);
14c4762a1bSJed Brown 
15*9371c9d4SSatish Balay static PetscErrorCode MatScaleUserImpl_MPIAIJ(Mat A, PetscScalar aa) {
16c4762a1bSJed Brown   Mat AA, AB;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A, &AA, &AB, NULL));
209566063dSJacob Faibussowitsch   PetscCall(MatScaleUserImpl(AA, aa));
219566063dSJacob Faibussowitsch   PetscCall(MatScaleUserImpl(AB, aa));
22c4762a1bSJed Brown   PetscFunctionReturn(0);
23c4762a1bSJed Brown }
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /* This routine registers MatScaleUserImpl_SeqAIJ() and
26c4762a1bSJed Brown    MatScaleUserImpl_MPIAIJ() as methods providing MatScaleUserImpl()
27c4762a1bSJed Brown    functionality for SeqAIJ and MPIAIJ matrix-types */
28*9371c9d4SSatish Balay PetscErrorCode RegisterMatScaleUserImpl(Mat mat) {
29c4762a1bSJed Brown   PetscMPIInt size;
30c4762a1bSJed Brown 
31362febeeSStefano Zampini   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
33c4762a1bSJed Brown   if (size == 1) { /* SeqAIJ Matrix */
349566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatScaleUserImpl_C", MatScaleUserImpl_SeqAIJ));
35c4762a1bSJed Brown   } else { /* MPIAIJ Matrix */
36c4762a1bSJed Brown     Mat AA, AB;
379566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJGetSeqAIJ(mat, &AA, &AB, NULL));
389566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatScaleUserImpl_C", MatScaleUserImpl_MPIAIJ));
399566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)AA, "MatScaleUserImpl_C", MatScaleUserImpl_SeqAIJ));
409566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)AB, "MatScaleUserImpl_C", MatScaleUserImpl_SeqAIJ));
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown   PetscFunctionReturn(0);
43c4762a1bSJed Brown }
44c4762a1bSJed Brown 
452e956fe4SStefano Zampini /* This routine deregisters MatScaleUserImpl_SeqAIJ() and
462e956fe4SStefano Zampini    MatScaleUserImpl_MPIAIJ() as methods providing MatScaleUserImpl()
472e956fe4SStefano Zampini    functionality for SeqAIJ and MPIAIJ matrix-types */
48*9371c9d4SSatish Balay PetscErrorCode DeRegisterMatScaleUserImpl(Mat mat) {
492e956fe4SStefano Zampini   PetscMPIInt size;
502e956fe4SStefano Zampini 
512e956fe4SStefano Zampini   PetscFunctionBegin;
522e956fe4SStefano Zampini   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
532e956fe4SStefano Zampini   if (size == 1) { /* SeqAIJ Matrix */
542e956fe4SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatScaleUserImpl_C", NULL));
552e956fe4SStefano Zampini   } else { /* MPIAIJ Matrix */
562e956fe4SStefano Zampini     Mat AA, AB;
572e956fe4SStefano Zampini     PetscCall(MatMPIAIJGetSeqAIJ(mat, &AA, &AB, NULL));
582e956fe4SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatScaleUserImpl_C", NULL));
592e956fe4SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)AA, "MatScaleUserImpl_C", NULL));
602e956fe4SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)AB, "MatScaleUserImpl_C", NULL));
612e956fe4SStefano Zampini   }
622e956fe4SStefano Zampini   PetscFunctionReturn(0);
632e956fe4SStefano Zampini }
642e956fe4SStefano Zampini 
65c4762a1bSJed Brown /* this routines queries the already registered MatScaleUserImp_XXX
66c4762a1bSJed Brown    implementations for the given matrix, and calls the correct
67c4762a1bSJed Brown    routine. i.e if MatType is SeqAIJ, MatScaleUserImpl_SeqAIJ() gets
68c4762a1bSJed Brown    called, and if MatType is MPIAIJ, MatScaleUserImpl_MPIAIJ() gets
69c4762a1bSJed Brown    called */
70*9371c9d4SSatish Balay PetscErrorCode MatScaleUserImpl(Mat mat, PetscScalar a) {
715f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(Mat, PetscScalar);
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   PetscFunctionBegin;
749566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)mat, "MatScaleUserImpl_C", &f));
759566063dSJacob Faibussowitsch   if (f) PetscCall((*f)(mat, a));
76c4762a1bSJed Brown   PetscFunctionReturn(0);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
79c4762a1bSJed Brown /* Main user code that uses MatScaleUserImpl() */
80c4762a1bSJed Brown 
81*9371c9d4SSatish Balay int main(int argc, char **args) {
82c4762a1bSJed Brown   Mat         mat;
83c4762a1bSJed Brown   PetscInt    i, j, m = 2, n, Ii, J;
84c4762a1bSJed Brown   PetscScalar v, none = -1.0;
85c4762a1bSJed Brown   PetscMPIInt rank, size;
86c4762a1bSJed Brown 
87327415f7SBarry Smith   PetscFunctionBeginUser;
889566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
91c4762a1bSJed Brown   n = 2 * size;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* create the matrix */
949566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &mat));
959566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
969566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, MATAIJ));
979566063dSJacob Faibussowitsch   PetscCall(MatSetUp(mat));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* register user defined MatScaleUser() operation for both SeqAIJ
100c4762a1bSJed Brown      and MPIAIJ types */
1019566063dSJacob Faibussowitsch   PetscCall(RegisterMatScaleUserImpl(mat));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* assemble the matrix */
104c4762a1bSJed Brown   for (i = 0; i < m; i++) {
105c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
106*9371c9d4SSatish Balay       v  = -1.0;
107*9371c9d4SSatish Balay       Ii = j + n * i;
108*9371c9d4SSatish Balay       if (i > 0) {
109*9371c9d4SSatish Balay         J = Ii - n;
110*9371c9d4SSatish Balay         PetscCall(MatSetValues(mat, 1, &Ii, 1, &J, &v, INSERT_VALUES));
111*9371c9d4SSatish Balay       }
112*9371c9d4SSatish Balay       if (i < m - 1) {
113*9371c9d4SSatish Balay         J = Ii + n;
114*9371c9d4SSatish Balay         PetscCall(MatSetValues(mat, 1, &Ii, 1, &J, &v, INSERT_VALUES));
115*9371c9d4SSatish Balay       }
116*9371c9d4SSatish Balay       if (j > 0) {
117*9371c9d4SSatish Balay         J = Ii - 1;
118*9371c9d4SSatish Balay         PetscCall(MatSetValues(mat, 1, &Ii, 1, &J, &v, INSERT_VALUES));
119*9371c9d4SSatish Balay       }
120*9371c9d4SSatish Balay       if (j < n - 1) {
121*9371c9d4SSatish Balay         J = Ii + 1;
122*9371c9d4SSatish Balay         PetscCall(MatSetValues(mat, 1, &Ii, 1, &J, &v, INSERT_VALUES));
123*9371c9d4SSatish Balay       }
124*9371c9d4SSatish Balay       v = 4.0;
125*9371c9d4SSatish Balay       PetscCall(MatSetValues(mat, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
126c4762a1bSJed Brown     }
127c4762a1bSJed Brown   }
1289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
1299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* check the matrix before and after scaling by -1.0 */
1329566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Matrix _before_ MatScaleUserImpl() operation\n"));
1339566063dSJacob Faibussowitsch   PetscCall(MatView(mat, PETSC_VIEWER_STDOUT_WORLD));
1349566063dSJacob Faibussowitsch   PetscCall(MatScaleUserImpl(mat, none));
1359566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Matrix _after_ MatScaleUserImpl() operation\n"));
1369566063dSJacob Faibussowitsch   PetscCall(MatView(mat, PETSC_VIEWER_STDOUT_WORLD));
137c4762a1bSJed Brown 
1382e956fe4SStefano Zampini   /* deregister user defined MatScaleUser() operation for both SeqAIJ
1392e956fe4SStefano Zampini      and MPIAIJ types */
1402e956fe4SStefano Zampini   PetscCall(DeRegisterMatScaleUserImpl(mat));
1419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
1429566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
143b122ec5aSJacob Faibussowitsch   return 0;
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown /*TEST
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    test:
149c4762a1bSJed Brown 
150c4762a1bSJed Brown    test:
151c4762a1bSJed Brown       suffix: 2
152c4762a1bSJed Brown       nsize: 2
153c4762a1bSJed Brown 
154c4762a1bSJed Brown TEST*/
155