awefem.py 2.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100
  1. """
  2. Solve the constant velocity scalar wave equation in an arbitrary number of dimensions.
  3. It injects a point source with a time-dependent source time function.
  4. """
  5. #Original script by Carlos da Costa:
  6. #https://github.com/cako/fenics-scripts/blob/master/awefem/awefem.py
  7. #
  8. from dolfin import *
  9. from vedo import ProgressBar, printc, download, settings
  10. from vedo.dolfin import plot
  11. import numpy as np
  12. set_log_level(30)
  13. def ricker_source(t, f=40):
  14. t -= 2 / f
  15. return (1 - 2 * (np.pi*f*t)**2) * np.exp(-(np.pi*f*t)**2)
  16. def sine_source(t, f=40):
  17. return np.sin(2 * np.pi*f*t)
  18. def awefem(mesh, t, source_loc=None):
  19. # Function space
  20. V = FunctionSpace(mesh, "Lagrange", 1)
  21. # Boundary condition
  22. bc = DirichletBC(V, Constant(0), "on_boundary")
  23. # Trial and test functions
  24. u = TrialFunction(V)
  25. v = TestFunction(V)
  26. # Discretization
  27. c = 6
  28. dt = t[1] - t[0]
  29. u0 = Function(V) # u0 = uN-1
  30. u1 = Function(V) # u1 = uN1
  31. # Variational formulation
  32. F = (u - 2 * u1 + u0) * v * dx + (dt * c) ** 2 * dot(
  33. grad(u + 2 * u1 + u0) / 4, grad(v) ) * dx
  34. a, L = lhs(F), rhs(F)
  35. # Solver
  36. A, b = assemble_system(a, L)
  37. solver = LUSolver(A, "mumps")
  38. solver.parameters["symmetric"] = True
  39. bc.apply(A, b)
  40. # Solution
  41. u = Function(V) # uN+1
  42. # Source
  43. if source_loc is None:
  44. mesh_center = np.mean(mesh.coordinates(), axis=0)
  45. source_loc = Point(mesh_center)
  46. else:
  47. source_loc = Point(source_loc)
  48. # Time stepping
  49. printc('\bomb Hit F1 to interrupt.', c='yellow')
  50. pb = ProgressBar(0, len(t))
  51. for i, t_ in enumerate(t[1:]):
  52. pb.print()
  53. b = assemble(L)
  54. delta = PointSource(V, source_loc, ricker_source(t_) * dt**2)
  55. delta.apply(b)
  56. solver.solve(u.vector(), b)
  57. u0.assign(u1)
  58. u1.assign(u)
  59. if t_>0.03:
  60. plt = plot(
  61. u,
  62. warp_zfactor=20, # set elevation along z
  63. vmin=.0, # sets a minimum to the color scale
  64. vmax=0.003,
  65. cmap='rainbow', # the color map style
  66. alpha=1, # transparency of the mesh
  67. lw=0.1, # linewidth of mesh
  68. scalarbar=None,
  69. #lighting='plastic',
  70. #elevation=-.3,
  71. interactive=False,
  72. ) # continue execution
  73. plt.clear()
  74. plt.interactive()
  75. if __name__ == "__main__":
  76. ot, dt, nt = 0, 1e-3, 150
  77. t = ot + np.arange(nt) * dt
  78. print("Computing wavefields over dolfin mesh")
  79. fpath = download("https://vedo.embl.es/examples/data/dolfin_fine.xml")
  80. mesh = Mesh(fpath)
  81. awefem(mesh, t, source_loc=(0.8, 0.8))