Actual source code: ex9.c
2: static char help[] = "Tests MPI parallel matrix creation.\n\n";
4: #include petscmat.h
8: int main(int argc,char **args)
9: {
10: Mat C;
11: MatInfo info;
12: PetscMPIInt rank,size;
13: PetscInt i,j,m = 3,n = 2,low,high,iglobal;
14: PetscInt Ii,J,ldim;
16: PetscTruth flg;
17: PetscScalar v,one = 1.0;
18: Vec u,b;
19: PetscInt bs,ndiag,diag[7]; bs = 1,ndiag = 5;
21: PetscInitialize(&argc,&args,(char *)0,help);
22: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
23: MPI_Comm_size(PETSC_COMM_WORLD,&size);
24: n = 2*size;
26: MatCreate(PETSC_COMM_WORLD,&C);
27: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
28: MatSetFromOptions(C);
30: diag[0] = n;
31: diag[1] = 1;
32: diag[2] = 0;
33: diag[3] = -1;
34: diag[4] = -n;
35: if (size>1) {ndiag = 7; diag[5] = 2; diag[6] = -2;}
36: MatMPIBDiagSetPreallocation(C,ndiag,bs,diag,PETSC_NULL);
37: MatMPIAIJSetPreallocation(C,5,PETSC_NULL,5,PETSC_NULL);
39: /* Create the matrix for the five point stencil, YET AGAIN */
40: for (i=0; i<m; i++) {
41: for (j=2*rank; j<2*rank+2; j++) {
42: v = -1.0; Ii = j + n*i;
43: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
44: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
45: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
46: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
47: v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
48: }
49: }
51: /* Add extra elements (to illustrate variants of MatGetInfo) */
52: Ii = n; J = n-2; v = 100.0;
53: MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
54: Ii = n-2; J = n; v = 100.0;
55: MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
57: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
58: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
60: /* Form vectors */
61: VecCreate(PETSC_COMM_WORLD,&u);
62: VecSetSizes(u,PETSC_DECIDE,m*n);
63: VecSetFromOptions(u);
64: VecDuplicate(u,&b);
65: VecGetLocalSize(u,&ldim);
66: VecGetOwnershipRange(u,&low,&high);
67: for (i=0; i<ldim; i++) {
68: iglobal = i + low;
69: v = one*((PetscReal)i) + 100.0*rank;
70: VecSetValues(u,1,&iglobal,&v,INSERT_VALUES);
71: }
72: VecAssemblyBegin(u);
73: VecAssemblyEnd(u);
75: MatMult(C,u,b);
77: VecView(u,PETSC_VIEWER_STDOUT_WORLD);
78: VecView(b,PETSC_VIEWER_STDOUT_WORLD);
80: VecDestroy(u);
81: VecDestroy(b);
83: PetscOptionsHasName(PETSC_NULL,"-view_info",&flg);
84: if (flg) {PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);}
85: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
87: MatGetInfo(C,MAT_GLOBAL_SUM,&info);
88: PetscPrintf(PETSC_COMM_WORLD,"matrix information (global sums):\n\
89: nonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
90: MatGetInfo (C,MAT_GLOBAL_MAX,&info);
91: PetscPrintf(PETSC_COMM_WORLD,"matrix information (global max):\n\
92: nonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
94: MatDestroy(C);
95: PetscFinalize();
96: return 0;
97: }