Actual source code: test15.c
slepc-3.19.0 2023-03-31
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 user interface for TRLANCZOS with 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"
16: " -s <s>, where <s> = scale parameter.\n\n";
18: #include <slepcsvd.h>
20: int main(int argc,char **argv)
21: {
22: Mat A,B;
23: SVD svd;
24: KSP ksp;
25: PC pc;
26: PetscInt m=15,n=20,p=21,i,j,d,Istart,Iend;
27: PetscReal keep,scale=1.0;
28: PetscBool flg,lock,oneside;
29: SVDTRLanczosGBidiag bidiag;
31: PetscFunctionBeginUser;
32: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
33: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
34: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
35: PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
36: PetscCall(PetscOptionsGetReal(NULL,NULL,"-s",&scale,NULL));
37: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%" PetscInt_FMT "+%" PetscInt_FMT ")x%" PetscInt_FMT "\n\n",m,p,n));
39: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: Generate the matrices
41: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
43: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
44: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n));
45: PetscCall(MatSetFromOptions(A));
46: PetscCall(MatSetUp(A));
48: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
49: for (i=Istart;i<Iend;i++) {
50: if (i>0 && i-1<n) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
51: if (i+1<n) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
52: if (i<n) PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
53: if (i>n) PetscCall(MatSetValue(A,i,n-1,1.0,INSERT_VALUES));
54: }
55: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
56: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
58: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
59: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n));
60: PetscCall(MatSetFromOptions(B));
61: PetscCall(MatSetUp(B));
63: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
64: d = PetscMax(0,n-p);
65: for (i=Istart;i<Iend;i++) {
66: for (j=0;j<=PetscMin(i,n-1);j++) PetscCall(MatSetValue(B,i,j+d,1.0,INSERT_VALUES));
67: }
68: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
69: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Create the singular value solver, set options and solve
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
76: PetscCall(SVDSetOperators(svd,A,B));
77: PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED));
78: PetscCall(SVDSetDimensions(svd,4,PETSC_DEFAULT,PETSC_DEFAULT));
79: PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_NORM));
81: PetscCall(SVDSetType(svd,SVDTRLANCZOS));
82: PetscCall(SVDTRLanczosSetGBidiag(svd,SVD_TRLANCZOS_GBIDIAG_UPPER));
83: PetscCall(SVDTRLanczosSetScale(svd,scale));
85: /* create a standalone KSP with appropriate settings */
86: PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp));
87: PetscCall(KSPSetType(ksp,KSPLSQR));
88: PetscCall(KSPGetPC(ksp,&pc));
89: PetscCall(PCSetType(pc,PCNONE));
90: PetscCall(KSPSetTolerances(ksp,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
91: PetscCall(KSPSetFromOptions(ksp));
92: PetscCall(SVDTRLanczosSetKSP(svd,ksp));
93: PetscCall(SVDTRLanczosSetRestart(svd,0.4));
94: PetscCall(SVDTRLanczosSetLocking(svd,PETSC_TRUE));
96: PetscCall(SVDSetFromOptions(svd));
98: PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDTRLANCZOS,&flg));
99: if (flg) {
100: PetscCall(SVDTRLanczosGetGBidiag(svd,&bidiag));
101: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: using %s bidiagonalization\n",SVDTRLanczosGBidiags[bidiag]));
102: PetscCall(SVDTRLanczosGetRestart(svd,&keep));
103: PetscCall(SVDTRLanczosGetLocking(svd,&lock));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: restarting parameter %.2f %s\n",(double)keep,lock?"(locking)":""));
105: PetscCall(SVDTRLanczosGetScale(svd,&scale));
106: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: scale parameter %g\n",(double)scale));
107: PetscCall(SVDTRLanczosGetOneSide(svd,&oneside));
108: if (oneside) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: using one-sided orthogonalization\n"));
109: }
111: PetscCall(SVDSolve(svd));
112: PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,NULL));
114: /* Free work space */
115: PetscCall(SVDDestroy(&svd));
116: PetscCall(KSPDestroy(&ksp));
117: PetscCall(MatDestroy(&A));
118: PetscCall(MatDestroy(&B));
119: PetscCall(SlepcFinalize());
120: return 0;
121: }
123: /*TEST
125: test:
126: suffix: 1
127: args: -svd_trlanczos_gbidiag {{single upper lower}}
128: filter: grep -v "TRLANCZOS: using"
129: requires: !single
131: test:
132: suffix: 2
133: args: -m 6 -n 12 -p 12 -svd_trlanczos_restart .7
134: requires: !single
136: test:
137: suffix: 3
138: args: -s 8 -svd_trlanczos_gbidiag lower
139: requires: !single
141: test:
142: suffix: 4
143: args: -s -5 -svd_trlanczos_gbidiag lower
144: requires: !single
146: TEST*/