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,&mat);
32: MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n);
33: MatSetFromOptions(mat);
34: MatGetOwnershipRange(mat,&rstart,&rend);
35: for (i=rstart; i<rend; i++) {
36: for (j=0; j<n; j++) {
37: v=10*i+j;
38: MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES);
39: }
40: }
41: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
42: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
44: /* Test whether matrix has been corrupted (just to demonstrate this
45: routine) not needed in most application codes. */
46: MatValid(mat,(PetscTruth*)&flg);
47: if (!flg) SETERRQ(1,"Corrupted matrix.");
49: /* Print info about original matrix */
50: MatGetInfo(mat,MAT_GLOBAL_SUM,&info);
51: PetscPrintf(PETSC_COMM_WORLD,"original matrix nonzeros = %D, allocated nonzeros = %D\n",
52: (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
53: MatNorm(mat,NORM_FROBENIUS,&normf);
54: MatNorm(mat,NORM_1,&norm1);
55: MatNorm(mat,NORM_INFINITY,&normi);
56: PetscPrintf(PETSC_COMM_WORLD,"original: Frobenious norm = %G, one norm = %G, infinity norm = %G\n",
57: normf,norm1,normi);
58: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
60: PetscTypeCompare((PetscObject)mat,MATSEQBDIAG,&isbdiag);
61: if (!isbdiag) {
62: PetscTypeCompare((PetscObject)mat,MATMPIBDIAG,&isbdiag);
63: }
64: if (isbdiag) {
65: MatBDiagGetData(mat,&nd,&bs,&diag,&bdlen,&diagv);
66: PetscPrintf(PETSC_COMM_WORLD,"original matrix: block diag format: %D diagonals, block size = %D\n",nd,bs);
67: for (i=0; i<nd; i++) {
68: PetscPrintf(PETSC_COMM_WORLD," diag=%D, bdlength=%D\n",diag[i],bdlen[i]);
69: }
70: }
72: /* Form matrix transpose */
73: PetscOptionsHasName(PETSC_NULL,"-in_place",&flg);
74: if (!rect && flg) {
75: MatTranspose(mat,0); /* in-place transpose */
76: tmat = mat; mat = 0;
77: } else { /* out-of-place transpose */
78: MatTranspose(mat,&tmat);
79: }
81: /* Print info about transpose matrix */
82: MatGetInfo(tmat,MAT_GLOBAL_SUM,&info);
83: PetscPrintf(PETSC_COMM_WORLD,"transpose matrix nonzeros = %D, allocated nonzeros = %D\n",
84: (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
85: MatNorm(tmat,NORM_FROBENIUS,&normf);
86: MatNorm(tmat,NORM_1,&norm1);
87: MatNorm(tmat,NORM_INFINITY,&normi);
88: PetscPrintf(PETSC_COMM_WORLD,"transpose: Frobenious norm = %G, one norm = %G, infinity norm = %G\n",
89: normf,norm1,normi);
90: MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
92: if (isbdiag) {
93: MatBDiagGetData(tmat,&nd,&bs,&diag,&bdlen,&diagv);
94: PetscPrintf(PETSC_COMM_WORLD,"transpose matrix: block diag format: %D diagonals, block size = %D\n",nd,bs);
95: for (i=0; i<nd; i++) {
96: PetscPrintf(PETSC_COMM_WORLD," diag=%D, bdlength=%D\n",diag[i],bdlen[i]);
97: }
98: }
100: /* Test MatAXPY */
101: if (mat && !rect) {
102: PetscScalar alpha = 1.0;
103: PetscOptionsGetScalar(PETSC_NULL,"-alpha",&alpha,PETSC_NULL);
104: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A\n");
105: MatAXPY(tmat,alpha,mat,DIFFERENT_NONZERO_PATTERN);
106: MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
107: }
109: /* Free data structures */
110: MatDestroy(tmat);
111: if (mat) {MatDestroy(mat);}
113: PetscFinalize();
114: return 0;
115: }
116: