Actual source code: test6.c
slepc-3.21.2 2024-09-25
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;
42: PetscFunctionBeginUser;
43: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
44: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
45: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg));
46: if (!flg) n=m+2;
47: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%" PetscInt_FMT " n=%" PetscInt_FMT "\n\n",m,n));
49: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50: Generate the matrix
51: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
54: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n));
55: PetscCall(MatSetFromOptions(A));
56: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
57: for (i=Istart;i<Iend;i++) {
58: col[0]=i; col[1]=i+1;
59: if (i<n-1) PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
60: else if (i==n-1) PetscCall(MatSetValue(A,i,col[0],value[0],INSERT_VALUES));
61: }
62: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
63: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Create a standalone EPS with appropriate settings
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
70: PetscCall(EPSGetST(eps,&st));
71: PetscCall(STSetType(st,STSINVERT));
72: PetscCall(STGetKSP(st,&ksp));
73: PetscCall(KSPSetType(ksp,KSPBCGS));
74: PetscCall(KSPGetPC(ksp,&pc));
75: PetscCall(PCSetType(pc,PCJACOBI));
76: PetscCall(EPSSetFromOptions(eps));
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Compute singular values
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
83: PetscCall(SVDSetOperators(svd,A,NULL));
84: PetscCall(SVDSetType(svd,SVDCROSS));
85: PetscCall(SVDCrossSetEPS(svd,eps));
86: PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
87: PetscCall(SVDSetFromOptions(svd));
88: PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDCROSS,&flg));
89: if (flg) {
90: PetscCall(SVDCrossGetExplicitMatrix(svd,&expmat));
91: if (expmat) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using explicit matrix with cross solver\n"));
92: }
93: PetscCall(SVDSolve(svd));
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Display solution and clean up
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,PETSC_VIEWER_STDOUT_WORLD));
99: PetscCall(SVDDestroy(&svd));
100: PetscCall(EPSDestroy(&eps));
101: PetscCall(MatDestroy(&A));
102: PetscCall(SlepcFinalize());
103: return 0;
104: }
106: /*TEST
108: testset:
109: output_file: output/test6_1.out
110: test:
111: suffix: 1_subspace
112: args: -eps_type subspace
113: test:
114: suffix: 1_lobpcg
115: args: -eps_type lobpcg -st_type precond
116: test:
117: suffix: 2_cuda
118: args: -eps_type subspace -mat_type aijcusparse
119: requires: cuda
121: TEST*/