Actual source code: ex15.c
2: static char help[] = "Tests MatNorm(), MatLUFactor(), MatSolve() and MatSolveAdd().\n\n";
4: #include petscmat.h
8: int main(int argc,char **args)
9: {
10: Mat C;
11: PetscInt i,j,m = 3,n = 3,Ii,J;
13: PetscTruth flg;
14: PetscScalar v;
15: IS perm,iperm;
16: Vec x,u,b,y;
17: PetscReal norm;
18: MatFactorInfo info;
19: PetscMPIInt size;
21: PetscInitialize(&argc,&args,(char *)0,help);
22: MPI_Comm_size(PETSC_COMM_WORLD,&size);
23: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
24: MatCreate(PETSC_COMM_WORLD,&C);
25: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
26: MatSetFromOptions(C);
27: PetscOptionsHasName(PETSC_NULL,"-symmetric",&flg);
28: if (flg) { /* Treat matrix as symmetric only if we set this flag */
29: MatSetOption(C,MAT_SYMMETRIC);
30: MatSetOption(C,MAT_SYMMETRY_ETERNAL);
31: }
33: /* Create the matrix for the five point stencil, YET AGAIN */
34: for (i=0; i<m; i++) {
35: for (j=0; j<n; j++) {
36: v = -1.0; Ii = j + n*i;
37: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
38: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
39: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
40: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
41: v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
42: }
43: }
44: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
45: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
46: MatGetOrdering(C,MATORDERING_RCM,&perm,&iperm);
47: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
48: ISView(perm,PETSC_VIEWER_STDOUT_SELF);
49: VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
50: VecSet(u,1.0);
51: VecDuplicate(u,&x);
52: VecDuplicate(u,&b);
53: VecDuplicate(u,&y);
54: MatMult(C,u,b);
55: VecCopy(b,y);
56: VecScale(y,2.0);
58: MatNorm(C,NORM_FROBENIUS,&norm);
59: PetscPrintf(PETSC_COMM_SELF,"Frobenius norm of matrix %G\n",norm);
60: MatNorm(C,NORM_1,&norm);
61: PetscPrintf(PETSC_COMM_SELF,"One norm of matrix %G\n",norm);
62: MatNorm(C,NORM_INFINITY,&norm);
63: PetscPrintf(PETSC_COMM_SELF,"Infinity norm of matrix %G\n",norm);
65: MatFactorInfoInitialize(&info);
66: info.fill = 2.0;
67: info.dtcol = 0.0;
68: info.shiftnz = 0.0;
69: info.zeropivot = 1.e-14;
70: info.pivotinblocks = 1.0;
71: MatLUFactor(C,perm,iperm,&info);
72: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
74: /* Test MatSolve */
75: MatSolve(C,b,x);
76: VecView(b,PETSC_VIEWER_STDOUT_SELF);
77: VecView(x,PETSC_VIEWER_STDOUT_SELF);
78: VecAXPY(x,-1.0,u);
79: VecNorm(x,NORM_2,&norm);
80: PetscPrintf(PETSC_COMM_SELF,"Norm of error %A\n",norm);
82: /* Test MatSolveAdd */
83: MatSolveAdd(C,b,y,x);
84: VecAXPY(x,-1.0,y);
85: VecAXPY(x,-1.0,u);
86: VecNorm(x,NORM_2,&norm);
88: PetscPrintf(PETSC_COMM_SELF,"Norm of error %A\n",norm);
90: ISDestroy(perm);
91: ISDestroy(iperm);
92: VecDestroy(u);
93: VecDestroy(y);
94: VecDestroy(b);
95: VecDestroy(x);
96: MatDestroy(C);
97: PetscFinalize();
98: return 0;
99: }