Actual source code: ex1.c
slepc-3.20.0 2023-09-29
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[] = "Standard symmetric eigenproblem corresponding to the Laplacian operator in 1 dimension.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";
15: #include <slepceps.h>
17: int main(int argc,char **argv)
18: {
19: Mat A; /* problem matrix */
20: EPS eps; /* eigenproblem solver context */
21: EPSType type;
22: PetscReal error,tol,re,im;
23: PetscScalar kr,ki;
24: Vec xr,xi;
25: PetscInt n=30,i,Istart,Iend,nev,maxit,its,nconv;
27: PetscFunctionBeginUser;
28: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
31: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
33: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: Compute the operator matrix that defines the eigensystem, Ax=kx
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
38: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
39: PetscCall(MatSetFromOptions(A));
40: PetscCall(MatSetUp(A));
42: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
43: for (i=Istart;i<Iend;i++) {
44: if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
45: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
46: PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
47: }
48: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
49: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
51: PetscCall(MatCreateVecs(A,NULL,&xr));
52: PetscCall(MatCreateVecs(A,NULL,&xi));
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Create the eigensolver and set various options
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: /*
58: Create eigensolver context
59: */
60: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
62: /*
63: Set operators. In this case, it is a standard eigenvalue problem
64: */
65: PetscCall(EPSSetOperators(eps,A,NULL));
66: PetscCall(EPSSetProblemType(eps,EPS_HEP));
68: /*
69: Set solver parameters at runtime
70: */
71: PetscCall(EPSSetFromOptions(eps));
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Solve the eigensystem
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: PetscCall(EPSSolve(eps));
78: /*
79: Optional: Get some information from the solver and display it
80: */
81: PetscCall(EPSGetIterationNumber(eps,&its));
82: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its));
83: PetscCall(EPSGetType(eps,&type));
84: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
85: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
86: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
87: PetscCall(EPSGetTolerances(eps,&tol,&maxit));
88: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Display solution and clean up
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: /*
94: Get number of converged approximate eigenpairs
95: */
96: PetscCall(EPSGetConverged(eps,&nconv));
97: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv));
99: if (nconv>0) {
100: /*
101: Display eigenvalues and relative errors
102: */
103: PetscCall(PetscPrintf(PETSC_COMM_WORLD,
104: " k ||Ax-kx||/||kx||\n"
105: " ----------------- ------------------\n"));
107: for (i=0;i<nconv;i++) {
108: /*
109: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
110: ki (imaginary part)
111: */
112: PetscCall(EPSGetEigenpair(eps,i,&kr,&ki,xr,xi));
113: /*
114: Compute the relative error associated to each eigenpair
115: */
116: PetscCall(EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error));
118: #if defined(PETSC_USE_COMPLEX)
119: re = PetscRealPart(kr);
120: im = PetscImaginaryPart(kr);
121: #else
122: re = kr;
123: im = ki;
124: #endif
125: if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error));
126: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error));
127: }
128: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
129: }
131: /*
132: Free work space
133: */
134: PetscCall(EPSDestroy(&eps));
135: PetscCall(MatDestroy(&A));
136: PetscCall(VecDestroy(&xr));
137: PetscCall(VecDestroy(&xi));
138: PetscCall(SlepcFinalize());
139: return 0;
140: }