Actual source code: test18.c
slepc-3.19.2 2023-09-05
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[] = "Test GSVD with user-provided initial vectors.\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 a GSVD problem for the bidiagonal matrices
22: | 1 2 | | 2 |
23: | 1 2 | | -1 2 |
24: | 1 2 | | -1 2 |
25: A = | . . | B = | . . |
26: | . . | | . . |
27: | 1 2 | | -1 2 |
28: | 1 2 | | -1 2 |
29: */
31: int main(int argc,char **argv)
32: {
33: Mat A,B;
34: SVD svd;
35: Vec v0,w0; /* initial vectors */
36: VecType vtype;
37: PetscInt m=22,n=20,p=22,Istart,Iend,i,col[2];
38: PetscScalar valsa[] = { 1, 2 }, valsb[] = { -1, 2 };
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,&A));
52: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n));
53: PetscCall(MatSetFromOptions(A));
54: PetscCall(MatSetUp(A));
55: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
56: for (i=Istart;i<Iend;i++) {
57: col[0]=i; col[1]=i+1;
58: if (i<n-1) PetscCall(MatSetValues(A,1,&i,2,col,valsa,INSERT_VALUES));
59: else if (i==n-1) PetscCall(MatSetValue(A,i,col[0],valsa[0],INSERT_VALUES));
60: }
61: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
62: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
64: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
65: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n));
66: PetscCall(MatSetFromOptions(B));
67: PetscCall(MatSetUp(B));
68: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
69: for (i=Istart;i<Iend;i++) {
70: col[0]=i-1; col[1]=i;
71: if (i==0) PetscCall(MatSetValue(B,i,col[1],valsb[1],INSERT_VALUES));
72: else if (i<n) PetscCall(MatSetValues(B,1,&i,2,col,valsb,INSERT_VALUES));
73: else if (i==n) PetscCall(MatSetValue(B,i,col[0],valsb[0],INSERT_VALUES));
74: }
75: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
76: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Create the singular value solver, set options and solve
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
83: PetscCall(SVDSetOperators(svd,A,B));
84: PetscCall(SVDSetFromOptions(svd));
86: /*
87: Set the initial vectors. This is optional, if not done the initial
88: vectors are set to random values
89: */
90: PetscCall(MatCreateVecs(A,&v0,NULL)); /* right initial vector, length n */
91: PetscCall(VecCreate(PETSC_COMM_WORLD,&w0)); /* left initial vector, length m+p */
92: PetscCall(VecSetSizes(w0,PETSC_DECIDE,m+p));
93: PetscCall(VecGetType(v0,&vtype));
94: PetscCall(VecSetType(w0,vtype));
95: PetscCall(VecSet(v0,1.0));
96: PetscCall(VecSet(w0,1.0));
97: PetscCall(SVDSetInitialSpaces(svd,1,&v0,1,&w0));
99: /*
100: Compute solution
101: */
102: PetscCall(SVDSolve(svd));
103: PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,NULL));
105: /* Free work space */
106: PetscCall(VecDestroy(&v0));
107: PetscCall(VecDestroy(&w0));
108: PetscCall(SVDDestroy(&svd));
109: PetscCall(MatDestroy(&A));
110: PetscCall(MatDestroy(&B));
111: PetscCall(SlepcFinalize());
112: return 0;
113: }
115: /*TEST
117: testset:
118: args: -svd_nsv 3
119: requires: !single
120: output_file: output/test18_1.out
121: test:
122: suffix: 1
123: args: -svd_type {{lapack cross cyclic}}
124: test:
125: suffix: 1_trlanczos
126: args: -svd_type trlanczos -svd_trlanczos_gbidiag {{single upper lower}}
127: requires: !__float128
129: test:
130: suffix: 2
131: args: -svd_nsv 3 -svd_type trlanczos -svd_monitor_conditioning
132: requires: double
134: TEST*/