1c4762a1bSJed Brown static char help[] = "Tests MatPtAP() for MPIMAIJ and MPIAIJ \n "; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmda.h> 4c4762a1bSJed Brown 5*9371c9d4SSatish Balay int main(int argc, char **argv) { 6c4762a1bSJed Brown DM coarsedm, finedm; 7c4762a1bSJed Brown PetscMPIInt size, rank; 8c4762a1bSJed Brown PetscInt M, N, Z, i, nrows; 9c4762a1bSJed Brown PetscScalar one = 1.0; 10c4762a1bSJed Brown PetscReal fill = 2.0; 11c4762a1bSJed Brown Mat A, P, C; 12c20d7725SJed Brown PetscScalar *array, alpha; 13c4762a1bSJed Brown PetscBool Test_3D = PETSC_FALSE, flg; 14c4762a1bSJed Brown const PetscInt *ia, *ja; 15c4762a1bSJed Brown PetscInt dof; 16c4762a1bSJed Brown MPI_Comm comm; 17c4762a1bSJed Brown 18327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 20c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 23*9371c9d4SSatish Balay M = 10; 24*9371c9d4SSatish Balay N = 10; 25*9371c9d4SSatish Balay Z = 10; 26c4762a1bSJed Brown dof = 10; 27c4762a1bSJed Brown 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_3D", &Test_3D, NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-Z", &Z, NULL)); 32c4762a1bSJed Brown /* Set up distributed array for fine grid */ 33c4762a1bSJed Brown if (!Test_3D) { 349566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, 1, NULL, NULL, &coarsedm)); 35c4762a1bSJed Brown } else { 369566063dSJacob Faibussowitsch PetscCall(DMDACreate3d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, Z, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof, 1, NULL, NULL, NULL, &coarsedm)); 37c4762a1bSJed Brown } 389566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(coarsedm)); 399566063dSJacob Faibussowitsch PetscCall(DMSetUp(coarsedm)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* This makes sure the coarse DMDA has the same partition as the fine DMDA */ 429566063dSJacob Faibussowitsch PetscCall(DMRefine(coarsedm, PetscObjectComm((PetscObject)coarsedm), &finedm)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /*------------------------------------------------------------*/ 459566063dSJacob Faibussowitsch PetscCall(DMSetMatType(finedm, MATAIJ)); 469566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(finedm, &A)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* set val=one to A */ 49c4762a1bSJed Brown if (size == 1) { 509566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 51c4762a1bSJed Brown if (flg) { 529566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(A, &array)); 53c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one; 549566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(A, &array)); 55c4762a1bSJed Brown } 569566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 57c4762a1bSJed Brown } else { 58c4762a1bSJed Brown Mat AA, AB; 599566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(A, &AA, &AB, NULL)); 609566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 61c4762a1bSJed Brown if (flg) { 629566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(AA, &array)); 63c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one; 649566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(AA, &array)); 65c4762a1bSJed Brown } 669566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 679566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 68c4762a1bSJed Brown if (flg) { 699566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(AB, &array)); 70c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one; 719566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(AB, &array)); 72c4762a1bSJed Brown } 739566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 769566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(coarsedm, finedm, &P, NULL)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Test P^T * A * P - MatPtAP() */ 79c4762a1bSJed Brown /*------------------------------*/ 80c20d7725SJed Brown /* (1) Developer API */ 819566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A, P, NULL, &C)); 829566063dSJacob Faibussowitsch PetscCall(MatProductSetType(C, MATPRODUCT_PtAP)); 839566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(C, "allatonce")); 849566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(C, PETSC_DEFAULT)); 859566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(C)); 869566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(C)); 879566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(C)); 889566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(C)); /* Test reuse of symbolic C */ 89c20d7725SJed Brown 90c2f6fc52SHong Zhang { /* Test MatProductView() */ 91c2f6fc52SHong Zhang PetscViewer viewer; 929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(comm, NULL, &viewer)); 939566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 949566063dSJacob Faibussowitsch PetscCall(MatProductView(C, viewer)); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 969566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 97c2f6fc52SHong Zhang } 98c2f6fc52SHong Zhang 999566063dSJacob Faibussowitsch PetscCall(MatPtAPMultEqual(A, P, C, 10, &flg)); 10028b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatProduct_PtAP"); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 102c20d7725SJed Brown 103c20d7725SJed Brown /* (2) User API */ 1049566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, fill, &C)); 105c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 106c4762a1bSJed Brown alpha = 1.0; 107c4762a1bSJed Brown for (i = 0; i < 1; i++) { 108c4762a1bSJed Brown alpha -= 0.1; 1099566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha)); 1109566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, fill, &C)); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */ 1149566063dSJacob Faibussowitsch PetscCall(MatProductClear(C)); 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(MatPtAPMultEqual(A, P, C, 10, &flg)); 11728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatPtAP"); 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 1229566063dSJacob Faibussowitsch PetscCall(DMDestroy(&finedm)); 1239566063dSJacob Faibussowitsch PetscCall(DMDestroy(&coarsedm)); 1249566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 125b122ec5aSJacob Faibussowitsch return 0; 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /*TEST 129c4762a1bSJed Brown 130c4762a1bSJed Brown test: 131c4762a1bSJed Brown args: -M 10 -N 10 -Z 10 132c4762a1bSJed Brown output_file: output/ex89_1.out 133c4762a1bSJed Brown 134c4762a1bSJed Brown test: 135c4762a1bSJed Brown suffix: allatonce 136c4762a1bSJed Brown nsize: 4 137c2f6fc52SHong Zhang args: -M 10 -N 10 -Z 10 138c2f6fc52SHong Zhang output_file: output/ex89_2.out 139c4762a1bSJed Brown 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: allatonce_merged 142c4762a1bSJed Brown nsize: 4 1433e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -mat_product_algorithm allatonce_merged 144c2f6fc52SHong Zhang output_file: output/ex89_3.out 145c4762a1bSJed Brown 146c4762a1bSJed Brown test: 147c2f6fc52SHong Zhang suffix: nonscalable_3D 148c4762a1bSJed Brown nsize: 4 1493e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm nonscalable 150c2f6fc52SHong Zhang output_file: output/ex89_4.out 151c4762a1bSJed Brown 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: allatonce_merged_3D 154c4762a1bSJed Brown nsize: 4 1553e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm allatonce_merged 156c2f6fc52SHong Zhang output_file: output/ex89_3.out 157c2f6fc52SHong Zhang 158c2f6fc52SHong Zhang test: 159c2f6fc52SHong Zhang suffix: nonscalable 160c2f6fc52SHong Zhang nsize: 4 1613e662e0bSHong Zhang args: -M 10 -N 10 -Z 10 -mat_product_algorithm nonscalable 162c2f6fc52SHong Zhang output_file: output/ex89_5.out 163c4762a1bSJed Brown 164c4762a1bSJed Brown TEST*/ 165