Actual source code: loaded_string.c
slepc-3.19.2 2023-09-05
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: PetscFunctionBeginUser;
44: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
46: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
47: PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
48: PetscCall(PetscOptionsGetReal(NULL,NULL,"-mass",&m,NULL));
49: sigma = kappa/m;
50: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string (QEP), n=%" PetscInt_FMT " kappa=%g m=%g\n\n",n,(double)kappa,(double)m));
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /* initialize matrices */
56: for (i=0;i<NMAT;i++) {
57: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
58: PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
59: PetscCall(MatSetFromOptions(A[i]));
60: PetscCall(MatSetUp(A[i]));
61: }
62: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
64: /* A0 */
65: for (i=Istart;i<Iend;i++) {
66: PetscCall(MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES));
67: if (i>0) PetscCall(MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES));
68: if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES));
69: }
71: /* A1 */
72: for (i=Istart;i<Iend;i++) {
73: PetscCall(MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES));
74: if (i>0) PetscCall(MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES));
75: if (i<n-1) PetscCall(MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES));
76: }
78: /* A2 */
79: if (Istart<=n-1 && n-1<Iend) PetscCall(MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES));
81: /* assemble matrices */
82: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
83: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));
85: /* build matrices for the QEP */
86: PetscCall(MatAXPY(A[2],1.0,A[0],DIFFERENT_NONZERO_PATTERN));
87: PetscCall(MatAXPY(A[2],sigma,A[1],SAME_NONZERO_PATTERN));
88: PetscCall(MatScale(A[2],-1.0));
89: PetscCall(MatScale(A[0],sigma));
90: M = A[1];
91: A[1] = A[2];
92: A[2] = M;
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create the eigensolver and solve the problem
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
99: PetscCall(PEPSetOperators(pep,3,A));
100: PetscCall(PEPSetFromOptions(pep));
101: PetscCall(PEPSolve(pep));
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Display solution and clean up
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: /* show detailed info unless -terse option is given by user */
108: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
109: if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
110: else {
111: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
112: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
113: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
114: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
115: }
116: PetscCall(PEPDestroy(&pep));
117: for (i=0;i<NMAT;i++) PetscCall(MatDestroy(&A[i]));
118: PetscCall(SlepcFinalize());
119: return 0;
120: }
122: /*TEST
124: test:
125: suffix: 1
126: args: -pep_hyperbolic -pep_interval 4,900 -pep_type stoar -st_type sinvert -st_pc_type cholesky -terse
127: requires: !single
129: TEST*/