Actual source code: test15.c
slepc-3.21.1 2024-04-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 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));
47: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
48: for (i=Istart;i<Iend;i++) {
49: if (i>0 && i-1<n) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
50: if (i+1<n) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
51: if (i<n) PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
52: if (i>n) PetscCall(MatSetValue(A,i,n-1,1.0,INSERT_VALUES));
53: }
54: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
55: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
57: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
58: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n));
59: PetscCall(MatSetFromOptions(B));
61: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
62: d = PetscMax(0,n-p);
63: for (i=Istart;i<Iend;i++) {
64: for (j=0;j<=PetscMin(i,n-1);j++) PetscCall(MatSetValue(B,i,j+d,1.0,INSERT_VALUES));
65: }
66: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
67: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
69: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: Create the singular value solver, set options and solve
71: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
74: PetscCall(SVDSetOperators(svd,A,B));
75: PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED));
76: PetscCall(SVDSetDimensions(svd,4,PETSC_DEFAULT,PETSC_DEFAULT));
77: PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_NORM));
79: PetscCall(SVDSetType(svd,SVDTRLANCZOS));
80: PetscCall(SVDTRLanczosSetGBidiag(svd,SVD_TRLANCZOS_GBIDIAG_UPPER));
81: PetscCall(SVDTRLanczosSetScale(svd,scale));
83: /* create a standalone KSP with appropriate settings */
84: PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp));
85: PetscCall(KSPSetType(ksp,KSPLSQR));
86: PetscCall(KSPGetPC(ksp,&pc));
87: PetscCall(PCSetType(pc,PCNONE));
88: PetscCall(KSPSetTolerances(ksp,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
89: PetscCall(KSPSetFromOptions(ksp));
90: PetscCall(SVDTRLanczosSetKSP(svd,ksp));
91: PetscCall(SVDTRLanczosSetRestart(svd,0.4));
92: PetscCall(SVDTRLanczosSetLocking(svd,PETSC_TRUE));
94: PetscCall(SVDSetFromOptions(svd));
96: PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDTRLANCZOS,&flg));
97: if (flg) {
98: PetscCall(SVDTRLanczosGetGBidiag(svd,&bidiag));
99: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: using %s bidiagonalization\n",SVDTRLanczosGBidiags[bidiag]));
100: PetscCall(SVDTRLanczosGetRestart(svd,&keep));
101: PetscCall(SVDTRLanczosGetLocking(svd,&lock));
102: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: restarting parameter %.2f %s\n",(double)keep,lock?"(locking)":""));
103: PetscCall(SVDTRLanczosGetScale(svd,&scale));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: scale parameter %g\n",(double)scale));
105: PetscCall(SVDTRLanczosGetOneSide(svd,&oneside));
106: if (oneside) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: using one-sided orthogonalization\n"));
107: }
109: PetscCall(SVDSolve(svd));
110: PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,NULL));
112: /* Free work space */
113: PetscCall(SVDDestroy(&svd));
114: PetscCall(KSPDestroy(&ksp));
115: PetscCall(MatDestroy(&A));
116: PetscCall(MatDestroy(&B));
117: PetscCall(SlepcFinalize());
118: return 0;
119: }
121: /*TEST
123: test:
124: suffix: 1
125: args: -svd_trlanczos_gbidiag {{single upper lower}}
126: filter: grep -v "TRLANCZOS: using"
127: requires: !single
129: test:
130: suffix: 2
131: args: -m 6 -n 12 -p 12 -svd_trlanczos_restart .7
132: requires: !single
134: test:
135: suffix: 3
136: args: -s 8 -svd_trlanczos_gbidiag lower
137: requires: !single
139: test:
140: suffix: 4
141: args: -s -5 -svd_trlanczos_gbidiag lower
142: requires: !single
144: TEST*/