Actual source code: test8.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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[] = "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";
15: #include <slepcsvd.h>
17: /*
18: This example computes the singular values of an nxn Grcar matrix,
19: which is a nonsymmetric Toeplitz matrix:
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 |
31: */
33: int main(int argc,char **argv)
34: {
35: Mat A;
36: SVD svd;
37: PetscInt N=30,Istart,Iend,i,col[5],nsv,ncv;
38: PetscScalar value[] = { -1, 1, 1, 1, 1 };
41: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
42: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
43: PetscPrintf(PETSC_COMM_WORLD,"\nSingular values of a Grcar matrix, n=%D",N);
44: PetscPrintf(PETSC_COMM_WORLD,"\n\n");
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Generate the matrix
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: MatCreate(PETSC_COMM_WORLD,&A);
51: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
52: MatSetFromOptions(A);
53: MatSetUp(A);
55: MatGetOwnershipRange(A,&Istart,&Iend);
56: for (i=Istart;i<Iend;i++) {
57: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
58: if (i==0) {
59: MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES);
60: } else {
61: MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
62: }
63: }
65: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
66: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Create the singular value solver and set the solution method
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: SVDCreate(PETSC_COMM_WORLD,&svd);
73: SVDSetOperator(svd,A);
74: SVDSetTolerances(svd,1e-6,1000);
75: SVDSetWhichSingularTriplets(svd,SVD_LARGEST);
76: SVDSetFromOptions(svd);
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Compute the singular values
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: /* First solve */
83: SVDSolve(svd);
84: PetscPrintf(PETSC_COMM_WORLD," - - - First solve, default subspace dimension - - -\n");
85: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
87: /* Second solve */
88: SVDGetDimensions(svd,&nsv,&ncv,NULL);
89: SVDSetDimensions(svd,nsv,ncv+2,PETSC_DEFAULT);
90: SVDSolve(svd);
91: PetscPrintf(PETSC_COMM_WORLD," - - - Second solve, subspace of increased size - - -\n");
92: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
94: /* Free work space */
95: SVDDestroy(&svd);
96: MatDestroy(&A);
97: SlepcFinalize();
98: return ierr;
99: }
101: /*TEST
103: test:
104: suffix: 1
105: args: -svd_type {{lanczos trlanczos cross cyclic lapack}} -svd_nsv 3 -svd_ncv 12
107: TEST*/