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