main.cpp 4.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142
  1. // This simple code shows how to setup a 3D mechanical problem with a periodic condition.
  2. // A central 1.3 um thick, 50 um diameter polysilicon micromembrane stands above a 300 nm deep cavity.
  3. // The central micromembrane is surrounded by 6 identical ones, one every 60 degrees.
  4. // Assuming a periodic vibration with 60 degrees periodicity allows to reduce the problem to only a 1/6 of the geometry.
  5. //
  6. // In the .nas mesh used the meshes on both faces of the periodic condition do not match.
  7. #include "sparselizardbase.h"
  8. using namespace mathop;
  9. void processmesh(void);
  10. void sparselizard(void)
  11. {
  12. // Name the regions for the inner and outer electrode, the clamp and the regions 'gamma' on which to apply the periodic condition:
  13. int electrodein = 1, electrodeout = 2, clamp = 3, gamma1 = 4, gamma2 = 5, cavity = 4007, solid = 4008, all = 4009;
  14. // Define all regions needed in the source .nas mesh and save it in .msh format.
  15. processmesh();
  16. mesh mymesh("cmutperiodic.msh");
  17. // The periodic condition is only applied to the solid region:
  18. gamma1 = regionintersection({gamma1, solid});
  19. gamma2 = regionintersection({gamma2, solid});
  20. wallclock clk;
  21. // Harmonic simulation at f0 = 1 MHz:
  22. setfundamentalfrequency(1e6);
  23. // Nodal shape functions 'h1' with 3 components for u, the membrane deflection.
  24. // Use harmonic 2 to have u = U*sin(2pi*f0*t).
  25. field u("h1xyz", {2});
  26. // Use interpolation order 2 everywhere:
  27. u.setorder(all, 2);
  28. // Clamp on surface 'clamp' (i.e. 0 valued-Dirichlet conditions):
  29. u.setconstraint(clamp);
  30. // E is Young's modulus. nu is Poisson's ratio.
  31. double E = 160e9, nu = 0.22, rho = 2320;
  32. formulation elasticity;
  33. // The linear elasticity formulation is classical and thus predefined:
  34. elasticity += integral(solid, predefinedelasticity(dof(u), tf(u), E, nu));
  35. // Add a pressure load at frequency f0 on both inner and outer electrodes:
  36. double p = 1e5;
  37. elasticity += integral(electrodein, -p*compz(tf(u.harmonic(2))));
  38. elasticity += integral(electrodeout, -p*compz(tf(u.harmonic(2))));
  39. // Add the inertia term:
  40. elasticity += integral(solid, -rho*dtdt(dof(u))*tf(u));
  41. // Add the periodic condition between gamma1 and gamma2.
  42. // Region gamma2 is obtained from gamma1 by a 60 degrees rotation around z (rotation center is the origin).
  43. elasticity += periodicitycondition(gamma1, gamma2, u, {0,0,0}, {0,0,60});
  44. // Generate, solve and store the solution to field u:
  45. solve(elasticity);
  46. // Write the deflection to ParaView .vtk format. Write with an order 2 interpolation.
  47. u.write(solid, "u.vtk", 2);
  48. // Confirm that the periodic condition is correct by comparing the inner and outer cavity deflection:
  49. double ucenterin = 1e9 * norm(u.harmonic(2)).interpolate(solid,{0,0,1.5e-6})[0];
  50. double ucenterout = 1e9 * norm(u.harmonic(2)).interpolate(solid,{60e-6,0,1.5e-6})[0];
  51. std::cout << "Deflection at center of inner/outer cavity is " << ucenterin << " / " << ucenterout << " nm" << std::endl;
  52. clk.print();
  53. // Code validation line. Can be removed.
  54. std::cout << (ucenterin < 26.5977 && ucenterin > 26.5975);
  55. }
  56. void processmesh(void)
  57. {
  58. // Define the central electrode, outer electrode and clamp regions as well as the regions to apply the periodic condition.
  59. int elecc = 1, eleco = 2, clamp = 3, gamma1 = 4, gamma2 = 5;
  60. setphysicalregionshift(1000);
  61. mesh mymesh1;
  62. mymesh1.load("cmutperiodic.nas", 0);
  63. int vac = regionunion({4001,4005});
  64. int solid = regionunion({4002,4003,4004,4006});
  65. int all = regionunion({vac,solid});
  66. // Rotate the mesh to easily select the bottom side for the periodic condition:
  67. mymesh1.rotate(0,0,30);
  68. mymesh1.write("cmutperiodic.msh", 0);
  69. setphysicalregionshift(0);
  70. mesh mymesh2;
  71. mymesh2.boxselection(elecc, 4001, 2, {-10,10,-10,10,0.3e-6-1e-10,0.3e-6+1e-10});
  72. mymesh2.boxselection(eleco, 4006, 2, {-10,10,-10,10,0.3e-6-1e-10,0.3e-6+1e-10});
  73. mymesh2.boxselection(clamp, all, 2, {-10,10,-10,10,-1e-10,1e-10});
  74. mymesh2.boxselection(gamma1, all, 2, {-10,10,-1e-10,1e-10,-10,10});
  75. mymesh2.load("cmutperiodic.msh", 0);
  76. // Rotate to the other direction to align the other region for the periodic condition:
  77. mymesh2.rotate(0,0,-60);
  78. mymesh2.write("cmutperiodic.msh", 0);
  79. mesh mymesh3;
  80. mymesh3.boxselection(gamma2, all, 2, {-10,10,-1e-10,1e-10,-10,10});
  81. mymesh3.load("cmutperiodic.msh", 0);
  82. // Bring the mesh back to its original angle:
  83. mymesh3.rotate(0,0,30);
  84. // Write the processed mesh:
  85. mymesh3.write("cmutperiodic.msh", 0);
  86. }
  87. int main(void)
  88. {
  89. SlepcInitialize(0,{},0,0);
  90. sparselizard();
  91. SlepcFinalize();
  92. return 0;
  93. }