Actual source code: ex15.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc. See the README file for conditions of use
7: and additional information.
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Singular value decomposition of the Lauchli matrix.\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = matrix dimension.\n"
14: " -mu <mu>, where <mu> = subdiagonal value.\n\n";
16: #include slepcsvd.h
20: int main( int argc, char **argv )
21: {
22: Mat A; /* operator matrix */
23: SVD svd; /* singular value problem solver context */
24: SVDType type;
25: PetscReal error, tol, sigma, mu=PETSC_SQRT_MACHINE_EPSILON;
27: PetscInt n=100, i, j, start, end;
28: int nsv, maxit, its, nconv;
30: SlepcInitialize(&argc,&argv,(char*)0,help);
32: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
33: PetscOptionsGetReal(PETSC_NULL,"-mu",&mu,PETSC_NULL);
34: PetscPrintf(PETSC_COMM_WORLD,"\nLauchli singular value decomposition, (%d x %d) mu=%g\n\n",n+1,n,mu);
36: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: Build the Lauchli matrix
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: MatCreate(PETSC_COMM_WORLD,&A);
41: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n);
42: MatSetFromOptions(A);
44: MatGetOwnershipRange(A,&start,&end);
45: for (i=start;i<end;i++) {
46: if (i == 0) {
47: for (j=0;j<n;j++) {
48: MatSetValue(A,0,j,1.0,INSERT_VALUES);
49: }
50: } else {
51: MatSetValue(A,i,i-1,mu,INSERT_VALUES);
52: }
53: }
54:
55: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
56: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
58: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: Create the singular value solver and set various options
60: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62: /*
63: Create singular value solver context
64: */
65: SVDCreate(PETSC_COMM_WORLD,&svd);
67: /*
68: Set operator
69: */
70: SVDSetOperator(svd,A);
71:
72: /*
73: Use thick-restart Lanczos as default solver
74: */
75: SVDSetType(svd,SVDTRLANCZOS);
77: /*
78: Set solver parameters at runtime
79: */
80: SVDSetFromOptions(svd);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Solve the singular value system
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: SVDSolve(svd);
87: SVDGetIterationNumber(svd, &its);
88: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
90: /*
91: Optional: Get some information from the solver and display it
92: */
93: SVDGetType(svd,&type);
94: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
95: SVDGetDimensions(svd,&nsv,PETSC_NULL);
96: PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %d\n",nsv);
97: SVDGetTolerances(svd,&tol,&maxit);
98: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Display solution and clean up
102: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: /*
105: Get number of converged singular triplets
106: */
107: SVDGetConverged(svd,&nconv);
108: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %d\n\n",nconv);
110: if (nconv>0) {
111: /*
112: Display singular values and relative errors
113: */
114: PetscPrintf(PETSC_COMM_WORLD,
115: " sigma residual norm\n"
116: " --------------------- ------------------\n" );
117: for( i=0; i<nconv; i++ ) {
118: /*
119: Get converged singular triplets: i-th singular value is stored in sigma
120: */
121: SVDGetSingularTriplet(svd,i,&sigma,PETSC_NULL,PETSC_NULL);
123: /*
124: Compute the error associated to each singular triplet
125: */
126: SVDComputeRelativeError(svd,i,&error);
128: PetscPrintf(PETSC_COMM_WORLD," % 6f ",sigma);
129: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);
130: }
131: PetscPrintf(PETSC_COMM_WORLD,"\n" );
132: }
133:
134: /*
135: Free work space
136: */
137: SVDDestroy(svd);
138: MatDestroy(A);
139: SlepcFinalize();
140: return 0;
141: }