main.cpp 4.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109
  1. // This code simulates the harmonic mechanical deflection of a 3D monocrystalline silicon - PZT microbilayer when
  2. // the PZT (a piezoelectric material) is sandwiched between two electrically actuated electrodes.
  3. //
  4. // The anisotropic behavior of both silicon and PZT is fully taken into account.
  5. // The silicon and PZT crystals are rotated to align with a selected orientation.
  6. #include "sparselizardbase.h"
  7. using namespace mathop;
  8. void sparselizard(void)
  9. {
  10. // The domain regions as defined in 'bilayer.geo':
  11. int pztlayer = 1, siliconlayer = 2, electrode = 3, ground = 4, clamp = 6, freeside = 7, elecbox = 8, mecabox = 9;
  12. mesh mymesh("bilayer.msh");
  13. int wholedomain = regionunion({pztlayer, siliconlayer});
  14. // Harmonic analysis. Set the fundamental frequency [Hz]:
  15. setfundamentalfrequency(1e4);
  16. // Nodal shape functions 'h1' for v (the electric potential) and u
  17. // (the mechanical displacement). Three components are used for u.
  18. // Use harmonic 2, i.e. u(x,t) = U(x)*sin(2pi*f0*t) for u and
  19. // v(x,t) = V(x)*sin(2pi*f0*t) for v.
  20. //
  21. field u("h1xyz",{2}), v("h1",{2});
  22. // Use interpolation order 2 on 'wholedomain' for u and 1 for v in the PZT layer:
  23. u.setorder(wholedomain, 2);
  24. v.setorder(pztlayer, 1);
  25. // Clamp on surface 'clamp' (i.e. 0 valued-Dirichlet conditions):
  26. u.setconstraint(clamp);
  27. v.harmonic(2).setconstraint(electrode, 10);
  28. v.setconstraint(ground, 0);
  29. // Silicon and PZT density [kg/m^3]:
  30. double rhosi = 2330, rhopzt = 7500;
  31. // Diagonal relative permittivity matrix for PZT:
  32. expression K(3,3,{1704,1704,1433});
  33. K = K * 8.854e-12;
  34. // Coupling matrix [C/m^2] for PZT (6 rows, 3 columns):
  35. expression C(6,3,{0,0,-6.62, 0,0,-6.62, 0,0,23.24, 0,17.03,0, 17.03,0,0, 0,0,0});
  36. // Anisotropic elasticity matrix [Pa] for silicon and PZT. Ordering is [exx,eyy,ezz,gyz,gxz,gxy] (Voigt form).
  37. // Only the lower triangular part (top to bottom and left to right) is provided since it is symmetric.
  38. expression Hsi(6,6, {194.5e9, 35.7e9,194.5e9, 64.1e9,64.1e9,165.7e9, 0,0,0,79.6e9, 0,0,0,0,79.6e9, 0,0,0,0,0,50.9e9});
  39. expression Hpzt(6,6, {1.27e11, 8.02e10,1.27e11, 8.46e10,8.46e10,1.17e11, 0,0,0,2.29e10, 0,0,0,0,2.29e10, 0,0,0,0,0,2.34e10});
  40. // Rotate the silicon and PZT crystals first by -30 degrees around y then by 45 degrees around z.
  41. // Refer to the documentation to make sure you understand how to use 'rotate'!
  42. Hsi.rotate(0,-30,45,"K","KT"); Hpzt.rotate(0,-30,45,"K","KT"); C.rotate(0,-30,45,"K","RT"); K.rotate(0,-30,45);
  43. formulation piezo;
  44. // Orthotropic elasticity in the silicon (not piezoelectric):
  45. piezo += integral(siliconlayer, predefinedelasticity(dof(u), tf(u), Hsi) );
  46. // Inertia term:
  47. piezo += integral(siliconlayer, -rhosi*dtdt(dof(u))*tf(u) );
  48. // The classical weak formulation for piezoelectricity below can be found e.g. in paper:
  49. //
  50. // "A new fnite-element formulation for electromechanical boundary value problems"
  51. // Define the mechanical equations of the problem in the piezo.
  52. // strain(u) returns the strain vector [exx,eyy,ezz,gyz,gxz,gxy] of u.
  53. piezo += integral(pztlayer, -( Hpzt*strain(dof(u)) )*strain(tf(u)) - ( C*grad(dof(v)) )*strain(tf(u)) );
  54. // Inertia term for PZT:
  55. piezo += integral(pztlayer, -rhopzt*dtdt(dof(u))*tf(u) );
  56. // Define the electrical equations of the problem in the piezo:
  57. piezo += integral(pztlayer, ( transpose(C)*strain(dof(u)) )*grad(tf(v)) - ( K*grad(dof(v)) )*grad(tf(v)) );
  58. solve(piezo);
  59. // Display the peak displacement:
  60. double umax = norm(u.harmonic(2)).max(wholedomain, 5)[0];
  61. std::cout << "Peak deflection is " << umax << " m" << std::endl;
  62. // Write the deflection and electric potential on the volume boundary.
  63. u.write(mecabox, "u.vtk", 2);
  64. v.write(elecbox, "v.vtk", 1);
  65. // Harmonic fields can be saved in time. Use 50 time steps:
  66. u.write(mecabox, "utime.vtk", 2, 50);
  67. v.write(elecbox, "vtime.vtk", 1, 50);
  68. // Code validation line. Can be removed.
  69. std::cout << (umax < 1.7386e-07 && umax > 1.7384e-07);
  70. }
  71. int main(void)
  72. {
  73. SlepcInitialize(0,{},0,0);
  74. sparselizard();
  75. SlepcFinalize();
  76. return 0;
  77. }