Actual source code: test28.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[] = "Tests multiple calls to EPSSolve with different matrix of different size.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
14: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
16: #include <slepceps.h>
18: int main(int argc,char **argv)
19: {
20: Mat A,B;
21: EPS eps;
22: PetscInt N,n=10,m=11,Istart,Iend,II,nev=3,i,j;
23: PetscBool flag,terse;
25: SlepcInitialize(&argc,&argv,(char*)0,help);
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
28: N = n*m;
29: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
31: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32: Create the 2-D Laplacian with coarse mesh
33: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
35: MatCreate(PETSC_COMM_WORLD,&A);
36: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
37: MatSetFromOptions(A);
38: MatSetUp(A);
39: MatGetOwnershipRange(A,&Istart,&Iend);
40: for (II=Istart;II<Iend;II++) {
41: i = II/n; j = II-i*n;
42: if (i>0) MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);
43: if (i<m-1) MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);
44: if (j>0) MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);
45: if (j<n-1) MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);
46: MatSetValue(A,II,II,4.0,INSERT_VALUES);
47: }
48: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
49: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Create the eigensolver, set options and solve the eigensystem
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: EPSCreate(PETSC_COMM_WORLD,&eps);
56: EPSSetOperators(eps,A,NULL);
57: EPSSetProblemType(eps,EPS_HEP);
58: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
59: EPSSetDimensions(eps,nev,PETSC_DEFAULT,PETSC_DEFAULT);
60: EPSSetFromOptions(eps);
62: EPSSolve(eps);
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Display solution of first solve
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
69: if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
70: else {
71: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
72: EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
73: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
74: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
75: }
77: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: Create the 2-D Laplacian with finer mesh
79: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: n *= 2;
82: m *= 2;
83: N = n*m;
84: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
86: MatCreate(PETSC_COMM_WORLD,&B);
87: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
88: MatSetFromOptions(B);
89: MatSetUp(B);
90: MatGetOwnershipRange(B,&Istart,&Iend);
91: for (II=Istart;II<Iend;II++) {
92: i = II/n; j = II-i*n;
93: if (i>0) MatSetValue(B,II,II-n,-1.0,INSERT_VALUES);
94: if (i<m-1) MatSetValue(B,II,II+n,-1.0,INSERT_VALUES);
95: if (j>0) MatSetValue(B,II,II-1,-1.0,INSERT_VALUES);
96: if (j<n-1) MatSetValue(B,II,II+1,-1.0,INSERT_VALUES);
97: MatSetValue(B,II,II,4.0,INSERT_VALUES);
98: }
99: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
100: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Solve again, calling EPSReset() since matrix size has changed
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: /*EPSReset(eps);*/ /* not required, will be called in EPSSetOperators() */
107: EPSSetOperators(eps,B,NULL);
108: EPSSolve(eps);
110: if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
111: else {
112: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
113: EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
114: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
115: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
116: }
118: EPSDestroy(&eps);
119: MatDestroy(&A);
120: MatDestroy(&B);
121: SlepcFinalize();
122: return 0;
123: }
125: /*TEST
127: testset:
128: requires: !single
129: output_file: output/test28_1.out
130: test:
131: suffix: 1
132: args: -eps_type {{krylovschur arnoldi gd jd rqcg lobpcg lapack}} -terse
133: test:
134: suffix: 1_lanczos
135: args: -eps_type lanczos -eps_lanczos_reorthog local -terse
137: test:
138: suffix: 2
139: args: -eps_type {{power subspace}} -eps_target 8 -st_type sinvert -terse
141: testset:
142: args: -eps_interval 0.5,0.67 -terse
143: output_file: output/test28_3.out
144: test:
145: suffix: 3
146: args: -st_type sinvert -st_pc_type cholesky
147: requires: !single
148: test:
149: suffix: 3_evsl
150: nsize: {{1 2}}
151: args: -eps_type evsl -eps_evsl_dos_method lanczos
152: requires: evsl
154: TEST*/