Actual source code: ex45.c
slepc-3.17.1 2022-04-11
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;
30: SlepcInitialize(&argc,&argv,(char*)0,help);
32: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
33: PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
34: if (!flg) n = m;
35: PetscOptionsGetInt(NULL,NULL,"-p",&p,&flg);
36: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%" PetscInt_FMT "+%" PetscInt_FMT ")x%" PetscInt_FMT "\n\n",m,p,n);
37: PetscOptionsGetBool(NULL,NULL,"-skiporth",&skiporth,NULL);
39: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: Build the matrices
41: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
43: MatCreate(PETSC_COMM_WORLD,&A);
44: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
45: MatSetFromOptions(A);
46: MatSetUp(A);
48: MatGetOwnershipRange(A,&Istart,&Iend);
49: for (i=Istart;i<Iend;i++) {
50: if (i>0 && i-1<n) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
51: if (i+1<n) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
52: if (i<n) MatSetValue(A,i,i,2.0,INSERT_VALUES);
53: if (i>n) MatSetValue(A,i,n-1,1.0,INSERT_VALUES);
54: }
55: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
56: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
58: MatCreate(PETSC_COMM_WORLD,&B);
59: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
60: MatSetFromOptions(B);
61: MatSetUp(B);
63: MatGetOwnershipRange(B,&Istart,&Iend);
64: d = PetscMax(0,n-p);
65: for (i=Istart;i<Iend;i++) {
66: for (j=0;j<=PetscMin(i,n-1);j++) MatSetValue(B,i,j+d,1.0,INSERT_VALUES);
67: }
68: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
69: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Create the singular value solver and set various options
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: /*
76: Create singular value solver context
77: */
78: SVDCreate(PETSC_COMM_WORLD,&svd);
80: /*
81: Set operators and problem type
82: */
83: SVDSetOperators(svd,A,B);
84: SVDSetProblemType(svd,SVD_GENERALIZED);
86: /*
87: Set solver parameters at runtime
88: */
89: SVDSetFromOptions(svd);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Solve the singular value system
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: SVDSolve(svd);
96: SVDGetIterationNumber(svd,&its);
97: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its);
99: /*
100: Optional: Get some information from the solver and display it
101: */
102: SVDGetType(svd,&type);
103: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
104: SVDGetDimensions(svd,&nsv,NULL,NULL);
105: PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %" PetscInt_FMT "\n",nsv);
106: SVDGetTolerances(svd,&tol,&maxit);
107: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit);
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Display solution and clean up
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: /*
114: Get number of converged singular triplets
115: */
116: SVDGetConverged(svd,&nconv);
117: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %" PetscInt_FMT "\n\n",nconv);
119: if (nconv>0) {
120: /*
121: Create vectors. The interface returns u and v as stacked on top of each other
122: [u;v] so need to create a special vector (VecNest) to extract them
123: */
124: MatCreateVecs(A,&x,&u);
125: MatCreateVecs(B,NULL,&v);
126: aux[0] = u;
127: aux[1] = v;
128: VecCreateNest(PETSC_COMM_WORLD,2,NULL,aux,&uv);
130: VecDuplicateVecs(u,nconv,&U);
131: VecDuplicateVecs(v,nconv,&V);
133: /*
134: Display singular values and errors relative to the norms
135: */
136: PetscCall(PetscPrintf(PETSC_COMM_WORLD,
137: " sigma ||r||/||[A;B]||\n"
138: " --------------------- ------------------\n"));
139: for (i=0;i<nconv;i++) {
140: /*
141: Get converged singular triplets: i-th singular value is stored in sigma
142: */
143: SVDGetSingularTriplet(svd,i,&sigma,uv,x);
145: /* at this point, u and v can be used normally as individual vectors */
146: VecCopy(u,U[i]);
147: VecCopy(v,V[i]);
149: /*
150: Compute the error associated to each singular triplet
151: */
152: SVDComputeError(svd,i,SVD_ERROR_NORM,&error);
154: PetscPrintf(PETSC_COMM_WORLD," % 6f ",(double)sigma);
155: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",(double)error);
156: }
157: PetscPrintf(PETSC_COMM_WORLD,"\n");
159: if (!skiporth) {
160: VecCheckOrthonormality(U,nconv,NULL,nconv,NULL,NULL,&lev1);
161: VecCheckOrthonormality(V,nconv,NULL,nconv,NULL,NULL,&lev2);
162: }
163: if (lev1+lev2<20*tol) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
164: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g (U) %g (V)\n",(double)lev1,(double)lev2);
165: VecDestroyVecs(nconv,&U);
166: VecDestroyVecs(nconv,&V);
167: VecDestroy(&x);
168: VecDestroy(&u);
169: VecDestroy(&v);
170: VecDestroy(&uv);
171: }
173: /*
174: Free work space
175: */
176: SVDDestroy(&svd);
177: MatDestroy(&A);
178: MatDestroy(&B);
179: SlepcFinalize();
180: return 0;
181: }
183: /*TEST
185: testset:
186: 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"
187: requires: double
188: test:
189: args: -svd_type lapack -m 20 -n 10 -p 6
190: suffix: 1
191: test:
192: args: -svd_type lapack -m 15 -n 20 -p 10 -svd_smallest
193: suffix: 2
194: test:
195: args: -svd_type lapack -m 15 -n 20 -p 21
196: suffix: 3
197: test:
198: args: -svd_type lapack -m 20 -n 15 -p 21
199: suffix: 4
201: testset:
202: args: -m 25 -n 20 -p 21 -svd_smallest -svd_nsv 4
203: 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"
204: requires: double
205: output_file: output/ex45_5.out
206: test:
207: args: -svd_type trlanczos -svd_ncv 8 -svd_trlanczos_gbidiag {{upper lower}}
208: suffix: 5
209: test:
210: args: -svd_type cross -svd_ncv 10 -svd_cross_explicitmatrix {{0 1}}
211: suffix: 5_cross
212: test:
213: args: -svd_type cyclic -svd_ncv 12 -svd_cyclic_explicitmatrix {{0 1}}
214: suffix: 5_cyclic
215: requires: !complex
217: testset:
218: args: -m 15 -n 20 -p 21 -svd_nsv 4 -svd_ncv 9
219: 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"
220: requires: double
221: output_file: output/ex45_6.out
222: test:
223: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}} -svd_trlanczos_locking {{0 1}}
224: suffix: 6
225: test:
226: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
227: suffix: 6_cross
228: test:
229: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
230: suffix: 6_cyclic
232: testset:
233: args: -m 20 -n 15 -p 21 -svd_nsv 4 -svd_ncv 9
234: 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"
235: requires: double
236: output_file: output/ex45_7.out
237: test:
238: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}} -svd_trlanczos_restart 0.4
239: suffix: 7
240: test:
241: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
242: suffix: 7_cross
243: test:
244: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
245: suffix: 7_cyclic
247: TEST*/