xref: /petsc/src/mat/tests/ex68.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatReorderForNonzeroDiagonal().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
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 
15*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   /* ------- Assemble matrix, --------- */
18c4762a1bSJed Brown 
195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&mat));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,4,4));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(mat));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(mat));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* set anti-diagonal of matrix */
25c4762a1bSJed Brown   v    = 1.0;
26c4762a1bSJed Brown   i    = 0; j = 3;
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
28c4762a1bSJed Brown   v    = 2.0;
29c4762a1bSJed Brown   i    = 1; j = 2;
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
31c4762a1bSJed Brown   v    = 3.0;
32c4762a1bSJed Brown   i    = 2; j = 1;
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
34c4762a1bSJed Brown   v    = 4.0;
35c4762a1bSJed Brown   i    = 3; j = 0;
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
39c4762a1bSJed Brown 
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix\n"));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(mat,viewer));
44c4762a1bSJed Brown 
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(mat,MATORDERINGNATURAL,&isrow,&iscol));
46c4762a1bSJed Brown 
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPermute(mat,isrow,iscol,&B));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by identity\n"));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,viewer));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
51c4762a1bSJed Brown 
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPermute(mat,isrow,iscol,&B));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by identity + NonzeroDiagonal()\n"));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,viewer));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Row permutation\n"));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(isrow,viewer));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Column permutation\n"));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(iscol,viewer));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscol));
64c4762a1bSJed Brown 
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(mat,MATORDERINGND,&isrow,&iscol));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPermute(mat,isrow,iscol,&B));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by ND\n"));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,viewer));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"ND row permutation\n"));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(isrow,viewer));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"ND column permutation\n"));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(iscol,viewer));
74c4762a1bSJed Brown 
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPermute(mat,isrow,iscol,&B));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by ND + NonzeroDiagonal()\n"));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,viewer));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"ND + NonzeroDiagonal() row permutation\n"));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(isrow,viewer));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"ND + NonzeroDiagonal() column permutation\n"));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(iscol,viewer));
84c4762a1bSJed Brown 
855f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscol));
87c4762a1bSJed Brown 
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(mat,MATORDERINGRCM,&isrow,&iscol));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPermute(mat,isrow,iscol,&B));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by RCM\n"));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,viewer));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"RCM row permutation\n"));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(isrow,viewer));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"RCM column permutation\n"));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(iscol,viewer));
97c4762a1bSJed Brown 
985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPermute(mat,isrow,iscol,&B));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by RCM + NonzeroDiagonal()\n"));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,viewer));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"RCM + NonzeroDiagonal() row permutation\n"));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(isrow,viewer));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"RCM + NonzeroDiagonal() column permutation\n"));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(iscol,viewer));
1065f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
107c4762a1bSJed Brown   if (size == 1) {
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE));
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,4,0,1,&identity));
1105f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPermute(B,identity,identity,&C));
1115f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(C,MATSEQSBAIJ,MAT_INPLACE_MATRIX,&C));
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&identity));
114c4762a1bSJed Brown   }
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
116c4762a1bSJed Brown   /* Test MatLUFactor(); set diagonal as zeros as requested by PETSc matrix factorization */
117c4762a1bSJed Brown   for (i=0; i<4; i++) {
118c4762a1bSJed Brown     v = 0.0;
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(mat,1,&i,1,&i,&v,INSERT_VALUES));
120c4762a1bSJed Brown   }
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactor(mat,isrow,iscol,NULL));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Free data structures */
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscol));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat));
129c4762a1bSJed Brown 
130*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
131*b122ec5aSJacob Faibussowitsch   return 0;
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /*TEST
135c4762a1bSJed Brown 
136c4762a1bSJed Brown    test:
137c4762a1bSJed Brown 
138c4762a1bSJed Brown TEST*/
139