xref: /petsc/src/mat/tests/ex194.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() with certain entire rows of matrix, modified from ex181.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;
9c4762a1bSJed Brown   PetscInt       i,j,m = 3,n = 2,rstart,rend,cstart,cend;
10c4762a1bSJed Brown   PetscMPIInt    size,rank;
11c4762a1bSJed Brown   PetscScalar    v;
12c4762a1bSJed Brown   IS             isrow,iscol;
13c4762a1bSJed Brown 
14*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
15*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
17c4762a1bSJed Brown   n    = 2*size;
18c4762a1bSJed Brown 
19*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
20*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
21*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
22*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /*
25c4762a1bSJed Brown         This is JUST to generate a nice test matrix, all processors fill up
26c4762a1bSJed Brown     the entire matrix. This is not something one would ever do in practice.
27c4762a1bSJed Brown   */
28*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C,&rstart,&rend));
29c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
30c4762a1bSJed Brown     for (j=0; j<m*n; j++) {
31c4762a1bSJed Brown       v    = i + j + 1;
32*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
33c4762a1bSJed Brown     }
34c4762a1bSJed Brown   }
35*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
36*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
37*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON));
38*9566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
39*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /*
42c4762a1bSJed Brown      Generate a new matrix consisting every column of the original matrix
43c4762a1bSJed Brown   */
44*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C,&rstart,&rend));
45*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(C,&cstart,&cend));
46c4762a1bSJed Brown 
47*9566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_WORLD,rend-rstart > 0 ? rend-rstart-1 : 0,rstart,1,&isrow));
48*9566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_WORLD,cend-cstart,cstart,1,&iscol));
49*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(C,isrow,NULL,MAT_INITIAL_MATRIX,&A));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Change C to test the case MAT_REUSE_MATRIX */
52dd400576SPatrick Sanan   if (rank == 0) {
53c4762a1bSJed Brown     i = 0; j = 0; v = 100;
54*9566063dSJacob Faibussowitsch     PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
55c4762a1bSJed Brown   }
56*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
57*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
58c4762a1bSJed Brown 
59*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(C,isrow,NULL,MAT_REUSE_MATRIX,&A));
60*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON));
61*9566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
62*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
63c4762a1bSJed Brown 
64*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
65*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
66*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
67*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
68*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
69b122ec5aSJacob Faibussowitsch   return 0;
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown /*TEST
73c4762a1bSJed Brown 
74c4762a1bSJed Brown    test:
75c4762a1bSJed Brown       nsize: 2
76c4762a1bSJed Brown 
77c4762a1bSJed Brown TEST*/
78