Actual source code: test14.c
slepc-3.18.2 2023-01-26
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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);
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=%" PetscInt_FMT " n=%" PetscInt_FMT "\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) MatSetValues(A,1,&i,2,col,valsa,INSERT_VALUES);
57: else if (i==n-1) MatSetValue(A,i,col[0],valsa[0],INSERT_VALUES);
58: }
59: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
60: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
62: MatCreate(PETSC_COMM_WORLD,&B);
63: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
64: MatSetFromOptions(B);
65: MatSetUp(B);
66: MatGetOwnershipRange(B,&Istart,&Iend);
67: for (i=Istart;i<Iend;i++) {
68: col[0]=i-1; col[1]=i;
69: if (i==0) MatSetValue(B,i,col[1],valsb[1],INSERT_VALUES);
70: else if (i<n) MatSetValues(B,1,&i,2,col,valsb,INSERT_VALUES);
71: else if (i==n) MatSetValue(B,i,col[0],valsb[0],INSERT_VALUES);
72: }
73: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Create the singular value solver, set options and solve
78: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: SVDCreate(PETSC_COMM_WORLD,&svd);
81: SVDSetOperators(svd,A,NULL);
82: SVDSetTolerances(svd,PETSC_DEFAULT,1000);
83: SVDSetFromOptions(svd);
84: SVDSolve(svd);
85: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Solve with second matrix
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: SVDSetOperators(svd,B,NULL);
92: SVDSolve(svd);
93: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
95: /* Free work space */
96: SVDDestroy(&svd);
97: MatDestroy(&A);
98: MatDestroy(&B);
99: SlepcFinalize();
100: return 0;
101: }
103: /*TEST
105: testset:
106: args: -svd_nsv 3
107: requires: !single
108: output_file: output/test14_1.out
109: test:
110: suffix: 1
111: args: -svd_type {{lanczos trlanczos lapack}}
112: test:
113: suffix: 1_cross
114: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
115: test:
116: suffix: 1_cyclic
117: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
119: testset:
120: args: -n 18 -svd_nsv 3
121: requires: !single
122: output_file: output/test14_2.out
123: test:
124: suffix: 2
125: args: -svd_type {{lanczos trlanczos lapack}}
126: test:
127: suffix: 2_cross
128: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
129: test:
130: suffix: 2_cyclic
131: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
133: TEST*/