Actual source code: test6.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[] = "SVD via the cross-product matrix with a user-provided EPS.\n\n"
12: "The command line options are:\n"
13: " -m <m>, where <m> = matrix rows.\n"
14: " -n <n>, where <n> = matrix columns (defaults to m+2).\n\n";
16: #include <slepcsvd.h>
18: /*
19: This example computes the singular values of a rectangular bidiagonal matrix
21: | 1 2 |
22: | 1 2 |
23: | 1 2 |
24: A = | . . |
25: | . . |
26: | 1 2 |
27: | 1 2 |
28: */
30: int main(int argc,char **argv)
31: {
32: Mat A;
33: SVD svd;
34: EPS eps;
35: ST st;
36: KSP ksp;
37: PC pc;
38: PetscInt m=20,n,Istart,Iend,i,col[2];
39: PetscScalar value[] = { 1, 2 };
40: PetscBool flg;
43: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
44: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
46: if (!flg) n=m+2;
47: PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%D n=%D\n\n",m,n);
49: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50: Generate the matrix
51: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: MatCreate(PETSC_COMM_WORLD,&A);
54: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
55: MatSetFromOptions(A);
56: MatSetUp(A);
57: MatGetOwnershipRange(A,&Istart,&Iend);
58: for (i=Istart;i<Iend;i++) {
59: col[0]=i; col[1]=i+1;
60: if (i<n-1) {
61: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
62: } else if (i==n-1) {
63: MatSetValue(A,i,col[0],value[0],INSERT_VALUES);
64: }
65: }
66: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
67: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
69: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: Create a standalone EPS with appropriate settings
71: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: EPSCreate(PETSC_COMM_WORLD,&eps);
74: EPSGetST(eps,&st);
75: STSetType(st,STSINVERT);
76: STGetKSP(st,&ksp);
77: KSPSetType(ksp,KSPBCGS);
78: KSPGetPC(ksp,&pc);
79: PCSetType(pc,PCJACOBI);
80: EPSSetFromOptions(eps);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Compute singular values
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: SVDCreate(PETSC_COMM_WORLD,&svd);
87: SVDSetOperator(svd,A);
88: SVDSetType(svd,SVDCROSS);
89: SVDCrossSetEPS(svd,eps);
90: SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
91: SVDSetFromOptions(svd);
92: SVDSolve(svd);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Display solution and clean up
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
98: SVDDestroy(&svd);
99: EPSDestroy(&eps);
100: MatDestroy(&A);
101: SlepcFinalize();
102: return ierr;
103: }
105: /*TEST
107: testset:
108: output_file: output/test6_1.out
109: test:
110: suffix: 1_subspace
111: args: -eps_type subspace
112: test:
113: suffix: 1_lobpcg
114: args: -eps_type lobpcg -st_type precond
115: test:
116: suffix: 2_cuda
117: args: -eps_type subspace -mat_type aijcusparse
118: requires: cuda
120: TEST*/