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,PETSC_DECIDE,PETSC_DECIDE,m,n,&mat);
38: MatSetFromOptions(mat);
39: MatGetOwnershipRange(mat,&rstart,&rend);
40: for (i=rstart; i<rend; i++) {
41: for (j=0; j<n; j++) {
42: v=10*i+j;
43: MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES);
44: }
45: }
46: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
47: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
49: /* Test whether matrix has been corrupted (just to demonstrate this
50: routine) not needed in most application codes. */
51: MatValid(mat,(PetscTruth*)&flg);
52: if (!flg) SETERRQ(1,"Corrupted matrix.");
54: /* ----------------- Test MatNorm() ----------------- */
55: #ifdef TestMatNorm
56: MatNorm(mat,NORM_FROBENIUS,&normf);
57: MatNorm(mat,NORM_1,&norm1);
58: MatNorm(mat,NORM_INFINITY,&normi);
59: PetscPrintf(PETSC_COMM_WORLD,"original: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",
60: normf,norm1,normi);
61: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
62: #endif
63: /* --------------- Test MatTranspose() -------------- */
64: #ifdef TestMatTranspose
65: PetscOptionsHasName(PETSC_NULL,"-in_place",&flg);
66: if (!rect && flg) {
67: MatTranspose(mat,0); /* in-place transpose */
68: tmat = mat; mat = 0;
69: } else { /* out-of-place transpose */
70: MatTranspose(mat,&tmat);
71: }
72: #endif
73: /* ----------------- Test MatNorm() ----------------- */
74: #ifdef TestMatNorm_tmat
75: /* Print info about transpose matrix */
76: MatNorm(tmat,NORM_FROBENIUS,&normf);
77: MatNorm(tmat,NORM_1,&norm1);
78: MatNorm(tmat,NORM_INFINITY,&normi);
79: PetscPrintf(PETSC_COMM_WORLD,"transpose: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",
80: normf,norm1,normi);
81: MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
82: #endif
83: /* ----------------- Test MatAXPY() ----------------- */
84: #ifdef TestMatAXPY
85: if (mat && !rect) {
86: alpha = 1.0;
87: PetscOptionsGetScalar(PETSC_NULL,"-alpha",&alpha,PETSC_NULL);
88: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A\n");
89: MatAXPY(&alpha,mat,tmat,DIFFERENT_NONZERO_PATTERN);
90: MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
91: }
92: MatDestroy(tmat);
93: #endif
94: #ifdef TestMatAXPY2
95: Mat matB;
96: alpha = 1.0;
97: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A, SAME_NONZERO_PATTERN\n");
98: MatDuplicate(mat,MAT_COPY_VALUES,&matB);
99: MatAXPY(&alpha,mat,matB,SAME_NONZERO_PATTERN);
100: MatView(matB,PETSC_VIEWER_STDOUT_WORLD);
101: MatDestroy(matB);
103: /* get matB that has nonzeros of mat in all even numbers of row and col */
104: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,&matB);
105: MatSetFromOptions(matB);
106: MatGetOwnershipRange(matB,&rstart,&rend);
107: for (i=rstart; i<rend; i += 2) {
108: for (j=0; j<n; j += 2) {
109: v=10*i+j;
110: MatSetValues(matB,1,&i,1,&j,&v,INSERT_VALUES);
111: }
112: }
113: MatAssemblyBegin(matB,MAT_FINAL_ASSEMBLY);
114: MatAssemblyEnd(matB,MAT_FINAL_ASSEMBLY);
115: PetscPrintf(PETSC_COMM_WORLD," B: original matrix:\n");
116: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
117: PetscPrintf(PETSC_COMM_WORLD," A(a subset of B):\n");
118: MatView(matB,PETSC_VIEWER_STDOUT_WORLD);
119: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A, SUBSET_NONZERO_PATTERN\n");
120: MatAXPY(&alpha,matB,mat,SUBSET_NONZERO_PATTERN);
121: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
122:
123: PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A, SUBSET_NONZERO_PATTERN\n");
124: MatAXPY(&alpha,matB,mat,SUBSET_NONZERO_PATTERN);
125: MatView(mat,PETSC_VIEWER_STDOUT_WORLD);
126:
127: MatDestroy(matB);
128: #endif
129:
131: /* Free data structures */
132: if (mat) {MatDestroy(mat);}
134: PetscFinalize();
135: return 0;
136: }
137: