Actual source code: ex13.c
slepc-3.22.1 2024-10-28
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[] = "Generalized Symmetric eigenproblem.\n\n"
12: "The problem is Ax = lambda Bx, with:\n"
13: " A = Laplacian operator in 2-D\n"
14: " B = diagonal matrix with all values equal to 4 except nulldim zeros\n\n"
15: "The command line options are:\n"
16: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
17: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
18: " -nulldim <k>, where <k> = dimension of the nullspace of B.\n\n";
20: #include <slepceps.h>
22: int main(int argc,char **argv)
23: {
24: Mat A,B; /* matrices */
25: EPS eps; /* eigenproblem solver context */
26: EPSType type;
27: PetscInt N,n=10,m,Istart,Iend,II,nev,i,j,nulldim=0;
28: PetscBool flag,terse;
30: PetscFunctionBeginUser;
31: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
33: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
34: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
35: if (!flag) m=n;
36: N = n*m;
37: PetscCall(PetscOptionsGetInt(NULL,NULL,"-nulldim",&nulldim,NULL));
38: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid), null(B)=%" PetscInt_FMT "\n\n",N,n,m,nulldim));
40: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41: Compute the matrices that define the eigensystem, Ax=kBx
42: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
44: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
45: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
46: PetscCall(MatSetFromOptions(A));
48: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
49: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
50: PetscCall(MatSetFromOptions(B));
52: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
53: for (II=Istart;II<Iend;II++) {
54: i = II/n; j = II-i*n;
55: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
56: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
57: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
58: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
59: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
60: if (II>=nulldim) PetscCall(MatSetValue(B,II,II,4.0,INSERT_VALUES));
61: }
63: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
64: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
65: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
66: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Create the eigensolver and set various options
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: /*
73: Create eigensolver context
74: */
75: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
77: /*
78: Set operators. In this case, it is a generalized eigenvalue problem
79: */
80: PetscCall(EPSSetOperators(eps,A,B));
81: PetscCall(EPSSetProblemType(eps,EPS_GHEP));
83: /*
84: Set solver parameters at runtime
85: */
86: PetscCall(EPSSetFromOptions(eps));
88: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89: Solve the eigensystem
90: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: PetscCall(EPSSolve(eps));
94: /*
95: Optional: Get some information from the solver and display it
96: */
97: PetscCall(EPSGetType(eps,&type));
98: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
99: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
100: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Display solution and clean up
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: /* show detailed info unless -terse option is given by user */
107: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
108: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
109: else {
110: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
111: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
112: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
113: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
114: }
115: PetscCall(EPSDestroy(&eps));
116: PetscCall(MatDestroy(&A));
117: PetscCall(MatDestroy(&B));
118: PetscCall(SlepcFinalize());
119: return 0;
120: }
122: /*TEST
124: test:
125: suffix: 1
126: args: -eps_nev 4 -eps_ncv 22 -eps_tol 1e-5 -st_type sinvert -terse
127: filter: grep -v Solution
129: test:
130: suffix: 2
131: args: -n 110 -nulldim 6 -eps_nev 4 -eps_ncv 18 -eps_tol 1e-5 -eps_purify 1 -st_type sinvert -st_matstructure {{different subset}} -terse
132: requires: !single
134: test:
135: suffix: 3
136: args: -eps_nev 3 -eps_tol 1e-5 -mat_type sbaij -st_type sinvert -terse
138: test:
139: suffix: 4
140: args: -eps_nev 4 -eps_tol 1e-4 -eps_smallest_real -eps_type {{gd lobpcg rqcg}} -terse
141: output_file: output/ex13_1.out
142: filter: grep -v Solution
144: test:
145: suffix: 5_primme
146: args: -n 10 -m 12 -eps_nev 4 -eps_target 0.9 -eps_max_it 15000 -eps_type primme -st_pc_type jacobi -terse
147: requires: primme defined(SLEPC_HAVE_PRIMME3) !single
149: test:
150: suffix: 6
151: nsize: 2
152: args: -eps_type ciss -rg_type ellipse -rg_ellipse_center 1.4 -rg_ellipse_radius 0.1 -eps_ciss_partitions 2 -terse
153: requires: !single
155: TEST*/