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
9: int main(int argc,char **args)
10: {
11: Mat C,A;
12: int i,j,m = 3,n = 2,rank,size,ierr,rstart,rend;
13: PetscScalar v;
14: IS isrow,iscol;
15: PetscTruth flg;
16: char type[256];
18: PetscInitialize(&argc,&args,(char *)0,help);
19: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
20: MPI_Comm_size(PETSC_COMM_WORLD,&size);
21: n = 2*size;
23: PetscStrcpy(type,MATSAME);
24: PetscOptionsGetString(PETSC_NULL,"-mat_type",type,256,PETSC_NULL);
26: PetscStrcmp(type,MATMPIDENSE,&flg);
27: if (flg) {
28: MatCreateMPIDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
29: m*n,m*n,PETSC_NULL,&C);
30: } else {
31: MatCreateMPIAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
32: m*n,m*n,5,PETSC_NULL,5,PETSC_NULL,&C);
33: }
35: /*
36: This is JUST to generate a nice test matrix, all processors fill up
37: the entire matrix. This is not something one would ever do in practice.
38: */
39: for (i=0; i<m*n; i++) {
40: for (j=0; j<m*n; j++) {
41: v = i + j + 1;
42: MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);
43: }
44: }
45: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
46: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
47: PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);
48: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
50: /*
51: Generate a new matrix consisting of every second row and column of
52: the original matrix
53: */
54: MatGetOwnershipRange(C,&rstart,&rend);
55: /* list the rows we want on THIS processor */
56: ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&isrow);
57: /* list ALL the columns we want */
58: ISCreateStride(PETSC_COMM_WORLD,(m*n)/2,0,2,&iscol);
59: MatGetSubMatrix(C,isrow,iscol,PETSC_DECIDE,MAT_INITIAL_MATRIX,&A);
60: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
62: MatGetSubMatrix(C,isrow,iscol,PETSC_DECIDE,MAT_REUSE_MATRIX,&A);
63: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
65: ISDestroy(isrow);
66: ISDestroy(iscol);
67: MatDestroy(A);
68: MatDestroy(C);
69: PetscFinalize();
70: return 0;
71: }