main.cpp 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138
  1. // This code simulates the laminar, incompressible water flow past a step.
  2. //
  3. // A parabolic normal flow velocity is forced at the inlet and
  4. // a zero pressure is imposed at the outlet.
  5. //
  6. // More information for this standard example can be found in:
  7. //
  8. // "Finite element methods for the incompressible Navier-Stokes equations", A. Segal
  9. //
  10. // The step geometry is as follows:
  11. //
  12. // l2
  13. // <--------------------------------->
  14. // __________________________________
  15. // | /|\
  16. // l1 | |
  17. // <----------->| |
  18. // _____________| h2 \|/
  19. // /|\ OUTLET
  20. // |
  21. // INLET |
  22. // \|/ h1
  23. // ________________________________________________
  24. #include "sparselizardbase.h"
  25. using namespace mathop;
  26. // First four arguments give the geometry dimension, last four arguments give the mesh size:
  27. mesh createmesh(double l1, double h1, double l2, double h2, int nl1, int nl2, int nh1, int nh2);
  28. void sparselizard(void)
  29. {
  30. // Region numbers used in this simulation:
  31. int fluid = 1, inlet = 2, outlet = 3, skin = 4;
  32. // Height of the inlet [m]:
  33. double h1 = 1e-3;
  34. mesh mymesh = createmesh(2e-3, h1, 12e-3, 1e-3, 30, 150, 20, 50);
  35. // Define the fluid wall (not including the inlet and outlet):
  36. int wall = regionexclusion(skin, regionunion({inlet,outlet}));
  37. // Dynamic viscosity of water [Pa.s] and density [kg/m3]:
  38. double mu = 8.9e-4, rho = 1000;
  39. // Field v is the flow velocity. It uses nodal shape functions "h1" with two components.
  40. // Field p is the relative pressure. Fields x and y are the space coordinates.
  41. field v("h1xy"), p("h1"), x("x"), y("y");
  42. // Force the flow velocity to 0 on the wall:
  43. v.setconstraint(wall);
  44. // Set a 0 relative pressure at the outlet:
  45. p.setconstraint(outlet);
  46. // Use an order 1 interpolation for p and 2 for v on the fluid region (satisfies the BB condition):
  47. p.setorder(fluid, 1);
  48. v.setorder(fluid, 2);
  49. // Define the weak formulation for incompressible laminar flow:
  50. formulation laminarflow;
  51. laminarflow += integral(fluid, predefinednavierstokes(dof(v), tf(v), v, dof(p), tf(p), mu, rho, 0, 0) );
  52. // This loop with the above formulation is a Newton iteration:
  53. int index = 0; double convergence = 1, velocity = 0.1;
  54. while (convergence > 1e-10)
  55. {
  56. // Slowly increase the velocity for a high Reynolds (1 m/s flow, 1000 Reynolds still converges).
  57. if (velocity < 0.299)
  58. velocity = velocity + 0.008;
  59. std::cout << "Flow velocity: " << velocity << " m/s" << std::endl;
  60. // Force the flow velocity at the inlet (quadratic profile w.r.t. the y axis):
  61. v.setconstraint(inlet, array2x1(velocity*y*(h1-y)/pow(h1*0.5,2), 0));
  62. // Get a measure of the solution for convergence evaluation:
  63. double measuresol = norm(v).integrate(fluid,2);
  64. // Generate and solve the laminar flow problem then save to the fields:
  65. solve(laminarflow);
  66. convergence = std::abs((norm(v).integrate(fluid,2) - measuresol)/norm(v).integrate(fluid,2));
  67. std::cout << "Relative solution change: " << convergence << std::endl;
  68. // Write the fields to file. Write v with an order 2 interpolation.
  69. p.write(fluid, "p.vtk");
  70. v.write(fluid, "v.vtk", 2);
  71. index++;
  72. }
  73. // Compute the flow velocity norm at position (5,1,0) mm in space:
  74. double vnorm = norm(v).interpolate(fluid, {5e-3, 1e-3, 0.0})[0];
  75. // Output the input and output flowrate for a unit width:
  76. double flowratein = (normal(inlet)*v).integrate(inlet, 4);
  77. double flowrateout = -(normal(outlet)*v).integrate(outlet, 4);
  78. std::cout << std::endl << "Flowrate in/out for a unit width: " << flowratein << " / " << flowrateout << " m^3/s" << std::endl;
  79. // Code validation line. Can be removed.
  80. std::cout << (vnorm*flowrateout < 2.64489e-05 && vnorm*flowrateout > 2.64485e-05);
  81. }
  82. mesh createmesh(double l1, double h1, double l2, double h2, int nl1, int nl2, int nh1, int nh2)
  83. {
  84. int fluid = 1, inlet = 2, outlet = 3, skin = 4;
  85. shape qthinleft("quadrangle", fluid, {0,0,0, l1,0,0, l1,h1,0, 0,h1,0}, {nl1,nh1,nl1,nh1});
  86. shape qthinright("quadrangle", fluid, {l1,0,0, l1+l2,0,0, l1+l2,h1,0, l1,h1,0}, {nl2,nh1,nl2,nh1});
  87. shape qthick("quadrangle", fluid, {l1,h1,0, l1+l2,h1,0, l1+l2,h1+h2,0, l1,h1+h2,0}, {nl2,nh2,nl2,nh2});
  88. shape linlet = qthinleft.getsons()[3];
  89. linlet.setphysicalregion(inlet);
  90. shape loutlet = shape("union", outlet, {qthick.getsons()[1],qthinright.getsons()[1]});
  91. mesh mymesh;
  92. mymesh.regionskin(skin, fluid);
  93. mymesh.load({qthinleft,qthinright,qthick,linlet,loutlet});
  94. mymesh.write("pipestep.msh");
  95. return mymesh;
  96. }
  97. int main(void)
  98. {
  99. SlepcInitialize(0,{},0,0);
  100. sparselizard();
  101. SlepcFinalize();
  102. return 0;
  103. }