xref: /petsc/src/mat/tests/ex181.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests MatCreateSubmatrix() with entire matrix, modified from ex59.c.";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            C,A,Adup;
9   PetscInt       i,j,m = 3,n = 2,rstart,rend;
10   PetscMPIInt    size,rank;
11   PetscErrorCode ierr;
12   PetscScalar    v;
13   IS             isrow;
14   PetscBool      detect_bug = PETSC_FALSE;
15 
16   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
17   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-detect_bug",&detect_bug));
18   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
19   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20   n    = 2*size;
21 
22   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
23   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
24   CHKERRQ(MatSetFromOptions(C));
25   CHKERRQ(MatSetUp(C));
26 
27   /*
28         This is JUST to generate a nice test matrix, all processors fill up
29     the entire matrix. This is not something one would ever do in practice.
30   */
31   CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend));
32   for (i=rstart; i<rend; i++) {
33     for (j=0; j<m*n; j++) {
34       v    = i + j + 1;
35       CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
36     }
37   }
38   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
39   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
40   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON));
41   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
42   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
43 
44   /*
45      Generate a new matrix consisting every row and column of the original matrix
46   */
47   CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend));
48 
49   /* Create parallel IS with the rows we want on THIS processor */
50   if (detect_bug && rank == 0) {
51     CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,1,rstart,1,&isrow));
52   } else {
53     CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,rend-rstart,rstart,1,&isrow));
54   }
55   CHKERRQ(MatCreateSubMatrix(C,isrow,NULL,MAT_INITIAL_MATRIX,&A));
56 
57   /* Change C to test the case MAT_REUSE_MATRIX */
58   if (rank == 0) {
59     i = 0; j = 0; v = 100;
60     CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
61   }
62   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
63   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
64 
65   CHKERRQ(MatCreateSubMatrix(C,isrow,NULL,MAT_REUSE_MATRIX,&A));
66   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON));
67   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
68   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
69 
70   /* Test MatDuplicate */
71   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&Adup));
72   CHKERRQ(MatDestroy(&Adup));
73 
74   CHKERRQ(ISDestroy(&isrow));
75   CHKERRQ(MatDestroy(&A));
76   CHKERRQ(MatDestroy(&C));
77   ierr = PetscFinalize();
78   return ierr;
79 }
80 
81 /*TEST
82 
83    test:
84       nsize: 2
85       filter: grep -v "Mat Object"
86       requires: !complex
87 
88    test:
89       suffix: 2
90       nsize: 3
91       args: -detect_bug
92       filter: grep -v "Mat Object"
93       requires: !complex
94 
95 TEST*/
96