Actual source code: ex15.c
slepc-3.22.2 2024-12-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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;
27: PetscFunctionBeginUser;
28: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
31: PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
32: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nLauchli singular value decomposition, (%" PetscInt_FMT " x %" PetscInt_FMT ") mu=%g\n\n",n+1,n,(double)mu));
34: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: Build the Lauchli matrix
36: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
38: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
39: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n));
40: PetscCall(MatSetFromOptions(A));
42: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
43: for (i=Istart;i<Iend;i++) {
44: if (i == 0) {
45: for (j=0;j<n;j++) PetscCall(MatSetValue(A,0,j,1.0,INSERT_VALUES));
46: } else PetscCall(MatSetValue(A,i,i-1,mu,INSERT_VALUES));
47: }
49: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
50: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
51: PetscCall(MatCreateVecs(A,&v,&u));
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Create the singular value solver and set various options
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: /*
58: Create singular value solver context
59: */
60: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
62: /*
63: Set operators and problem type
64: */
65: PetscCall(SVDSetOperators(svd,A,NULL));
66: PetscCall(SVDSetProblemType(svd,SVD_STANDARD));
68: /*
69: Use thick-restart Lanczos as default solver
70: */
71: PetscCall(SVDSetType(svd,SVDTRLANCZOS));
73: /*
74: Set solver parameters at runtime
75: */
76: PetscCall(SVDSetFromOptions(svd));
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Solve the singular value system
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: PetscCall(SVDSolve(svd));
83: PetscCall(SVDGetIterationNumber(svd,&its));
84: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its));
86: /*
87: Optional: Get some information from the solver and display it
88: */
89: PetscCall(SVDGetType(svd,&type));
90: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
91: PetscCall(SVDGetDimensions(svd,&nsv,NULL,NULL));
92: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %" PetscInt_FMT "\n",nsv));
93: PetscCall(SVDGetTolerances(svd,&tol,&maxit));
94: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Display solution and clean up
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: /*
101: Get number of converged singular triplets
102: */
103: PetscCall(SVDGetConverged(svd,&nconv));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %" PetscInt_FMT "\n\n",nconv));
106: if (nconv>0) {
107: /*
108: Display singular values and relative errors
109: */
110: PetscCall(PetscPrintf(PETSC_COMM_WORLD,
111: " sigma relative error\n"
112: " --------------------- ------------------\n"));
113: for (i=0;i<nconv;i++) {
114: /*
115: Get converged singular triplets: i-th singular value is stored in sigma
116: */
117: PetscCall(SVDGetSingularTriplet(svd,i,&sigma,u,v));
119: /*
120: Compute the error associated to each singular triplet
121: */
122: PetscCall(SVDComputeError(svd,i,SVD_ERROR_RELATIVE,&error));
124: PetscCall(PetscPrintf(PETSC_COMM_WORLD," % 6f ",(double)sigma));
125: PetscCall(PetscPrintf(PETSC_COMM_WORLD," % 12g\n",(double)error));
126: }
127: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
128: }
130: /*
131: Free work space
132: */
133: PetscCall(SVDDestroy(&svd));
134: PetscCall(MatDestroy(&A));
135: PetscCall(VecDestroy(&u));
136: PetscCall(VecDestroy(&v));
137: PetscCall(SlepcFinalize());
138: return 0;
139: }
141: /*TEST
143: testset:
144: filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
145: requires: double
146: test:
147: suffix: 1
148: test:
149: suffix: 1_scalapack
150: nsize: {{1 2}}
151: args: -svd_type scalapack
152: requires: scalapack
154: TEST*/