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