Actual source code: ex45.c
slepc-3.18.1 2022-11-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[] = "Computes a partial generalized singular value decomposition (GSVD).\n"
12: "The command line options are:\n"
13: " -m <m>, where <m> = number of rows of A.\n"
14: " -n <n>, where <n> = number of columns of A.\n"
15: " -p <p>, where <p> = number of rows of B.\n\n";
17: #include <slepcsvd.h>
19: int main(int argc,char **argv)
20: {
21: Mat A,B; /* operator matrices */
22: Vec u,v,x; /* singular vectors */
23: SVD svd; /* singular value problem solver context */
24: SVDType type;
25: Vec uv,aux[2],*U,*V;
26: PetscReal error,tol,sigma,lev1=0.0,lev2=0.0;
27: PetscInt m=100,n,p=14,i,j,d,Istart,Iend,nsv,maxit,its,nconv;
28: PetscBool flg,skiporth=PETSC_FALSE;
31: SlepcInitialize(&argc,&argv,(char*)0,help);
33: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
34: PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
35: if (!flg) n = m;
36: PetscOptionsGetInt(NULL,NULL,"-p",&p,&flg);
37: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%" PetscInt_FMT "+%" PetscInt_FMT ")x%" PetscInt_FMT "\n\n",m,p,n);
38: PetscOptionsGetBool(NULL,NULL,"-skiporth",&skiporth,NULL);
40: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41: Build the matrices
42: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
44: MatCreate(PETSC_COMM_WORLD,&A);
45: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
46: MatSetFromOptions(A);
47: MatSetUp(A);
49: MatGetOwnershipRange(A,&Istart,&Iend);
50: for (i=Istart;i<Iend;i++) {
51: if (i>0 && i-1<n) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
52: if (i+1<n) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
53: if (i<n) MatSetValue(A,i,i,2.0,INSERT_VALUES);
54: if (i>n) MatSetValue(A,i,n-1,1.0,INSERT_VALUES);
55: }
56: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
59: MatCreate(PETSC_COMM_WORLD,&B);
60: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
61: MatSetFromOptions(B);
62: MatSetUp(B);
64: MatGetOwnershipRange(B,&Istart,&Iend);
65: d = PetscMax(0,n-p);
66: for (i=Istart;i<Iend;i++) {
67: for (j=0;j<=PetscMin(i,n-1);j++) MatSetValue(B,i,j+d,1.0,INSERT_VALUES);
68: }
69: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Create the singular value solver and set various options
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: /*
77: Create singular value solver context
78: */
79: SVDCreate(PETSC_COMM_WORLD,&svd);
81: /*
82: Set operators and problem type
83: */
84: SVDSetOperators(svd,A,B);
85: SVDSetProblemType(svd,SVD_GENERALIZED);
87: /*
88: Set solver parameters at runtime
89: */
90: SVDSetFromOptions(svd);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Solve the singular value system
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: SVDSolve(svd);
97: SVDGetIterationNumber(svd,&its);
98: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its);
100: /*
101: Optional: Get some information from the solver and display it
102: */
103: SVDGetType(svd,&type);
104: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
105: SVDGetDimensions(svd,&nsv,NULL,NULL);
106: PetscPrintf(PETSC_COMM_WORLD," Number of requested generalized singular values: %" PetscInt_FMT "\n",nsv);
107: SVDGetTolerances(svd,&tol,&maxit);
108: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit);
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Display solution and clean up
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: /*
115: Get number of converged singular triplets
116: */
117: SVDGetConverged(svd,&nconv);
118: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %" PetscInt_FMT "\n\n",nconv);
120: if (nconv>0) {
121: /*
122: Create vectors. The interface returns u and v as stacked on top of each other
123: [u;v] so need to create a special vector (VecNest) to extract them
124: */
125: MatCreateVecs(A,&x,&u);
126: MatCreateVecs(B,NULL,&v);
127: aux[0] = u;
128: aux[1] = v;
129: VecCreateNest(PETSC_COMM_WORLD,2,NULL,aux,&uv);
131: VecDuplicateVecs(u,nconv,&U);
132: VecDuplicateVecs(v,nconv,&V);
134: /*
135: Display singular values and errors relative to the norms
136: */
137: PetscCall(PetscPrintf(PETSC_COMM_WORLD,
138: " sigma ||r||/||[A;B]||\n"
139: " --------------------- ------------------\n"));
140: for (i=0;i<nconv;i++) {
141: /*
142: Get converged singular triplets: i-th singular value is stored in sigma
143: */
144: SVDGetSingularTriplet(svd,i,&sigma,uv,x);
146: /* at this point, u and v can be used normally as individual vectors */
147: VecCopy(u,U[i]);
148: VecCopy(v,V[i]);
150: /*
151: Compute the error associated to each singular triplet
152: */
153: SVDComputeError(svd,i,SVD_ERROR_NORM,&error);
155: PetscPrintf(PETSC_COMM_WORLD," % 6f ",(double)sigma);
156: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",(double)error);
157: }
158: PetscPrintf(PETSC_COMM_WORLD,"\n");
160: if (!skiporth) {
161: VecCheckOrthonormality(U,nconv,NULL,nconv,NULL,NULL,&lev1);
162: VecCheckOrthonormality(V,nconv,NULL,nconv,NULL,NULL,&lev2);
163: }
164: if (lev1+lev2<20*tol) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
165: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g (U) %g (V)\n",(double)lev1,(double)lev2);
166: VecDestroyVecs(nconv,&U);
167: VecDestroyVecs(nconv,&V);
168: VecDestroy(&x);
169: VecDestroy(&u);
170: VecDestroy(&v);
171: VecDestroy(&uv);
172: }
174: /*
175: Free work space
176: */
177: SVDDestroy(&svd);
178: MatDestroy(&A);
179: MatDestroy(&B);
180: SlepcFinalize();
181: return 0;
182: }
184: /*TEST
186: testset:
187: filter: grep -v "Solution method" | grep -v "Number of iterations" | sed -e "s/, maxit=1[0]*$//" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
188: requires: double
189: test:
190: args: -svd_type lapack -m 20 -n 10 -p 6
191: suffix: 1
192: test:
193: args: -svd_type lapack -m 15 -n 20 -p 10 -svd_smallest
194: suffix: 2
195: test:
196: args: -svd_type lapack -m 15 -n 20 -p 21
197: suffix: 3
198: test:
199: args: -svd_type lapack -m 20 -n 15 -p 21
200: suffix: 4
202: testset:
203: args: -m 25 -n 20 -p 21 -svd_smallest -svd_nsv 2
204: filter: grep -v "Solution method" | grep -v "Number of iterations" | sed -e "s/, maxit=1[0]*$//" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
205: requires: double
206: output_file: output/ex45_5.out
207: test:
208: args: -svd_type trlanczos -svd_ncv 8 -svd_trlanczos_gbidiag {{upper lower}} -svd_trlanczos_oneside {{0 1}}
209: suffix: 5
210: test:
211: args: -svd_type cross -svd_ncv 10 -svd_cross_explicitmatrix
212: suffix: 5_cross
213: test:
214: args: -svd_type cross -svd_ncv 10 -svd_cross_eps_type krylovschur -svd_cross_st_type sinvert -svd_cross_eps_target 0 -svd_cross_st_ksp_type gmres -svd_cross_st_pc_type jacobi
215: suffix: 5_cross_implicit
216: test:
217: args: -svd_type cyclic -svd_ncv 12 -svd_cyclic_explicitmatrix {{0 1}}
218: suffix: 5_cyclic
219: requires: !complex
221: testset:
222: args: -m 15 -n 20 -p 21 -svd_nsv 4 -svd_ncv 9
223: filter: grep -v "Solution method" | grep -v "Number of iterations" | sed -e "s/7.884967/7.884968/" | sed -e "s/, maxit=1[0]*$//" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
224: requires: double
225: output_file: output/ex45_6.out
226: test:
227: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}} -svd_trlanczos_locking {{0 1}} -svd_trlanczos_oneside {{0 1}}
228: suffix: 6
229: test:
230: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
231: suffix: 6_cross
233: test:
234: args: -m 15 -n 20 -p 21 -svd_nsv 4 -svd_ncv 9 -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
235: filter: grep -v "Number of iterations" | sed -e "s/7.884967/7.884968/" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
236: requires: double
237: suffix: 6_cyclic
238: output_file: output/ex45_6_cyclic.out
240: testset:
241: args: -m 20 -n 15 -p 21 -svd_nsv 4 -svd_ncv 9
242: filter: grep -v "Solution method" | grep -v "Number of iterations" | sed -e "s/, maxit=1[0]*$//" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
243: requires: double
244: output_file: output/ex45_7.out
245: test:
246: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}} -svd_trlanczos_restart 0.4 -svd_trlanczos_oneside {{0 1}}
247: suffix: 7
248: test:
249: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
250: suffix: 7_cross
252: test:
253: args: -m 20 -n 15 -p 21 -svd_nsv 4 -svd_ncv 16 -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
254: filter: grep -v "Number of iterations" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
255: requires: double
256: suffix: 7_cyclic
257: output_file: output/ex45_7_cyclic.out
259: test:
260: args: -m 25 -n 20 -p 21 -svd_smallest -svd_nsv 2 -svd_ncv 5 -svd_type trlanczos -svd_trlanczos_gbidiag {{upper lower}} -svd_trlanczos_scale {{0.1 -20}}
261: filter: grep -v "Solution method" | grep -v "Number of iterations" | grep -v "Stopping condition" | sed -e "s/, maxit=1[0]*$//" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
262: suffix: 8
264: TEST*/