Actual source code: test27.c
slepc-3.17.1 2022-04-11
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[] = "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;
25: Vec v,*X;
27: SlepcInitialize(&argc,&argv,(char*)0,help);
28: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
29: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
30: if (!flag) m=n;
31: N = n*m;
32: PetscOptionsGetInt(NULL,NULL,"-neigs",&neigs,NULL);
33: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid), neigs=%" PetscInt_FMT "\n\n",N,n,m,neigs);
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Create the 2-D Laplacian
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: MatCreate(PETSC_COMM_WORLD,&A);
40: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
41: MatSetFromOptions(A);
42: MatSetUp(A);
43: MatGetOwnershipRange(A,&Istart,&Iend);
44: for (II=Istart;II<Iend;II++) {
45: i = II/n; j = II-i*n;
46: if (i>0) MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);
47: if (i<m-1) MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);
48: if (j>0) MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);
49: if (j<n-1) MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);
50: MatSetValue(A,II,II,4.0,INSERT_VALUES);
51: }
52: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
53: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
54: MatCreateVecs(A,&v,NULL);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Create the eigensolver and set various options
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: EPSCreate(PETSC_COMM_WORLD,&eps);
61: EPSSetOperators(eps,A,NULL);
62: EPSSetProblemType(eps,EPS_HEP);
63: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
64: EPSSetDimensions(eps,neigs,PETSC_DEFAULT,PETSC_DEFAULT);
65: EPSSetFromOptions(eps);
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Solve the eigensystem for the first time
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: EPSSolve(eps);
72: EPSGetDimensions(eps,&nev,NULL,NULL);
73: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
75: EPSGetConverged(eps,&nconv);
77: VecDuplicateVecs(v,neigs,&X);
78: for (i=0;i<neigs;i++) EPSGetEigenvector(eps,i,X[i],NULL);
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Display solution of first solve
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
85: if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
86: else {
87: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
88: EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
89: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
90: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
91: }
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Solve the eigensystem again, feeding the initial vectors
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: PetscPrintf(PETSC_COMM_WORLD," Solving again with eigenvectors as initial guesses\n");
98: EPSSetInitialSpace(eps,neigs,X);
99: EPSSolve(eps);
101: if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
102: else {
103: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
104: EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
105: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
106: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
107: }
109: VecDestroy(&v);
110: VecDestroyVecs(neigs,&X);
111: EPSDestroy(&eps);
112: MatDestroy(&A);
113: SlepcFinalize();
114: return 0;
115: }
117: /*TEST
119: test:
120: suffix: 1
121: args: -eps_type {{gd jd rqcg lobpcg}} -terse
123: testset:
124: args: -eps_interval .17,1.3 -terse
125: filter: grep -v "requested"
126: output_file: output/test27_2.out
127: test:
128: suffix: 2
129: args: -st_type filter -st_filter_degree 150 -eps_nev 1
130: requires: !single
131: test:
132: suffix: 2_evsl
133: nsize: {{1 2}}
134: args: -eps_type evsl
135: requires: evsl
137: TEST*/