Actual source code: damped_beam.c

slepc-3.20.1 2023-11-27
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:   PetscFunctionBeginUser;
 39:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 41:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 42:   nele = n/2;
 43:   n    = 2*nele;
 44:   PetscCall(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:   PetscCall(MatCreate(PETSC_COMM_WORLD,&Ko));
 66:   PetscCall(MatSetSizes(Ko,PETSC_DECIDE,PETSC_DECIDE,n+2,n+2));
 67:   PetscCall(MatSetBlockSize(Ko,2));
 68:   PetscCall(MatSetFromOptions(Ko));
 69:   PetscCall(MatSetUp(Ko));

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

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

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

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

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

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

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

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

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

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