Actual source code: test34.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[] = "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: PetscFunctionBeginUser;
30: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
32: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
33: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
34: if (!flag) m=n;
35: N = n*m;
36: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nStandard eigenproblem with PRIMME, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
38: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39: Compute the matrices that define the eigensystem, Ax=kBx
40: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
42: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
43: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
44: PetscCall(MatSetFromOptions(A));
46: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
47: for (II=Istart;II<Iend;II++) {
48: i = II/n; j = II-i*n;
49: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
50: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
51: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
52: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
53: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
54: }
56: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
57: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create the eigensolver and set various options
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
64: PetscCall(EPSSetOperators(eps,A,NULL));
65: PetscCall(EPSSetProblemType(eps,EPS_HEP));
66: PetscCall(EPSSetType(eps,EPSPRIMME));
68: /*
69: Set several options
70: */
71: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
72: PetscCall(EPSGetST(eps,&st));
73: PetscCall(STSetType(st,STPRECOND));
74: PetscCall(STGetKSP(st,&ksp));
75: PetscCall(KSPGetPC(ksp,&pc));
76: PetscCall(KSPSetType(ksp,KSPPREONLY));
77: PetscCall(PCSetType(pc,PCICC));
79: PetscCall(EPSPRIMMESetBlockSize(eps,4));
80: PetscCall(EPSPRIMMESetMethod(eps,EPS_PRIMME_GD_OLSEN_PLUSK));
81: PetscCall(EPSSetFromOptions(eps));
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Compute eigenvalues and display info
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: PetscCall(EPSSolve(eps));
88: PetscCall(EPSPRIMMEGetBlockSize(eps,&bs));
89: PetscCall(EPSPRIMMEGetMethod(eps,&meth));
90: PetscCall(PetscPrintf(PETSC_COMM_WORLD," PRIMME: using block size %" PetscInt_FMT ", method %s\n",bs,EPSPRIMMEMethods[meth]));
92: PetscCall(EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL));
94: PetscCall(EPSDestroy(&eps));
95: PetscCall(MatDestroy(&A));
96: PetscCall(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*/