main.cpp 2.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384
  1. // This code gives the mechanical vibration eigenfrequencies and eigenmodes of a
  2. // 3D disk that is clamped at its outer face.
  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 3 on 'vol', the whole domain:
  15. u.setorder(vol, 3);
  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. formulation elasticity;
  22. // The linear elasticity formulation is classical and thus predefined:
  23. elasticity += integral(vol, predefinedelasticity(dof(u), tf(u), E, nu));
  24. // Add the inertia terms:
  25. elasticity += integral(vol, -rho*dtdt(dof(u))*tf(u));
  26. elasticity.generate();
  27. // Get the stiffness and mass matrix:
  28. mat K = elasticity.K();
  29. mat M = elasticity.M();
  30. // Remove the rows and columns corresponding to the 0 constraints:
  31. K.removeconstraints();
  32. M.removeconstraints();
  33. // Create the object to solve the generalized eigenvalue problem K*x = lambda*M*x :
  34. eigenvalue eig(K, M);
  35. // Compute the 10 eigenvalues closest to the target magnitude 0.0 (i.e. the 10 first ones):
  36. eig.compute(10, 0.0);
  37. // Print the eigenfrequencies:
  38. eig.printeigenfrequencies();
  39. // The eigenvectors are real thus we only need the real part:
  40. std::vector<vec> myeigenvectors = eig.geteigenvectorrealpart();
  41. // Loop on all eigenvectors found:
  42. for (int i = 0; i < myeigenvectors.size(); i++)
  43. {
  44. // Transfer the data from the ith eigenvector to field u:
  45. u.setdata(top, myeigenvectors[i]);
  46. // Write the deflection on the top surface of the membrane with an order 3 interpolation:
  47. u.write(top, "u"+std::to_string(i)+".pos", 3);
  48. }
  49. // Code validation line. Can be removed.
  50. std::cout << (eig.geteigenvaluerealpart()[0] < 6.25240e+06 && eig.geteigenvaluerealpart()[0] > 6.25235e+06);
  51. }
  52. int main(void)
  53. {
  54. SlepcInitialize(0,{},0,0);
  55. sparselizard();
  56. SlepcFinalize();
  57. return 0;
  58. }