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 I,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,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
27: MatSetFromOptions(C);
29: diag[0] = n;
30: diag[1] = 1;
31: diag[2] = 0;
32: diag[3] = -1;
33: diag[4] = -n;
34: if (size>1) {ndiag = 7; diag[5] = 2; diag[6] = -2;}
35: MatMPIBDiagSetPreallocation(C,ndiag,bs,diag,PETSC_NULL);
36: MatMPIAIJSetPreallocation(C,5,PETSC_NULL,5,PETSC_NULL);
38: /* Create the matrix for the five point stencil, YET AGAIN */
39: for (i=0; i<m; i++) {
40: for (j=2*rank; j<2*rank+2; j++) {
41: v = -1.0; I = j + n*i;
42: if (i>0) {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
43: if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
44: if (j>0) {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
45: if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
46: v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
47: }
48: }
50: /* Add extra elements (to illustrate variants of MatGetInfo) */
51: I = n; J = n-2; v = 100.0;
52: MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);
53: I = n-2; J = n; v = 100.0;
54: MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);
56: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
59: /* Form vectors */
60: VecCreate(PETSC_COMM_WORLD,&u);
61: VecSetSizes(u,PETSC_DECIDE,m*n);
62: VecSetFromOptions(u);
63: VecDuplicate(u,&b);
64: VecGetLocalSize(u,&ldim);
65: VecGetOwnershipRange(u,&low,&high);
66: for (i=0; i<ldim; i++) {
67: iglobal = i + low;
68: v = one*((PetscReal)i) + 100.0*rank;
69: VecSetValues(u,1,&iglobal,&v,INSERT_VALUES);
70: }
71: VecAssemblyBegin(u);
72: VecAssemblyEnd(u);
74: MatMult(C,u,b);
76: VecView(u,PETSC_VIEWER_STDOUT_WORLD);
77: VecView(b,PETSC_VIEWER_STDOUT_WORLD);
79: VecDestroy(u);
80: VecDestroy(b);
82: PetscOptionsHasName(PETSC_NULL,"-view_info",&flg);
83: if (flg) {PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);}
84: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
86: MatGetInfo(C,MAT_GLOBAL_SUM,&info);
87: PetscPrintf(PETSC_COMM_WORLD,"matrix information (global sums):\n\
88: nonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
89: MatGetInfo (C,MAT_GLOBAL_MAX,&info);
90: PetscPrintf(PETSC_COMM_WORLD,"matrix information (global max):\n\
91: nonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
93: MatDestroy(C);
94: PetscFinalize();
95: return 0;
96: }