Actual source code: ex45.c
slepc-3.16.0 2021-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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);if (ierr) return ierr;
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, (%D+%D)x%D\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++) {
68: MatSetValue(B,i,j+d,1.0,INSERT_VALUES);
69: }
70: }
71: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
72: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
74: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: Create the singular value solver and set various options
76: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: /*
79: Create singular value solver context
80: */
81: SVDCreate(PETSC_COMM_WORLD,&svd);
83: /*
84: Set operators and problem type
85: */
86: SVDSetOperators(svd,A,B);
87: SVDSetProblemType(svd,SVD_GENERALIZED);
89: /*
90: Set solver parameters at runtime
91: */
92: SVDSetFromOptions(svd);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Solve the singular value system
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: SVDSolve(svd);
99: SVDGetIterationNumber(svd,&its);
100: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
102: /*
103: Optional: Get some information from the solver and display it
104: */
105: SVDGetType(svd,&type);
106: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
107: SVDGetDimensions(svd,&nsv,NULL,NULL);
108: PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %D\n",nsv);
109: SVDGetTolerances(svd,&tol,&maxit);
110: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%D\n",(double)tol,maxit);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Display solution and clean up
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: /*
117: Get number of converged singular triplets
118: */
119: SVDGetConverged(svd,&nconv);
120: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %D\n\n",nconv);
122: if (nconv>0) {
123: /*
124: Create vectors. The interface returns u and v as stacked on top of each other
125: [u;v] so need to create a special vector (VecNest) to extract them
126: */
127: MatCreateVecs(A,&x,&u);
128: MatCreateVecs(B,NULL,&v);
129: aux[0] = u;
130: aux[1] = v;
131: VecCreateNest(PETSC_COMM_WORLD,2,NULL,aux,&uv);
133: VecDuplicateVecs(u,nconv,&U);
134: VecDuplicateVecs(v,nconv,&V);
136: /*
137: Display singular values and relative errors
138: */
139: PetscPrintf(PETSC_COMM_WORLD,
140: " sigma relative error\n"
141: " --------------------- ------------------\n");
142: for (i=0;i<nconv;i++) {
143: /*
144: Get converged singular triplets: i-th singular value is stored in sigma
145: */
146: SVDGetSingularTriplet(svd,i,&sigma,uv,x);
148: /* at this point, u and v can be used normally as individual vectors */
149: VecCopy(u,U[i]);
150: VecCopy(v,V[i]);
152: /*
153: Compute the error associated to each singular triplet
154: */
155: SVDComputeError(svd,i,SVD_ERROR_RELATIVE,&error);
157: PetscPrintf(PETSC_COMM_WORLD," % 6f ",(double)sigma);
158: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",(double)error);
159: }
160: PetscPrintf(PETSC_COMM_WORLD,"\n");
162: if (!skiporth) {
163: VecCheckOrthonormality(U,nconv,NULL,nconv,NULL,NULL,&lev1);
164: VecCheckOrthonormality(V,nconv,NULL,nconv,NULL,NULL,&lev2);
165: }
166: if (lev1+lev2<20*tol) {
167: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
168: } else {
169: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g (U) %g (V)\n",(double)lev1,(double)lev2);
170: }
171: VecDestroyVecs(nconv,&U);
172: VecDestroyVecs(nconv,&V);
173: VecDestroy(&x);
174: VecDestroy(&u);
175: VecDestroy(&v);
176: VecDestroy(&uv);
177: }
179: /*
180: Free work space
181: */
182: SVDDestroy(&svd);
183: MatDestroy(&A);
184: MatDestroy(&B);
185: SlepcFinalize();
186: return ierr;
187: }
189: /*TEST
191: testset:
192: 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"
193: requires: double
194: test:
195: args: -svd_type lapack -m 20 -n 10 -p 6
196: suffix: 1
197: test:
198: args: -svd_type lapack -m 15 -n 20 -p 10 -svd_smallest
199: suffix: 2
200: test:
201: args: -svd_type lapack -m 15 -n 20 -p 21
202: suffix: 3
203: test:
204: args: -svd_type lapack -m 20 -n 15 -p 21
205: suffix: 4
207: testset:
208: args: -m 25 -n 20 -p 21 -svd_smallest -svd_nsv 4
209: 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"
210: requires: double
211: output_file: output/ex45_5.out
212: test:
213: args: -svd_type trlanczos -svd_ncv 8 -svd_trlanczos_gbidiag {{upper lower}}
214: suffix: 5
215: test:
216: args: -svd_type cross -svd_ncv 10 -svd_cross_explicitmatrix {{0 1}}
217: suffix: 5_cross
218: test:
219: args: -svd_type cyclic -svd_ncv 12 -svd_cyclic_explicitmatrix {{0 1}}
220: suffix: 5_cyclic
221: requires: !complex
223: testset:
224: args: -m 15 -n 20 -p 21 -svd_nsv 4 -svd_ncv 9
225: 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"
226: requires: double
227: output_file: output/ex45_6.out
228: test:
229: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}} -svd_trlanczos_locking {{0 1}}
230: suffix: 6
231: test:
232: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
233: suffix: 6_cross
234: test:
235: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
236: suffix: 6_cyclic
238: testset:
239: args: -m 20 -n 15 -p 21 -svd_nsv 4 -svd_ncv 9
240: 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"
241: requires: double
242: output_file: output/ex45_7.out
243: test:
244: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}} -svd_trlanczos_restart 0.4
245: suffix: 7
246: test:
247: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
248: suffix: 7_cross
249: test:
250: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
251: suffix: 7_cyclic
253: TEST*/