Actual source code: test15.c
slepc-3.17.2 2022-08-09
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\n";
17: #include <slepcsvd.h>
19: int main(int argc,char **argv)
20: {
21: Mat A,B;
22: SVD svd;
23: KSP ksp;
24: PC pc;
25: PetscInt m=15,n=20,p=21,i,j,d,Istart,Iend;
26: PetscReal keep;
27: PetscBool flg,lock;
28: SVDTRLanczosGBidiag bidiag;
30: SlepcInitialize(&argc,&argv,(char*)0,help);
31: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
32: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
33: PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);
34: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%" PetscInt_FMT "+%" PetscInt_FMT ")x%" PetscInt_FMT "\n\n",m,p,n);
36: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: Generate the matrices
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: MatCreate(PETSC_COMM_WORLD,&A);
41: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
42: MatSetFromOptions(A);
43: MatSetUp(A);
45: MatGetOwnershipRange(A,&Istart,&Iend);
46: for (i=Istart;i<Iend;i++) {
47: if (i>0 && i-1<n) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
48: if (i+1<n) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
49: if (i<n) MatSetValue(A,i,i,2.0,INSERT_VALUES);
50: if (i>n) MatSetValue(A,i,n-1,1.0,INSERT_VALUES);
51: }
52: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
53: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
55: MatCreate(PETSC_COMM_WORLD,&B);
56: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
57: MatSetFromOptions(B);
58: MatSetUp(B);
60: MatGetOwnershipRange(B,&Istart,&Iend);
61: d = PetscMax(0,n-p);
62: for (i=Istart;i<Iend;i++) {
63: for (j=0;j<=PetscMin(i,n-1);j++) MatSetValue(B,i,j+d,1.0,INSERT_VALUES);
64: }
65: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
66: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Create the singular value solver, set options and solve
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: SVDCreate(PETSC_COMM_WORLD,&svd);
73: SVDSetOperators(svd,A,B);
74: SVDSetProblemType(svd,SVD_GENERALIZED);
75: SVDSetDimensions(svd,4,PETSC_DEFAULT,PETSC_DEFAULT);
76: SVDSetConvergenceTest(svd,SVD_CONV_NORM);
78: SVDSetType(svd,SVDTRLANCZOS);
79: SVDTRLanczosSetGBidiag(svd,SVD_TRLANCZOS_GBIDIAG_UPPER);
81: /* create a standalone KSP with appropriate settings */
82: KSPCreate(PETSC_COMM_WORLD,&ksp);
83: KSPSetType(ksp,KSPLSQR);
84: KSPGetPC(ksp,&pc);
85: PCSetType(pc,PCNONE);
86: KSPSetTolerances(ksp,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
87: KSPSetFromOptions(ksp);
88: SVDTRLanczosSetKSP(svd,ksp);
89: SVDTRLanczosSetRestart(svd,0.4);
90: SVDTRLanczosSetLocking(svd,PETSC_TRUE);
92: SVDSetFromOptions(svd);
94: PetscObjectTypeCompare((PetscObject)svd,SVDTRLANCZOS,&flg);
95: if (flg) {
96: SVDTRLanczosGetGBidiag(svd,&bidiag);
97: PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: using %s bidiagonalization\n",SVDTRLanczosGBidiags[bidiag]);
98: SVDTRLanczosGetRestart(svd,&keep);
99: SVDTRLanczosGetLocking(svd,&lock);
100: PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: restarting parameter %.2f %s\n",(double)keep,lock?"(locking)":"");
101: }
103: SVDSolve(svd);
104: SVDErrorView(svd,SVD_ERROR_NORM,NULL);
106: /* Free work space */
107: SVDDestroy(&svd);
108: KSPDestroy(&ksp);
109: MatDestroy(&A);
110: MatDestroy(&B);
111: SlepcFinalize();
112: return 0;
113: }
115: /*TEST
117: test:
118: suffix: 1
119: args: -svd_trlanczos_gbidiag {{single upper lower}}
120: filter: grep -v "TRLANCZOS: using"
121: requires: !single
123: test:
124: suffix: 2
125: args: -m 6 -n 12 -p 12 -svd_trlanczos_restart .7
126: requires: !single
128: TEST*/