Actual source code: test19.c
slepc-3.21.1 2024-04-26
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[] = "Test RSVD on a low-rank matrix.\n\n";
13: #include <slepcsvd.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,Ur,Vr;
18: SVD svd;
19: PetscInt m=80,n=40,rank=3,Istart,Iend,i,j;
20: PetscScalar *u;
21: PetscReal tol=PETSC_SMALL;
23: PetscFunctionBeginUser;
24: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
26: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
27: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
28: PetscCall(PetscOptionsGetInt(NULL,NULL,"-rank",&rank,NULL));
29: PetscCheck(rank>0,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"The rank must be >=1");
30: PetscCheck(rank<PetscMin(m,n),PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"The rank must be <min(m,n)");
31: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSVD of low-rank matrix, size %" PetscInt_FMT "x%" PetscInt_FMT ", rank %" PetscInt_FMT "\n\n",m,n,rank));
33: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: Create a low-rank matrix A = Ur*Vr'
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37: PetscCall(MatCreate(PETSC_COMM_WORLD,&Ur));
38: PetscCall(MatSetSizes(Ur,PETSC_DECIDE,PETSC_DECIDE,m,rank));
39: PetscCall(MatSetType(Ur,MATDENSE));
40: PetscCall(MatGetOwnershipRange(Ur,&Istart,&Iend));
41: PetscCall(MatDenseGetArray(Ur,&u));
42: for (i=Istart;i<Iend;i++) {
43: if (i<m/2) u[i-Istart] = 1.0;
44: if (i==0) u[i+Iend-2*Istart] = -2.0;
45: if (i==1) u[i+Iend-2*Istart] = -1.0;
46: if (i==2) u[i+Iend-2*Istart] = -1.0;
47: for (j=2;j<rank;j++) if (i==j) u[j+j*Iend-(j+1)*Istart] = j;
48: }
49: PetscCall(MatDenseRestoreArray(Ur,&u));
50: PetscCall(MatAssemblyBegin(Ur,MAT_FINAL_ASSEMBLY));
51: PetscCall(MatAssemblyEnd(Ur,MAT_FINAL_ASSEMBLY));
52: PetscCall(MatCreate(PETSC_COMM_WORLD,&Vr));
53: PetscCall(MatSetSizes(Vr,PETSC_DECIDE,PETSC_DECIDE,n,rank));
54: PetscCall(MatSetType(Vr,MATDENSE));
55: PetscCall(MatGetOwnershipRange(Vr,&Istart,&Iend));
56: PetscCall(MatDenseGetArray(Vr,&u));
57: for (i=Istart;i<Iend;i++) {
58: if (i>n/2) u[i-Istart] = 1.0;
59: if (i==0) u[i+Iend-2*Istart] = 1.0;
60: if (i==1) u[i+Iend-2*Istart] = 2.0;
61: if (i==2) u[i+Iend-2*Istart] = 2.0;
62: for (j=2;j<rank;j++) if (i==j) u[j+j*Iend-(j+1)*Istart] = j;
63: }
64: PetscCall(MatDenseRestoreArray(Vr,&u));
65: PetscCall(MatAssemblyBegin(Vr,MAT_FINAL_ASSEMBLY));
66: PetscCall(MatAssemblyEnd(Vr,MAT_FINAL_ASSEMBLY));
67: PetscCall(MatCreateLRC(NULL,Ur,NULL,Vr,&A));
68: PetscCall(MatDestroy(&Ur));
69: PetscCall(MatDestroy(&Vr));
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Create the SVD solver
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
76: PetscCall(SVDSetOperators(svd,A,NULL));
77: PetscCall(SVDSetType(svd,SVDRANDOMIZED));
78: PetscCall(SVDSetDimensions(svd,2,PETSC_DEFAULT,PETSC_DEFAULT));
79: PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_MAXIT));
80: PetscCall(SVDSetTolerances(svd,tol,1)); /* maxit=1 to disable outer iteration */
81: PetscCall(SVDSetFromOptions(svd));
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Compute the singular triplets and display solution
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: PetscCall(SVDSolve(svd));
88: PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));
89: PetscCall(SVDDestroy(&svd));
90: PetscCall(MatDestroy(&A));
91: PetscCall(SlepcFinalize());
92: return 0;
93: }
95: /*TEST
97: test:
98: args: -bv_orthog_block tsqr # currently fail with other block orthogonalization methods
99: requires: !single
101: TEST*/