Actual source code: ex59.c
petsc-3.4.5 2014-06-29
2: static char help[] = "Tests MatGetSubmatrix() in parallel.";
4: #include <petscmat.h>
8: int main(int argc,char **args)
9: {
10: Mat C,A;
11: PetscInt i,j,m = 3,n = 2,rstart,rend;
12: PetscMPIInt size,rank;
14: PetscScalar v;
15: IS isrow,iscol;
16: PetscBool flg;
17: char type[256];
19: PetscInitialize(&argc,&args,(char*)0,help);
20: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
21: MPI_Comm_size(PETSC_COMM_WORLD,&size);
22: n = 2*size;
24: PetscStrcpy(type,MATSAME);
25: PetscOptionsGetString(NULL,"-mat_type",type,256,NULL);
27: PetscStrcmp(type,MATMPIDENSE,&flg);
28: if (flg) {
29: MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE, m*n,m*n,NULL,&C);
30: } else {
31: MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE, m*n,m*n,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&C);
32: }
34: /*
35: This is JUST to generate a nice test matrix, all processors fill up
36: the entire matrix. This is not something one would ever do in practice.
37: */
38: for (i=0; i<m*n; i++) {
39: for (j=0; j<m*n; j++) {
40: v = i + j + 1;
41: MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);
42: }
43: }
44: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
45: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
46: PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);
47: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
49: /*
50: Generate a new matrix consisting of every second row and column of
51: the original matrix
52: */
53: MatGetOwnershipRange(C,&rstart,&rend);
54: /* Create parallel IS with the rows we want on THIS processor */
55: ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&isrow);
56: /* Create parallel IS with the rows we want on THIS processor (same as rows for now) */
57: ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&iscol);
58: MatGetSubMatrix(C,isrow,iscol,MAT_INITIAL_MATRIX,&A);
59: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
61: MatGetSubMatrix(C,isrow,iscol,MAT_REUSE_MATRIX,&A);
62: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
64: ISDestroy(&isrow);
65: ISDestroy(&iscol);
66: MatDestroy(&A);
67: MatDestroy(&C);
68: PetscFinalize();
69: return 0;
70: }