Actual source code: wiresaw.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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 two 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: WIRESAW1 is a gyroscopic QEP from vibration analysis of a wiresaw,
19: where the parameter V represents the speed of the wire. When the
20: parameter eta is nonzero, then it turns into the WIRESAW2 problem
21: (with added viscous damping, e.g. eta=0.8).
23: [2] S. Wei and I. Kao, "Vibration analysis of wire and frequency
24: response in the modern wiresaw manufacturing process", J. Sound
25: Vib. 213(5):1383-1395, 2000.
26: */
28: static char help[] = "Vibration analysis of a wiresaw.\n\n"
29: "The command line options are:\n"
30: " -n <n> ... dimension of the matrices (default 10).\n"
31: " -v <value> ... velocity of the wire (default 0.01).\n"
32: " -eta <value> ... viscous damping (default 0.0).\n\n";
34: #include <slepcpep.h>
36: int main(int argc,char **argv)
37: {
38: Mat M,D,K,A[3]; /* problem matrices */
39: PEP pep; /* polynomial eigenproblem solver context */
40: PetscInt n=10,Istart,Iend,j,k;
41: PetscReal v=0.01,eta=0.0;
42: PetscBool terse;
45: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
47: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
48: PetscOptionsGetReal(NULL,NULL,"-v",&v,NULL);
49: PetscOptionsGetReal(NULL,NULL,"-eta",&eta,NULL);
50: PetscPrintf(PETSC_COMM_WORLD,"\nVibration analysis of a wiresaw, n=%D v=%g eta=%g\n\n",n,(double)v,(double)eta);
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Compute the matrices that define the eigensystem, (k^2*M+k*D+K)x=0
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: /* K is a diagonal matrix */
57: MatCreate(PETSC_COMM_WORLD,&K);
58: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
59: MatSetFromOptions(K);
60: MatSetUp(K);
62: MatGetOwnershipRange(K,&Istart,&Iend);
63: for (j=Istart;j<Iend;j++) {
64: MatSetValue(K,j,j,(j+1)*(j+1)*PETSC_PI*PETSC_PI*(1.0-v*v),INSERT_VALUES);
65: }
67: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
68: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
69: MatScale(K,0.5);
71: /* D is a tridiagonal */
72: MatCreate(PETSC_COMM_WORLD,&D);
73: MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n);
74: MatSetFromOptions(D);
75: MatSetUp(D);
77: MatGetOwnershipRange(D,&Istart,&Iend);
78: for (j=Istart;j<Iend;j++) {
79: for (k=0;k<n;k++) {
80: if ((j+k)%2) {
81: MatSetValue(D,j,k,8.0*(j+1)*(k+1)*v/((j+1)*(j+1)-(k+1)*(k+1)),INSERT_VALUES);
82: }
83: }
84: }
86: MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY);
87: MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY);
88: MatScale(D,0.5);
90: /* M is a diagonal matrix */
91: MatCreate(PETSC_COMM_WORLD,&M);
92: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
93: MatSetFromOptions(M);
94: MatSetUp(M);
95: MatGetOwnershipRange(M,&Istart,&Iend);
96: for (j=Istart;j<Iend;j++) {
97: MatSetValue(M,j,j,1.0,INSERT_VALUES);
98: }
99: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
100: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
101: MatScale(M,0.5);
103: /* add damping */
104: if (eta>0.0) {
105: MatAXPY(K,eta,D,DIFFERENT_NONZERO_PATTERN); /* K = K + eta*D */
106: MatShift(D,eta); /* D = D + eta*eye(n) */
107: }
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Create the eigensolver and solve the problem
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: PEPCreate(PETSC_COMM_WORLD,&pep);
114: A[0] = K; A[1] = D; A[2] = M;
115: PEPSetOperators(pep,3,A);
116: if (eta==0.0) {
117: PEPSetProblemType(pep,PEP_GYROSCOPIC);
118: } else {
119: PEPSetProblemType(pep,PEP_GENERAL);
120: }
121: PEPSetFromOptions(pep);
122: PEPSolve(pep);
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Display solution and clean up
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128: /* show detailed info unless -terse option is given by user */
129: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
130: if (terse) {
131: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
132: } else {
133: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
134: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
135: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
136: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
137: }
138: PEPDestroy(&pep);
139: MatDestroy(&M);
140: MatDestroy(&D);
141: MatDestroy(&K);
142: SlepcFinalize();
143: return ierr;
144: }
146: /*TEST
148: testset:
149: args: -pep_nev 4 -terse
150: requires: !single !complex
151: output_file: output/wiresaw_1.out
152: test:
153: suffix: 1
154: args: -pep_type {{toar qarnoldi}}
155: test:
156: suffix: 1_linear_h1
157: args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 1,0 -pep_linear_st_ksp_type bcgs -pep_linear_st_pc_type kaczmarz
158: test:
159: suffix: 1_linear_h2
160: args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 0,1
162: TEST*/