Actual source code: spring.c
slepc-3.7.3 2016-09-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, 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 spring problem is a QEP from the finite element model of a damped
30: mass-spring system. This implementation supports only scalar parameters,
31: that is all masses, dampers and springs have the same constants.
32: Furthermore, this implementation does not consider different constants
33: for dampers and springs connecting adjacent masses or masses to the ground.
34: */
36: static char help[] = "FEM model of a damped mass-spring system.\n\n"
37: "The command line options are:\n"
38: " -n <n> ... dimension of the matrices.\n"
39: " -mu <value> ... mass (default 1).\n"
40: " -tau <value> ... damping constant of the dampers (default 10).\n"
41: " -kappa <value> ... damping constant of the springs (default 5).\n\n";
43: #include <slepcpep.h>
47: int main(int argc,char **argv)
48: {
49: Mat M,C,K,A[3]; /* problem matrices */
50: PEP pep; /* polynomial eigenproblem solver context */
51: PetscInt n=5,Istart,Iend,i;
52: PetscReal mu=1,tau=10,kappa=5;
53: PetscBool terse;
56: SlepcInitialize(&argc,&argv,(char*)0,help);
58: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
59: PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL);
60: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
61: PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
62: PetscPrintf(PETSC_COMM_WORLD,"\nDamped mass-spring system, n=%D mu=%g tau=%g kappa=%g\n\n",n,(double)mu,(double)tau,(double)kappa);
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: /* K is a tridiagonal */
69: MatCreate(PETSC_COMM_WORLD,&K);
70: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
71: MatSetFromOptions(K);
72: MatSetUp(K);
73:
74: MatGetOwnershipRange(K,&Istart,&Iend);
75: for (i=Istart;i<Iend;i++) {
76: if (i>0) {
77: MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
78: }
79: MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
80: if (i<n-1) {
81: MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
82: }
83: }
85: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
86: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
88: /* C is a tridiagonal */
89: MatCreate(PETSC_COMM_WORLD,&C);
90: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
91: MatSetFromOptions(C);
92: MatSetUp(C);
93:
94: MatGetOwnershipRange(C,&Istart,&Iend);
95: for (i=Istart;i<Iend;i++) {
96: if (i>0) {
97: MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
98: }
99: MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
100: if (i<n-1) {
101: MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
102: }
103: }
105: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
106: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
107:
108: /* M is a diagonal matrix */
109: MatCreate(PETSC_COMM_WORLD,&M);
110: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
111: MatSetFromOptions(M);
112: MatSetUp(M);
113: MatGetOwnershipRange(M,&Istart,&Iend);
114: for (i=Istart;i<Iend;i++) {
115: MatSetValue(M,i,i,mu,INSERT_VALUES);
116: }
117: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
118: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Create the eigensolver and solve the problem
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: PEPCreate(PETSC_COMM_WORLD,&pep);
125: A[0] = K; A[1] = C; A[2] = M;
126: PEPSetOperators(pep,3,A);
127: PEPSetFromOptions(pep);
128: PEPSolve(pep);
130: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: Display solution and clean up
132: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133:
134: /* show detailed info unless -terse option is given by user */
135: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
136: if (terse) {
137: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
138: } else {
139: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
140: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
141: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
142: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
143: }
144: PEPDestroy(&pep);
145: MatDestroy(&M);
146: MatDestroy(&C);
147: MatDestroy(&K);
148: SlepcFinalize();
149: return ierr;
150: }