Actual source code: ex7.c
1: /*$Id: ex7.c,v 1.20 2001/08/07 21:30:08 bsmith Exp $*/
3: static char help[] = "Tests matrix factorization. Note that most users should\n\
4: employ the SLES interface to the linear solvers instead of using the factorization\n\
5: routines directly.\n\n";
7: #include petscmat.h
11: int main(int argc,char **args)
12: {
13: Mat C,LU;
14: MatInfo info;
15: int i,j,m = 3,n = 3,I,J,ierr;
16: PetscScalar v,mone = -1.0,one = 1.0;
17: IS perm,iperm;
18: Vec x,u,b;
19: PetscReal norm;
20: MatFactorInfo luinfo;
22: PetscInitialize(&argc,&args,(char *)0,help);
24: /* Create the matrix for the five point stencil, YET AGAIN */
25: MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
26: MatSetFromOptions(C);
27: for (i=0; i<m; i++) {
28: for (j=0; j<n; j++) {
29: v = -1.0; I = j + n*i;
30: if (i>0) {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
31: if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
32: if (j>0) {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
33: if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
34: v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
35: }
36: }
37: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
38: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
39: MatGetOrdering(C,MATORDERING_RCM,&perm,&iperm);
40: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
41: ISView(perm,PETSC_VIEWER_STDOUT_SELF);
43: luinfo.fill = 2.0;
44: luinfo.dtcol = 0.0;
45: luinfo.damping = 0.0;
46: luinfo.zeropivot = 1.e-14;
47: luinfo.pivotinblocks = 1.0;
48: MatLUFactorSymbolic(C,perm,iperm,&luinfo,&LU);
49: MatLUFactorNumeric(C,&LU);
50: MatView(LU,PETSC_VIEWER_STDOUT_WORLD);
52: VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
53: VecSet(&one,u);
54: VecDuplicate(u,&x);
55: VecDuplicate(u,&b);
57: MatMult(C,u,b);
58: MatSolve(LU,b,x);
60: VecView(b,PETSC_VIEWER_STDOUT_SELF);
61: VecView(x,PETSC_VIEWER_STDOUT_SELF);
63: VecAXPY(&mone,u,x);
64: VecNorm(x,NORM_2,&norm);
65: PetscPrintf(PETSC_COMM_SELF,"Norm of error %A\n",norm);
67: MatGetInfo(C,MAT_LOCAL,&info);
68: PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %d\n",(int)info.nz_used);
69: MatGetInfo(LU,MAT_LOCAL,&info);
70: PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %d\n",(int)info.nz_used);
72: VecDestroy(u);
73: VecDestroy(b);
74: VecDestroy(x);
75: ISDestroy(perm);
76: ISDestroy(iperm);
77: MatDestroy(C);
78: MatDestroy(LU);
80: PetscFinalize();
81: return 0;
82: }