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 with different matrix size.\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,B;
36 6 : SVD svd;
37 6 : PetscInt N=30,Istart,Iend,i,col[5];
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\n",N));
44 :
45 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46 : Generate the matrix of size N
47 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48 :
49 6 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
50 6 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
51 6 : PetscCall(MatSetFromOptions(A));
52 6 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
53 186 : for (i=Istart;i<Iend;i++) {
54 180 : col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
55 180 : if (i==0) PetscCall(MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES));
56 180 : else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
57 : }
58 6 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
59 6 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
60 :
61 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62 : Create the singular value solver, set options and solve
63 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64 :
65 6 : PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
66 6 : PetscCall(SVDSetOperators(svd,A,NULL));
67 6 : PetscCall(SVDSetTolerances(svd,1e-6,1000));
68 6 : PetscCall(SVDSetFromOptions(svd));
69 6 : PetscCall(SVDSolve(svd));
70 6 : PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));
71 :
72 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73 : Generate the matrix of size 2*N
74 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75 :
76 6 : N *= 2;
77 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSingular values of a Grcar matrix, n=%" PetscInt_FMT "\n\n",N));
78 :
79 6 : PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
80 6 : PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
81 6 : PetscCall(MatSetFromOptions(B));
82 6 : PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
83 366 : for (i=Istart;i<Iend;i++) {
84 360 : col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
85 360 : if (i==0) PetscCall(MatSetValues(B,1,&i,4,col+1,value+1,INSERT_VALUES));
86 360 : else PetscCall(MatSetValues(B,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
87 : }
88 6 : PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
89 6 : PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
90 :
91 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92 : Solve again, calling SVDReset() since matrix size has changed
93 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94 :
95 6 : PetscCall(SVDReset(svd)); /* if this is omitted, it will be called in SVDSetOperators() */
96 6 : PetscCall(SVDSetOperators(svd,B,NULL));
97 6 : PetscCall(SVDSolve(svd));
98 6 : PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL));
99 :
100 : /* Free work space */
101 6 : PetscCall(SVDDestroy(&svd));
102 6 : PetscCall(MatDestroy(&A));
103 6 : PetscCall(MatDestroy(&B));
104 6 : PetscCall(SlepcFinalize());
105 : return 0;
106 : }
107 :
108 : /*TEST
109 :
110 : test:
111 : suffix: 1
112 : args: -svd_type {{lanczos trlanczos cross cyclic lapack randomized}} -svd_nsv 3
113 : requires: double
114 :
115 : TEST*/
|