Actual source code: test2.c
slepc-3.9.0 2018-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, Universitat Politecnica de Valencia, Spain
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].
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: */
22: static char help[] = "Test the solution of a PEP from a finite element model of "
23: "damped mass-spring system (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";
31: #include <slepcpep.h>
33: int main(int argc,char **argv)
34: {
35: Mat M,C,K,A[3]; /* problem matrices */
36: PEP pep; /* polynomial eigenproblem solver context */
38: PetscInt n=30,Istart,Iend,i,nev;
39: PetscScalar mu=1.0,tau=10.0,kappa=5.0;
40: PetscBool initv=PETSC_FALSE;
41: Vec v0;
43: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
46: PetscOptionsGetScalar(NULL,NULL,"-mu",&mu,NULL);
47: PetscOptionsGetScalar(NULL,NULL,"-tau",&tau,NULL);
48: PetscOptionsGetScalar(NULL,NULL,"-kappa",&kappa,NULL);
49: PetscOptionsGetBool(NULL,NULL,"-initv",&initv,NULL);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /* K is a tridiagonal */
56: MatCreate(PETSC_COMM_WORLD,&K);
57: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
58: MatSetFromOptions(K);
59: MatSetUp(K);
61: MatGetOwnershipRange(K,&Istart,&Iend);
62: for (i=Istart;i<Iend;i++) {
63: if (i>0) {
64: MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
65: }
66: MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
67: if (i<n-1) {
68: MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
69: }
70: }
72: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
73: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
75: /* C is a tridiagonal */
76: MatCreate(PETSC_COMM_WORLD,&C);
77: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
78: MatSetFromOptions(C);
79: MatSetUp(C);
81: MatGetOwnershipRange(C,&Istart,&Iend);
82: for (i=Istart;i<Iend;i++) {
83: if (i>0) {
84: MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
85: }
86: MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
87: if (i<n-1) {
88: MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
89: }
90: }
92: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
93: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
95: /* M is a diagonal matrix */
96: MatCreate(PETSC_COMM_WORLD,&M);
97: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
98: MatSetFromOptions(M);
99: MatSetUp(M);
100: MatGetOwnershipRange(M,&Istart,&Iend);
101: for (i=Istart;i<Iend;i++) {
102: MatSetValue(M,i,i,mu,INSERT_VALUES);
103: }
104: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
105: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
107: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: Create the eigensolver and set various options
109: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: PEPCreate(PETSC_COMM_WORLD,&pep);
112: A[0] = K; A[1] = C; A[2] = M;
113: PEPSetOperators(pep,3,A);
114: PEPSetProblemType(pep,PEP_GENERAL);
115: PEPSetTolerances(pep,PETSC_SMALL,PETSC_DEFAULT);
116: if (initv) { /* initial vector */
117: MatCreateVecs(K,&v0,NULL);
118: VecSetValue(v0,0,-1.0,INSERT_VALUES);
119: VecSetValue(v0,1,0.5,INSERT_VALUES);
120: VecAssemblyBegin(v0);
121: VecAssemblyEnd(v0);
122: PEPSetInitialSpace(pep,1,&v0);
123: VecDestroy(&v0);
124: }
125: PEPSetFromOptions(pep);
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Solve the eigensystem
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: PEPSolve(pep);
132: PEPGetDimensions(pep,&nev,NULL,NULL);
133: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
135: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136: Display solution and clean up
137: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
140: PEPDestroy(&pep);
141: MatDestroy(&M);
142: MatDestroy(&C);
143: MatDestroy(&K);
144: SlepcFinalize();
145: return ierr;
146: }