Actual source code: test16.c
slepc-3.16.3 2022-04-11
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 (GSVD).\n\n"
12: "The command line options are:\n"
13: " -m <m>, where <m> = number of rows of A.\n"
14: " -n <n>, where <n> = number of columns of A.\n"
15: " -p <p>, where <p> = number of rows of B.\n\n";
17: #include <slepcsvd.h>
19: /*
20: This example solves two GSVD problems for the bidiagonal matrices
22: | 1 2 | | 1 |
23: | 1 2 | | 2 1 |
24: | 1 2 | | 2 1 |
25: A1 = | . . | A2 = | . . |
26: | . . | | . . |
27: | 1 2 | | 2 1 |
28: | 1 2 | | 2 1 |
30: with B = tril(ones(p,n))
31: */
33: int main(int argc,char **argv)
34: {
35: Mat A1,A2,B;
36: SVD svd;
37: PetscInt m=15,n=20,p=21,Istart,Iend,i,j,d,col[2];
38: PetscScalar valsa[] = { 1, 2 }, valsb[] = { 2, 1 };
41: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
42: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
43: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
44: PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);
45: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%D+%D)x%D\n\n",m,p,n);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Generate the matrices
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: MatCreate(PETSC_COMM_WORLD,&A1);
52: MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m,n);
53: MatSetFromOptions(A1);
54: MatSetUp(A1);
55: MatGetOwnershipRange(A1,&Istart,&Iend);
56: for (i=Istart;i<Iend;i++) {
57: col[0]=i; col[1]=i+1;
58: if (i<n-1) {
59: MatSetValues(A1,1,&i,2,col,valsa,INSERT_VALUES);
60: } else if (i==n-1) {
61: MatSetValue(A1,i,col[0],valsa[0],INSERT_VALUES);
62: }
63: }
64: MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);
65: MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);
67: MatCreate(PETSC_COMM_WORLD,&A2);
68: MatSetSizes(A2,PETSC_DECIDE,PETSC_DECIDE,m,n);
69: MatSetFromOptions(A2);
70: MatSetUp(A2);
71: MatGetOwnershipRange(A2,&Istart,&Iend);
72: for (i=Istart;i<Iend;i++) {
73: col[0]=i-1; col[1]=i;
74: if (i==0) {
75: MatSetValue(A2,i,col[1],valsb[1],INSERT_VALUES);
76: } else if (i<n) {
77: MatSetValues(A2,1,&i,2,col,valsb,INSERT_VALUES);
78: } else if (i==n) {
79: MatSetValue(A2,i,col[0],valsb[0],INSERT_VALUES);
80: }
81: }
82: MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);
85: MatCreate(PETSC_COMM_WORLD,&B);
86: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
87: MatSetFromOptions(B);
88: MatSetUp(B);
89: MatGetOwnershipRange(B,&Istart,&Iend);
90: d = PetscMax(0,n-p);
91: for (i=Istart;i<Iend;i++) {
92: for (j=0;j<=PetscMin(i,n-1);j++) {
93: MatSetValue(B,i,j+d,1.0,INSERT_VALUES);
94: }
95: }
96: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
97: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Create the singular value solver, set options and solve
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: SVDCreate(PETSC_COMM_WORLD,&svd);
104: SVDSetOperators(svd,A1,B);
105: SVDSetFromOptions(svd);
106: SVDSolve(svd);
107: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Solve second problem
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: SVDSetOperators(svd,A2,B);
114: SVDSolve(svd);
115: SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
117: /* Free work space */
118: SVDDestroy(&svd);
119: MatDestroy(&A1);
120: MatDestroy(&A2);
121: MatDestroy(&B);
122: SlepcFinalize();
123: return ierr;
124: }
126: /*TEST
128: testset:
129: args: -svd_nsv 3
130: requires: !single
131: output_file: output/test16_1.out
132: test:
133: suffix: 1_lapack
134: args: -svd_type lapack
135: test:
136: suffix: 1_cross
137: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
138: test:
139: suffix: 1_cyclic
140: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
141: test:
142: suffix: 1_trlanczos
143: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}}
145: TEST*/