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: int main(int argc,char **argv)
8: {
9: Mat mat,fact;
10: MatInfo info;
11: int m = 10,n = 10,i = 4,ierr,rstart,rend;
12: PetscScalar value = 1.0;
13: Vec x,y,b;
14: PetscReal norm;
15: IS perm;
17: PetscInitialize(&argc,&argv,(char*) 0,help);
19: VecCreate(PETSC_COMM_WORLD,&y);
20: VecSetSizes(y,PETSC_DECIDE,m);
21: VecSetFromOptions(y);
22: VecDuplicate(y,&x);
23: VecSet(&value,x);
24: VecCreate(PETSC_COMM_WORLD,&b);
25: VecSetSizes(b,PETSC_DECIDE,n);
26: VecSetFromOptions(b);
27: ISCreateStride(PETSC_COMM_WORLD,m,0,1,&perm);
29: MatCreateSeqDense(PETSC_COMM_WORLD,m,n,PETSC_NULL,&mat);
31: MatGetOwnershipRange(mat,&rstart,&rend);
32: for (i=rstart; i<rend; i++) {
33: value = (PetscReal)i+1;
34: MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES);
35: }
36: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
37: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
39: MatGetInfo(mat,MAT_LOCAL,&info);
40: PetscPrintf(PETSC_COMM_WORLD,"matrix nonzeros = %d, allocated nonzeros = %dn",
41: (int)info.nz_used,(int)info.nz_allocated);
43: /* Cholesky factorization is not yet in place for this matrix format */
44: MatMult(mat,x,b);
45: MatConvert(mat,MATSAME,&fact);
46: MatCholeskyFactor(fact,perm,1.0);
47: MatSolve(fact,b,y);
48: MatDestroy(fact);
49: value = -1.0; VecAXPY(&value,x,y);
50: VecNorm(y,NORM_2,&norm);
51: PetscPrintf(PETSC_COMM_WORLD,"Norm of error for Cholesky %An",norm);
52: MatCholeskyFactorSymbolic(mat,perm,1.0,&fact);
53: MatCholeskyFactorNumeric(mat,&fact);
54: MatSolve(fact,b,y);
55: value = -1.0; VecAXPY(&value,x,y);
56: VecNorm(y,NORM_2,&norm);
57: PetscPrintf(PETSC_COMM_WORLD,"Norm of error for Cholesky %An",norm);
58: MatDestroy(fact);
60: i = m-1; value = 1.0;
61: MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES);
62: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
63: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
64: MatMult(mat,x,b);
65: MatConvert(mat,MATSAME,&fact);
66: MatLUFactor(fact,perm,perm,PETSC_NULL);
67: MatSolve(fact,b,y);
68: value = -1.0; VecAXPY(&value,x,y);
69: VecNorm(y,NORM_2,&norm);
70: PetscPrintf(PETSC_COMM_WORLD,"Norm of error for LU %An",norm);
71: MatDestroy(fact);
73: MatLUFactorSymbolic(mat,perm,perm,PETSC_NULL,&fact);
74: MatLUFactorNumeric(mat,&fact);
75: MatSolve(fact,b,y);
76: value = -1.0; VecAXPY(&value,x,y);
77: VecNorm(y,NORM_2,&norm);
78: PetscPrintf(PETSC_COMM_WORLD,"Norm of error for LU %An",norm);
79: MatDestroy(fact);
80: MatDestroy(mat);
81: VecDestroy(x);
82: VecDestroy(b);
83: VecDestroy(y);
84: ISDestroy(perm);
86: PetscFinalize();
87: return 0;
88: }
89: