1c4762a1bSJed Brown static char help[] = "Shows how to add a new MatOperation to AIJ MatType\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown #include <petscblaslapack.h> 5c4762a1bSJed Brown 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScaleUserImpl_SeqAIJ(Mat inA, PetscScalar alpha) 7d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown PetscFunctionBegin; 99566063dSJacob Faibussowitsch PetscCall(MatScale(inA, alpha)); 103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11c4762a1bSJed Brown } 12c4762a1bSJed Brown 13c4762a1bSJed Brown extern PetscErrorCode MatScaleUserImpl(Mat, PetscScalar); 14c4762a1bSJed Brown 15d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScaleUserImpl_MPIAIJ(Mat A, PetscScalar aa) 16d71ae5a4SJacob Faibussowitsch { 17c4762a1bSJed Brown Mat AA, AB; 18c4762a1bSJed Brown 19c4762a1bSJed Brown PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(A, &AA, &AB, NULL)); 219566063dSJacob Faibussowitsch PetscCall(MatScaleUserImpl(AA, aa)); 229566063dSJacob Faibussowitsch PetscCall(MatScaleUserImpl(AB, aa)); 233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* This routine registers MatScaleUserImpl_SeqAIJ() and 27c4762a1bSJed Brown MatScaleUserImpl_MPIAIJ() as methods providing MatScaleUserImpl() 28c4762a1bSJed Brown functionality for SeqAIJ and MPIAIJ matrix-types */ 29d71ae5a4SJacob Faibussowitsch PetscErrorCode RegisterMatScaleUserImpl(Mat mat) 30d71ae5a4SJacob Faibussowitsch { 31c4762a1bSJed Brown PetscMPIInt size; 32c4762a1bSJed Brown 33362febeeSStefano Zampini PetscFunctionBegin; 349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 35c4762a1bSJed Brown if (size == 1) { /* SeqAIJ Matrix */ 369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatScaleUserImpl_C", MatScaleUserImpl_SeqAIJ)); 37c4762a1bSJed Brown } else { /* MPIAIJ Matrix */ 38c4762a1bSJed Brown Mat AA, AB; 399566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(mat, &AA, &AB, NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatScaleUserImpl_C", MatScaleUserImpl_MPIAIJ)); 419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)AA, "MatScaleUserImpl_C", MatScaleUserImpl_SeqAIJ)); 429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)AB, "MatScaleUserImpl_C", MatScaleUserImpl_SeqAIJ)); 43c4762a1bSJed Brown } 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 472e956fe4SStefano Zampini /* This routine deregisters MatScaleUserImpl_SeqAIJ() and 482e956fe4SStefano Zampini MatScaleUserImpl_MPIAIJ() as methods providing MatScaleUserImpl() 492e956fe4SStefano Zampini functionality for SeqAIJ and MPIAIJ matrix-types */ 50d71ae5a4SJacob Faibussowitsch PetscErrorCode DeRegisterMatScaleUserImpl(Mat mat) 51d71ae5a4SJacob Faibussowitsch { 522e956fe4SStefano Zampini PetscMPIInt size; 532e956fe4SStefano Zampini 542e956fe4SStefano Zampini PetscFunctionBegin; 552e956fe4SStefano Zampini PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 562e956fe4SStefano Zampini if (size == 1) { /* SeqAIJ Matrix */ 572e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatScaleUserImpl_C", NULL)); 582e956fe4SStefano Zampini } else { /* MPIAIJ Matrix */ 592e956fe4SStefano Zampini Mat AA, AB; 602e956fe4SStefano Zampini PetscCall(MatMPIAIJGetSeqAIJ(mat, &AA, &AB, NULL)); 612e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatScaleUserImpl_C", NULL)); 622e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)AA, "MatScaleUserImpl_C", NULL)); 632e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)AB, "MatScaleUserImpl_C", NULL)); 642e956fe4SStefano Zampini } 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 662e956fe4SStefano Zampini } 672e956fe4SStefano Zampini 68c4762a1bSJed Brown /* this routines queries the already registered MatScaleUserImp_XXX 69c4762a1bSJed Brown implementations for the given matrix, and calls the correct 70c4762a1bSJed Brown routine. i.e if MatType is SeqAIJ, MatScaleUserImpl_SeqAIJ() gets 71c4762a1bSJed Brown called, and if MatType is MPIAIJ, MatScaleUserImpl_MPIAIJ() gets 72c4762a1bSJed Brown called */ 73d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScaleUserImpl(Mat mat, PetscScalar a) 74d71ae5a4SJacob Faibussowitsch { 755f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(Mat, PetscScalar); 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscFunctionBegin; 789566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)mat, "MatScaleUserImpl_C", &f)); 799566063dSJacob Faibussowitsch if (f) PetscCall((*f)(mat, a)); 803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Main user code that uses MatScaleUserImpl() */ 84c4762a1bSJed Brown 85d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 86d71ae5a4SJacob Faibussowitsch { 87c4762a1bSJed Brown Mat mat; 88c4762a1bSJed Brown PetscInt i, j, m = 2, n, Ii, J; 89c4762a1bSJed Brown PetscScalar v, none = -1.0; 90c4762a1bSJed Brown PetscMPIInt rank, size; 91c4762a1bSJed Brown 92327415f7SBarry Smith PetscFunctionBeginUser; 93*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, 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++) { 1119371c9d4SSatish Balay v = -1.0; 1129371c9d4SSatish Balay Ii = j + n * i; 1139371c9d4SSatish Balay if (i > 0) { 1149371c9d4SSatish Balay J = Ii - n; 1159371c9d4SSatish Balay PetscCall(MatSetValues(mat, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 1169371c9d4SSatish Balay } 1179371c9d4SSatish Balay if (i < m - 1) { 1189371c9d4SSatish Balay J = Ii + n; 1199371c9d4SSatish Balay PetscCall(MatSetValues(mat, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 1209371c9d4SSatish Balay } 1219371c9d4SSatish Balay if (j > 0) { 1229371c9d4SSatish Balay J = Ii - 1; 1239371c9d4SSatish Balay PetscCall(MatSetValues(mat, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 1249371c9d4SSatish Balay } 1259371c9d4SSatish Balay if (j < n - 1) { 1269371c9d4SSatish Balay J = Ii + 1; 1279371c9d4SSatish Balay PetscCall(MatSetValues(mat, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 1289371c9d4SSatish Balay } 1299371c9d4SSatish Balay v = 4.0; 1309371c9d4SSatish Balay PetscCall(MatSetValues(mat, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown } 1339566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 1349566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* check the matrix before and after scaling by -1.0 */ 1379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Matrix _before_ MatScaleUserImpl() operation\n")); 1389566063dSJacob Faibussowitsch PetscCall(MatView(mat, PETSC_VIEWER_STDOUT_WORLD)); 1399566063dSJacob Faibussowitsch PetscCall(MatScaleUserImpl(mat, none)); 1409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Matrix _after_ MatScaleUserImpl() operation\n")); 1419566063dSJacob Faibussowitsch PetscCall(MatView(mat, PETSC_VIEWER_STDOUT_WORLD)); 142c4762a1bSJed Brown 1432e956fe4SStefano Zampini /* deregister user defined MatScaleUser() operation for both SeqAIJ 1442e956fe4SStefano Zampini and MPIAIJ types */ 1452e956fe4SStefano Zampini PetscCall(DeRegisterMatScaleUserImpl(mat)); 1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 1479566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 148b122ec5aSJacob Faibussowitsch return 0; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown /*TEST 152c4762a1bSJed Brown 153c4762a1bSJed Brown test: 154c4762a1bSJed Brown 155c4762a1bSJed Brown test: 156c4762a1bSJed Brown suffix: 2 157c4762a1bSJed Brown nsize: 2 158c4762a1bSJed Brown 159c4762a1bSJed Brown TEST*/ 160