elasticity3.py 3.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104
  1. #!/usr/bin/env python3
  2. """An initial circle is stretched by means of a variable force into its final shape.
  3. Colored lines are the trajectories of a few initial points."""
  4. from dolfin import *
  5. from mshr import *
  6. import numpy as np
  7. import vedo
  8. from vedo.dolfin import plot
  9. set_log_level(30)
  10. class AllBoundaries(SubDomain):
  11. def inside(self, x, on_boundary):
  12. return on_boundary and x[0]<-10
  13. def solve_problem(mesh, mfunc, force):
  14. V = VectorFunctionSpace(mesh, "CG", 1)
  15. u = TrialFunction(V)
  16. v = TestFunction(V)
  17. displacement = Function(V)
  18. bc = [DirichletBC(V, Constant((0,0)), mfunc, 1)]
  19. F = Constant(force)
  20. E = Constant(5000)
  21. nu = Constant(0.3)
  22. mu = E / (2.0*(1+nu))
  23. lmbda = E * nu / ((1.0+nu) * (1.0-2*nu))
  24. sigma = 2.0 * mu * sym(grad(u)) + lmbda * tr( sym(grad(u)) ) * Identity(2)
  25. solve(inner(sigma, grad(v)) * dx == inner(F, v) * dx, displacement, bc)
  26. displacement.set_allow_extrapolation(True)
  27. return displacement
  28. def update(mesh, displacement):
  29. new_mesh = Mesh(mesh)
  30. ALE.move(new_mesh, displacement)
  31. return new_mesh
  32. def remesh(mesh, res=50):
  33. if isinstance(mesh, vedo.Mesh):
  34. vmesh = mesh
  35. else:
  36. vmesh = vedo.Mesh([mesh.coordinates(), mesh.cells()])
  37. bpts = vmesh.computeNormals(cells=True).boundaries().join(reset=1) #extract boundary
  38. vz = vmesh.celldata["Normals"][0][2] # check z component of normals at first point
  39. bpts.generate_mesh(invert=vz<0).smooth().write('tmpmesh.xml') #vedo REMESHING + smoothing
  40. return Mesh("tmpmesh.xml")
  41. #################################################################################
  42. N = 40 # number of iterations of stretching
  43. res = 15 # resolution of meshes
  44. do_remesh = False # grab the boundary and remesh the interior at each iteration
  45. vedo.settings.use_parallel_projection = True # avoid perspective parallax
  46. circle = Circle(Point(0, 0), 50)
  47. mesh = generate_mesh(circle, res)
  48. meshes = [mesh]
  49. displacements = []
  50. for i in range(N):
  51. mfunc = MeshFunction('size_t', mesh, 1, mesh.domains())
  52. mfunc.set_all(0)
  53. allb = AllBoundaries()
  54. allb.mark(mfunc, 1)
  55. F = np.array([2, (i-N/2)/N]) # some variable force
  56. displacement = solve_problem(mesh, mfunc, F)
  57. new_mesh = update(mesh, displacement)
  58. if do_remesh:
  59. mesh = remesh(new_mesh)
  60. else:
  61. mesh = new_mesh
  62. meshes.append(mesh)
  63. displacements.append(displacement)
  64. dmesh_i = meshes[0] # initial mesh
  65. dmesh_f = meshes[-1] # final mesh
  66. vmesh_i = vedo.Mesh([dmesh_i.coordinates(), dmesh_i.cells()], c='grey5').z(-1)
  67. vmesh_f = vedo.Mesh([dmesh_f.coordinates(), dmesh_f.cells()], c='grey3').wireframe()
  68. plt = vedo.Plotter()
  69. # move a few points along the deformation of the circle
  70. seeds = vedo.Circle(r=50, res=res).vertices[:,(0,1)] # make points 2d with [:,(0,1)]
  71. endpoints = []
  72. for i, p in enumerate(seeds):
  73. line = [p]
  74. for u in displacements:
  75. p = p + u(p)
  76. line.append(p)
  77. plt += vedo.Line(line, c=i, lw=4).z(1)
  78. endpoints.append(p)
  79. plt += [vmesh_i, vmesh_f, __doc__]
  80. plt.show(axes=1)
  81. # to invert everything and move the end points back in place, check out discussion:
  82. #https://fenicsproject.discourse.group/t/precision-in-hyperelastic-model/6824/3