Actual source code: test13.c
slepc-3.18.1 2022-11-02
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 EPSSetArbitrarySelection.\n\n";
13: #include <slepceps.h>
15: PetscErrorCode MyArbitrarySelection(PetscScalar eigr,PetscScalar eigi,Vec xr,Vec xi,PetscScalar *rr,PetscScalar *ri,void *ctx)
16: {
17: Vec xref = *(Vec*)ctx;
20: VecDot(xr,xref,rr);
21: *rr = PetscAbsScalar(*rr);
22: if (ri) *ri = 0.0;
23: return 0;
24: }
26: int main(int argc,char **argv)
27: {
28: Mat A; /* problem matrices */
29: EPS eps; /* eigenproblem solver context */
30: PetscScalar seigr,seigi;
31: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
32: Vec sxr,sxi;
33: PetscInt n=30,i,Istart,Iend,nconv;
36: SlepcInitialize(&argc,&argv,(char*)0,help);
38: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
39: PetscPrintf(PETSC_COMM_WORLD,"\nTridiagonal with zero diagonal, n=%" PetscInt_FMT "\n\n",n);
41: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42: Create matrix tridiag([-1 0 -1])
43: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
44: MatCreate(PETSC_COMM_WORLD,&A);
45: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
46: MatSetFromOptions(A);
47: MatSetUp(A);
49: MatGetOwnershipRange(A,&Istart,&Iend);
50: for (i=Istart;i<Iend;i++) {
51: if (i>0) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
52: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
53: }
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create the eigensolver
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: EPSCreate(PETSC_COMM_WORLD,&eps);
61: EPSSetProblemType(eps,EPS_HEP);
62: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
63: EPSSetOperators(eps,A,NULL);
64: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
65: EPSSetFromOptions(eps);
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Solve eigenproblem and store some solution
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: EPSSolve(eps);
71: MatCreateVecs(A,&sxr,NULL);
72: MatCreateVecs(A,&sxi,NULL);
73: EPSGetConverged(eps,&nconv);
74: if (nconv>0) {
75: EPSGetEigenpair(eps,0,&seigr,&seigi,sxr,sxi);
76: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Solve eigenproblem using an arbitrary selection
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: EPSSetArbitrarySelection(eps,MyArbitrarySelection,&sxr);
82: EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE);
83: EPSSolve(eps);
84: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
85: } else PetscPrintf(PETSC_COMM_WORLD,"Problem: no eigenpairs converged.\n");
87: EPSDestroy(&eps);
88: VecDestroy(&sxr);
89: VecDestroy(&sxi);
90: MatDestroy(&A);
91: SlepcFinalize();
92: return 0;
93: }
95: /*TEST
97: testset:
98: args: -eps_max_it 5000 -st_pc_type jacobi
99: output_file: output/test13_1.out
100: filter: sed -e "s/-1.98975/-1.98974/"
101: test:
102: suffix: 1
103: args: -eps_type {{krylovschur gd jd}}
104: test:
105: suffix: 1_gd2
106: args: -eps_type gd -eps_gd_double_expansion
107: test:
108: suffix: 2
109: args: -eps_non_hermitian -eps_type {{krylovschur gd jd}}
110: test:
111: suffix: 2_gd2
112: args: -eps_non_hermitian -eps_type gd -eps_gd_double_expansion
113: timeoutfactor: 2
115: TEST*/