Actual source code: ex1.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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;
28: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
30: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
31: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%D\n\n",n);
33: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: Compute the operator matrix that defines the eigensystem, Ax=kx
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37: MatCreate(PETSC_COMM_WORLD,&A);
38: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
39: MatSetFromOptions(A);
40: MatSetUp(A);
42: MatGetOwnershipRange(A,&Istart,&Iend);
43: for (i=Istart;i<Iend;i++) {
44: if (i>0) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
45: if (i<n-1) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
46: MatSetValue(A,i,i,2.0,INSERT_VALUES);
47: }
48: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
49: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
51: MatCreateVecs(A,NULL,&xr);
52: MatCreateVecs(A,NULL,&xi);
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Create the eigensolver and set various options
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: /*
58: Create eigensolver context
59: */
60: EPSCreate(PETSC_COMM_WORLD,&eps);
62: /*
63: Set operators. In this case, it is a standard eigenvalue problem
64: */
65: EPSSetOperators(eps,A,NULL);
66: EPSSetProblemType(eps,EPS_HEP);
68: /*
69: Set solver parameters at runtime
70: */
71: EPSSetFromOptions(eps);
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Solve the eigensystem
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: EPSSolve(eps);
78: /*
79: Optional: Get some information from the solver and display it
80: */
81: EPSGetIterationNumber(eps,&its);
82: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
83: EPSGetType(eps,&type);
84: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
85: EPSGetDimensions(eps,&nev,NULL,NULL);
86: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
87: EPSGetTolerances(eps,&tol,&maxit);
88: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%D\n",(double)tol,maxit);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Display solution and clean up
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: /*
94: Get number of converged approximate eigenpairs
95: */
96: EPSGetConverged(eps,&nconv);
97: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %D\n\n",nconv);
99: if (nconv>0) {
100: /*
101: Display eigenvalues and relative errors
102: */
103: 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: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
113: /*
114: Compute the relative error associated to each eigenpair
115: */
116: 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) {
126: PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error);
127: } else {
128: PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error);
129: }
130: }
131: PetscPrintf(PETSC_COMM_WORLD,"\n");
132: }
134: /*
135: Free work space
136: */
137: EPSDestroy(&eps);
138: MatDestroy(&A);
139: VecDestroy(&xr);
140: VecDestroy(&xi);
141: SlepcFinalize();
142: return ierr;
143: }