Actual source code: ex2.c
2: static char help[] = "Tests MatTranspose(), MatNorm(), MatValid(), and MatAXPY().\n\n";
4: #include petscmat.h
6: #define TestMatNorm
7: #define TestMatTranspose
8: #define TestMatNorm_tmat
9: #define TestMatAXPY
10: #undef TestMatAXPY2
14: int main(int argc,char **argv)
15: {
16: Mat mat,tmat = 0;
17: PetscInt m = 7,n,i,j,rstart,rend,rect = 0;
18: PetscErrorCode ierr;
19: PetscMPIInt size,rank;
20: PetscTruth flg;
21: PetscScalar v, alpha;
22: PetscReal normf,normi,norm1;
24: PetscInitialize(&argc,&argv,(char*)0,help);
25: PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);
26: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
27: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
28: MPI_Comm_size(PETSC_COMM_WORLD,&size);
29: n = m;
30: PetscOptionsHasName(PETSC_NULL,"-rectA",&flg);
31: if (flg) {n += 2; rect = 1;}
32: PetscOptionsHasName(PETSC_NULL,"-rectB",&flg);
33: if (flg) {n -= 2; rect = 1;}
35: /* ------- Assemble matrix, test MatValid() --------- */
37: MatCreate(PETSC_COMM_WORLD,&mat);
38: MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n);
39: MatSetFromOptions(mat);
40: MatGetOwnershipRange(mat,&rstart,&rend);
41: for (i=rstart; i<rend; i++) {
42: for (j=0; j<n; j++) {
43: v=10*i+j;
44: MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES);
45: }
46: }
47: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
48: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
50: /* Test whether matrix has been corrupted (just to demonstrate this
51: routine) not needed in most application codes. */
52: MatValid(mat,(PetscTruth*)&flg);
53: if (!flg) SETERRQ(1,"Corrupted matrix.");
55: /* ----------------- Test MatNorm() ----------------- */
56: #ifdef TestMatNorm
57: MatNorm(mat,NORM_FROBENIUS,&normf);
58: MatNorm(mat,NORM_1,&norm1);
59: MatNorm(mat,NORM_INFINITY,&normi);
60: PetscPrintf(PETSC_COMM_WORLD,"original: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",
61: normf,norm1,normi);
62: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
63: #endif
64: /* --------------- Test MatTranspose() -------------- */
65: #ifdef TestMatTranspose
66: PetscOptionsHasName(PETSC_NULL,"-in_place",&flg);
67: if (!rect && flg) {
68: MatTranspose(mat,0); /* in-place transpose */
69: tmat = mat; mat = 0;
70: } else { /* out-of-place transpose */
71: MatTranspose(mat,&tmat);
72: }
73: #endif
74: /* ----------------- Test MatNorm() ----------------- */
75: #ifdef TestMatNorm_tmat
76: /* Print info about transpose matrix */
77: MatNorm(tmat,NORM_FROBENIUS,&normf);
78: MatNorm(tmat,NORM_1,&norm1);
79: MatNorm(tmat,NORM_INFINITY,&normi);
80: PetscPrintf(PETSC_COMM_WORLD,"transpose: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",
81: normf,norm1,normi);
82: MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
83: #endif
84: /* ----------------- Test MatAXPY() ----------------- */
85: #ifdef TestMatAXPY
86: if (mat && !rect) {
87: alpha = 1.0;
88: PetscOptionsGetScalar(PETSC_NULL,"-alpha",&alpha,PETSC_NULL);
89: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A\n");
90: MatAXPY(tmat,alpha,mat,DIFFERENT_NONZERO_PATTERN);
91: MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
92: }
93: MatDestroy(tmat);
94: #endif
95: #ifdef TestMatAXPY2
96: Mat matB;
97: alpha = 1.0;
98: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A, SAME_NONZERO_PATTERN\n");
99: MatDuplicate(mat,MAT_COPY_VALUES,&matB);
100: MatAXPY(matB,alpha,mat,SAME_NONZERO_PATTERN);
101: MatView(matB,PETSC_VIEWER_STDOUT_WORLD);
102: MatDestroy(matB);
104: /* get matB that has nonzeros of mat in all even numbers of row and col */
105: MatCreate(PETSC_COMM_WORLD,&matB);
106: MatSetSizes(matB,PETSC_DECIDE,PETSC_DECIDE,m,n);
107: MatSetFromOptions(matB);
108: MatGetOwnershipRange(matB,&rstart,&rend);
109: for (i=rstart; i<rend; i += 2) {
110: for (j=0; j<n; j += 2) {
111: v=10*i+j;
112: MatSetValues(matB,1,&i,1,&j,&v,INSERT_VALUES);
113: }
114: }
115: MatAssemblyBegin(matB,MAT_FINAL_ASSEMBLY);
116: MatAssemblyEnd(matB,MAT_FINAL_ASSEMBLY);
117: PetscPrintf(PETSC_COMM_WORLD," B: original matrix:\n");
118: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
119: PetscPrintf(PETSC_COMM_WORLD," A(a subset of B):\n");
120: MatView(matB,PETSC_VIEWER_STDOUT_WORLD);
121: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A, SUBSET_NONZERO_PATTERN\n");
122: MatAXPY(mat,alpha,matB,SUBSET_NONZERO_PATTERN);
123: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
124:
125: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A, SUBSET_NONZERO_PATTERN\n");
126: MatAXPY(mat,alpha,matB,SUBSET_NONZERO_PATTERN);
127: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
128:
129: MatDestroy(matB);
130: #endif
131:
133: /* Free data structures */
134: if (mat) {MatDestroy(mat);}
136: PetscFinalize();
137: return 0;
138: }
139: