Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test RSVD on a low-rank matrix.\n\n";
12 :
13 : #include <slepcsvd.h>
14 :
15 1 : int main(int argc,char **argv)
16 : {
17 1 : Mat A,Ur,Vr;
18 1 : SVD svd;
19 1 : PetscInt m=80,n=40,rank=3,Istart,Iend,i,j;
20 1 : PetscScalar *u;
21 1 : PetscReal tol=PETSC_SMALL;
22 :
23 1 : PetscFunctionBeginUser;
24 1 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
25 :
26 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
27 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
28 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-rank",&rank,NULL));
29 1 : PetscCheck(rank>0,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"The rank must be >=1");
30 1 : PetscCheck(rank<PetscMin(m,n),PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"The rank must be <min(m,n)");
31 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSVD of low-rank matrix, size %" PetscInt_FMT "x%" PetscInt_FMT ", rank %" PetscInt_FMT "\n\n",m,n,rank));
32 :
33 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34 : Create a low-rank matrix A = Ur*Vr'
35 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36 :
37 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&Ur));
38 1 : PetscCall(MatSetSizes(Ur,PETSC_DECIDE,PETSC_DECIDE,m,rank));
39 1 : PetscCall(MatSetType(Ur,MATDENSE));
40 1 : PetscCall(MatGetOwnershipRange(Ur,&Istart,&Iend));
41 1 : PetscCall(MatDenseGetArray(Ur,&u));
42 81 : for (i=Istart;i<Iend;i++) {
43 80 : if (i<m/2) u[i-Istart] = 1.0;
44 80 : if (i==0) u[i+Iend-2*Istart] = -2.0;
45 80 : if (i==1) u[i+Iend-2*Istart] = -1.0;
46 80 : if (i==2) u[i+Iend-2*Istart] = -1.0;
47 160 : for (j=2;j<rank;j++) if (i==j) u[j+j*Iend-(j+1)*Istart] = j;
48 : }
49 1 : PetscCall(MatDenseRestoreArray(Ur,&u));
50 1 : PetscCall(MatAssemblyBegin(Ur,MAT_FINAL_ASSEMBLY));
51 1 : PetscCall(MatAssemblyEnd(Ur,MAT_FINAL_ASSEMBLY));
52 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&Vr));
53 1 : PetscCall(MatSetSizes(Vr,PETSC_DECIDE,PETSC_DECIDE,n,rank));
54 1 : PetscCall(MatSetType(Vr,MATDENSE));
55 1 : PetscCall(MatGetOwnershipRange(Vr,&Istart,&Iend));
56 1 : PetscCall(MatDenseGetArray(Vr,&u));
57 41 : for (i=Istart;i<Iend;i++) {
58 40 : if (i>n/2) u[i-Istart] = 1.0;
59 40 : if (i==0) u[i+Iend-2*Istart] = 1.0;
60 40 : if (i==1) u[i+Iend-2*Istart] = 2.0;
61 40 : if (i==2) u[i+Iend-2*Istart] = 2.0;
62 80 : for (j=2;j<rank;j++) if (i==j) u[j+j*Iend-(j+1)*Istart] = j;
63 : }
64 1 : PetscCall(MatDenseRestoreArray(Vr,&u));
65 1 : PetscCall(MatAssemblyBegin(Vr,MAT_FINAL_ASSEMBLY));
66 1 : PetscCall(MatAssemblyEnd(Vr,MAT_FINAL_ASSEMBLY));
67 1 : PetscCall(MatCreateLRC(NULL,Ur,NULL,Vr,&A));
68 1 : PetscCall(MatDestroy(&Ur));
69 1 : PetscCall(MatDestroy(&Vr));
70 :
71 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72 : Create the SVD solver
73 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74 :
75 1 : PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
76 1 : PetscCall(SVDSetOperators(svd,A,NULL));
77 1 : PetscCall(SVDSetType(svd,SVDRANDOMIZED));
78 1 : PetscCall(SVDSetDimensions(svd,2,PETSC_DETERMINE,PETSC_DETERMINE));
79 1 : PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_MAXIT));
80 1 : PetscCall(SVDSetTolerances(svd,tol,1)); /* maxit=1 to disable outer iteration */
81 1 : PetscCall(SVDSetFromOptions(svd));
82 :
83 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84 : Compute the singular triplets and display solution
85 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86 :
87 1 : PetscCall(SVDSolve(svd));
88 1 : PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));
89 1 : PetscCall(SVDDestroy(&svd));
90 1 : PetscCall(MatDestroy(&A));
91 1 : PetscCall(SlepcFinalize());
92 : return 0;
93 : }
94 :
95 : /*TEST
96 :
97 : test:
98 : args: -bv_orthog_block tsqr # currently fail with other block orthogonalization methods
99 : requires: !single
100 :
101 : TEST*/
|