Actual source code: damped_beam.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  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;

 38:   SlepcInitialize(&argc,&argv,(char*)0,help);

 40:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 41:   nele = n/2;
 42:   n    = 2*nele;
 43:   PetscPrintf(PETSC_COMM_WORLD,"\nSimply supported beam damped in the middle, n=%" PetscInt_FMT " (nele=%" PetscInt_FMT ")\n\n",n,nele);

 45:   dlen = glength/nele;
 46:   EI   = 7e10*width*height*height*height/12.0;
 47:   area = width*height;
 48:   rho  = 0.674/(area*glength);

 50:   K1[0]  =  12;  K1[1]  =   6*dlen;  K1[2]  =   6*dlen;  K1[3]  =  4*dlen*dlen;
 51:   K2[0]  = -12;  K2[1]  =   6*dlen;  K2[2]  =  -6*dlen;  K2[3]  =  2*dlen*dlen;
 52:   K2t[0] = -12;  K2t[1] =  -6*dlen;  K2t[2] =   6*dlen;  K2t[3] =  2*dlen*dlen;
 53:   K3[0]  =  12;  K3[1]  =  -6*dlen;  K3[2]  =  -6*dlen;  K3[3]  =  4*dlen*dlen;
 54:   M1[0]  = 156;  M1[1]  =  22*dlen;  M1[2]  =  22*dlen;  M1[3]  =  4*dlen*dlen;
 55:   M2[0]  =  54;  M2[1]  = -13*dlen;  M2[2]  =  13*dlen;  M2[3]  = -3*dlen*dlen;
 56:   M2t[0] =  54;  M2t[1] =  13*dlen;  M2t[2] = -13*dlen;  M2t[3] = -3*dlen*dlen;
 57:   M3[0]  = 156;  M3[1]  = -22*dlen;  M3[2]  = -22*dlen;  M3[3]  =  4*dlen*dlen;

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 63:   /* K is block-tridiagonal */
 64:   MatCreate(PETSC_COMM_WORLD,&Ko);
 65:   MatSetSizes(Ko,PETSC_DECIDE,PETSC_DECIDE,n+2,n+2);
 66:   MatSetBlockSize(Ko,2);
 67:   MatSetFromOptions(Ko);
 68:   MatSetUp(Ko);

 70:   MatGetOwnershipRange(Ko,&Istart,&Iend);
 71:   for (i=Istart/2;i<Iend/2;i++) {
 72:     if (i>0) {
 73:       j = i-1;
 74:       MatSetValuesBlocked(Ko,1,&i,1,&j,K2t,ADD_VALUES);
 75:       MatSetValuesBlocked(Ko,1,&i,1,&i,K3,ADD_VALUES);
 76:     }
 77:     if (i<nele) {
 78:       j = i+1;
 79:       MatSetValuesBlocked(Ko,1,&i,1,&j,K2,ADD_VALUES);
 80:       MatSetValuesBlocked(Ko,1,&i,1,&i,K1,ADD_VALUES);
 81:     }
 82:   }
 83:   MatAssemblyBegin(Ko,MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(Ko,MAT_FINAL_ASSEMBLY);
 85:   MatScale(Ko,EI/(dlen*dlen*dlen));

 87:   /* M is block-tridiagonal */
 88:   MatCreate(PETSC_COMM_WORLD,&Mo);
 89:   MatSetSizes(Mo,PETSC_DECIDE,PETSC_DECIDE,n+2,n+2);
 90:   MatSetBlockSize(Mo,2);
 91:   MatSetFromOptions(Mo);
 92:   MatSetUp(Mo);

 94:   MatGetOwnershipRange(Mo,&Istart,&Iend);
 95:   for (i=Istart/2;i<Iend/2;i++) {
 96:     if (i>0) {
 97:       j = i-1;
 98:       MatSetValuesBlocked(Mo,1,&i,1,&j,M2t,ADD_VALUES);
 99:       MatSetValuesBlocked(Mo,1,&i,1,&i,M3,ADD_VALUES);
100:     }
101:     if (i<nele) {
102:       j = i+1;
103:       MatSetValuesBlocked(Mo,1,&i,1,&j,M2,ADD_VALUES);
104:       MatSetValuesBlocked(Mo,1,&i,1,&i,M1,ADD_VALUES);
105:     }
106:   }
107:   MatAssemblyBegin(Mo,MAT_FINAL_ASSEMBLY);
108:   MatAssemblyEnd(Mo,MAT_FINAL_ASSEMBLY);
109:   MatScale(Mo,rho*area*dlen/420);

111:   /* remove rows/columns from K and M corresponding to boundary conditions */
112:   ISCreateStride(PETSC_COMM_WORLD,Iend-Istart,Istart,1,&isf);
113:   bc[0] = 0; bc[1] = n;
114:   ISCreateGeneral(PETSC_COMM_SELF,2,bc,PETSC_USE_POINTER,&isbc);
115:   ISDifference(isf,isbc,&is);
116:   MatCreateSubMatrix(Ko,is,is,MAT_INITIAL_MATRIX,&K);
117:   MatCreateSubMatrix(Mo,is,is,MAT_INITIAL_MATRIX,&M);
118:   MatGetLocalSize(M,&mloc,&nloc);

120:   /* C is zero except for the (nele,nele)-entry */
121:   MatCreate(PETSC_COMM_WORLD,&C);
122:   MatSetSizes(C,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE);
123:   MatSetFromOptions(C);
124:   MatSetUp(C);

126:   MatGetOwnershipRange(C,&Istart,&Iend);
127:   if (nele-1>=Istart && nele-1<Iend) MatSetValue(C,nele-1,nele-1,damp,INSERT_VALUES);
128:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
129:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:                 Create the eigensolver and solve the problem
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

135:   PEPCreate(PETSC_COMM_WORLD,&pep);
136:   A[0] = K; A[1] = C; A[2] = M;
137:   PEPSetOperators(pep,3,A);
138:   PEPSetFromOptions(pep);
139:   PEPSolve(pep);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:                     Display solution and clean up
143:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

145:   /* show detailed info unless -terse option is given by user */
146:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
147:   if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
148:   else {
149:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
150:     PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
151:     PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
152:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
153:   }
154:   PEPDestroy(&pep);
155:   ISDestroy(&isf);
156:   ISDestroy(&isbc);
157:   ISDestroy(&is);
158:   MatDestroy(&M);
159:   MatDestroy(&C);
160:   MatDestroy(&K);
161:   MatDestroy(&Ko);
162:   MatDestroy(&Mo);
163:   SlepcFinalize();
164:   return 0;
165: }

167: /*TEST

169:    testset:
170:       args: -pep_nev 2 -pep_ncv 12 -pep_target 0 -terse
171:       requires: !single
172:       output_file: output/damped_beam_1.out
173:       test:
174:          suffix: 1
175:          args: -pep_type {{toar linear}} -st_type sinvert
176:       test:
177:          suffix: 1_qarnoldi
178:          args: -pep_type qarnoldi -pep_qarnoldi_locking 0 -st_type sinvert
179:       test:
180:          suffix: 1_jd
181:          args: -pep_type jd

183:    testset:
184:       args: -pep_nev 2 -pep_ncv 12 -pep_target 1i -terse
185:       requires: complex !single
186:       output_file: output/damped_beam_1.out
187:       test:
188:          suffix: 1_complex
189:          args: -pep_type {{toar linear}} -st_type sinvert
190:       test:
191:          suffix: 1_qarnoldi_complex
192:          args: -pep_type qarnoldi -pep_qarnoldi_locking 0 -st_type sinvert
193:       test:
194:          suffix: 1_jd_complex
195:          args: -pep_type jd

197: TEST*/