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