Actual source code: ex15.c
petsc-3.7.3 2016-08-01
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: PetscBool flg;
14: PetscScalar v;
15: IS perm,iperm;
16: Vec x,u,b,y;
17: PetscReal norm,tol=PETSC_SMALL;
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(PETSC_COMM_WORLD,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: MatSetUp(C);
28: PetscOptionsHasName(NULL,NULL,"-symmetric",&flg);
29: if (flg) { /* Treat matrix as symmetric only if we set this flag */
30: MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE);
31: MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
32: }
34: /* Create the matrix for the five point stencil, YET AGAIN */
35: for (i=0; i<m; i++) {
36: for (j=0; j<n; j++) {
37: v = -1.0; Ii = j + n*i;
38: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
39: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
40: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
41: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
42: v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
43: }
44: }
45: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
46: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
47: MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm);
48: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
49: ISView(perm,PETSC_VIEWER_STDOUT_SELF);
50: VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
51: VecSet(u,1.0);
52: VecDuplicate(u,&x);
53: VecDuplicate(u,&b);
54: VecDuplicate(u,&y);
55: MatMult(C,u,b);
56: VecCopy(b,y);
57: VecScale(y,2.0);
59: MatNorm(C,NORM_FROBENIUS,&norm);
60: PetscPrintf(PETSC_COMM_SELF,"Frobenius norm of matrix %g\n",(double)norm);
61: MatNorm(C,NORM_1,&norm);
62: PetscPrintf(PETSC_COMM_SELF,"One norm of matrix %g\n",(double)norm);
63: MatNorm(C,NORM_INFINITY,&norm);
64: PetscPrintf(PETSC_COMM_SELF,"Infinity norm of matrix %g\n",(double)norm);
66: MatFactorInfoInitialize(&info);
67: info.fill = 2.0;
68: info.dtcol = 0.0;
69: info.zeropivot = 1.e-14;
70: info.pivotinblocks = 1.0;
72: MatLUFactor(C,perm,iperm,&info);
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: if (norm > tol) {
81: PetscPrintf(PETSC_COMM_SELF,"MatSolve: Norm of error %g\n",(double)norm);
82: }
84: /* Test MatSolveAdd */
85: MatSolveAdd(C,b,y,x);
86: VecAXPY(x,-1.0,y);
87: VecAXPY(x,-1.0,u);
88: VecNorm(x,NORM_2,&norm);
89: if (norm > tol) {
90: PetscPrintf(PETSC_COMM_SELF,"MatSolveAdd(): Norm of error %g\n",(double)norm);
91: }
93: ISDestroy(&perm);
94: ISDestroy(&iperm);
95: VecDestroy(&u);
96: VecDestroy(&y);
97: VecDestroy(&b);
98: VecDestroy(&x);
99: MatDestroy(&C);
100: PetscFinalize();
101: return 0;
102: }