Actual source code: test7.c
slepc-3.19.2 2023-09-05
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 cyclic 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));
45: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
46: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg));
47: if (!flg) n=m+2;
48: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%" PetscInt_FMT " n=%" PetscInt_FMT "\n\n",m,n));
50: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: Generate the matrix
52: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
55: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n));
56: PetscCall(MatSetFromOptions(A));
57: PetscCall(MatSetUp(A));
58: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
59: for (i=Istart;i<Iend;i++) {
60: col[0]=i; col[1]=i+1;
61: if (i<n-1) PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
62: else if (i==n-1) PetscCall(MatSetValue(A,i,col[0],value[0],INSERT_VALUES));
63: }
64: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
65: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Create a standalone EPS with appropriate settings
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
72: PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
73: PetscCall(EPSSetTarget(eps,1.0));
74: PetscCall(EPSGetST(eps,&st));
75: PetscCall(STSetType(st,STSINVERT));
76: PetscCall(STSetShift(st,1.01));
77: PetscCall(STGetKSP(st,&ksp));
78: PetscCall(KSPSetType(ksp,KSPPREONLY));
79: PetscCall(KSPGetPC(ksp,&pc));
80: PetscCall(PCSetType(pc,PCLU));
81: PetscCall(EPSSetFromOptions(eps));
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Compute singular values
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
88: PetscCall(SVDSetOperators(svd,A,NULL));
89: PetscCall(SVDSetType(svd,SVDCYCLIC));
90: PetscCall(SVDCyclicSetEPS(svd,eps));
91: PetscCall(SVDCyclicSetExplicitMatrix(svd,PETSC_TRUE));
92: PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
93: PetscCall(SVDSetFromOptions(svd));
94: PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDCYCLIC,&flg));
95: if (flg) {
96: PetscCall(SVDCyclicGetExplicitMatrix(svd,&expmat));
97: if (expmat) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using explicit matrix with cyclic solver\n"));
98: }
99: PetscCall(SVDSolve(svd));
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Display solution and clean up
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
105: PetscCall(SVDDestroy(&svd));
106: PetscCall(EPSDestroy(&eps));
107: PetscCall(MatDestroy(&A));
108: PetscCall(SlepcFinalize());
109: return 0;
110: }
112: /*TEST
114: test:
115: suffix: 1
116: args: -log_exclude svd
118: test:
119: suffix: 2_cuda
120: args: -log_exclude svd -mat_type aijcusparse
121: requires: cuda
122: output_file: output/test7_1.out
124: TEST*/