Actual source code: test6.c
slepc-3.18.0 2022-10-01
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[] = "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,expmat;
43: SlepcInitialize(&argc,&argv,(char*)0,help);
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=%" PetscInt_FMT " n=%" PetscInt_FMT "\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) MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
61: else if (i==n-1) MatSetValue(A,i,col[0],value[0],INSERT_VALUES);
62: }
63: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
64: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Create a standalone EPS with appropriate settings
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: EPSCreate(PETSC_COMM_WORLD,&eps);
71: EPSGetST(eps,&st);
72: STSetType(st,STSINVERT);
73: STGetKSP(st,&ksp);
74: KSPSetType(ksp,KSPBCGS);
75: KSPGetPC(ksp,&pc);
76: PCSetType(pc,PCJACOBI);
77: EPSSetFromOptions(eps);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Compute singular values
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: SVDCreate(PETSC_COMM_WORLD,&svd);
84: SVDSetOperators(svd,A,NULL);
85: SVDSetType(svd,SVDCROSS);
86: SVDCrossSetEPS(svd,eps);
87: SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
88: SVDSetFromOptions(svd);
89: PetscObjectTypeCompare((PetscObject)svd,SVDCROSS,&flg);
90: if (flg) {
91: SVDCrossGetExplicitMatrix(svd,&expmat);
92: if (expmat) PetscPrintf(PETSC_COMM_WORLD," Using explicit matrix with cross solver\n");
93: }
94: SVDSolve(svd);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Display solution and clean up
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
100: SVDDestroy(&svd);
101: EPSDestroy(&eps);
102: MatDestroy(&A);
103: SlepcFinalize();
104: return 0;
105: }
107: /*TEST
109: testset:
110: output_file: output/test6_1.out
111: test:
112: suffix: 1_subspace
113: args: -eps_type subspace
114: test:
115: suffix: 1_lobpcg
116: args: -eps_type lobpcg -st_type precond
117: test:
118: suffix: 2_cuda
119: args: -eps_type subspace -mat_type aijcusparse
120: requires: cuda
122: TEST*/