Actual source code: ex1.c
1: /*$Id: ex1.c,v 1.23 2001/08/07 21:30:08 bsmith Exp $*/
3: static char help[] = "Tests LU and Cholesky factorization for a dense matrix.nn";
5: #include petscmat.h
7: #undef __FUNCT__
9: int main(int argc,char **argv)
10: {
11: Mat mat,fact;
12: MatInfo info;
13: int m = 10,n = 10,i = 4,ierr,rstart,rend;
14: PetscScalar value = 1.0;
15: Vec x,y,b;
16: PetscReal norm;
17: IS perm;
18: MatLUInfo luinfo;
20: PetscInitialize(&argc,&argv,(char*) 0,help);
22: VecCreate(PETSC_COMM_WORLD,&y);
23: VecSetSizes(y,PETSC_DECIDE,m);
24: VecSetFromOptions(y);
25: VecDuplicate(y,&x);
26: VecSet(&value,x);
27: VecCreate(PETSC_COMM_WORLD,&b);
28: VecSetSizes(b,PETSC_DECIDE,n);
29: VecSetFromOptions(b);
30: ISCreateStride(PETSC_COMM_WORLD,m,0,1,&perm);
32: MatCreateSeqDense(PETSC_COMM_WORLD,m,n,PETSC_NULL,&mat);
34: MatGetOwnershipRange(mat,&rstart,&rend);
35: for (i=rstart; i<rend; i++) {
36: value = (PetscReal)i+1;
37: MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES);
38: }
39: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
40: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
42: MatGetInfo(mat,MAT_LOCAL,&info);
43: PetscPrintf(PETSC_COMM_WORLD,"matrix nonzeros = %d, allocated nonzeros = %dn",
44: (int)info.nz_used,(int)info.nz_allocated);
46: /* Cholesky factorization is not yet in place for this matrix format */
47: MatMult(mat,x,b);
48: MatConvert(mat,MATSAME,&fact);
49: MatCholeskyFactor(fact,perm,1.0);
50: MatSolve(fact,b,y);
51: MatDestroy(fact);
52: value = -1.0; VecAXPY(&value,x,y);
53: VecNorm(y,NORM_2,&norm);
54: PetscPrintf(PETSC_COMM_WORLD,"Norm of error for Cholesky %An",norm);
55: MatCholeskyFactorSymbolic(mat,perm,1.0,&fact);
56: MatCholeskyFactorNumeric(mat,&fact);
57: MatSolve(fact,b,y);
58: value = -1.0; VecAXPY(&value,x,y);
59: VecNorm(y,NORM_2,&norm);
60: PetscPrintf(PETSC_COMM_WORLD,"Norm of error for Cholesky %An",norm);
61: MatDestroy(fact);
63: i = m-1; value = 1.0;
64: MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES);
65: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
66: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
67: MatMult(mat,x,b);
68: MatConvert(mat,MATSAME,&fact);
69: MatLUFactor(fact,perm,perm,PETSC_NULL);
70: MatSolve(fact,b,y);
71: value = -1.0; VecAXPY(&value,x,y);
72: VecNorm(y,NORM_2,&norm);
73: PetscPrintf(PETSC_COMM_WORLD,"Norm of error for LU %An",norm);
74: MatDestroy(fact);
76: luinfo.fill = 2.0;
77: luinfo.dtcol = 0.0;
78: luinfo.damping = 0.0;
79: luinfo.zeropivot = 1.e-14;
80: luinfo.pivotinblocks = 1.0;
81: MatLUFactorSymbolic(mat,perm,perm,&luinfo,&fact);
82: MatLUFactorNumeric(mat,&fact);
83: MatSolve(fact,b,y);
84: value = -1.0; VecAXPY(&value,x,y);
85: VecNorm(y,NORM_2,&norm);
86: PetscPrintf(PETSC_COMM_WORLD,"Norm of error for LU %An",norm);
87: MatDestroy(fact);
88: MatDestroy(mat);
89: VecDestroy(x);
90: VecDestroy(b);
91: VecDestroy(y);
92: ISDestroy(perm);
94: PetscFinalize();
95: return 0;
96: }
97: