Actual source code: test10.c
slepc-3.20.2 2024-03-15
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[] = "Lanczos SVD. Also illustrates the use of SVDSetBV().\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: PetscInt m=20,n,Istart,Iend,i,k=6,col[2];
35: PetscScalar value[] = { 1, 2 };
36: PetscBool flg,oneside=PETSC_FALSE;
37: const char *prefix;
38: BV U,V;
39: Vec u,v;
41: PetscFunctionBeginUser;
42: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
43: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
44: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg));
45: if (!flg) n=m+2;
46: PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
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(MatSetUp(A));
57: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
58: for (i=Istart;i<Iend;i++) {
59: col[0]=i; col[1]=i+1;
60: if (i<n-1) PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
61: else if (i==n-1) PetscCall(MatSetValue(A,i,col[0],value[0],INSERT_VALUES));
62: }
63: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
64: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
65: PetscCall(MatCreateVecs(A,&v,&u));
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Create standalone BV objects to illustrate use of SVDSetBV()
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: PetscCall(BVCreate(PETSC_COMM_WORLD,&U));
72: PetscCall(PetscObjectSetName((PetscObject)U,"U"));
73: PetscCall(BVSetSizesFromVec(U,u,k));
74: PetscCall(BVSetFromOptions(U));
75: PetscCall(BVCreate(PETSC_COMM_WORLD,&V));
76: PetscCall(PetscObjectSetName((PetscObject)V,"V"));
77: PetscCall(BVSetSizesFromVec(V,v,k));
78: PetscCall(BVSetFromOptions(V));
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Compute singular values
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
85: PetscCall(SVDSetBV(svd,V,U));
86: PetscCall(SVDSetOptionsPrefix(svd,"check_"));
87: PetscCall(SVDAppendOptionsPrefix(svd,"myprefix_"));
88: PetscCall(SVDGetOptionsPrefix(svd,&prefix));
89: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"SVD prefix is currently: %s\n\n",prefix));
90: PetscCall(PetscObjectSetName((PetscObject)svd,"SVD_solver"));
92: PetscCall(SVDSetOperators(svd,A,NULL));
93: PetscCall(SVDSetType(svd,SVDLANCZOS));
94: PetscCall(SVDSetFromOptions(svd));
96: PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDLANCZOS,&flg));
97: if (flg) {
98: PetscCall(SVDLanczosGetOneSide(svd,&oneside));
99: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Running Lanczos %s\n\n",oneside?"(onesided)":""));
100: }
101: PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDTRLANCZOS,&flg));
102: if (flg) {
103: PetscCall(SVDTRLanczosGetOneSide(svd,&oneside));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Running thick-restart Lanczos %s\n\n",oneside?"(onesided)":""));
105: }
107: PetscCall(SVDSolve(svd));
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Display solution and clean up
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
113: PetscCall(BVDestroy(&U));
114: PetscCall(BVDestroy(&V));
115: PetscCall(VecDestroy(&u));
116: PetscCall(VecDestroy(&v));
117: PetscCall(SVDDestroy(&svd));
118: PetscCall(MatDestroy(&A));
119: PetscCall(SlepcFinalize());
120: return 0;
121: }
123: /*TEST
125: testset:
126: args: -check_myprefix_svd_nsv 3
127: requires: double
128: test:
129: suffix: 1
130: args: -check_myprefix_svd_view_vectors ::ascii_info
131: test:
132: suffix: 2
133: args: -check_myprefix_svd_type trlanczos -check_myprefix_svd_monitor -check_myprefix_svd_view_values ::ascii_matlab
134: filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
135: test:
136: suffix: 3
137: args: -m 22 -n 20
139: TEST*/