Actual source code: ex9.c
1: /*$Id: ex9.c,v 1.25 2001/09/11 16:32:50 bsmith Exp $*/
3: static char help[] = "Tests MPI parallel matrix creation.nn";
5: #include petscmat.h
7: int main(int argc,char **args)
8: {
9: Mat C;
10: MatInfo info;
11: int i,j,m = 3,n = 2,rank,size,low,high,iglobal;
12: int I,J,ierr,ldim;
13: PetscTruth flg;
14: PetscScalar v,one = 1.0;
15: Vec u,b;
16: int bs,ndiag,diag[7]; bs = 1,ndiag = 5;
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: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
24: MatSetFromOptions(C);
26: diag[0] = n;
27: diag[1] = 1;
28: diag[2] = 0;
29: diag[3] = -1;
30: diag[4] = -n;
31: if (size>1) {ndiag = 7; diag[5] = 2; diag[6] = -2;}
32: MatMPIBDiagSetPreallocation(C,ndiag,bs,diag,PETSC_NULL);
33: MatMPIAIJSetPreallocation(C,5,PETSC_NULL,5,PETSC_NULL);
35: /* Create the matrix for the five point stencil, YET AGAIN */
36: for (i=0; i<m; i++) {
37: for (j=2*rank; j<2*rank+2; j++) {
38: v = -1.0; I = j + n*i;
39: if (i>0) {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
40: if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
41: if (j>0) {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
42: if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
43: v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
44: }
45: }
47: /* Add extra elements (to illustrate variants of MatGetInfo) */
48: I = n; J = n-2; v = 100.0;
49: MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);
50: I = n-2; J = n; v = 100.0;
51: MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);
53: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
56: /* Form vectors */
57: VecCreate(PETSC_COMM_WORLD,&u);
58: VecSetSizes(u,PETSC_DECIDE,m*n);
59: VecSetFromOptions(u);
60: VecDuplicate(u,&b);
61: VecGetLocalSize(u,&ldim);
62: VecGetOwnershipRange(u,&low,&high);
63: for (i=0; i<ldim; i++) {
64: iglobal = i + low;
65: v = one*((PetscReal)i) + 100.0*rank;
66: VecSetValues(u,1,&iglobal,&v,INSERT_VALUES);
67: }
68: VecAssemblyBegin(u);
69: VecAssemblyEnd(u);
71: MatMult(C,u,b);
73: VecView(u,PETSC_VIEWER_STDOUT_WORLD);
74: VecView(b,PETSC_VIEWER_STDOUT_WORLD);
76: VecDestroy(u);
77: VecDestroy(b);
79: PetscOptionsHasName(PETSC_NULL,"-view_info",&flg);
80: if (flg) {PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);}
81: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
83: MatGetInfo(C,MAT_GLOBAL_SUM,&info);
84: PetscPrintf(PETSC_COMM_WORLD,"matrix information (global sums):n
85: nonzeros = %d, allocated nonzeros = %dn",(int)info.nz_used,(int)info.nz_allocated);
86: MatGetInfo (C,MAT_GLOBAL_MAX,&info);
87: PetscPrintf(PETSC_COMM_WORLD,"matrix information (global max):n
88: nonzeros = %d, allocated nonzeros = %dn",(int)info.nz_used,(int)info.nz_allocated);
90: MatDestroy(C);
91: PetscFinalize();
92: return 0;
93: }