Actual source code: test34.c
slepc-3.17.0 2022-03-31
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[] = "Test interface to external package PRIMME.\n\n"
12: "This is based on ex3.c. 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; /* matrix */
21: EPS eps; /* eigenproblem solver context */
22: ST st; /* spectral transformation context */
23: KSP ksp;
24: PC pc;
25: PetscInt N,n=35,m,Istart,Iend,II,i,j,bs;
26: PetscBool flag;
27: EPSPRIMMEMethod meth;
29: SlepcInitialize(&argc,&argv,(char*)0,help);
31: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
32: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
33: if (!flag) m=n;
34: N = n*m;
35: PetscPrintf(PETSC_COMM_WORLD,"\nStandard eigenproblem with PRIMME, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
37: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: Compute the matrices that define the eigensystem, Ax=kBx
39: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
41: MatCreate(PETSC_COMM_WORLD,&A);
42: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
43: MatSetFromOptions(A);
44: MatSetUp(A);
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: }
56: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create the eigensolver and set various options
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: EPSCreate(PETSC_COMM_WORLD,&eps);
64: EPSSetOperators(eps,A,NULL);
65: EPSSetProblemType(eps,EPS_HEP);
66: EPSSetType(eps,EPSPRIMME);
68: /*
69: Set several options
70: */
71: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
72: EPSGetST(eps,&st);
73: STSetType(st,STPRECOND);
74: STGetKSP(st,&ksp);
75: KSPGetPC(ksp,&pc);
76: KSPSetType(ksp,KSPPREONLY);
77: PCSetType(pc,PCICC);
79: EPSPRIMMESetBlockSize(eps,4);
80: EPSPRIMMESetMethod(eps,EPS_PRIMME_GD_OLSEN_PLUSK);
81: EPSSetFromOptions(eps);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Compute eigenvalues and display info
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: EPSSolve(eps);
88: EPSPRIMMEGetBlockSize(eps,&bs);
89: EPSPRIMMEGetMethod(eps,&meth);
90: PetscPrintf(PETSC_COMM_WORLD," PRIMME: using block size %" PetscInt_FMT ", method %s\n",bs,EPSPRIMMEMethods[meth]);
92: EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL);
94: EPSDestroy(&eps);
95: MatDestroy(&A);
96: SlepcFinalize();
97: return 0;
98: }
100: /*TEST
102: build:
103: requires: primme
105: testset:
106: args: -eps_nev 4
107: requires: primme
108: output_file: output/test34_1.out
109: test:
110: suffix: 1
111: test:
112: suffix: 2
113: args: -st_pc_type bjacobi -eps_target 0.01 -eps_target_real -eps_refined
114: nsize: 2
115: test:
116: suffix: 3
117: args: -eps_smallest_magnitude -eps_harmonic
119: TEST*/