main.cpp 3.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899
  1. // This code gives the mechanical vibration eigenfrequencies and eigenmodes of a
  2. // 3D disk that is clamped at its outer face and subject to proportional damping.
  3. #include "sparselizardbase.h"
  4. using namespace mathop;
  5. void sparselizard(void)
  6. {
  7. // The domain regions as defined in 'disk.geo':
  8. int vol = 1, sur = 2, top = 3;
  9. // The mesh can be curved!
  10. mesh mymesh("disk.msh");
  11. // Nodal shape functions 'h1' with 3 components.
  12. // Field u is the membrane deflection.
  13. field u("h1xyz");
  14. // Use interpolation order 2 on 'vol', the whole domain:
  15. u.setorder(vol, 2);
  16. // Clamp on surface 'sur' (i.e. 0 valued-Dirichlet conditions):
  17. u.setconstraint(sur);
  18. // E is Young's modulus. nu is Poisson's ratio. rho is the volumic mass.
  19. parameter E, nu, rho;
  20. E|vol = 150e9; nu|vol = 0.3; rho|vol = 2330;
  21. // Proportional damping C = alpha*M + beta*K:
  22. double alpha = 100, beta = 0.00006;
  23. formulation elasticity;
  24. // The linear elasticity formulation is classical and thus predefined:
  25. elasticity += integral(vol, predefinedelasticity(dof(u), tf(u), E, nu));
  26. // Add the inertia terms:
  27. elasticity += integral(vol, -rho*dtdt(dof(u))*tf(u));
  28. // Add the proportional damping terms to get damping matrix C = alpha*M + beta*K.
  29. // The mechanical equation can be written as M*dtdt(u) + C*dt(u) + K*u = 0.
  30. // All terms with a first order time derivative on the dof are added to C.
  31. elasticity += integral(vol, alpha * -rho*dt(dof(u))*tf(u));
  32. elasticity += integral(vol, beta * predefinedelasticity(dt(dof(u)), tf(u), E, nu));
  33. elasticity.generate();
  34. // Get the stiffness, damping and mass matrix:
  35. mat K = elasticity.K();
  36. mat C = elasticity.C();
  37. mat M = elasticity.M();
  38. // Remove the rows and columns corresponding to the 0 constraints:
  39. K.removeconstraints();
  40. C.removeconstraints();
  41. M.removeconstraints();
  42. // In case of proportional damping this yields the same results:
  43. // C = alpha*M + beta*K;
  44. // Create the object to solve the second order polynomial eigenvalue problem (M*lambda^2 + C*lambda + K)u = 0 :
  45. eigenvalue eig(K, C, M); // Replace by eig(K, M) for an undamped simulation
  46. // Compute the 10 eigenvalues closest to the target magnitude 0.0 (i.e. the 10 first ones):
  47. eig.compute(10, 0.0);
  48. // Print the eigenfrequencies:
  49. eig.printeigenfrequencies();
  50. // The eigenvectors are real only in the undamped case:
  51. std::vector<vec> myrealeigenvectors = eig.geteigenvectorrealpart();
  52. std::vector<vec> myimageigenvectors = eig.geteigenvectorimaginarypart();
  53. // Loop on all eigenvectors found:
  54. for (int i = 0; i < myrealeigenvectors.size(); i++)
  55. {
  56. // Transfer the data from the ith eigenvector to field u:
  57. u.setdata(vol, myrealeigenvectors[i]);
  58. // Write the deflection on the membrane with an order 2 interpolation:
  59. u.write(vol, "u"+std::to_string(i)+".vtk", 2);
  60. }
  61. // Code validation line. Can be removed.
  62. double checkval = eig.geteigenvaluerealpart()[1]*myimageigenvectors[3].norm();
  63. std::cout << (checkval > -0.0029259 && checkval < -0.00292588);
  64. }
  65. int main(void)
  66. {
  67. SlepcInitialize(0,{},0,0);
  68. sparselizard();
  69. SlepcFinalize();
  70. return 0;
  71. }