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