Actual source code: test27.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[] = "Illustrates feeding exact eigenvectors as initial vectors of a second solve.\n\n"
12: "Based on ex2.\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
15: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
17: #include <slepceps.h>
19: int main(int argc,char **argv)
20: {
21: Mat A;
22: EPS eps;
23: PetscInt N,n=10,m,Istart,Iend,II,nev,nconv,i,j,neigs=5;
24: PetscBool flag,terse;
26: Vec v,*X;
28: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
29: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
30: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
31: if (!flag) m=n;
32: N = n*m;
33: PetscOptionsGetInt(NULL,NULL,"-neigs",&neigs,NULL);
34: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%D (%Dx%D grid), neigs=%D\n\n",N,n,m,neigs);
36: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: Create the 2-D Laplacian
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: MatCreate(PETSC_COMM_WORLD,&A);
41: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
42: MatSetFromOptions(A);
43: MatSetUp(A);
44: MatGetOwnershipRange(A,&Istart,&Iend);
45: for (II=Istart;II<Iend;II++) {
46: i = II/n; j = II-i*n;
47: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
48: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
49: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
50: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
51: MatSetValue(A,II,II,4.0,INSERT_VALUES);
52: }
53: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
55: MatCreateVecs(A,&v,NULL);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create the eigensolver and set various options
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: EPSCreate(PETSC_COMM_WORLD,&eps);
62: EPSSetOperators(eps,A,NULL);
63: EPSSetProblemType(eps,EPS_HEP);
64: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
65: EPSSetDimensions(eps,neigs,PETSC_DEFAULT,PETSC_DEFAULT);
66: EPSSetFromOptions(eps);
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Solve the eigensystem for the first time
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: EPSSolve(eps);
73: EPSGetDimensions(eps,&nev,NULL,NULL);
74: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
76: EPSGetConverged(eps,&nconv);
77: if (nconv<neigs) SETERRQ2(PETSC_COMM_WORLD,1,"Only %D eigenvalues have converged, %D requested",nconv,neigs);
78: VecDuplicateVecs(v,neigs,&X);
79: for (i=0;i<neigs;i++) {
80: EPSGetEigenvector(eps,i,X[i],NULL);
81: }
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Display solution of first solve
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
88: if (terse) {
89: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
90: } else {
91: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
92: EPSReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
93: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
94: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
95: }
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Solve the eigensystem again, feeding the initial vectors
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101: PetscPrintf(PETSC_COMM_WORLD," Solving again with eigenvectors as initial guesses\n");
102: EPSSetInitialSpace(eps,neigs,X);
103: EPSSolve(eps);
105: if (terse) {
106: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
107: } else {
108: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
109: EPSReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
110: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
111: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
112: }
114: VecDestroy(&v);
115: VecDestroyVecs(neigs,&X);
116: EPSDestroy(&eps);
117: MatDestroy(&A);
118: SlepcFinalize();
119: return ierr;
120: }
122: /*TEST
124: test:
125: suffix: 1
126: args: -eps_type {{gd jd rqcg lobpcg}} -terse
127: requires: !single
128: output_file: output/test27_1.out
130: test:
131: suffix: 2
132: args: -eps_interval .17,1.3 -st_type filter -eps_nev 1 -terse
133: filter: grep -v "requested"
134: requires: !single
136: TEST*/