main.cpp 6.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167
  1. // This code simulates the eigenmodes of two nearby rectangular SiN dielectric photonic
  2. // waveguides buried in a SiO2 clad. To do so the 2D cross sectional geometry is
  3. // considered and the propagation is assumed to be along the z axis.
  4. // To solve this eigenvalue problem in 2D while taking into account the z component
  5. // the electric field is split into a transverse and a longitudinal part.
  6. // Since the searched modes are confined in and around the waveguide, a perfect
  7. // conductor boundary condition is valid, as long as the domain is large compared
  8. // to the waveguide dimensions.
  9. //
  10. // More details can be found in 'NASA technical paper 3485'.
  11. //
  12. //
  13. // The waveguides are 500 nm wide and 250 nm high. Their spacing is 100 nm.
  14. // _____________________________________________________
  15. // | |
  16. // | SiO2 ___________ ___________ |
  17. // | | | | | |
  18. // | | SiN | | SiN | |
  19. // | |___________| |___________| |
  20. // | |
  21. // |_____________________________________________________|
  22. //
  23. // Credits: R. Haouari
  24. #include "sparselizardbase.h"
  25. using namespace mathop;
  26. void sparselizard(void)
  27. {
  28. // Give names to the physical region numbers :
  29. int wavg1 = 1, wavg2 = 2, clad = 3, all = 4, bound = 5;
  30. mesh mymesh;
  31. mymesh.regionskin(bound, all);
  32. mymesh.load("optical_waveguide.msh");
  33. int waveguides = regionunion({wavg1, wavg2});
  34. // Waveguide boundary plotting:
  35. expression(1).write(regionintersection({waveguides, clad}), "skin.vtu");
  36. wallclock clk;
  37. // Edge shape functions 'hcurl' for the tranverse electric field Et.
  38. // Nodal shape functions 'h1' for the longitudinal electric field El.
  39. field Et("hcurl"), El("h1");
  40. // Use interpolation order 2 on the whole domain:
  41. Et.setorder(all, 2);
  42. El.setorder(all, 2);
  43. // Material properties definition
  44. parameter n, epsr, mur;
  45. n|clad = 1.4;
  46. n|wavg1 = 2.0;
  47. n|wavg2 = 2.0;
  48. epsr|all = n*n;
  49. mur|all = 1.0;
  50. // Light property
  51. double lambda = 680e-9, c = 299792458, f0 = c/lambda, k0 = 2.0*getpi()/lambda;
  52. std::cout << std::endl << "lambda = " << lambda*1e9 <<" nm, f0 = " << f0 << " Hz, k0 = " << k0 << " 1/m" << std::endl << std::endl;
  53. // Perfect conductor boundary condition:
  54. El.setconstraint(bound);
  55. Et.setconstraint(bound);
  56. // Weak formulation for the eigenvalue problem for confined propagation of an EM wave in a waveguide:
  57. //
  58. // 1/mur*curl(Et)*curl(Et') - k0^2*epsr*Et*Et' = -bt^2/mur*( grad(El)*Et') + Et*Et' )
  59. //
  60. // bt^2/mur*( grad(El)*grad(El') + Et*grad(El') ) = k0^2*bt^2*epsr* El*El'
  61. //
  62. // where curl() and grad() have a special definition in the transverse (xy) plane.
  63. //
  64. // The bt constant coming from the space derivation of exp(i*bt*z) is replaced by a time derivative dt().
  65. // Operators grad() and curl() in the transverse plane:
  66. expression dtdtgradEl(3,1,{dtdt(dx(dof(El))), dtdt(dy(dof(El))), 0});
  67. expression gradtfEl(3,1,{dx(tf(El)), dy(tf(El)), 0});
  68. formulation mode;
  69. mode += integral(all, curl(dof(Et))*curl(tf(Et))- k0*k0*mur*epsr*(dof(Et))*tf(Et));
  70. mode += integral(all, dtdtgradEl*tf(Et)+ dtdt(dof(Et))*tf(Et));
  71. mode += integral(all, dtdtgradEl*gradtfEl+dtdt(dof(Et))*gradtfEl);
  72. mode += integral(all, -k0*k0*mur*epsr*dtdt(dof(El))*tf(El));
  73. mode.generate();
  74. // Get the stiffness matrix K and mass matrix M:
  75. mat K = mode.K();
  76. mat M = mode.M();
  77. // Remove the rows and columns corresponding to 0 constraints:
  78. K.removeconstraints();
  79. M.removeconstraints();
  80. // Create the object to solve the eigenvalue problem:
  81. eigenvalue eig(K, M);
  82. // Compute the 5 eigenvalues closest to the target magnitude bt^2 (propagation constant^2).
  83. // We are looking for modes around neff_target, the target effective refractive index.
  84. double neff_target = 1.6, bt = k0*neff_target, bt2 = std::pow(bt,2.0);
  85. std::cout << "Target is neff = " << neff_target << ", beta = " << bt << " rad/m" << std::endl << std::endl;
  86. // Propagation mode eigenvalue is purely imaginary (i*bt), we thus target -bt^2:
  87. eig.compute(5, -bt2);
  88. // Get all eigenvectors and eigenvalues found:
  89. std::vector<vec> myrealeigenvectors = eig.geteigenvectorrealpart();
  90. std::vector<double> myrealeigenvalues = eig.geteigenvaluerealpart();
  91. std::vector<double> neffcs;
  92. expression Etotal = array3x1(compx(Et), compy(Et), El);
  93. // Loop on all eigenvalues found:
  94. int index = 1;
  95. for (int i = 0; i < myrealeigenvalues.size(); i++)
  96. {
  97. if (myrealeigenvalues[i] < 0)
  98. {
  99. // Transfer the data from the ith eigenvector to fields Et and El:
  100. Et.setdata(all, myrealeigenvectors[i]);
  101. El.setdata(all, myrealeigenvectors[i]);
  102. // Compute the propagation constant and the effective refractive index:
  103. double btc = std::sqrt(std::abs(myrealeigenvalues[i]));
  104. double neffc = btc/k0;
  105. // We need to write separately on the clad and waveguide regions to visualize the discontinuity:
  106. Etotal.write(clad, "Eclad_"+ std::to_string(index) +".vtu", 2);
  107. Etotal.write(regionunion({wavg1,wavg2}), "Ewav_"+ std::to_string(index) +".vtu", 2);
  108. // Display mode information:
  109. std::cout << "Mode " << index << ": beta = " << btc << " rad/m, neff = " << neffc << std::endl;
  110. neffcs.push_back(neffc);
  111. index++;
  112. }
  113. }
  114. // Write a .pvd file to visualize all the modes together:
  115. grouptimesteps("Eclad.pvd", "Eclad_", 1, neffcs);
  116. grouptimesteps("Ewav.pvd", "Ewav_", 1, neffcs);
  117. std::cout << std::endl;
  118. clk.print("Total time elapsed:");
  119. // Code validation line. Can be removed.
  120. std::cout << (neffcs[0] < 1.65780 && neffcs[0] > 1.65778);
  121. }
  122. int main(void)
  123. {
  124. SlepcInitialize(0,{},0,0);
  125. sparselizard();
  126. SlepcFinalize();
  127. return 0;
  128. }