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 : Example based on spring problem in NLEVP collection [1]. See the parameters
12 : meaning at Example 2 in [2].
13 :
14 : [1] T. Betcke, N. J. Higham, V. Mehrmann, C. Schroder, and F. Tisseur,
15 : NLEVP: A Collection of Nonlinear Eigenvalue Problems, MIMS EPrint
16 : 2010.98, November 2010.
17 : [2] F. Tisseur, Backward error and condition of polynomial eigenvalue
18 : problems, Linear Algebra and its Applications, 309 (2000), pp. 339--361,
19 : April 2000.
20 : */
21 :
22 : static char help[] = "Tests multiple calls to PEPSolve with different matrix of different size.\n\n"
23 : "This is based on the spring problem from NLEVP collection.\n\n"
24 : "The command line options are:\n"
25 : " -n <n> ... number of grid subdivisions.\n"
26 : " -mu <value> ... mass (default 1).\n"
27 : " -tau <value> ... damping constant of the dampers (default 10).\n"
28 : " -kappa <value> ... damping constant of the springs (default 5).\n"
29 : " -initv ... set an initial vector.\n\n";
30 :
31 : #include <slepcpep.h>
32 :
33 4 : int main(int argc,char **argv)
34 : {
35 4 : Mat M,C,K,A[3]; /* problem matrices */
36 4 : PEP pep; /* polynomial eigenproblem solver context */
37 4 : PetscInt n=30,Istart,Iend,i,nev;
38 4 : PetscReal mu=1.0,tau=10.0,kappa=5.0;
39 4 : PetscBool terse=PETSC_FALSE;
40 :
41 4 : PetscFunctionBeginUser;
42 4 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
43 :
44 4 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
45 4 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
46 4 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
47 4 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
48 :
49 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50 : Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
51 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52 :
53 : /* K is a tridiagonal */
54 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
55 4 : PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
56 4 : PetscCall(MatSetFromOptions(K));
57 :
58 4 : PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
59 124 : for (i=Istart;i<Iend;i++) {
60 120 : if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
61 120 : PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
62 120 : if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
63 : }
64 :
65 4 : PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
66 4 : PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
67 :
68 : /* C is a tridiagonal */
69 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
70 4 : PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
71 4 : PetscCall(MatSetFromOptions(C));
72 :
73 4 : PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
74 124 : for (i=Istart;i<Iend;i++) {
75 120 : if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
76 120 : PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
77 120 : if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
78 : }
79 :
80 4 : PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
81 4 : PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
82 :
83 : /* M is a diagonal matrix */
84 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
85 4 : PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
86 4 : PetscCall(MatSetFromOptions(M));
87 4 : PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
88 124 : for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
89 4 : PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
90 4 : PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
91 :
92 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93 : Create the eigensolver and set various options
94 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95 :
96 4 : PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
97 4 : A[0] = K; A[1] = C; A[2] = M;
98 4 : PetscCall(PEPSetOperators(pep,3,A));
99 4 : PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
100 4 : PetscCall(PEPSetTolerances(pep,PETSC_SMALL,PETSC_CURRENT));
101 4 : PetscCall(PEPSetFromOptions(pep));
102 :
103 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104 : Solve the eigensystem
105 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106 :
107 4 : PetscCall(PEPSolve(pep));
108 4 : PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
109 4 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
110 :
111 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112 : Display solution of first solve
113 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114 4 : PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
115 4 : if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
116 : else {
117 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
118 0 : PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
119 0 : PetscCall(PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
120 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
121 : }
122 4 : PetscCall(MatDestroy(&M));
123 4 : PetscCall(MatDestroy(&C));
124 4 : PetscCall(MatDestroy(&K));
125 :
126 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127 : Compute the eigensystem, (k^2*M+k*C+K)x=0 for bigger n
128 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129 :
130 4 : n *= 2;
131 : /* K is a tridiagonal */
132 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
133 4 : PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
134 4 : PetscCall(MatSetFromOptions(K));
135 :
136 4 : PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
137 244 : for (i=Istart;i<Iend;i++) {
138 240 : if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
139 240 : PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
140 240 : if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
141 : }
142 :
143 4 : PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
144 4 : PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
145 :
146 : /* C is a tridiagonal */
147 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
148 4 : PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
149 4 : PetscCall(MatSetFromOptions(C));
150 :
151 4 : PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
152 244 : for (i=Istart;i<Iend;i++) {
153 240 : if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
154 240 : PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
155 240 : if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
156 : }
157 :
158 4 : PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
159 4 : PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
160 :
161 : /* M is a diagonal matrix */
162 4 : PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
163 4 : PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
164 4 : PetscCall(MatSetFromOptions(M));
165 4 : PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
166 244 : for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
167 4 : PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
168 4 : PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
169 :
170 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171 : Solve again, calling PEPReset() since matrix size has changed
172 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173 : /* PetscCall(PEPReset(pep)); */ /* not required, will be called in PEPSetOperators() */
174 4 : A[0] = K; A[1] = C; A[2] = M;
175 4 : PetscCall(PEPSetOperators(pep,3,A));
176 4 : PetscCall(PEPSolve(pep));
177 :
178 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179 : Display solution and clean up
180 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181 4 : if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
182 : else {
183 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
184 0 : PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
185 0 : PetscCall(PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
186 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
187 : }
188 4 : PetscCall(PEPDestroy(&pep));
189 4 : PetscCall(MatDestroy(&M));
190 4 : PetscCall(MatDestroy(&C));
191 4 : PetscCall(MatDestroy(&K));
192 4 : PetscCall(SlepcFinalize());
193 : return 0;
194 : }
195 :
196 : /*TEST
197 :
198 : test:
199 : suffix: 1
200 : args: -pep_type {{toar qarnoldi linear}} -pep_nev 4 -terse
201 : requires: double
202 :
203 : test:
204 : suffix: 2
205 : args: -pep_type stoar -pep_hermitian -pep_nev 4 -terse
206 : requires: !single
207 :
208 : TEST*/
|