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[] = "Test the solution of a PEP without calling PEPSetFromOptions (based on ex16.c).\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"
15 : " -type <pep_type> = pep type to test.\n"
16 : " -epstype <eps_type> = eps type to test (for linear).\n\n";
17 :
18 : #include <slepcpep.h>
19 :
20 4 : int main(int argc,char **argv)
21 : {
22 4 : Mat M,C,K,A[3]; /* problem matrices */
23 4 : PEP pep; /* polynomial eigenproblem solver context */
24 4 : PetscInt N,n=10,m,Istart,Iend,II,nev,i,j;
25 4 : PetscReal keep;
26 4 : PetscBool flag,isgd2,epsgiven,lock;
27 4 : char peptype[30] = "linear",epstype[30] = "";
28 4 : EPS eps;
29 4 : ST st;
30 4 : KSP ksp;
31 4 : PC pc;
32 :
33 4 : PetscFunctionBeginUser;
34 4 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
35 :
36 4 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
37 4 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
38 4 : if (!flag) m=n;
39 4 : N = n*m;
40 4 : PetscCall(PetscOptionsGetString(NULL,NULL,"-type",peptype,sizeof(peptype),NULL));
41 4 : PetscCall(PetscOptionsGetString(NULL,NULL,"-epstype",epstype,sizeof(epstype),&epsgiven));
42 4 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
43 :
44 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45 : Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
46 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47 :
48 : /* K is the 2-D Laplacian */
49 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
50 4 : PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
51 4 : PetscCall(MatSetFromOptions(K));
52 4 : PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
53 444 : for (II=Istart;II<Iend;II++) {
54 440 : i = II/n; j = II-i*n;
55 440 : if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
56 440 : if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
57 440 : if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
58 440 : if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
59 440 : PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
60 : }
61 4 : PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
62 4 : PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
63 :
64 : /* C is the 1-D Laplacian on horizontal lines */
65 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
66 4 : PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
67 4 : PetscCall(MatSetFromOptions(C));
68 4 : PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
69 444 : for (II=Istart;II<Iend;II++) {
70 440 : i = II/n; j = II-i*n;
71 440 : if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
72 440 : if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
73 440 : PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
74 : }
75 4 : PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
76 4 : PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
77 :
78 : /* M is a diagonal matrix */
79 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
80 4 : PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
81 4 : PetscCall(MatSetFromOptions(M));
82 4 : PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
83 444 : for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
84 4 : PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
85 4 : PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
86 :
87 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88 : Create the eigensolver and set various options
89 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90 :
91 4 : PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
92 4 : A[0] = K; A[1] = C; A[2] = M;
93 4 : PetscCall(PEPSetOperators(pep,3,A));
94 4 : PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
95 4 : PetscCall(PEPSetDimensions(pep,4,20,PETSC_DETERMINE));
96 4 : PetscCall(PEPSetTolerances(pep,PETSC_SMALL,PETSC_CURRENT));
97 :
98 : /*
99 : Set solver type at runtime
100 : */
101 4 : PetscCall(PEPSetType(pep,peptype));
102 4 : if (epsgiven) {
103 1 : PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&flag));
104 1 : if (flag) {
105 1 : PetscCall(PEPLinearGetEPS(pep,&eps));
106 1 : PetscCall(PetscStrcmp(epstype,"gd2",&isgd2));
107 1 : if (isgd2) {
108 0 : PetscCall(EPSSetType(eps,EPSGD));
109 0 : PetscCall(EPSGDSetDoubleExpansion(eps,PETSC_TRUE));
110 1 : } else PetscCall(EPSSetType(eps,epstype));
111 1 : PetscCall(EPSGetST(eps,&st));
112 1 : PetscCall(STGetKSP(st,&ksp));
113 1 : PetscCall(KSPGetPC(ksp,&pc));
114 1 : PetscCall(PCSetType(pc,PCJACOBI));
115 1 : PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSGD,&flag));
116 : }
117 1 : PetscCall(PEPLinearSetExplicitMatrix(pep,PETSC_TRUE));
118 : }
119 4 : PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPQARNOLDI,&flag));
120 4 : if (flag) {
121 1 : PetscCall(STCreate(PETSC_COMM_WORLD,&st));
122 1 : PetscCall(STSetTransform(st,PETSC_TRUE));
123 1 : PetscCall(PEPSetST(pep,st));
124 1 : PetscCall(STDestroy(&st));
125 1 : PetscCall(PEPQArnoldiGetRestart(pep,&keep));
126 1 : PetscCall(PEPQArnoldiGetLocking(pep,&lock));
127 1 : if (!lock && keep<0.6) PetscCall(PEPQArnoldiSetRestart(pep,0.6));
128 : }
129 4 : PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPTOAR,&flag));
130 4 : if (flag) {
131 1 : PetscCall(PEPTOARGetRestart(pep,&keep));
132 1 : PetscCall(PEPTOARGetLocking(pep,&lock));
133 1 : if (!lock && keep<0.6) PetscCall(PEPTOARSetRestart(pep,0.6));
134 : }
135 :
136 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137 : Solve the eigensystem
138 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139 :
140 4 : PetscCall(PEPSolve(pep));
141 4 : PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
142 4 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
143 :
144 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145 : Display solution and clean up
146 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147 :
148 4 : PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
149 4 : PetscCall(PEPDestroy(&pep));
150 4 : PetscCall(MatDestroy(&M));
151 4 : PetscCall(MatDestroy(&C));
152 4 : PetscCall(MatDestroy(&K));
153 4 : PetscCall(SlepcFinalize());
154 : return 0;
155 : }
156 :
157 : /*TEST
158 :
159 : testset:
160 : args: -m 11
161 : output_file: output/test1_1.out
162 : filter: sed -e "s/1.16403/1.16404/g" | sed -e "s/1.65362i/1.65363i/g" | sed -e "s/-1.16404-1.65363i, -1.16404+1.65363i/-1.16404+1.65363i, -1.16404-1.65363i/" | sed -e "s/-0.51784-1.31039i, -0.51784+1.31039i/-0.51784+1.31039i, -0.51784-1.31039i/"
163 : requires: !single
164 : test:
165 : suffix: 1
166 : args: -type {{toar qarnoldi linear}}
167 : test:
168 : suffix: 1_linear_gd
169 : args: -type linear -epstype gd
170 : requires: !__float128
171 :
172 : TEST*/
|