Actual source code: sleeper.c
slepc-3.6.3 2016-03-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
21: /*
22: This example implements one of the problems found at
23: NLEVP: A Collection of Nonlinear Eigenvalue Problems,
24: The University of Manchester.
25: The details of the collection can be found at:
26: [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
27: Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
29: The sleeper problem is a proportionally damped QEP describing the
30: oscillations of a rail track resting on sleepers.
31: */
33: static char help[] = "NLEVP problem: sleeper.\n\n"
34: "The command line options are:\n"
35: " -n <n>, where <n> = dimension of the matrices.\n\n";
37: #include <slepcpep.h>
41: int main(int argc,char **argv)
42: {
43: Mat M,C,K,A[3]; /* problem matrices */
44: PEP pep; /* polynomial eigenproblem solver context */
45: PetscInt n=10,Istart,Iend,i;
46: PetscBool terse;
49: SlepcInitialize(&argc,&argv,(char*)0,help);
51: PetscOptionsGetInt(NULL,"-n",&n,NULL);
52: PetscPrintf(PETSC_COMM_WORLD,"\nRailtrack resting on sleepers, n=%D\n\n",n);
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: /* K is a pentadiagonal */
59: MatCreate(PETSC_COMM_WORLD,&K);
60: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
61: MatSetFromOptions(K);
62: MatSetUp(K);
63:
64: MatGetOwnershipRange(K,&Istart,&Iend);
65: for (i=Istart;i<Iend;i++) {
66: if (i==0) {
67: MatSetValue(K,i,n-1,-3.0,INSERT_VALUES);
68: MatSetValue(K,i,n-2,1.0,INSERT_VALUES);
69: }
70: if (i==1) { MatSetValue(K,i,n-1,1.0,INSERT_VALUES); }
71: if (i>0) { MatSetValue(K,i,i-1,-3.0,INSERT_VALUES); }
72: if (i>1) { MatSetValue(K,i,i-2,1.0,INSERT_VALUES); }
73: MatSetValue(K,i,i,5.0,INSERT_VALUES);
74: if (i==n-1) {
75: MatSetValue(K,i,0,-3.0,INSERT_VALUES);
76: MatSetValue(K,i,1,1.0,INSERT_VALUES);
77: }
78: if (i==n-2) { MatSetValue(K,i,0,1.0,INSERT_VALUES); }
79: if (i<n-1) { MatSetValue(K,i,i+1,-3.0,INSERT_VALUES); }
80: if (i<n-2) { MatSetValue(K,i,i+2,1.0,INSERT_VALUES); }
81: }
83: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
84: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
86: /* C is a circulant matrix */
87: MatCreate(PETSC_COMM_WORLD,&C);
88: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
89: MatSetFromOptions(C);
90: MatSetUp(C);
91:
92: MatGetOwnershipRange(C,&Istart,&Iend);
93: for (i=Istart;i<Iend;i++) {
94: if (i==0) {
95: MatSetValue(C,i,n-1,-4.0,INSERT_VALUES);
96: MatSetValue(C,i,n-2,1.0,INSERT_VALUES);
97: }
98: if (i==1) { MatSetValue(C,i,n-1,1.0,INSERT_VALUES); }
99: if (i>0) { MatSetValue(C,i,i-1,-4.0,INSERT_VALUES); }
100: if (i>1) { MatSetValue(C,i,i-2,1.0,INSERT_VALUES); }
101: MatSetValue(C,i,i,7.0,INSERT_VALUES);
102: if (i==n-1) {
103: MatSetValue(C,i,0,-4.0,INSERT_VALUES);
104: MatSetValue(C,i,1,1.0,INSERT_VALUES);
105: }
106: if (i==n-2) { MatSetValue(C,i,0,1.0,INSERT_VALUES); }
107: if (i<n-1) { MatSetValue(C,i,i+1,-4.0,INSERT_VALUES); }
108: if (i<n-2) { MatSetValue(C,i,i+2,1.0,INSERT_VALUES); }
109: }
111: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
112: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
113:
114: /* M is the identity matrix */
115: MatCreate(PETSC_COMM_WORLD,&M);
116: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
117: MatSetFromOptions(M);
118: MatSetUp(M);
119: MatGetOwnershipRange(M,&Istart,&Iend);
120: for (i=Istart;i<Iend;i++) {
121: MatSetValue(M,i,i,1.0,INSERT_VALUES);
122: }
123: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
124: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Create the eigensolver and solve the problem
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: PEPCreate(PETSC_COMM_WORLD,&pep);
131: A[0] = K; A[1] = C; A[2] = M;
132: PEPSetOperators(pep,3,A);
133: PEPSetFromOptions(pep);
134: PEPSolve(pep);
136: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137: Display solution and clean up
138: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139:
140: /* show detailed info unless -terse option is given by user */
141: PetscOptionsHasName(NULL,"-terse",&terse);
142: if (terse) {
143: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
144: } else {
145: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
146: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
147: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
148: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
149: }
150: PEPDestroy(&pep);
151: MatDestroy(&M);
152: MatDestroy(&C);
153: MatDestroy(&K);
154: SlepcFinalize();
155: return 0;
156: }