Actual source code: test30.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 changing EPS type.\n\n";
13: #include <slepceps.h>
15: int main(int argc,char **argv)
16: {
17: Mat A; /* problem matrix */
18: EPS eps; /* eigenproblem solver context */
19: ST st;
20: KSP ksp;
21: PetscInt n=20,i,Istart,Iend,nrest;
23: PetscFunctionBeginUser;
24: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
25: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
27: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
28: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
29: PetscCall(MatSetFromOptions(A));
30: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
31: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A,i,i,i+1,INSERT_VALUES));
32: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
33: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Create eigensolver and view options
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
38: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
39: PetscCall(EPSSetOperators(eps,A,NULL));
40: PetscCall(EPSSetProblemType(eps,EPS_HEP));
41: PetscCall(EPSSetTarget(eps,4.8));
42: PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
43: PetscCall(EPSSetDimensions(eps,4,PETSC_DETERMINE,PETSC_DETERMINE));
45: PetscCall(EPSSetType(eps,EPSRQCG));
46: PetscCall(EPSRQCGSetReset(eps,10));
47: PetscCall(EPSSetFromOptions(eps));
48: PetscCall(EPSRQCGGetReset(eps,&nrest));
49: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"RQCG reset parameter = %" PetscInt_FMT "\n\n",nrest));
50: PetscCall(EPSView(eps,NULL));
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Change eigensolver type and solve problem
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: PetscCall(EPSSetType(eps,EPSGD));
56: PetscCall(EPSGetST(eps,&st));
57: PetscCall(STGetKSP(st,&ksp));
58: PetscCall(KSPSetType(ksp,KSPPREONLY));
59: PetscCall(EPSSolve(eps));
60: PetscCall(EPSView(eps,NULL));
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Display solution and clean up
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
66: PetscCall(EPSDestroy(&eps));
67: PetscCall(MatDestroy(&A));
68: PetscCall(SlepcFinalize());
69: return 0;
70: }
72: /*TEST
74: test:
75: suffix: 1
76: args: -eps_harmonic -eps_conv_norm
77: filter: grep -v tolerance | sed -e "s/symmetric/hermitian/"
79: TEST*/