Actual source code: damped_beam.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */
 21: /*
 22:    This example implements one of the problems found at
 23:        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
 24:        The University of Manchester.
 25:    The details of the collection can be found at:
 26:        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
 27:            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.

 29:    The damped_beam problem is a QEP from the vibrarion analysis of a beam
 30:    simply supported at both ends and damped in the middle.
 31: */

 33: static char help[] = "NLEVP problem: damped_beam.\n\n"
 34:   "The command line options are:\n"
 35:   "  -n <n> ... dimension of the matrices.\n\n";

 37: #include <slepcpep.h>

 41: int main(int argc,char **argv)
 42: {
 43:   Mat            M,Mo,C,K,Ko,A[3]; /* problem matrices */
 44:   PEP            pep;              /* polynomial eigenproblem solver context */
 45:   IS             isf,isbc,is;
 46:   PetscInt       n=200,nele,Istart,Iend,i,j,mloc,nloc,bc[2];
 47:   PetscReal      width=0.05,height=0.005,glength=1.0,dlen,EI,area,rho;
 48:   PetscScalar    K1[4],K2[4],K2t[4],K3[4],M1[4],M2[4],M2t[4],M3[4],damp=5.0;

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

 53:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 54:   nele = n/2;
 55:   n    = 2*nele;
 56:   PetscPrintf(PETSC_COMM_WORLD,"\nSimply supported beam damped in the middle, n=%D (nele=%D)\n\n",n,nele);

 58:   dlen = glength/nele;
 59:   EI   = 7e10*width*height*height*height/12.0;
 60:   area = width*height;
 61:   rho  = 0.674/(area*glength);

 63:   K1[0]  =  12;  K1[1]  =   6*dlen;  K1[2]  =   6*dlen;  K1[3]  =  4*dlen*dlen;
 64:   K2[0]  = -12;  K2[1]  =   6*dlen;  K2[2]  =  -6*dlen;  K2[3]  =  2*dlen*dlen;
 65:   K2t[0] = -12;  K2t[1] =  -6*dlen;  K2t[2] =   6*dlen;  K2t[3] =  2*dlen*dlen;
 66:   K3[0]  =  12;  K3[1]  =  -6*dlen;  K3[2]  =  -6*dlen;  K3[3]  =  4*dlen*dlen;
 67:   M1[0]  = 156;  M1[1]  =  22*dlen;  M1[2]  =  22*dlen;  M1[3]  =  4*dlen*dlen;
 68:   M2[0]  =  54;  M2[1]  = -13*dlen;  M2[2]  =  13*dlen;  M2[3]  = -3*dlen*dlen;
 69:   M2t[0] =  54;  M2t[1] =  13*dlen;  M2t[2] = -13*dlen;  M2t[3] = -3*dlen*dlen;
 70:   M3[0]  = 156;  M3[1]  = -22*dlen;  M3[2]  = -22*dlen;  M3[3]  =  4*dlen*dlen;

 72:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 73:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 74:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 76:   /* K is block-tridiagonal */
 77:   MatCreate(PETSC_COMM_WORLD,&Ko);
 78:   MatSetSizes(Ko,PETSC_DECIDE,PETSC_DECIDE,n+2,n+2);
 79:   MatSetBlockSize(Ko,2);
 80:   MatSetFromOptions(Ko);
 81:   MatSetUp(Ko);
 82:   
 83:   MatGetOwnershipRange(Ko,&Istart,&Iend);
 84:   for (i=Istart/2;i<Iend/2;i++) {
 85:     if (i>0) {
 86:       j = i-1;
 87:       MatSetValuesBlocked(Ko,1,&i,1,&j,K2t,ADD_VALUES);
 88:       MatSetValuesBlocked(Ko,1,&i,1,&i,K3,ADD_VALUES);
 89:     }
 90:     if (i<nele) {
 91:       j = i+1;
 92:       MatSetValuesBlocked(Ko,1,&i,1,&j,K2,ADD_VALUES);
 93:       MatSetValuesBlocked(Ko,1,&i,1,&i,K1,ADD_VALUES);
 94:     }
 95:   }
 96:   MatAssemblyBegin(Ko,MAT_FINAL_ASSEMBLY);
 97:   MatAssemblyEnd(Ko,MAT_FINAL_ASSEMBLY);
 98:   MatScale(Ko,EI/(dlen*dlen*dlen));

100:   /* M is block-tridiagonal */
101:   MatCreate(PETSC_COMM_WORLD,&Mo);
102:   MatSetSizes(Mo,PETSC_DECIDE,PETSC_DECIDE,n+2,n+2);
103:   MatSetBlockSize(Mo,2);
104:   MatSetFromOptions(Mo);
105:   MatSetUp(Mo);

107:   MatGetOwnershipRange(Mo,&Istart,&Iend);
108:   for (i=Istart/2;i<Iend/2;i++) {
109:     if (i>0) {
110:       j = i-1;
111:       MatSetValuesBlocked(Mo,1,&i,1,&j,M2t,ADD_VALUES);
112:       MatSetValuesBlocked(Mo,1,&i,1,&i,M3,ADD_VALUES);
113:     }
114:     if (i<nele) {
115:       j = i+1;
116:       MatSetValuesBlocked(Mo,1,&i,1,&j,M2,ADD_VALUES);
117:       MatSetValuesBlocked(Mo,1,&i,1,&i,M1,ADD_VALUES);
118:     }
119:   }
120:   MatAssemblyBegin(Mo,MAT_FINAL_ASSEMBLY);
121:   MatAssemblyEnd(Mo,MAT_FINAL_ASSEMBLY);
122:   MatScale(Mo,rho*area*dlen/420);

124:   /* remove rows/columns from K and M corresponding to boundary conditions */
125:   ISCreateStride(PETSC_COMM_WORLD,Iend-Istart,Istart,1,&isf);
126:   bc[0] = 0; bc[1] = n;
127:   ISCreateGeneral(PETSC_COMM_SELF,2,bc,PETSC_USE_POINTER,&isbc);
128:   ISDifference(isf,isbc,&is);
129:   MatGetSubMatrix(Ko,is,is,MAT_INITIAL_MATRIX,&K);
130:   MatGetSubMatrix(Mo,is,is,MAT_INITIAL_MATRIX,&M);
131:   MatGetLocalSize(M,&mloc,&nloc);

133:   /* C is zero except for the (nele,nele)-entry */
134:   MatCreate(PETSC_COMM_WORLD,&C);
135:   MatSetSizes(C,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE);
136:   MatSetFromOptions(C);
137:   MatSetUp(C);
138:   
139:   MatGetOwnershipRange(C,&Istart,&Iend);
140:   if (nele-1>=Istart && nele-1<Iend) { 
141:     MatSetValue(C,nele-1,nele-1,damp,INSERT_VALUES);
142:   }
143:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
144:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
145:   
146:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
147:                 Create the eigensolver and solve the problem
148:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

150:   PEPCreate(PETSC_COMM_WORLD,&pep);
151:   A[0] = K; A[1] = C; A[2] = M;
152:   PEPSetOperators(pep,3,A);
153:   PEPSetFromOptions(pep);
154:   PEPSolve(pep);

156:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157:                     Display solution and clean up
158:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159:   
160:   PEPPrintSolution(pep,NULL);
161:   PEPDestroy(&pep);
162:   ISDestroy(&isf);
163:   ISDestroy(&isbc);
164:   ISDestroy(&is);
165:   MatDestroy(&M);
166:   MatDestroy(&C);
167:   MatDestroy(&K);
168:   MatDestroy(&Ko);
169:   MatDestroy(&Mo);
170:   SlepcFinalize();
171:   return 0;
172: }