Actual source code: ex1.c
slepc-main 2024-11-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[] = "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,NULL,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));
41: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
42: for (i=Istart;i<Iend;i++) {
43: if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
44: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
45: PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
46: }
47: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
48: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
50: PetscCall(MatCreateVecs(A,NULL,&xr));
51: PetscCall(MatCreateVecs(A,NULL,&xi));
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Create the eigensolver and set various options
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: /*
57: Create eigensolver context
58: */
59: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
61: /*
62: Set operators. In this case, it is a standard eigenvalue problem
63: */
64: PetscCall(EPSSetOperators(eps,A,NULL));
65: PetscCall(EPSSetProblemType(eps,EPS_HEP));
67: /*
68: Set solver parameters at runtime
69: */
70: PetscCall(EPSSetFromOptions(eps));
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Solve the eigensystem
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: PetscCall(EPSSolve(eps));
77: /*
78: Optional: Get some information from the solver and display it
79: */
80: PetscCall(EPSGetIterationNumber(eps,&its));
81: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its));
82: PetscCall(EPSGetType(eps,&type));
83: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
84: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
85: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
86: PetscCall(EPSGetTolerances(eps,&tol,&maxit));
87: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Display solution and clean up
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: /*
93: Get number of converged approximate eigenpairs
94: */
95: PetscCall(EPSGetConverged(eps,&nconv));
96: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv));
98: if (nconv>0) {
99: /*
100: Display eigenvalues and relative errors
101: */
102: PetscCall(PetscPrintf(PETSC_COMM_WORLD,
103: " k ||Ax-kx||/||kx||\n"
104: " ----------------- ------------------\n"));
106: for (i=0;i<nconv;i++) {
107: /*
108: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
109: ki (imaginary part)
110: */
111: PetscCall(EPSGetEigenpair(eps,i,&kr,&ki,xr,xi));
112: /*
113: Compute the relative error associated to each eigenpair
114: */
115: PetscCall(EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error));
117: #if defined(PETSC_USE_COMPLEX)
118: re = PetscRealPart(kr);
119: im = PetscImaginaryPart(kr);
120: #else
121: re = kr;
122: im = ki;
123: #endif
124: if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error));
125: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error));
126: }
127: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
128: }
130: /*
131: Free work space
132: */
133: PetscCall(EPSDestroy(&eps));
134: PetscCall(MatDestroy(&A));
135: PetscCall(VecDestroy(&xr));
136: PetscCall(VecDestroy(&xi));
137: PetscCall(SlepcFinalize());
138: return 0;
139: }