Actual source code: test16.c
slepc-3.20.1 2023-11-27
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 (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 };
40: PetscFunctionBeginUser;
41: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
42: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
43: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
44: PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
45: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%" PetscInt_FMT "+%" PetscInt_FMT ")x%" PetscInt_FMT "\n\n",m,p,n));
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Generate the matrices
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: PetscCall(MatCreate(PETSC_COMM_WORLD,&A1));
52: PetscCall(MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m,n));
53: PetscCall(MatSetFromOptions(A1));
54: PetscCall(MatSetUp(A1));
55: PetscCall(MatGetOwnershipRange(A1,&Istart,&Iend));
56: for (i=Istart;i<Iend;i++) {
57: col[0]=i; col[1]=i+1;
58: if (i<n-1) PetscCall(MatSetValues(A1,1,&i,2,col,valsa,INSERT_VALUES));
59: else if (i==n-1) PetscCall(MatSetValue(A1,i,col[0],valsa[0],INSERT_VALUES));
60: }
61: PetscCall(MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY));
62: PetscCall(MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY));
64: PetscCall(MatCreate(PETSC_COMM_WORLD,&A2));
65: PetscCall(MatSetSizes(A2,PETSC_DECIDE,PETSC_DECIDE,m,n));
66: PetscCall(MatSetFromOptions(A2));
67: PetscCall(MatSetUp(A2));
68: PetscCall(MatGetOwnershipRange(A2,&Istart,&Iend));
69: for (i=Istart;i<Iend;i++) {
70: col[0]=i-1; col[1]=i;
71: if (i==0) PetscCall(MatSetValue(A2,i,col[1],valsb[1],INSERT_VALUES));
72: else if (i<n) PetscCall(MatSetValues(A2,1,&i,2,col,valsb,INSERT_VALUES));
73: else if (i==n) PetscCall(MatSetValue(A2,i,col[0],valsb[0],INSERT_VALUES));
74: }
75: PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY));
76: PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY));
78: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
79: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n));
80: PetscCall(MatSetFromOptions(B));
81: PetscCall(MatSetUp(B));
82: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
83: d = PetscMax(0,n-p);
84: for (i=Istart;i<Iend;i++) {
85: for (j=PetscMax(0,i-5);j<=PetscMin(i,n-1);j++) PetscCall(MatSetValue(B,i,j+d,1.0,INSERT_VALUES));
86: }
87: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
88: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Create the singular value solver, set options and solve
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
95: PetscCall(SVDSetOperators(svd,A1,B));
96: PetscCall(SVDSetFromOptions(svd));
97: PetscCall(SVDSolve(svd));
98: PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,NULL));
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Solve second problem
102: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: PetscCall(SVDSetOperators(svd,A2,B));
105: PetscCall(SVDSolve(svd));
106: PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,NULL));
108: /* Free work space */
109: PetscCall(SVDDestroy(&svd));
110: PetscCall(MatDestroy(&A1));
111: PetscCall(MatDestroy(&A2));
112: PetscCall(MatDestroy(&B));
113: PetscCall(SlepcFinalize());
114: return 0;
115: }
117: /*TEST
119: testset:
120: args: -svd_nsv 3
121: requires: !single
122: output_file: output/test16_1.out
123: test:
124: suffix: 1_lapack
125: args: -svd_type lapack
126: test:
127: suffix: 1_cross
128: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
129: test:
130: suffix: 1_cyclic
131: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
132: test:
133: suffix: 1_trlanczos
134: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single lower}} -svd_trlanczos_ksp_rtol 1e-10
135: requires: double
136: test:
137: suffix: 1_trlanczos_par
138: nsize: 2
139: args: -svd_type trlanczos -ds_parallel {{redundant synchronized}}
141: testset:
142: args: -svd_nsv 3 -mat_type aijcusparse
143: requires: cuda !single
144: output_file: output/test16_1.out
145: test:
146: suffix: 2_cross
147: args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
148: test:
149: suffix: 2_cyclic
150: args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}}
151: test:
152: suffix: 2_trlanczos
153: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single lower}} -svd_trlanczos_ksp_rtol 1e-10
154: requires: double
156: TEST*/