Actual source code: ex15.c
1: /*$Id: ex15.c,v 1.22 2001/08/07 21:30:08 bsmith Exp $*/
3: static char help[] = "Tests MatNorm(), MatLUFactor(), MatSolve() and MatSolveAdd().\n\n";
5: #include petscmat.h
9: int main(int argc,char **args)
10: {
11: Mat C;
12: int i,j,m = 3,n = 3,I,J,ierr;
13: PetscTruth flg;
14: PetscScalar v,mone = -1.0,one = 1.0,alpha = 2.0;
15: IS perm,iperm;
16: Vec x,u,b,y;
17: PetscReal norm;
18: MatFactorInfo info;
20: PetscInitialize(&argc,&args,(char *)0,help);
22: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
23: MatSetFromOptions(C);
24: PetscOptionsHasName(PETSC_NULL,"-symmetric",&flg);
25: if (flg) { /* Treat matrix as symmetric only if we set this flag */
26: MatSetOption(C,MAT_SYMMETRIC);
27: }
29: /* Create the matrix for the five point stencil, YET AGAIN */
30: for (i=0; i<m; i++) {
31: for (j=0; j<n; j++) {
32: v = -1.0; I = j + n*i;
33: if (i>0) {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
34: if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
35: if (j>0) {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
36: if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
37: v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
38: }
39: }
40: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
41: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
42: MatGetOrdering(C,MATORDERING_RCM,&perm,&iperm);
43: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
44: ISView(perm,PETSC_VIEWER_STDOUT_SELF);
45: VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
46: VecSet(&one,u);
47: VecDuplicate(u,&x);
48: VecDuplicate(u,&b);
49: VecDuplicate(u,&y);
50: MatMult(C,u,b);
51: VecCopy(b,y);
52: VecScale(&alpha,y);
54: MatNorm(C,NORM_FROBENIUS,&norm);
55: PetscPrintf(PETSC_COMM_SELF,"Frobenius norm of matrix %g\n",norm);
56: MatNorm(C,NORM_1,&norm);
57: PetscPrintf(PETSC_COMM_SELF,"One norm of matrix %g\n",norm);
58: MatNorm(C,NORM_INFINITY,&norm);
59: PetscPrintf(PETSC_COMM_SELF,"Infinity norm of matrix %g\n",norm);
61: info.fill = 2.0;
62: info.dtcol = 0.0;
63: info.damping = 0.0;
64: info.zeropivot = 1.e-14;
65: info.pivotinblocks = 1.0;
66: MatLUFactor(C,perm,iperm,&info);
67: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
69: /* Test MatSolve */
70: MatSolve(C,b,x);
71: VecView(b,PETSC_VIEWER_STDOUT_SELF);
72: VecView(x,PETSC_VIEWER_STDOUT_SELF);
73: VecAXPY(&mone,u,x);
74: VecNorm(x,NORM_2,&norm);
75: PetscPrintf(PETSC_COMM_SELF,"Norm of error %A\n",norm);
77: /* Test MatSolveAdd */
78: MatSolveAdd(C,b,y,x);
79: VecAXPY(&mone,y,x);
80: VecAXPY(&mone,u,x);
81: VecNorm(x,NORM_2,&norm);
83: PetscPrintf(PETSC_COMM_SELF,"Norm of error %A\n",norm);
85: ISDestroy(perm);
86: ISDestroy(iperm);
87: VecDestroy(u);
88: VecDestroy(y);
89: VecDestroy(b);
90: VecDestroy(x);
91: MatDestroy(C);
92: PetscFinalize();
93: return 0;
94: }