Actual source code: test14.c
slepc-3.15.1 2021-05-28
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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 equal matrix size.\n\n"
12: "The command line options are:\n"
13: " -m <m>, where <m> = matrix rows.\n"
14: " -n <n>, where <n> = matrix columns (defaults to m+2).\n\n";
16: #include <slepcsvd.h>
18: /*
19: This example computes the singular values of two rectangular bidiagonal matrices
21: | 1 2 | | 1 |
22: | 1 2 | | 2 1 |
23: | 1 2 | | 2 1 |
24: A = | . . | B = | . . |
25: | . . | | . . |
26: | 1 2 | | 2 1 |
27: | 1 2 | | 2 1 |
28: */
30: int main(int argc,char **argv)
31: {
32: Mat A,B;
33: SVD svd;
34: PetscInt m=20,n,Istart,Iend,i,col[2];
35: PetscScalar valsa[] = { 1, 2 }, valsb[] = { 2, 1 };
36: PetscBool flg;
39: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
40: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
41: PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
42: if (!flg) n=m+2;
43: PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%D n=%D\n\n",m,n);
45: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: Generate the matrices
47: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49: MatCreate(PETSC_COMM_WORLD,&A);
50: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
51: MatSetFromOptions(A);
52: MatSetUp(A);
53: MatGetOwnershipRange(A,&Istart,&Iend);
54: for (i=Istart;i<Iend;i++) {
55: col[0]=i; col[1]=i+1;
56: if (i<n-1) {
57: MatSetValues(A,1,&i,2,col,valsa,INSERT_VALUES);
58: } else if (i==n-1) {
59: MatSetValue(A,i,col[0],valsa[0],INSERT_VALUES);
60: }
61: }
62: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
63: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
65: MatCreate(PETSC_COMM_WORLD,&B);
66: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
67: MatSetFromOptions(B);
68: MatSetUp(B);
69: MatGetOwnershipRange(B,&Istart,&Iend);
70: for (i=Istart;i<Iend;i++) {
71: col[0]=i; col[1]=i+1;
72: if (i==0) {
73: MatSetValue(B,i,col[0],valsb[1],INSERT_VALUES);
74: } else if (i<n-1) {
75: MatSetValues(B,1,&i,2,col,valsb,INSERT_VALUES);
76: } else if (i==n-1) {
77: MatSetValue(B,i,col[0],valsb[0],INSERT_VALUES);
78: }
79: }
80: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Create the singular value solver, set options and solve
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: SVDCreate(PETSC_COMM_WORLD,&svd);
88: SVDSetOperators(svd,A,NULL);
89: SVDSetTolerances(svd,PETSC_DEFAULT,1000);
90: SVDSetFromOptions(svd);
91: SVDSolve(svd);
92: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Solve with second matrix
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: SVDSetOperators(svd,B,NULL);
99: SVDSolve(svd);
100: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
102: /* Free work space */
103: SVDDestroy(&svd);
104: MatDestroy(&A);
105: MatDestroy(&B);
106: SlepcFinalize();
107: return ierr;
108: }
110: /*TEST
112: testset:
113: args: -svd_nsv 3
114: requires: !single
115: output_file: output/test14_1.out
116: test:
117: suffix: 1
118: args: -svd_type {{lanczos trlanczos lapack}}
119: test:
120: suffix: 1_cross
121: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
122: test:
123: suffix: 1_cyclic
124: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
126: testset:
127: args: -n 18 -svd_nsv 3
128: requires: !single
129: output_file: output/test14_2.out
130: test:
131: suffix: 2
132: args: -svd_type {{lanczos trlanczos lapack}}
133: test:
134: suffix: 2_cross
135: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
136: test:
137: suffix: 2_cyclic
138: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
140: TEST*/