Actual source code: ex59.c
1: /*$Id: ex59.c,v 1.20 2001/08/07 03:03:07 balay Exp $*/
3: static char help[] = "Tests MatGetSubmatrix() in parallel.";
5: #include petscmat.h
7: int main(int argc,char **args)
8: {
9: Mat C,A;
10: int i,j,m = 3,n = 2,rank,size,ierr,rstart,rend;
11: PetscScalar v;
12: IS isrow,iscol;
13: PetscTruth flg;
14: char type[256];
16: PetscInitialize(&argc,&args,(char *)0,help);
17: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
18: MPI_Comm_size(PETSC_COMM_WORLD,&size);
19: n = 2*size;
21: PetscStrcpy(type,MATSAME);
22: PetscOptionsGetString(PETSC_NULL,"-mat_type",type,256,PETSC_NULL);
24: PetscStrcmp(type,MATMPIDENSE,&flg);
25: if (flg) {
26: MatCreateMPIDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
27: m*n,m*n,PETSC_NULL,&C);
28: } else {
29: MatCreateMPIAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
30: m*n,m*n,5,PETSC_NULL,5,PETSC_NULL,&C);
31: }
33: /*
34: This is JUST to generate a nice test matrix, all processors fill up
35: the entire matrix. This is not something one would ever do in practice.
36: */
37: for (i=0; i<m*n; i++) {
38: for (j=0; j<m*n; j++) {
39: v = i + j + 1;
40: MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);
41: }
42: }
43: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
44: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
45: PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);
46: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
48: /*
49: Generate a new matrix consisting of every second row and column of
50: the original matrix
51: */
52: MatGetOwnershipRange(C,&rstart,&rend);
53: /* list the rows we want on THIS processor */
54: ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&isrow);
55: /* list ALL the columns we want */
56: ISCreateStride(PETSC_COMM_WORLD,(m*n)/2,0,2,&iscol);
57: MatGetSubMatrix(C,isrow,iscol,PETSC_DECIDE,MAT_INITIAL_MATRIX,&A);
58: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
60: MatGetSubMatrix(C,isrow,iscol,PETSC_DECIDE,MAT_REUSE_MATRIX,&A);
61: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
63: ISDestroy(isrow);
64: ISDestroy(iscol);
65: MatDestroy(A);
66: MatDestroy(C);
67: PetscFinalize();
68: return 0;
69: }