Actual source code: ex49.c
2: static char help[] = "Tests MatTranspose(), MatNorm(), MatValid(), and MatAXPY().\n\n";
4: #include petscmat.h
8: int main(int argc,char **argv)
9: {
10: Mat mat,tmat = 0;
11: PetscInt m = 4,n,i,j;
13: PetscMPIInt size,rank;
14: PetscInt rstart,rend,rect = 0,nd,bs,*diag,*bdlen;
15: PetscTruth flg,isbdiag;
16: PetscScalar v,**diagv;
17: PetscReal normf,normi,norm1;
18: MatInfo info;
19:
20: PetscInitialize(&argc,&argv,(char*)0,help);
21: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
22: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
23: MPI_Comm_size(PETSC_COMM_WORLD,&size);
24: n = m;
25: PetscOptionsHasName(PETSC_NULL,"-rect1",&flg);
26: if (flg) {n += 2; rect = 1;}
27: PetscOptionsHasName(PETSC_NULL,"-rect2",&flg);
28: if (flg) {n -= 2; rect = 1;}
30: /* Create and assemble matrix */
31: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,&mat);
32: MatSetFromOptions(mat);
33: MatGetOwnershipRange(mat,&rstart,&rend);
34: for (i=rstart; i<rend; i++) {
35: for (j=0; j<n; j++) {
36: v=10*i+j;
37: MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES);
38: }
39: }
40: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
41: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
43: /* Test whether matrix has been corrupted (just to demonstrate this
44: routine) not needed in most application codes. */
45: MatValid(mat,(PetscTruth*)&flg);
46: if (!flg) SETERRQ(1,"Corrupted matrix.");
48: /* Print info about original matrix */
49: MatGetInfo(mat,MAT_GLOBAL_SUM,&info);
50: PetscPrintf(PETSC_COMM_WORLD,"original matrix nonzeros = %D, allocated nonzeros = %D\n",
51: (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
52: MatNorm(mat,NORM_FROBENIUS,&normf);
53: MatNorm(mat,NORM_1,&norm1);
54: MatNorm(mat,NORM_INFINITY,&normi);
55: PetscPrintf(PETSC_COMM_WORLD,"original: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",
56: normf,norm1,normi);
57: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
59: PetscTypeCompare((PetscObject)mat,MATSEQBDIAG,&isbdiag);
60: if (!isbdiag) {
61: PetscTypeCompare((PetscObject)mat,MATMPIBDIAG,&isbdiag);
62: }
63: if (isbdiag) {
64: MatBDiagGetData(mat,&nd,&bs,&diag,&bdlen,&diagv);
65: PetscPrintf(PETSC_COMM_WORLD,"original matrix: block diag format: %D diagonals, block size = %D\n",nd,bs);
66: for (i=0; i<nd; i++) {
67: PetscPrintf(PETSC_COMM_WORLD," diag=%D, bdlength=%D\n",diag[i],bdlen[i]);
68: }
69: }
71: /* Form matrix transpose */
72: PetscOptionsHasName(PETSC_NULL,"-in_place",&flg);
73: if (!rect && flg) {
74: MatTranspose(mat,0); /* in-place transpose */
75: tmat = mat; mat = 0;
76: } else { /* out-of-place transpose */
77: MatTranspose(mat,&tmat);
78: }
80: /* Print info about transpose matrix */
81: MatGetInfo(tmat,MAT_GLOBAL_SUM,&info);
82: PetscPrintf(PETSC_COMM_WORLD,"transpose matrix nonzeros = %D, allocated nonzeros = %D\n",
83: (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
84: MatNorm(tmat,NORM_FROBENIUS,&normf);
85: MatNorm(tmat,NORM_1,&norm1);
86: MatNorm(tmat,NORM_INFINITY,&normi);
87: PetscPrintf(PETSC_COMM_WORLD,"transpose: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",
88: normf,norm1,normi);
89: MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
91: if (isbdiag) {
92: MatBDiagGetData(tmat,&nd,&bs,&diag,&bdlen,&diagv);
93: PetscPrintf(PETSC_COMM_WORLD,"transpose matrix: block diag format: %D diagonals, block size = %D\n",nd,bs);
94: for (i=0; i<nd; i++) {
95: PetscPrintf(PETSC_COMM_WORLD," diag=%D, bdlength=%D\n",diag[i],bdlen[i]);
96: }
97: }
99: /* Test MatAXPY */
100: if (mat && !rect) {
101: PetscScalar alpha = 1.0;
102: PetscOptionsGetScalar(PETSC_NULL,"-alpha",&alpha,PETSC_NULL);
103: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A\n");
104: MatAXPY(&alpha,mat,tmat,DIFFERENT_NONZERO_PATTERN);
105: MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
106: }
108: /* Free data structures */
109: MatDestroy(tmat);
110: if (mat) {MatDestroy(mat);}
112: PetscFinalize();
113: return 0;
114: }
115: