Actual source code: test28.c
slepc-3.20.2 2024-03-15
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: PetscFunctionBeginUser;
26: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
27: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
28: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
29: N = n*m;
30: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
32: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: Create the 2-D Laplacian with coarse mesh
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
37: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
38: PetscCall(MatSetFromOptions(A));
39: PetscCall(MatSetUp(A));
40: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
41: for (II=Istart;II<Iend;II++) {
42: i = II/n; j = II-i*n;
43: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
44: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
45: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
46: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
47: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
48: }
49: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
50: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Create the eigensolver, set options and solve the eigensystem
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
57: PetscCall(EPSSetOperators(eps,A,NULL));
58: PetscCall(EPSSetProblemType(eps,EPS_HEP));
59: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
60: PetscCall(EPSSetDimensions(eps,nev,PETSC_DEFAULT,PETSC_DEFAULT));
61: PetscCall(EPSSetFromOptions(eps));
63: PetscCall(EPSSolve(eps));
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Display solution of first solve
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
70: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
71: else {
72: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
73: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
74: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
75: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
76: }
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Create the 2-D Laplacian with finer mesh
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: n *= 2;
83: m *= 2;
84: N = n*m;
85: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
87: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
88: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
89: PetscCall(MatSetFromOptions(B));
90: PetscCall(MatSetUp(B));
91: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
92: for (II=Istart;II<Iend;II++) {
93: i = II/n; j = II-i*n;
94: if (i>0) PetscCall(MatSetValue(B,II,II-n,-1.0,INSERT_VALUES));
95: if (i<m-1) PetscCall(MatSetValue(B,II,II+n,-1.0,INSERT_VALUES));
96: if (j>0) PetscCall(MatSetValue(B,II,II-1,-1.0,INSERT_VALUES));
97: if (j<n-1) PetscCall(MatSetValue(B,II,II+1,-1.0,INSERT_VALUES));
98: PetscCall(MatSetValue(B,II,II,4.0,INSERT_VALUES));
99: }
100: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
101: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Solve again, calling EPSReset() since matrix size has changed
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: /* PetscCall(EPSReset(eps)); */ /* not required, will be called in EPSSetOperators() */
108: PetscCall(EPSSetOperators(eps,B,NULL));
109: PetscCall(EPSSolve(eps));
111: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
112: else {
113: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
114: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
115: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
116: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
117: }
119: PetscCall(EPSDestroy(&eps));
120: PetscCall(MatDestroy(&A));
121: PetscCall(MatDestroy(&B));
122: PetscCall(SlepcFinalize());
123: return 0;
124: }
126: /*TEST
128: testset:
129: requires: !single
130: output_file: output/test28_1.out
131: test:
132: suffix: 1
133: args: -eps_type {{krylovschur arnoldi gd jd rqcg lobpcg lapack}} -terse
134: test:
135: suffix: 1_lanczos
136: args: -eps_type lanczos -eps_lanczos_reorthog local -terse
138: test:
139: suffix: 2
140: args: -eps_type {{power subspace}} -eps_target 8 -st_type sinvert -terse
142: testset:
143: args: -eps_interval 0.5,0.67 -terse
144: output_file: output/test28_3.out
145: test:
146: suffix: 3
147: args: -st_type sinvert -st_pc_type cholesky
148: requires: !single
149: test:
150: suffix: 3_evsl
151: nsize: {{1 2}}
152: args: -eps_type evsl -eps_evsl_dos_method lanczos
153: requires: evsl
155: TEST*/