Actual source code: ex15.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
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>
18: int main(int argc,char **argv)
19: {
20: Mat A; /* operator matrix */
21: Vec u,v; /* left and right singular vectors */
22: SVD svd; /* singular value problem solver context */
23: SVDType type;
24: PetscReal error,tol,sigma,mu=PETSC_SQRT_MACHINE_EPSILON;
25: PetscInt n=100,i,j,Istart,Iend,nsv,maxit,its,nconv;
28: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
30: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
31: PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
32: PetscPrintf(PETSC_COMM_WORLD,"\nLauchli singular value decomposition, (%D x %D) mu=%g\n\n",n+1,n,(double)mu);
34: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: Build the Lauchli matrix
36: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
38: MatCreate(PETSC_COMM_WORLD,&A);
39: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n);
40: MatSetFromOptions(A);
41: MatSetUp(A);
43: MatGetOwnershipRange(A,&Istart,&Iend);
44: for (i=Istart;i<Iend;i++) {
45: if (i == 0) {
46: for (j=0;j<n;j++) {
47: MatSetValue(A,0,j,1.0,INSERT_VALUES);
48: }
49: } else {
50: MatSetValue(A,i,i-1,mu,INSERT_VALUES);
51: }
52: }
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
56: MatCreateVecs(A,&v,&u);
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);
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,NULL,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",(double)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 relative error\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,u,v);
123: /*
124: Compute the error associated to each singular triplet
125: */
126: SVDComputeError(svd,i,SVD_ERROR_RELATIVE,&error);
128: PetscPrintf(PETSC_COMM_WORLD," % 6f ",(double)sigma);
129: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",(double)error);
130: }
131: PetscPrintf(PETSC_COMM_WORLD,"\n");
132: }
134: /*
135: Free work space
136: */
137: SVDDestroy(&svd);
138: MatDestroy(&A);
139: VecDestroy(&u);
140: VecDestroy(&v);
141: SlepcFinalize();
142: return ierr;
143: }
145: /*TEST
147: test:
148: suffix: 1
149: filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
150: requires: double
152: TEST*/