Actual source code: test18.c
slepc-3.17.2 2022-08-09
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[] = "Symmetric-indefinite eigenproblem.\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; /* problem matrices */
21: EPS eps; /* eigenproblem solver context */
22: PetscInt N,n=10,m,Istart,Iend,II,nev,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: if (!flag) m=n;
29: N = n*m;
30: PetscPrintf(PETSC_COMM_WORLD,"\nSymmetric-indefinite eigenproblem, N=%" PetscInt_FMT "\n\n",N);
32: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: Compute the matrices that define the eigensystem, Ax=kBx
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36: MatCreate(PETSC_COMM_WORLD,&A);
37: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
38: MatSetFromOptions(A);
39: MatSetUp(A);
41: MatCreate(PETSC_COMM_WORLD,&B);
42: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
43: MatSetFromOptions(B);
44: MatSetUp(B);
46: MatGetOwnershipRange(A,&Istart,&Iend);
47: for (II=Istart;II<Iend;II++) {
48: i = II/n; j = II-i*n;
49: if (i>0) MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);
50: if (i<m-1) MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);
51: if (j>0) MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);
52: if (j<n-1) MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);
53: MatSetValue(A,II,II,4.0,INSERT_VALUES);
54: MatSetValue(B,II,N-II-1,1.0,INSERT_VALUES);
55: }
57: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
58: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
60: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Create the eigensolver and set various options
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: EPSCreate(PETSC_COMM_WORLD,&eps);
67: EPSSetOperators(eps,A,B);
68: EPSSetProblemType(eps,EPS_GHIEP);
69: EPSSetFromOptions(eps);
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Solve the eigensystem
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: EPSSolve(eps);
76: EPSGetDimensions(eps,&nev,NULL,NULL);
77: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Display solution and clean up
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: /* show detailed info unless -terse option is given by user */
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: EPSDestroy(&eps);
94: MatDestroy(&A);
95: MatDestroy(&B);
96: SlepcFinalize();
97: return 0;
98: }
100: /*TEST
102: testset:
103: args: -eps_nev 4 -eps_ncv 12 -terse -st_type sinvert -eps_krylovschur_restart .3
104: requires: !single
105: output_file: output/test18_1.out
106: test:
107: suffix: 1_ks
108: test:
109: suffix: 1_ks_gnhep
110: args: -eps_gen_non_hermitian
111: requires: !__float128
112: test:
113: suffix: 2_cuda_ks
114: args: -mat_type aijcusparse
115: requires: cuda
116: test:
117: suffix: 2_cuda_ks_gnhep
118: args: -eps_gen_non_hermitian -mat_type aijcusparse
119: requires: cuda
121: test:
122: suffix: 2
123: args: -n 10 -m 11 -eps_type {{gd jd}} -eps_target 0.2 -eps_harmonic -eps_nev 2 -eps_ncv 11 -terse
124: requires: !single
125: filter: sed -e "s/[+-]0\.0*i//g"
127: TEST*/