Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Tests a HEP problem with Hermitian matrix.\n\n";
12 :
13 : #include <slepceps.h>
14 :
15 9 : int main(int argc,char **argv)
16 : {
17 9 : Mat A; /* matrix */
18 9 : EPS eps; /* eigenproblem solver context */
19 9 : PetscInt N,n=20,m,Istart,Iend,II,i,j;
20 9 : PetscBool flag;
21 :
22 9 : PetscFunctionBeginUser;
23 9 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
24 9 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
25 9 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
26 9 : if (!flag) m=n;
27 9 : N = n*m;
28 9 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nHermitian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
29 : #if !defined(PETSC_USE_COMPLEX)
30 : SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example requires complex scalars!");
31 : #endif
32 :
33 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34 : Compute the matrix that defines the eigensystem, Ax=kx
35 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36 :
37 9 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
38 9 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
39 9 : PetscCall(MatSetFromOptions(A));
40 :
41 9 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
42 3261 : for (II=Istart;II<Iend;II++) {
43 3252 : i = II/n; j = II-i*n;
44 3252 : if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0-0.1*PETSC_i,INSERT_VALUES));
45 3252 : if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0+0.1*PETSC_i,INSERT_VALUES));
46 3252 : if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0-0.1*PETSC_i,INSERT_VALUES));
47 3252 : if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0+0.1*PETSC_i,INSERT_VALUES));
48 3252 : PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
49 : }
50 9 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
51 9 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
52 :
53 9 : PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
54 :
55 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56 : Create the eigensolver and solve the problem
57 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58 :
59 9 : PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
60 9 : PetscCall(EPSSetOperators(eps,A,NULL));
61 9 : PetscCall(EPSSetProblemType(eps,EPS_HEP));
62 9 : PetscCall(EPSSetFromOptions(eps));
63 9 : PetscCall(EPSSolve(eps));
64 9 : PetscCall(EPSErrorView(eps,EPS_ERROR_BACKWARD,NULL));
65 :
66 9 : PetscCall(EPSDestroy(&eps));
67 9 : PetscCall(MatDestroy(&A));
68 9 : PetscCall(SlepcFinalize());
69 : return 0;
70 : }
71 :
72 : /*TEST
73 :
74 : build:
75 : requires: complex
76 :
77 : testset:
78 : args: -m 18 -n 19 -eps_nev 4 -eps_max_it 1000
79 : requires: !single complex
80 : output_file: output/test36_1.out
81 : test:
82 : suffix: 1
83 : args: -eps_type {{krylovschur subspace arnoldi gd jd lapack}}
84 : test:
85 : suffix: 1_elemental
86 : args: -eps_type elemental
87 : requires: elemental
88 :
89 : test:
90 : suffix: 2
91 : args: -eps_nev 4 -eps_smallest_real -eps_type {{lobpcg rqcg lapack}}
92 : requires: !single complex
93 :
94 : TEST*/
|