Actual source code: damped_beam.c
slepc-3.18.2 2023-01-26
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 damped_beam problem is a QEP from the vibrarion analysis of a beam
19: simply supported at both ends and damped in the middle.
20: */
22: static char help[] = "Quadratic eigenproblem from the vibrarion analysis of a beam.\n\n"
23: "The command line options are:\n"
24: " -n <n> ... dimension of the matrices.\n\n";
26: #include <slepcpep.h>
28: int main(int argc,char **argv)
29: {
30: Mat M,Mo,C,K,Ko,A[3]; /* problem matrices */
31: PEP pep; /* polynomial eigenproblem solver context */
32: IS isf,isbc,is;
33: PetscInt n=200,nele,Istart,Iend,i,j,mloc,nloc,bc[2];
34: PetscReal width=0.05,height=0.005,glength=1.0,dlen,EI,area,rho;
35: PetscScalar K1[4],K2[4],K2t[4],K3[4],M1[4],M2[4],M2t[4],M3[4],damp=5.0;
36: PetscBool terse;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
41: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
42: nele = n/2;
43: n = 2*nele;
44: PetscPrintf(PETSC_COMM_WORLD,"\nSimply supported beam damped in the middle, n=%" PetscInt_FMT " (nele=%" PetscInt_FMT ")\n\n",n,nele);
46: dlen = glength/nele;
47: EI = 7e10*width*height*height*height/12.0;
48: area = width*height;
49: rho = 0.674/(area*glength);
51: K1[0] = 12; K1[1] = 6*dlen; K1[2] = 6*dlen; K1[3] = 4*dlen*dlen;
52: K2[0] = -12; K2[1] = 6*dlen; K2[2] = -6*dlen; K2[3] = 2*dlen*dlen;
53: K2t[0] = -12; K2t[1] = -6*dlen; K2t[2] = 6*dlen; K2t[3] = 2*dlen*dlen;
54: K3[0] = 12; K3[1] = -6*dlen; K3[2] = -6*dlen; K3[3] = 4*dlen*dlen;
55: M1[0] = 156; M1[1] = 22*dlen; M1[2] = 22*dlen; M1[3] = 4*dlen*dlen;
56: M2[0] = 54; M2[1] = -13*dlen; M2[2] = 13*dlen; M2[3] = -3*dlen*dlen;
57: M2t[0] = 54; M2t[1] = 13*dlen; M2t[2] = -13*dlen; M2t[3] = -3*dlen*dlen;
58: M3[0] = 156; M3[1] = -22*dlen; M3[2] = -22*dlen; M3[3] = 4*dlen*dlen;
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: /* K is block-tridiagonal */
65: MatCreate(PETSC_COMM_WORLD,&Ko);
66: MatSetSizes(Ko,PETSC_DECIDE,PETSC_DECIDE,n+2,n+2);
67: MatSetBlockSize(Ko,2);
68: MatSetFromOptions(Ko);
69: MatSetUp(Ko);
71: MatGetOwnershipRange(Ko,&Istart,&Iend);
72: for (i=Istart/2;i<Iend/2;i++) {
73: if (i>0) {
74: j = i-1;
75: MatSetValuesBlocked(Ko,1,&i,1,&j,K2t,ADD_VALUES);
76: MatSetValuesBlocked(Ko,1,&i,1,&i,K3,ADD_VALUES);
77: }
78: if (i<nele) {
79: j = i+1;
80: MatSetValuesBlocked(Ko,1,&i,1,&j,K2,ADD_VALUES);
81: MatSetValuesBlocked(Ko,1,&i,1,&i,K1,ADD_VALUES);
82: }
83: }
84: MatAssemblyBegin(Ko,MAT_FINAL_ASSEMBLY);
85: MatAssemblyEnd(Ko,MAT_FINAL_ASSEMBLY);
86: MatScale(Ko,EI/(dlen*dlen*dlen));
88: /* M is block-tridiagonal */
89: MatCreate(PETSC_COMM_WORLD,&Mo);
90: MatSetSizes(Mo,PETSC_DECIDE,PETSC_DECIDE,n+2,n+2);
91: MatSetBlockSize(Mo,2);
92: MatSetFromOptions(Mo);
93: MatSetUp(Mo);
95: MatGetOwnershipRange(Mo,&Istart,&Iend);
96: for (i=Istart/2;i<Iend/2;i++) {
97: if (i>0) {
98: j = i-1;
99: MatSetValuesBlocked(Mo,1,&i,1,&j,M2t,ADD_VALUES);
100: MatSetValuesBlocked(Mo,1,&i,1,&i,M3,ADD_VALUES);
101: }
102: if (i<nele) {
103: j = i+1;
104: MatSetValuesBlocked(Mo,1,&i,1,&j,M2,ADD_VALUES);
105: MatSetValuesBlocked(Mo,1,&i,1,&i,M1,ADD_VALUES);
106: }
107: }
108: MatAssemblyBegin(Mo,MAT_FINAL_ASSEMBLY);
109: MatAssemblyEnd(Mo,MAT_FINAL_ASSEMBLY);
110: MatScale(Mo,rho*area*dlen/420);
112: /* remove rows/columns from K and M corresponding to boundary conditions */
113: ISCreateStride(PETSC_COMM_WORLD,Iend-Istart,Istart,1,&isf);
114: bc[0] = 0; bc[1] = n;
115: ISCreateGeneral(PETSC_COMM_SELF,2,bc,PETSC_USE_POINTER,&isbc);
116: ISDifference(isf,isbc,&is);
117: MatCreateSubMatrix(Ko,is,is,MAT_INITIAL_MATRIX,&K);
118: MatCreateSubMatrix(Mo,is,is,MAT_INITIAL_MATRIX,&M);
119: MatGetLocalSize(M,&mloc,&nloc);
121: /* C is zero except for the (nele,nele)-entry */
122: MatCreate(PETSC_COMM_WORLD,&C);
123: MatSetSizes(C,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE);
124: MatSetFromOptions(C);
125: MatSetUp(C);
127: MatGetOwnershipRange(C,&Istart,&Iend);
128: if (nele-1>=Istart && nele-1<Iend) MatSetValue(C,nele-1,nele-1,damp,INSERT_VALUES);
129: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
130: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
132: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133: Create the eigensolver and solve the problem
134: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: PEPCreate(PETSC_COMM_WORLD,&pep);
137: A[0] = K; A[1] = C; A[2] = M;
138: PEPSetOperators(pep,3,A);
139: PEPSetFromOptions(pep);
140: PEPSolve(pep);
142: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143: Display solution and clean up
144: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146: /* show detailed info unless -terse option is given by user */
147: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
148: if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
149: else {
150: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
151: PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
152: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
153: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
154: }
155: PEPDestroy(&pep);
156: ISDestroy(&isf);
157: ISDestroy(&isbc);
158: ISDestroy(&is);
159: MatDestroy(&M);
160: MatDestroy(&C);
161: MatDestroy(&K);
162: MatDestroy(&Ko);
163: MatDestroy(&Mo);
164: SlepcFinalize();
165: return 0;
166: }
168: /*TEST
170: testset:
171: args: -pep_nev 2 -pep_ncv 12 -pep_target 0 -terse
172: requires: !single
173: output_file: output/damped_beam_1.out
174: test:
175: suffix: 1
176: args: -pep_type {{toar linear}} -st_type sinvert
177: test:
178: suffix: 1_qarnoldi
179: args: -pep_type qarnoldi -pep_qarnoldi_locking 0 -st_type sinvert
180: test:
181: suffix: 1_jd
182: args: -pep_type jd
184: testset:
185: args: -pep_nev 2 -pep_ncv 12 -pep_target 1i -terse
186: requires: complex !single
187: output_file: output/damped_beam_1.out
188: test:
189: suffix: 1_complex
190: args: -pep_type {{toar linear}} -st_type sinvert
191: test:
192: suffix: 1_qarnoldi_complex
193: args: -pep_type qarnoldi -pep_qarnoldi_locking 0 -st_type sinvert
194: test:
195: suffix: 1_jd_complex
196: args: -pep_type jd
198: TEST*/