1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatReorderForNonzeroDiagonal().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 7*d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Mat mat, B, C; 9c4762a1bSJed Brown PetscInt i, j; 10c4762a1bSJed Brown PetscMPIInt size; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown IS isrow, iscol, identity; 13c4762a1bSJed Brown PetscViewer viewer; 14c4762a1bSJed Brown 15327415f7SBarry Smith PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* ------- Assemble matrix, --------- */ 19c4762a1bSJed Brown 209566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &mat)); 219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, PETSC_DECIDE, PETSC_DECIDE, 4, 4)); 229566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 239566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* set anti-diagonal of matrix */ 26c4762a1bSJed Brown v = 1.0; 279371c9d4SSatish Balay i = 0; 289371c9d4SSatish Balay j = 3; 299566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &j, &v, INSERT_VALUES)); 30c4762a1bSJed Brown v = 2.0; 319371c9d4SSatish Balay i = 1; 329371c9d4SSatish Balay j = 2; 339566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &j, &v, INSERT_VALUES)); 34c4762a1bSJed Brown v = 3.0; 359371c9d4SSatish Balay i = 2; 369371c9d4SSatish Balay j = 1; 379566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &j, &v, INSERT_VALUES)); 38c4762a1bSJed Brown v = 4.0; 399371c9d4SSatish Balay i = 3; 409371c9d4SSatish Balay j = 0; 419566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &j, &v, INSERT_VALUES)); 429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD, &viewer)); 469566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_DENSE)); 479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Original matrix\n")); 489566063dSJacob Faibussowitsch PetscCall(MatView(mat, viewer)); 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(mat, MATORDERINGNATURAL, &isrow, &iscol)); 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(MatPermute(mat, isrow, iscol, &B)); 539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Original matrix permuted by identity\n")); 549566063dSJacob Faibussowitsch PetscCall(MatView(B, viewer)); 559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(mat, 1.e-8, isrow, iscol)); 589566063dSJacob Faibussowitsch PetscCall(MatPermute(mat, isrow, iscol, &B)); 599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Original matrix permuted by identity + NonzeroDiagonal()\n")); 609566063dSJacob Faibussowitsch PetscCall(MatView(B, viewer)); 619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Row permutation\n")); 629566063dSJacob Faibussowitsch PetscCall(ISView(isrow, viewer)); 639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Column permutation\n")); 649566063dSJacob Faibussowitsch PetscCall(ISView(iscol, viewer)); 659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(mat, MATORDERINGND, &isrow, &iscol)); 719566063dSJacob Faibussowitsch PetscCall(MatPermute(mat, isrow, iscol, &B)); 729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Original matrix permuted by ND\n")); 739566063dSJacob Faibussowitsch PetscCall(MatView(B, viewer)); 749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ND row permutation\n")); 769566063dSJacob Faibussowitsch PetscCall(ISView(isrow, viewer)); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ND column permutation\n")); 789566063dSJacob Faibussowitsch PetscCall(ISView(iscol, viewer)); 79c4762a1bSJed Brown 809566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(mat, 1.e-8, isrow, iscol)); 819566063dSJacob Faibussowitsch PetscCall(MatPermute(mat, isrow, iscol, &B)); 829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Original matrix permuted by ND + NonzeroDiagonal()\n")); 839566063dSJacob Faibussowitsch PetscCall(MatView(B, viewer)); 849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ND + NonzeroDiagonal() row permutation\n")); 869566063dSJacob Faibussowitsch PetscCall(ISView(isrow, viewer)); 879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ND + NonzeroDiagonal() column permutation\n")); 889566063dSJacob Faibussowitsch PetscCall(ISView(iscol, viewer)); 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(mat, MATORDERINGRCM, &isrow, &iscol)); 949566063dSJacob Faibussowitsch PetscCall(MatPermute(mat, isrow, iscol, &B)); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Original matrix permuted by RCM\n")); 969566063dSJacob Faibussowitsch PetscCall(MatView(B, viewer)); 979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "RCM row permutation\n")); 999566063dSJacob Faibussowitsch PetscCall(ISView(isrow, viewer)); 1009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "RCM column permutation\n")); 1019566063dSJacob Faibussowitsch PetscCall(ISView(iscol, viewer)); 102c4762a1bSJed Brown 1039566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(mat, 1.e-8, isrow, iscol)); 1049566063dSJacob Faibussowitsch PetscCall(MatPermute(mat, isrow, iscol, &B)); 1059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Original matrix permuted by RCM + NonzeroDiagonal()\n")); 1069566063dSJacob Faibussowitsch PetscCall(MatView(B, viewer)); 1079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "RCM + NonzeroDiagonal() row permutation\n")); 1089566063dSJacob Faibussowitsch PetscCall(ISView(isrow, viewer)); 1099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "RCM + NonzeroDiagonal() column permutation\n")); 1109566063dSJacob Faibussowitsch PetscCall(ISView(iscol, viewer)); 1119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 112c4762a1bSJed Brown if (size == 1) { 1139566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_SYMMETRIC, PETSC_TRUE)); 1149566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, 4, 0, 1, &identity)); 1159566063dSJacob Faibussowitsch PetscCall(MatPermute(B, identity, identity, &C)); 1169566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSEQSBAIJ, MAT_INPLACE_MATRIX, &C)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1189566063dSJacob Faibussowitsch PetscCall(ISDestroy(&identity)); 119c4762a1bSJed Brown } 1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 121c4762a1bSJed Brown /* Test MatLUFactor(); set diagonal as zeros as requested by PETSc matrix factorization */ 122c4762a1bSJed Brown for (i = 0; i < 4; i++) { 123c4762a1bSJed Brown v = 0.0; 1249566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &i, &v, INSERT_VALUES)); 125c4762a1bSJed Brown } 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 1289566063dSJacob Faibussowitsch PetscCall(MatLUFactor(mat, isrow, iscol, NULL)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Free data structures */ 1319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 1329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 1339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 134c4762a1bSJed Brown 1359566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 136b122ec5aSJacob Faibussowitsch return 0; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /*TEST 140c4762a1bSJed Brown 141c4762a1bSJed Brown test: 142c4762a1bSJed Brown 143c4762a1bSJed Brown TEST*/ 144