Actual source code: test13.c
slepc-3.22.1 2024-10-28
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;
19: PetscFunctionBeginUser;
20: PetscCall(VecDot(xr,xref,rr));
21: *rr = PetscAbsScalar(*rr);
22: if (ri) *ri = 0.0;
23: PetscFunctionReturn(PETSC_SUCCESS);
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;
35: PetscFunctionBeginUser;
36: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
38: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
39: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTridiagonal with zero diagonal, n=%" PetscInt_FMT "\n\n",n));
41: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42: Create matrix tridiag([-1 0 -1])
43: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
44: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
45: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
46: PetscCall(MatSetFromOptions(A));
48: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
49: for (i=Istart;i<Iend;i++) {
50: if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
51: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
52: }
53: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
54: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Create the eigensolver
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
60: PetscCall(EPSSetProblemType(eps,EPS_HEP));
61: PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT));
62: PetscCall(EPSSetOperators(eps,A,NULL));
63: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
64: PetscCall(EPSSetFromOptions(eps));
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Solve eigenproblem and store some solution
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: PetscCall(EPSSolve(eps));
70: PetscCall(MatCreateVecs(A,&sxr,NULL));
71: PetscCall(MatCreateVecs(A,&sxi,NULL));
72: PetscCall(EPSGetConverged(eps,&nconv));
73: if (nconv>0) {
74: PetscCall(EPSGetEigenpair(eps,0,&seigr,&seigi,sxr,sxi));
75: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
77: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: Solve eigenproblem using an arbitrary selection
79: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: PetscCall(EPSSetArbitrarySelection(eps,MyArbitrarySelection,&sxr));
81: PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE));
82: PetscCall(EPSSolve(eps));
83: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
84: } else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Problem: no eigenpairs converged.\n"));
86: PetscCall(EPSDestroy(&eps));
87: PetscCall(VecDestroy(&sxr));
88: PetscCall(VecDestroy(&sxi));
89: PetscCall(MatDestroy(&A));
90: PetscCall(SlepcFinalize());
91: return 0;
92: }
94: /*TEST
96: testset:
97: args: -eps_max_it 5000 -st_pc_type jacobi
98: output_file: output/test13_1.out
99: filter: sed -e "s/-1.98975/-1.98974/"
100: test:
101: suffix: 1
102: args: -eps_type {{krylovschur gd}}
103: test:
104: suffix: 1_jd
105: args: -eps_type jd -st_ksp_type gmres
106: test:
107: suffix: 1_gd2
108: args: -eps_type gd -eps_gd_double_expansion
109: test:
110: suffix: 2
111: args: -eps_non_hermitian -eps_type {{krylovschur gd jd}}
112: test:
113: suffix: 2_gd2
114: args: -eps_non_hermitian -eps_type gd -eps_gd_double_expansion
115: timeoutfactor: 2
117: TEST*/