Actual source code: loaded_string.c
slepc-3.17.2 2022-08-09
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: This example implements one of the problems found at
12: NLEVP: A Collection of Nonlinear Eigenvalue Problems,
13: The University of Manchester.
14: The details of the collection can be found at:
15: [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
16: Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
18: The loaded_string problem is a rational eigenvalue problem for the
19: finite element model of a loaded vibrating string.
20: This example solves the loaded_string problem by first transforming
21: it to a quadratic eigenvalue problem.
22: */
24: static char help[] = "Finite element model of a loaded vibrating string.\n\n"
25: "The command line options are:\n"
26: " -n <n>, dimension of the matrices.\n"
27: " -kappa <kappa>, stiffness of elastic spring.\n"
28: " -mass <m>, mass of the attached load.\n\n";
30: #include <slepcpep.h>
32: #define NMAT 3
34: int main(int argc,char **argv)
35: {
36: Mat A[3],M; /* problem matrices */
37: PEP pep; /* polynomial eigenproblem solver context */
38: PetscInt n=100,Istart,Iend,i;
39: PetscBool terse;
40: PetscReal kappa=1.0,m=1.0;
41: PetscScalar sigma;
43: SlepcInitialize(&argc,&argv,(char*)0,help);
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
46: PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
47: PetscOptionsGetReal(NULL,NULL,"-mass",&m,NULL);
48: sigma = kappa/m;
49: PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string (QEP), n=%" PetscInt_FMT " kappa=%g m=%g\n\n",n,(double)kappa,(double)m);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54: /* initialize matrices */
55: for (i=0;i<NMAT;i++) {
56: MatCreate(PETSC_COMM_WORLD,&A[i]);
57: MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
58: MatSetFromOptions(A[i]);
59: MatSetUp(A[i]);
60: }
61: MatGetOwnershipRange(A[0],&Istart,&Iend);
63: /* A0 */
64: for (i=Istart;i<Iend;i++) {
65: MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES);
66: if (i>0) MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES);
67: if (i<n-1) MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES);
68: }
70: /* A1 */
71: for (i=Istart;i<Iend;i++) {
72: MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES);
73: if (i>0) MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES);
74: if (i<n-1) MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES);
75: }
77: /* A2 */
78: if (Istart<=n-1 && n-1<Iend) MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES);
80: /* assemble matrices */
81: for (i=0;i<NMAT;i++) MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
82: for (i=0;i<NMAT;i++) MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
84: /* build matrices for the QEP */
85: MatAXPY(A[2],1.0,A[0],DIFFERENT_NONZERO_PATTERN);
86: MatAXPY(A[2],sigma,A[1],SAME_NONZERO_PATTERN);
87: MatScale(A[2],-1.0);
88: MatScale(A[0],sigma);
89: M = A[1];
90: A[1] = A[2];
91: A[2] = M;
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Create the eigensolver and solve the problem
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: PEPCreate(PETSC_COMM_WORLD,&pep);
98: PEPSetOperators(pep,3,A);
99: PEPSetFromOptions(pep);
100: PEPSolve(pep);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Display solution and clean up
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: /* show detailed info unless -terse option is given by user */
107: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
108: if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
109: else {
110: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
111: PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
112: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
113: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
114: }
115: PEPDestroy(&pep);
116: for (i=0;i<NMAT;i++) MatDestroy(&A[i]);
117: SlepcFinalize();
118: return 0;
119: }
121: /*TEST
123: test:
124: suffix: 1
125: args: -pep_hyperbolic -pep_interval 4,900 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
126: requires: !single
128: TEST*/