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[] = "Tests multiple calls to SVDSolve changing ncv.\n\n"
12 : "The command line options are:\n"
13 : " -n <n>, where <n> = matrix dimension.\n\n";
14 :
15 : #include <slepcsvd.h>
16 :
17 : /*
18 : This example computes the singular values of an nxn Grcar matrix,
19 : which is a nonsymmetric Toeplitz matrix:
20 :
21 : | 1 1 1 1 |
22 : | -1 1 1 1 1 |
23 : | -1 1 1 1 1 |
24 : | . . . . . |
25 : A = | . . . . . |
26 : | -1 1 1 1 1 |
27 : | -1 1 1 1 |
28 : | -1 1 1 |
29 : | -1 1 |
30 :
31 : */
32 :
33 6 : int main(int argc,char **argv)
34 : {
35 6 : Mat A;
36 6 : SVD svd;
37 6 : PetscInt N=30,Istart,Iend,i,col[5],nsv,ncv;
38 6 : PetscScalar value[] = { -1, 1, 1, 1, 1 };
39 :
40 6 : PetscFunctionBeginUser;
41 6 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
42 6 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
43 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSingular values of a Grcar matrix, n=%" PetscInt_FMT,N));
44 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n\n"));
45 :
46 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47 : Generate the matrix
48 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49 :
50 6 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
51 6 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
52 6 : PetscCall(MatSetFromOptions(A));
53 :
54 6 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
55 186 : for (i=Istart;i<Iend;i++) {
56 180 : col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
57 180 : if (i==0) PetscCall(MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES));
58 180 : else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
59 : }
60 :
61 6 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
62 6 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
63 :
64 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65 : Create the singular value solver and set the solution method
66 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67 :
68 6 : PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
69 6 : PetscCall(SVDSetOperators(svd,A,NULL));
70 6 : PetscCall(SVDSetTolerances(svd,1e-6,1000));
71 6 : PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
72 6 : PetscCall(SVDSetFromOptions(svd));
73 :
74 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75 : Compute the singular values
76 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77 :
78 : /* First solve */
79 6 : PetscCall(SVDSolve(svd));
80 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," - - - First solve, default subspace dimension - - -\n"));
81 6 : PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));
82 :
83 : /* Second solve */
84 6 : PetscCall(SVDGetDimensions(svd,&nsv,&ncv,NULL));
85 6 : PetscCall(SVDSetDimensions(svd,nsv,ncv+2,PETSC_DETERMINE));
86 6 : PetscCall(SVDSolve(svd));
87 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," - - - Second solve, subspace of increased size - - -\n"));
88 6 : PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));
89 :
90 : /* Free work space */
91 6 : PetscCall(SVDDestroy(&svd));
92 6 : PetscCall(MatDestroy(&A));
93 6 : PetscCall(SlepcFinalize());
94 : return 0;
95 : }
96 :
97 : /*TEST
98 :
99 : test:
100 : suffix: 1
101 : args: -svd_type {{lanczos trlanczos cross cyclic lapack randomized}} -svd_nsv 3 -svd_ncv 12
102 : requires: !single
103 :
104 : TEST*/
|