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