simulation.html 76 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562
  1. <!doctype html>
  2. <html lang="en">
  3. <head>
  4. <meta charset="utf-8">
  5. <meta name="viewport" content="width=device-width, initial-scale=1, minimum-scale=1" />
  6. <meta name="generator" content="pdoc 0.5.2" />
  7. <title>simulation API documentation</title>
  8. <meta name="description" content="Definition of the Simulation class and the Galaxy constructor." />
  9. <link href='https://cdnjs.cloudflare.com/ajax/libs/normalize/8.0.0/normalize.min.css' rel='stylesheet'>
  10. <link href='https://cdnjs.cloudflare.com/ajax/libs/10up-sanitize.css/8.0.0/sanitize.min.css' rel='stylesheet'>
  11. <link href="https://cdnjs.cloudflare.com/ajax/libs/highlight.js/9.12.0/styles/github.min.css" rel="stylesheet">
  12. <style>.flex{display:flex !important}body{line-height:1.5em}#content{padding:20px}#sidebar{padding:30px;overflow:hidden}.http-server-breadcrumbs{font-size:130%;margin:0 0 15px 0}#footer{font-size:.75em;padding:5px 30px;border-top:1px solid #ddd;text-align:right}#footer p{margin:0 0 0 1em;display:inline-block}#footer p:last-child{margin-right:30px}h1,h2,h3,h4,h5{font-weight:300}h1{font-size:2.5em;line-height:1.1em}h2{font-size:1.75em;margin:1em 0 .50em 0}h3{font-size:1.4em;margin:25px 0 10px 0}h4{margin:0;font-size:105%}a{color:#058;text-decoration:none;transition:color .3s ease-in-out}a:hover{color:#e82}.title code{font-weight:bold}h2[id^="header-"]{margin-top:2em}.ident{color:#900}pre code{background:#f8f8f8;font-size:.8em;line-height:1.4em}code{background:#f2f2f1;padding:1px 4px;overflow-wrap:break-word}h1 code{background:transparent}pre{background:#f8f8f8;border:0;border-top:1px solid #ccc;border-bottom:1px solid #ccc;margin:1em 0;padding:1ex}#http-server-module-list{display:flex;flex-flow:column}#http-server-module-list div{display:flex}#http-server-module-list dt{min-width:10%}#http-server-module-list p{margin-top:0}.toc ul,#index{list-style-type:none;margin:0;padding:0}#index code{background:transparent}#index h3{border-bottom:1px solid #ddd}#index ul{padding:0}#index h4{font-weight:bold}#index h4 + ul{margin-bottom:.6em}#index .two-column{column-count:2}dl{margin-bottom:2em}dl dl:last-child{margin-bottom:4em}dd{margin:0 0 1em 3em}#header-classes + dl > dd{margin-bottom:3em}dd dd{margin-left:2em}dd p{margin:10px 0}.name{background:#eee;font-weight:bold;font-size:.85em;padding:5px 10px;display:inline-block;min-width:40%}.name:hover{background:#e0e0e0}.name > span:first-child{white-space:nowrap}.name.class > span:nth-child(2){margin-left:.4em}.name small{font-weight:normal}.inherited{color:#999;border-left:5px solid #eee;padding-left:1em}.inheritance em{font-style:normal;font-weight:bold}.desc h2{font-weight:400;font-size:1.25em}.desc h3{font-size:1em}.desc dt code{background:inherit}.source summary{color:#666;text-align:right;font-weight:400;font-size:.8em;text-transform:uppercase;cursor:pointer}.source pre{max-height:500px;overflow:auto;margin:0}.source pre code{overflow:visible}.hlist{list-style:none}.hlist li{display:inline}.hlist li:after{content:',\2002'}.hlist li:last-child:after{content:none}.hlist .hlist{display:inline;padding-left:1em}img{max-width:100%}.admonition{padding:.1em .5em}.admonition-title{font-weight:bold}.admonition.note,.admonition.info,.admonition.important{background:#aef}.admonition.todo,.admonition.versionadded,.admonition.tip,.admonition.hint{background:#dfd}.admonition.warning,.admonition.versionchanged,.admonition.deprecated{background:#fd4}.admonition.error,.admonition.danger,.admonition.caution{background:lightpink}</style>
  13. <style media="screen and (min-width: 700px)">@media screen and (min-width:700px){#sidebar{width:30%}#content{width:70%;max-width:100ch;padding:3em 4em;border-left:1px solid #ddd}pre code{font-size:1em}.item .name{font-size:1em}main{display:flex;flex-direction:row-reverse;justify-content:flex-end}.toc ul ul,#index ul{padding-left:1.5em}.toc > ul > li{margin-top:.5em}}</style>
  14. <style media="print">@media print{#sidebar h1{page-break-before:always}.source{display:none}}@media print{*{background:transparent !important;color:#000 !important;box-shadow:none !important;text-shadow:none !important}a[href]:after{content:" (" attr(href) ")";font-size:90%}a[href][title]:after{content:none}abbr[title]:after{content:" (" attr(title) ")"}.ir a:after,a[href^="javascript:"]:after,a[href^="#"]:after{content:""}pre,blockquote{border:1px solid #999;page-break-inside:avoid}thead{display:table-header-group}tr,img{page-break-inside:avoid}img{max-width:100% !important}@page{margin:0.5cm}p,h2,h3{orphans:3;widows:3}h1,h2,h3,h4,h5,h6{page-break-after:avoid}}</style>
  15. </head>
  16. <body>
  17. <main>
  18. <article id="content">
  19. <header>
  20. <h1 class="title"><code>simulation</code> module</h1>
  21. </header>
  22. <section id="section-intro">
  23. <p>Definition of the Simulation class and the Galaxy constructor.</p>
  24. <details class="source">
  25. <summary>Source code</summary>
  26. <pre><code class="python">&#34;&#34;&#34;Definition of the Simulation class and the Galaxy constructor.&#34;&#34;&#34;
  27. import os
  28. import pickle
  29. import numpy as np
  30. import matplotlib.pyplot as plt
  31. from utils import random_unit_vectors, cascade_round
  32. from distributions import PLUMMER, HERNQUIST, UNIFORM, EXP, NFW
  33. import acceleration
  34. ##############################################################################
  35. ##############################################################################
  36. class Simulation:
  37. &#34;&#34;&#34;&#34;Main class for the gravitational simulation.
  38. Attributes:
  39. r_vec (array): position of the particles in the current timestep.
  40. Shape: (number of particles, 3)
  41. rprev_vec (array): position of the particles in the previous timestep.
  42. Shape: (number of particles, 3)
  43. v_vec (array): velocity in the current timestep.
  44. Shape: (number of particles, 3)
  45. a_vec (array): acceleration in the current timestep.
  46. Shape: (number of particles, 3)
  47. mass (array): mass of each particle in the simulation.
  48. Shape: (number of particles,)
  49. type (array): non-unique identifier for each particle.
  50. Shape: (number of particles,)
  51. tracks (array): list of positions through the simulation for central
  52. masses. Shape: (tracked particles, n+1, 3).
  53. CONFIG (array): configuration used to create the simulation.
  54. It will be saved along the state of the simulation.
  55. dt (float): timestep of the simulation
  56. n (int): current timestep. Initialized as n=0.
  57. soft (float): softening length used by the simulation.
  58. verbose (boolean): When True progress statements will be printed.
  59. &#34;&#34;&#34;
  60. def __init__(self, dt, soft, verbose, CONFIG, method, **kwargs):
  61. &#34;&#34;&#34;Constructor for the Simulation class.
  62. Arguments:
  63. dt (float): timestep of the simulation
  64. n (int): current timestep. Initialized as n=0.
  65. soft (float): softening length used by the simulation.
  66. verbose (bool): When True progress statements will be printed.
  67. CONFIG (dict): configuration file used to create the simulation.
  68. method (string): Optional. Algorithm to use when computing the
  69. gravitational forces. One of &#39;bruteForce&#39;, &#39;bruteForce_numba&#39;,
  70. &#39;bruteForce_numbaopt&#39;, &#39;bruteForce_CPP&#39;, &#39;barnesHut_CPP&#39;.
  71. &#34;&#34;&#34;
  72. self.n = 0
  73. self.t = 0
  74. self.dt = dt
  75. self.soft = soft
  76. self.verbose = verbose
  77. self.CONFIG = CONFIG
  78. # Initialize empty arrays for all necessary properties
  79. self.r_vec = np.empty((0, 3))
  80. self.v_vec = np.empty((0, 3))
  81. self.a_vec = np.empty((0, 3))
  82. self.mass = np.empty((0,))
  83. self.type = np.empty((0, 2))
  84. algorithms = {
  85. &#39;bruteForce&#39;: acceleration.bruteForce,
  86. &#39;bruteForceNumba&#39;: acceleration.bruteForceNumba,
  87. &#39;bruteForceNumbaOptimized&#39;: acceleration.bruteForceNumbaOptimized,
  88. &#39;bruteForceCPP&#39;: acceleration.bruteForceCPP,
  89. &#39;barnesHutCPP&#39;: acceleration.barnesHutCPP
  90. }
  91. try:
  92. self.acceleration = algorithms[method]
  93. except: raise Exception(&#34;Method &#39;{}&#39; unknown&#34;.format(method))
  94. def add(self, body):
  95. &#34;&#34;&#34;Add a body to the simulation. It must expose the public attributes
  96. body.r_vec, body.v_vec, body.a_vec, body.type, body.mass.
  97. Arguments:
  98. body: Object to be added to the simulation (e.g. a Galaxy object)
  99. &#34;&#34;&#34;
  100. # Extend all relevant attributes by concatenating the body
  101. for name in [&#39;r_vec&#39;, &#39;v_vec&#39;, &#39;a_vec&#39;, &#39;type&#39;, &#39;mass&#39;]:
  102. simattr, bodyattr = getattr(self, name), getattr(body, name)
  103. setattr(self, name, np.concatenate([simattr, bodyattr], axis=0))
  104. # Order based on mass
  105. order = np.argsort(-self.mass)
  106. for name in [&#39;r_vec&#39;, &#39;v_vec&#39;, &#39;a_vec&#39;, &#39;type&#39;, &#39;mass&#39;]:
  107. setattr(self, name, getattr(self, name)[order])
  108. # Update the list of objects to keep track of
  109. self.tracks = np.empty((np.sum(self.type[:,0]==&#39;center&#39;), 0, 3))
  110. def step(self):
  111. &#34;&#34;&#34;Perform a single step of the simulation.
  112. Makes use of a 4th order Verlet integrator.
  113. &#34;&#34;&#34;
  114. # Calculate the acceleration
  115. self.a_vec = self.acceleration(self.r_vec, self.mass, soft=self.soft)
  116. # Update the state using the Verlet algorithm
  117. # (A custom algorithm is written mainly for learning purposes)
  118. self.r_vec, self.rprev_vec = (2*self.r_vec - self.rprev_vec
  119. + self.a_vec * self.dt**2, self.r_vec)
  120. self.n += 1
  121. # Update tracks
  122. self.tracks = np.concatenate([self.tracks,
  123. self.r_vec[self.type[:,0]==&#39;center&#39;][:,np.newaxis]], axis=1)
  124. def run(self, tmax, saveEvery, outputFolder, **kwargs):
  125. &#34;&#34;&#34;Run the galactic simulation.
  126. Attributes:
  127. tmax (float): Time to which the simulation will run to.
  128. This is measured here since the start of the simulation,
  129. not since pericenter.
  130. saveEvery (int): The state is saved every saveEvery steps.
  131. outputFolder (string): It will be saved to /data/outputFolder/
  132. &#34;&#34;&#34;
  133. # When the simulation starts, intialize self.rprev_vec
  134. self.rprev_vec = self.r_vec - self.v_vec * self.dt
  135. if self.verbose: print(&#39;Simulation starting. Bon voyage!&#39;)
  136. while(self.t &lt; tmax):
  137. self.step()
  138. if(self.n % saveEvery == 0):
  139. self.save(&#39;data/{}&#39;.format(outputFolder))
  140. print(&#39;Simulation complete.&#39;)
  141. def save(self, outputFolder):
  142. &#34;&#34;&#34;Save the state of the simulation to the outputFolder.
  143. Two files are saved:
  144. sim{self.n}.pickle: serializing the state.
  145. sim{self.n}.png: a simplified 2D plot of x, y.
  146. &#34;&#34;&#34;
  147. # Create the output folder if it doesn&#39;t exist
  148. if not os.path.exists(outputFolder): os.makedirs(outputFolder)
  149. # Compute some useful quantities
  150. # v_vec is not required by the integrator, but useful
  151. self.v_vec = (self.r_vec - self.rprev_vec)/self.dt
  152. self.t = self.n * self.dt # prevents numerical rounding errors
  153. # Serialize state
  154. file = open(outputFolder+&#39;/data{}.pickle&#39;.format(self.n), &#34;wb&#34;)
  155. pickle.dump({&#39;r_vec&#39;: self.r_vec, &#39;v_vec&#39;: self.v_vec,
  156. &#39;type&#39;: self.type, &#39;mass&#39;: self.mass,
  157. &#39;CONFIG&#39;: self.CONFIG, &#39;t&#39;: self.t,
  158. &#39;tracks&#39;: self.tracks}, file)
  159. # Save simplified plot of the current state.
  160. # Its main use is to detect ill-behaved situations early on.
  161. fig = plt.figure()
  162. plt.xlim(-5, 5); plt.ylim(-5, 5); plt.axis(&#39;equal&#39;)
  163. # Dark halo is plotted in red, disk in blue, bulge in green
  164. PLTCON = [(&#39;dark&#39;, &#39;r&#39;, 0.3), (&#39;disk&#39;, &#39;b&#39;, 1.0), (&#39;bulge&#39;, &#39;g&#39;, 0.5)]
  165. for type_, c, a in PLTCON:
  166. plt.scatter(self.r_vec[self.type[:,0]==type_][:,0],
  167. self.r_vec[self.type[:,0]==type_][:,1], s=0.1, c=c, alpha=a)
  168. # Central mass as a magenta star
  169. plt.scatter(self.r_vec[self.type[:,0]==&#39;center&#39;][:,0],
  170. self.r_vec[self.type[:,0]==&#39;center&#39;][:,1], s=100, marker=&#34;*&#34;, c=&#39;m&#39;)
  171. # Save to png file
  172. fig.savefig(outputFolder+&#39;/sim{}.png&#39;.format(self.n), dpi=150)
  173. plt.close(fig)
  174. def project(self, theta, phi, view=0):
  175. &#34;&#34;&#34;Projects the 3D simulation onto a plane as viewed from the
  176. direction described by the (theta, phi, view). Angles in radians.
  177. (This is used by the simulated annealing algorithm)
  178. Parameters:
  179. theta (float): polar angle.
  180. phi (float): azimuthal angle.
  181. view (float): rotation along line of sight.
  182. &#34;&#34;&#34;
  183. M1 = np.array([[np.cos(phi), np.sin(phi), 0],
  184. [-np.sin(phi), np.cos(phi), 0],
  185. [0, 0, 1]])
  186. M2 = np.array([[1, 0, 0],
  187. [0, np.cos(theta), np.sin(theta)],
  188. [0, -np.sin(theta), np.cos(theta)]])
  189. M3 = np.array([[np.cos(view), np.sin(view), 0],
  190. [-np.sin(view), np.cos(view), 0],
  191. [0, 0, 1]])
  192. M = np.matmul(M1, np.matmul(M2, M3)) # combine rotations
  193. r = np.tensordot(self.r_vec, M, axes=[1, 0])
  194. return r[:,0:2]
  195. def setOrbit(self, galaxy1, galaxy2, e, rmin, R0):
  196. &#34;&#34;&#34;Sets the two galaxies galaxy1, galaxy2 in an orbit.
  197. Parameters:
  198. galaxy1 (Galaxy): 1st galaxy (main)
  199. galaxy2 (Galaxy): 2nd galaxy (companion)
  200. e: eccentricity of the orbit
  201. rmin: distance of closest approach
  202. R0: initial separation
  203. &#34;&#34;&#34;
  204. m1, m2 = np.sum(galaxy1.mass), np.sum(galaxy2.mass)
  205. # Relevant formulae:
  206. # $r_0 = r (1 + e) \cos(\phi)$, where $r_0$ ($\neq R_0$) is the semi-latus rectum
  207. # $r_0 = r_\textup{min} (1 + e)$
  208. # $v^2 = G M (2/r - 1/a)$, where a is the semimajor axis
  209. # Solve the reduced two-body problem with reduced mass $\mu$ (mu)
  210. M = m1 + m2
  211. r0 = rmin * (1 + e)
  212. try:
  213. phi = np.arccos((r0/R0 - 1) / e) # inverting the orbit equation
  214. phi = -np.abs(phi) # Choose negative (incoming) angle
  215. ainv = (1 - e) / rmin # ainv = $1/a$, as a may be infinite
  216. v = np.sqrt(M * (2/R0 - ainv))
  217. # Finally, calculate the angle of motion. angle = tan(dy/dx)
  218. # $dy/dx = ((dr/d\phi) sin(\phi) + r \cos(\phi))/((dr/d\phi) cos(\phi) - r \sin(\phi))$
  219. dy = R0/r0 * e * np.sin(phi)**2 + np.cos(phi)
  220. dx = R0/r0 * e * np.sin(phi) * np.cos(phi) - np.sin(phi)
  221. vangle = np.arctan2(dy, dx)
  222. except: raise Exception(&#39;&#39;&#39;The orbital parameters cannot be satisfied.
  223. For elliptical orbits check that R0 is posible (&lt;rmax).&#39;&#39;&#39;)
  224. # We now need the actual motion of each body
  225. R_vec = np.array([[R0*np.cos(phi), R0*np.sin(phi), 0.]])
  226. V_vec = np.array([[v*np.cos(vangle), v*np.sin(vangle), 0.]])
  227. galaxy1.r_vec += m2/M * R_vec
  228. galaxy1.v_vec += m2/M * V_vec
  229. galaxy2.r_vec += -m1/M * R_vec
  230. galaxy2.v_vec += -m1/M * V_vec
  231. # Explicitely add the galaxies to the simulation
  232. self.add(galaxy1)
  233. self.add(galaxy2)
  234. if self.verbose: print(&#39;Galaxies set in orbit.&#39;)
  235. ##############################################################################
  236. ##############################################################################
  237. class Galaxy():
  238. &#34;&#34;&#34;&#34;Helper class for creating galaxies.
  239. Attributes:
  240. r_vec (array): position of the particles in the current timestep.
  241. Shape: (number of particles, 3)
  242. v_vec (array): velocity in the current timestep.
  243. Shape: (number of particles, 3)
  244. a_vec (array): acceleration in the current timestep.
  245. Shape: (number of particles, 3)
  246. mass (array): mass of each particle in the simulation.
  247. Shape: (number of particles,)
  248. type (array): non-unique identifier for each particle.
  249. Shape: (number of particles,) &#34;&#34;&#34;
  250. def __init__(self, orientation, centralMass, bulge, disk, halo, sim):
  251. &#34;&#34;&#34;Constructor for the Galaxy class.
  252. Parameters:
  253. orientation (tupple): (inclination i, argument of pericenter w)
  254. centralMass (float): mass at the center of the galaxy
  255. bulge (dict): passed to the addBulge method.
  256. disk (dict): passed to the addDisk method.
  257. halo (dict): passed to the addHalo method.
  258. sim (Simulation): simulation object
  259. &#34;&#34;&#34;
  260. if sim.verbose: print(&#39;Initializing galaxy&#39;)
  261. # Build the central mass
  262. self.r_vec = np.zeros((1, 3))
  263. self.v_vec = np.zeros((1, 3))
  264. self.a_vec = np.zeros((1, 3))
  265. self.mass = np.array([centralMass])
  266. self.type = np.array([[&#39;center&#39;, 0]])
  267. # Build the other components
  268. self.addBulge(**bulge)
  269. if sim.verbose: print(&#39;Bulge created.&#39;)
  270. self.addDisk(**disk)
  271. if sim.verbose: print(&#39;Disk created.&#39;)
  272. self.addHalo(**halo)
  273. if sim.verbose: print(&#39;Halo created.&#39;)
  274. # Correct particles velocities to generate circular orbits
  275. # $a_\textup{centripetal} = v^2/r$
  276. a_vec = sim.acceleration(self.r_vec, self.mass, soft=sim.soft)
  277. a = np.linalg.norm(a_vec, axis=-1, keepdims=True)
  278. r = np.linalg.norm(self.r_vec, axis=-1, keepdims=True)
  279. v = np.linalg.norm(self.v_vec[1:], axis=-1, keepdims=True)
  280. direction_unit = self.v_vec[1:]/v
  281. # Set orbital velocities (except for central mass)
  282. self.v_vec[1:] = np.sqrt(a[1:] * r[1:]) * direction_unit
  283. self.a_vec = np.zeros_like(self.r_vec)
  284. # Rotate the galaxy into its correct orientation
  285. self.rotate(*(np.array(orientation)/360 * 2*np.pi))
  286. if sim.verbose: print(&#39;Galaxy set in rotation and oriented.&#39;)
  287. def addBulge(self, model, totalMass, particles, l):
  288. &#34;&#34;&#34;Adds a bulge to the galaxy.
  289. Parameters:
  290. model (string): parametrization of the bulge.
  291. &#39;plummer&#39; and &#39;hernquist&#39; are supported.
  292. totalMass (float): total mass of the bulge
  293. particles (int): number of particles in the bulge
  294. l (float): characteristic length scale of the model.
  295. &#34;&#34;&#34;
  296. if particles == 0: return None
  297. # Divide the mass equally among all particles
  298. mass = np.ones(particles) * totalMass/particles
  299. self.mass = np.concatenate([self.mass, mass], axis=0)
  300. # Create particles according to the radial distribution from model
  301. if model == &#39;plummer&#39;:
  302. r = PLUMMER.ppf(np.random.rand(particles), scale=l)
  303. elif model == &#39;hernquist&#39;:
  304. r = HERNQUIST.ppf(np.random.rand(particles), scale=l)
  305. else: raise Exception(&#34;&#34;&#34;Bulge distribution not allowed.
  306. &#39;plummer&#39; and &#39;hernquist&#39; are the supported values&#34;&#34;&#34;)
  307. r_vec = r[:,np.newaxis] * random_unit_vectors(size=particles)
  308. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  309. # Set them orbitting along random directions normal to r_vec
  310. v_vec = np.cross(r_vec, random_unit_vectors(size=particles))
  311. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)
  312. # Label the particles
  313. type_ = [[&#39;bulge&#39;, 0]]*particles
  314. self.type = np.concatenate([self.type, type_], axis=0)
  315. def addDisk(self, model, totalMass, particles, l):
  316. &#34;&#34;&#34;Adds a disk to the galaxy.
  317. Parameters:
  318. model (string): parametrization of the disk.
  319. &#39;rings&#39;, &#39;uniform&#39; and &#39;exp&#39; are supported.
  320. totalMass (float): total mass of the bulge
  321. particles (int): number of particles in the bulge
  322. l: fot &#39;exp&#39; and &#39;uniform&#39; characteristic length of the
  323. model. For &#39;rings&#39; tupple of the form (inner radius,
  324. outer radius, number of rings)
  325. &#34;&#34;&#34;
  326. if particles == 0: return None
  327. # Create particles according to the radial distribution from model
  328. if model == &#39;uniform&#39;:
  329. r = UNIFORM.ppf(np.random.rand(particles), scale=l)
  330. type_ = [[&#39;disk&#39;, 0]]*particles
  331. r_vec = r[:,np.newaxis] * random_unit_vectors(particles, &#39;2D&#39;)
  332. self.type = np.concatenate([self.type, type_], axis=0)
  333. elif model == &#39;rings&#39;:
  334. # l = [inner radius, outter radius, number of rings]
  335. distances = np.linspace(*l)
  336. # Aim for roughly constant areal density
  337. # Cascade rounding preserves the total number of particles
  338. perRing = cascade_round(particles * distances / np.sum(distances))
  339. particles = int(np.sum(perRing)) # prevents numerical errors
  340. r_vec = np.empty((0, 3))
  341. for d, n, i in zip(distances, perRing, range(l[2])):
  342. type_ = [[&#39;disk&#39;, i+1]]*int(n)
  343. self.type = np.concatenate([self.type, type_], axis=0)
  344. phi = np.linspace(0, 2 * np.pi, n, endpoint=False)
  345. ringr = d * np.array([[np.cos(i), np.sin(i), 0] for i in phi])
  346. r_vec = np.concatenate([r_vec, ringr], axis=0)
  347. elif model == &#39;exp&#39;:
  348. r = EXP.ppf(np.random.rand(particles), scale=l)
  349. r_vec = r[:,np.newaxis] * random_unit_vectors(particles, &#39;2D&#39;)
  350. type_ = [[&#39;disk&#39;, 0]]*particles
  351. self.type = np.concatenate([self.type, type_], axis=0)
  352. else:
  353. raise Exception(&#34;&#34;&#34;Disk distribution not allowed.
  354. &#39;uniform&#39;, &#39;rings&#39; and &#39;exp&#39; are the supported values&#34;&#34;&#34;)
  355. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  356. # Divide the mass equally among all particles
  357. mass = np.ones(particles) * totalMass/particles
  358. self.mass = np.concatenate([self.mass, mass], axis=0)
  359. # Set them orbitting along the spin plane
  360. v_vec = np.cross(r_vec, [0, 0, 1])
  361. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)
  362. def addHalo(self, model, totalMass, particles, rs):
  363. &#34;&#34;&#34;Adds a halo to the galaxy.
  364. Parameters:
  365. model (string): parametrization of the halo.
  366. Only &#39;NFW&#39; is supported.
  367. totalMass (float): total mass of the halo
  368. particles (int): number of particles in the halo
  369. rs (float): characteristic length scale of the NFW profile.
  370. &#34;&#34;&#34;
  371. if particles == 0: return None
  372. # Divide the mass equally among all particles
  373. mass = np.ones(particles)*totalMass/particles
  374. self.mass = np.concatenate([self.mass, mass], axis=0)
  375. # Create particles according to the radial distribution from model
  376. if model == &#39;NFW&#39;:
  377. r = NFW.ppf(np.random.rand(particles), scale=rs)
  378. else: raise Exception(&#34;&#34;&#34;Bulge distribution not allowed.
  379. &#39;plummer&#39; and &#39;hernquist&#39; are the supported values&#34;&#34;&#34;)
  380. r_vec = r[:,np.newaxis] * random_unit_vectors(size=particles)
  381. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  382. # Orbit along random directions normal to the radial vector
  383. v_vec = np.cross(r_vec, random_unit_vectors(size=particles))
  384. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)
  385. # Label the particles
  386. type_ = [[&#39;dark&#39;], 0]*particles
  387. self.type = np.concatenate([self.type, type_], axis=0)
  388. def rotate(self, theta, phi):
  389. &#34;&#34;&#34;Rotates the galaxy so that its spin is along the (theta, phi)
  390. direction.
  391. Parameters:
  392. theta (float): polar angle.
  393. phi (float): azimuthal angle.
  394. &#34;&#34;&#34;
  395. M1 = np.array([[1, 0, 0],
  396. [0, np.cos(theta), np.sin(theta)],
  397. [0, -np.sin(theta), np.cos(theta)]])
  398. M2 = np.array([[np.cos(phi), np.sin(phi), 0],
  399. [-np.sin(phi), np.cos(phi), 0],
  400. [0, 0, 1]])
  401. M = np.matmul(M1, M2) # combine rotations
  402. self.r_vec = np.tensordot(self.r_vec, M, axes=[1, 0])
  403. self.v_vec = np.tensordot(self.v_vec, M, axes=[1, 0])</code></pre>
  404. </details>
  405. </section>
  406. <section>
  407. </section>
  408. <section>
  409. </section>
  410. <section>
  411. </section>
  412. <section>
  413. <h2 class="section-title" id="header-classes">Classes</h2>
  414. <dl>
  415. <dt id="simulation.Galaxy"><code class="flex name class">
  416. <span>class <span class="ident">Galaxy</span></span>
  417. </code></dt>
  418. <dd>
  419. <section class="desc"><p>"Helper class for creating galaxies.</p>
  420. <h2 id="attributes">Attributes</h2>
  421. <dl>
  422. <dt><strong><code>r_vec</code></strong> :&ensp;<code>array</code></dt>
  423. <dd>position of the particles in the current timestep.
  424. Shape: (number of particles, 3)</dd>
  425. <dt><strong><code>v_vec</code></strong> :&ensp;<code>array</code></dt>
  426. <dd>velocity in the current timestep.
  427. Shape: (number of particles, 3)</dd>
  428. <dt><strong><code>a_vec</code></strong> :&ensp;<code>array</code></dt>
  429. <dd>acceleration in the current timestep.
  430. Shape: (number of particles, 3)</dd>
  431. <dt><strong><code>mass</code></strong> :&ensp;<code>array</code></dt>
  432. <dd>mass of each particle in the simulation.
  433. Shape: (number of particles,)</dd>
  434. <dt><strong><code>type</code></strong> :&ensp;<code>array</code></dt>
  435. <dd>non-unique identifier for each particle.
  436. Shape: (number of particles,)</dd>
  437. </dl></section>
  438. <details class="source">
  439. <summary>Source code</summary>
  440. <pre><code class="python">class Galaxy():
  441. &#34;&#34;&#34;&#34;Helper class for creating galaxies.
  442. Attributes:
  443. r_vec (array): position of the particles in the current timestep.
  444. Shape: (number of particles, 3)
  445. v_vec (array): velocity in the current timestep.
  446. Shape: (number of particles, 3)
  447. a_vec (array): acceleration in the current timestep.
  448. Shape: (number of particles, 3)
  449. mass (array): mass of each particle in the simulation.
  450. Shape: (number of particles,)
  451. type (array): non-unique identifier for each particle.
  452. Shape: (number of particles,) &#34;&#34;&#34;
  453. def __init__(self, orientation, centralMass, bulge, disk, halo, sim):
  454. &#34;&#34;&#34;Constructor for the Galaxy class.
  455. Parameters:
  456. orientation (tupple): (inclination i, argument of pericenter w)
  457. centralMass (float): mass at the center of the galaxy
  458. bulge (dict): passed to the addBulge method.
  459. disk (dict): passed to the addDisk method.
  460. halo (dict): passed to the addHalo method.
  461. sim (Simulation): simulation object
  462. &#34;&#34;&#34;
  463. if sim.verbose: print(&#39;Initializing galaxy&#39;)
  464. # Build the central mass
  465. self.r_vec = np.zeros((1, 3))
  466. self.v_vec = np.zeros((1, 3))
  467. self.a_vec = np.zeros((1, 3))
  468. self.mass = np.array([centralMass])
  469. self.type = np.array([[&#39;center&#39;, 0]])
  470. # Build the other components
  471. self.addBulge(**bulge)
  472. if sim.verbose: print(&#39;Bulge created.&#39;)
  473. self.addDisk(**disk)
  474. if sim.verbose: print(&#39;Disk created.&#39;)
  475. self.addHalo(**halo)
  476. if sim.verbose: print(&#39;Halo created.&#39;)
  477. # Correct particles velocities to generate circular orbits
  478. # $a_\textup{centripetal} = v^2/r$
  479. a_vec = sim.acceleration(self.r_vec, self.mass, soft=sim.soft)
  480. a = np.linalg.norm(a_vec, axis=-1, keepdims=True)
  481. r = np.linalg.norm(self.r_vec, axis=-1, keepdims=True)
  482. v = np.linalg.norm(self.v_vec[1:], axis=-1, keepdims=True)
  483. direction_unit = self.v_vec[1:]/v
  484. # Set orbital velocities (except for central mass)
  485. self.v_vec[1:] = np.sqrt(a[1:] * r[1:]) * direction_unit
  486. self.a_vec = np.zeros_like(self.r_vec)
  487. # Rotate the galaxy into its correct orientation
  488. self.rotate(*(np.array(orientation)/360 * 2*np.pi))
  489. if sim.verbose: print(&#39;Galaxy set in rotation and oriented.&#39;)
  490. def addBulge(self, model, totalMass, particles, l):
  491. &#34;&#34;&#34;Adds a bulge to the galaxy.
  492. Parameters:
  493. model (string): parametrization of the bulge.
  494. &#39;plummer&#39; and &#39;hernquist&#39; are supported.
  495. totalMass (float): total mass of the bulge
  496. particles (int): number of particles in the bulge
  497. l (float): characteristic length scale of the model.
  498. &#34;&#34;&#34;
  499. if particles == 0: return None
  500. # Divide the mass equally among all particles
  501. mass = np.ones(particles) * totalMass/particles
  502. self.mass = np.concatenate([self.mass, mass], axis=0)
  503. # Create particles according to the radial distribution from model
  504. if model == &#39;plummer&#39;:
  505. r = PLUMMER.ppf(np.random.rand(particles), scale=l)
  506. elif model == &#39;hernquist&#39;:
  507. r = HERNQUIST.ppf(np.random.rand(particles), scale=l)
  508. else: raise Exception(&#34;&#34;&#34;Bulge distribution not allowed.
  509. &#39;plummer&#39; and &#39;hernquist&#39; are the supported values&#34;&#34;&#34;)
  510. r_vec = r[:,np.newaxis] * random_unit_vectors(size=particles)
  511. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  512. # Set them orbitting along random directions normal to r_vec
  513. v_vec = np.cross(r_vec, random_unit_vectors(size=particles))
  514. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)
  515. # Label the particles
  516. type_ = [[&#39;bulge&#39;, 0]]*particles
  517. self.type = np.concatenate([self.type, type_], axis=0)
  518. def addDisk(self, model, totalMass, particles, l):
  519. &#34;&#34;&#34;Adds a disk to the galaxy.
  520. Parameters:
  521. model (string): parametrization of the disk.
  522. &#39;rings&#39;, &#39;uniform&#39; and &#39;exp&#39; are supported.
  523. totalMass (float): total mass of the bulge
  524. particles (int): number of particles in the bulge
  525. l: fot &#39;exp&#39; and &#39;uniform&#39; characteristic length of the
  526. model. For &#39;rings&#39; tupple of the form (inner radius,
  527. outer radius, number of rings)
  528. &#34;&#34;&#34;
  529. if particles == 0: return None
  530. # Create particles according to the radial distribution from model
  531. if model == &#39;uniform&#39;:
  532. r = UNIFORM.ppf(np.random.rand(particles), scale=l)
  533. type_ = [[&#39;disk&#39;, 0]]*particles
  534. r_vec = r[:,np.newaxis] * random_unit_vectors(particles, &#39;2D&#39;)
  535. self.type = np.concatenate([self.type, type_], axis=0)
  536. elif model == &#39;rings&#39;:
  537. # l = [inner radius, outter radius, number of rings]
  538. distances = np.linspace(*l)
  539. # Aim for roughly constant areal density
  540. # Cascade rounding preserves the total number of particles
  541. perRing = cascade_round(particles * distances / np.sum(distances))
  542. particles = int(np.sum(perRing)) # prevents numerical errors
  543. r_vec = np.empty((0, 3))
  544. for d, n, i in zip(distances, perRing, range(l[2])):
  545. type_ = [[&#39;disk&#39;, i+1]]*int(n)
  546. self.type = np.concatenate([self.type, type_], axis=0)
  547. phi = np.linspace(0, 2 * np.pi, n, endpoint=False)
  548. ringr = d * np.array([[np.cos(i), np.sin(i), 0] for i in phi])
  549. r_vec = np.concatenate([r_vec, ringr], axis=0)
  550. elif model == &#39;exp&#39;:
  551. r = EXP.ppf(np.random.rand(particles), scale=l)
  552. r_vec = r[:,np.newaxis] * random_unit_vectors(particles, &#39;2D&#39;)
  553. type_ = [[&#39;disk&#39;, 0]]*particles
  554. self.type = np.concatenate([self.type, type_], axis=0)
  555. else:
  556. raise Exception(&#34;&#34;&#34;Disk distribution not allowed.
  557. &#39;uniform&#39;, &#39;rings&#39; and &#39;exp&#39; are the supported values&#34;&#34;&#34;)
  558. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  559. # Divide the mass equally among all particles
  560. mass = np.ones(particles) * totalMass/particles
  561. self.mass = np.concatenate([self.mass, mass], axis=0)
  562. # Set them orbitting along the spin plane
  563. v_vec = np.cross(r_vec, [0, 0, 1])
  564. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)
  565. def addHalo(self, model, totalMass, particles, rs):
  566. &#34;&#34;&#34;Adds a halo to the galaxy.
  567. Parameters:
  568. model (string): parametrization of the halo.
  569. Only &#39;NFW&#39; is supported.
  570. totalMass (float): total mass of the halo
  571. particles (int): number of particles in the halo
  572. rs (float): characteristic length scale of the NFW profile.
  573. &#34;&#34;&#34;
  574. if particles == 0: return None
  575. # Divide the mass equally among all particles
  576. mass = np.ones(particles)*totalMass/particles
  577. self.mass = np.concatenate([self.mass, mass], axis=0)
  578. # Create particles according to the radial distribution from model
  579. if model == &#39;NFW&#39;:
  580. r = NFW.ppf(np.random.rand(particles), scale=rs)
  581. else: raise Exception(&#34;&#34;&#34;Bulge distribution not allowed.
  582. &#39;plummer&#39; and &#39;hernquist&#39; are the supported values&#34;&#34;&#34;)
  583. r_vec = r[:,np.newaxis] * random_unit_vectors(size=particles)
  584. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  585. # Orbit along random directions normal to the radial vector
  586. v_vec = np.cross(r_vec, random_unit_vectors(size=particles))
  587. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)
  588. # Label the particles
  589. type_ = [[&#39;dark&#39;], 0]*particles
  590. self.type = np.concatenate([self.type, type_], axis=0)
  591. def rotate(self, theta, phi):
  592. &#34;&#34;&#34;Rotates the galaxy so that its spin is along the (theta, phi)
  593. direction.
  594. Parameters:
  595. theta (float): polar angle.
  596. phi (float): azimuthal angle.
  597. &#34;&#34;&#34;
  598. M1 = np.array([[1, 0, 0],
  599. [0, np.cos(theta), np.sin(theta)],
  600. [0, -np.sin(theta), np.cos(theta)]])
  601. M2 = np.array([[np.cos(phi), np.sin(phi), 0],
  602. [-np.sin(phi), np.cos(phi), 0],
  603. [0, 0, 1]])
  604. M = np.matmul(M1, M2) # combine rotations
  605. self.r_vec = np.tensordot(self.r_vec, M, axes=[1, 0])
  606. self.v_vec = np.tensordot(self.v_vec, M, axes=[1, 0])</code></pre>
  607. </details>
  608. <h3>Methods</h3>
  609. <dl>
  610. <dt id="simulation.Galaxy.__init__"><code class="name flex">
  611. <span>def <span class="ident">__init__</span></span>(<span>self, orientation, centralMass, bulge, disk, halo, sim)</span>
  612. </code></dt>
  613. <dd>
  614. <section class="desc"><p>Constructor for the Galaxy class.</p>
  615. <h2 id="parameters">Parameters</h2>
  616. <dl>
  617. <dt><strong><code>orientation</code></strong> :&ensp;<code>tupple</code></dt>
  618. <dd>(inclination i, argument of pericenter w)</dd>
  619. <dt><strong><code>centralMass</code></strong> :&ensp;<code>float</code></dt>
  620. <dd>mass at the center of the galaxy</dd>
  621. <dt><strong><code>bulge</code></strong> :&ensp;<code>dict</code></dt>
  622. <dd>passed to the addBulge method.</dd>
  623. <dt><strong><code>disk</code></strong> :&ensp;<code>dict</code></dt>
  624. <dd>passed to the addDisk method.</dd>
  625. <dt><strong><code>halo</code></strong> :&ensp;<code>dict</code></dt>
  626. <dd>passed to the addHalo method.</dd>
  627. <dt><strong><code>sim</code></strong> :&ensp;<a title="simulation.Simulation" href="#simulation.Simulation"><code>Simulation</code></a></dt>
  628. <dd>simulation object</dd>
  629. </dl></section>
  630. <details class="source">
  631. <summary>Source code</summary>
  632. <pre><code class="python">def __init__(self, orientation, centralMass, bulge, disk, halo, sim):
  633. &#34;&#34;&#34;Constructor for the Galaxy class.
  634. Parameters:
  635. orientation (tupple): (inclination i, argument of pericenter w)
  636. centralMass (float): mass at the center of the galaxy
  637. bulge (dict): passed to the addBulge method.
  638. disk (dict): passed to the addDisk method.
  639. halo (dict): passed to the addHalo method.
  640. sim (Simulation): simulation object
  641. &#34;&#34;&#34;
  642. if sim.verbose: print(&#39;Initializing galaxy&#39;)
  643. # Build the central mass
  644. self.r_vec = np.zeros((1, 3))
  645. self.v_vec = np.zeros((1, 3))
  646. self.a_vec = np.zeros((1, 3))
  647. self.mass = np.array([centralMass])
  648. self.type = np.array([[&#39;center&#39;, 0]])
  649. # Build the other components
  650. self.addBulge(**bulge)
  651. if sim.verbose: print(&#39;Bulge created.&#39;)
  652. self.addDisk(**disk)
  653. if sim.verbose: print(&#39;Disk created.&#39;)
  654. self.addHalo(**halo)
  655. if sim.verbose: print(&#39;Halo created.&#39;)
  656. # Correct particles velocities to generate circular orbits
  657. # $a_\textup{centripetal} = v^2/r$
  658. a_vec = sim.acceleration(self.r_vec, self.mass, soft=sim.soft)
  659. a = np.linalg.norm(a_vec, axis=-1, keepdims=True)
  660. r = np.linalg.norm(self.r_vec, axis=-1, keepdims=True)
  661. v = np.linalg.norm(self.v_vec[1:], axis=-1, keepdims=True)
  662. direction_unit = self.v_vec[1:]/v
  663. # Set orbital velocities (except for central mass)
  664. self.v_vec[1:] = np.sqrt(a[1:] * r[1:]) * direction_unit
  665. self.a_vec = np.zeros_like(self.r_vec)
  666. # Rotate the galaxy into its correct orientation
  667. self.rotate(*(np.array(orientation)/360 * 2*np.pi))
  668. if sim.verbose: print(&#39;Galaxy set in rotation and oriented.&#39;)</code></pre>
  669. </details>
  670. </dd>
  671. <dt id="simulation.Galaxy.addBulge"><code class="name flex">
  672. <span>def <span class="ident">addBulge</span></span>(<span>self, model, totalMass, particles, l)</span>
  673. </code></dt>
  674. <dd>
  675. <section class="desc"><p>Adds a bulge to the galaxy.</p>
  676. <h2 id="parameters">Parameters</h2>
  677. <dl>
  678. <dt><strong><code>model</code></strong> :&ensp;<code>string</code></dt>
  679. <dd>parametrization of the bulge.
  680. 'plummer' and 'hernquist' are supported.</dd>
  681. <dt><strong><code>totalMass</code></strong> :&ensp;<code>float</code></dt>
  682. <dd>total mass of the bulge</dd>
  683. <dt><strong><code>particles</code></strong> :&ensp;<code>int</code></dt>
  684. <dd>number of particles in the bulge</dd>
  685. <dt><strong><code>l</code></strong> :&ensp;<code>float</code></dt>
  686. <dd>characteristic length scale of the model.</dd>
  687. </dl></section>
  688. <details class="source">
  689. <summary>Source code</summary>
  690. <pre><code class="python">def addBulge(self, model, totalMass, particles, l):
  691. &#34;&#34;&#34;Adds a bulge to the galaxy.
  692. Parameters:
  693. model (string): parametrization of the bulge.
  694. &#39;plummer&#39; and &#39;hernquist&#39; are supported.
  695. totalMass (float): total mass of the bulge
  696. particles (int): number of particles in the bulge
  697. l (float): characteristic length scale of the model.
  698. &#34;&#34;&#34;
  699. if particles == 0: return None
  700. # Divide the mass equally among all particles
  701. mass = np.ones(particles) * totalMass/particles
  702. self.mass = np.concatenate([self.mass, mass], axis=0)
  703. # Create particles according to the radial distribution from model
  704. if model == &#39;plummer&#39;:
  705. r = PLUMMER.ppf(np.random.rand(particles), scale=l)
  706. elif model == &#39;hernquist&#39;:
  707. r = HERNQUIST.ppf(np.random.rand(particles), scale=l)
  708. else: raise Exception(&#34;&#34;&#34;Bulge distribution not allowed.
  709. &#39;plummer&#39; and &#39;hernquist&#39; are the supported values&#34;&#34;&#34;)
  710. r_vec = r[:,np.newaxis] * random_unit_vectors(size=particles)
  711. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  712. # Set them orbitting along random directions normal to r_vec
  713. v_vec = np.cross(r_vec, random_unit_vectors(size=particles))
  714. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)
  715. # Label the particles
  716. type_ = [[&#39;bulge&#39;, 0]]*particles
  717. self.type = np.concatenate([self.type, type_], axis=0)</code></pre>
  718. </details>
  719. </dd>
  720. <dt id="simulation.Galaxy.addDisk"><code class="name flex">
  721. <span>def <span class="ident">addDisk</span></span>(<span>self, model, totalMass, particles, l)</span>
  722. </code></dt>
  723. <dd>
  724. <section class="desc"><p>Adds a disk to the galaxy.</p>
  725. <h2 id="parameters">Parameters</h2>
  726. <dl>
  727. <dt><strong><code>model</code></strong> :&ensp;<code>string</code></dt>
  728. <dd>parametrization of the disk.
  729. 'rings', 'uniform' and 'exp' are supported.</dd>
  730. <dt><strong><code>totalMass</code></strong> :&ensp;<code>float</code></dt>
  731. <dd>total mass of the bulge</dd>
  732. <dt><strong><code>particles</code></strong> :&ensp;<code>int</code></dt>
  733. <dd>number of particles in the bulge</dd>
  734. <dt><strong><code>l</code></strong></dt>
  735. <dd>fot 'exp' and 'uniform' characteristic length of the
  736. model. For 'rings' tupple of the form (inner radius,
  737. outer radius, number of rings)</dd>
  738. </dl></section>
  739. <details class="source">
  740. <summary>Source code</summary>
  741. <pre><code class="python">def addDisk(self, model, totalMass, particles, l):
  742. &#34;&#34;&#34;Adds a disk to the galaxy.
  743. Parameters:
  744. model (string): parametrization of the disk.
  745. &#39;rings&#39;, &#39;uniform&#39; and &#39;exp&#39; are supported.
  746. totalMass (float): total mass of the bulge
  747. particles (int): number of particles in the bulge
  748. l: fot &#39;exp&#39; and &#39;uniform&#39; characteristic length of the
  749. model. For &#39;rings&#39; tupple of the form (inner radius,
  750. outer radius, number of rings)
  751. &#34;&#34;&#34;
  752. if particles == 0: return None
  753. # Create particles according to the radial distribution from model
  754. if model == &#39;uniform&#39;:
  755. r = UNIFORM.ppf(np.random.rand(particles), scale=l)
  756. type_ = [[&#39;disk&#39;, 0]]*particles
  757. r_vec = r[:,np.newaxis] * random_unit_vectors(particles, &#39;2D&#39;)
  758. self.type = np.concatenate([self.type, type_], axis=0)
  759. elif model == &#39;rings&#39;:
  760. # l = [inner radius, outter radius, number of rings]
  761. distances = np.linspace(*l)
  762. # Aim for roughly constant areal density
  763. # Cascade rounding preserves the total number of particles
  764. perRing = cascade_round(particles * distances / np.sum(distances))
  765. particles = int(np.sum(perRing)) # prevents numerical errors
  766. r_vec = np.empty((0, 3))
  767. for d, n, i in zip(distances, perRing, range(l[2])):
  768. type_ = [[&#39;disk&#39;, i+1]]*int(n)
  769. self.type = np.concatenate([self.type, type_], axis=0)
  770. phi = np.linspace(0, 2 * np.pi, n, endpoint=False)
  771. ringr = d * np.array([[np.cos(i), np.sin(i), 0] for i in phi])
  772. r_vec = np.concatenate([r_vec, ringr], axis=0)
  773. elif model == &#39;exp&#39;:
  774. r = EXP.ppf(np.random.rand(particles), scale=l)
  775. r_vec = r[:,np.newaxis] * random_unit_vectors(particles, &#39;2D&#39;)
  776. type_ = [[&#39;disk&#39;, 0]]*particles
  777. self.type = np.concatenate([self.type, type_], axis=0)
  778. else:
  779. raise Exception(&#34;&#34;&#34;Disk distribution not allowed.
  780. &#39;uniform&#39;, &#39;rings&#39; and &#39;exp&#39; are the supported values&#34;&#34;&#34;)
  781. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  782. # Divide the mass equally among all particles
  783. mass = np.ones(particles) * totalMass/particles
  784. self.mass = np.concatenate([self.mass, mass], axis=0)
  785. # Set them orbitting along the spin plane
  786. v_vec = np.cross(r_vec, [0, 0, 1])
  787. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)</code></pre>
  788. </details>
  789. </dd>
  790. <dt id="simulation.Galaxy.addHalo"><code class="name flex">
  791. <span>def <span class="ident">addHalo</span></span>(<span>self, model, totalMass, particles, rs)</span>
  792. </code></dt>
  793. <dd>
  794. <section class="desc"><p>Adds a halo to the galaxy.</p>
  795. <h2 id="parameters">Parameters</h2>
  796. <dl>
  797. <dt><strong><code>model</code></strong> :&ensp;<code>string</code></dt>
  798. <dd>parametrization of the halo.
  799. Only 'NFW' is supported.</dd>
  800. <dt><strong><code>totalMass</code></strong> :&ensp;<code>float</code></dt>
  801. <dd>total mass of the halo</dd>
  802. <dt><strong><code>particles</code></strong> :&ensp;<code>int</code></dt>
  803. <dd>number of particles in the halo</dd>
  804. <dt><strong><code>rs</code></strong> :&ensp;<code>float</code></dt>
  805. <dd>characteristic length scale of the NFW profile.</dd>
  806. </dl></section>
  807. <details class="source">
  808. <summary>Source code</summary>
  809. <pre><code class="python">def addHalo(self, model, totalMass, particles, rs):
  810. &#34;&#34;&#34;Adds a halo to the galaxy.
  811. Parameters:
  812. model (string): parametrization of the halo.
  813. Only &#39;NFW&#39; is supported.
  814. totalMass (float): total mass of the halo
  815. particles (int): number of particles in the halo
  816. rs (float): characteristic length scale of the NFW profile.
  817. &#34;&#34;&#34;
  818. if particles == 0: return None
  819. # Divide the mass equally among all particles
  820. mass = np.ones(particles)*totalMass/particles
  821. self.mass = np.concatenate([self.mass, mass], axis=0)
  822. # Create particles according to the radial distribution from model
  823. if model == &#39;NFW&#39;:
  824. r = NFW.ppf(np.random.rand(particles), scale=rs)
  825. else: raise Exception(&#34;&#34;&#34;Bulge distribution not allowed.
  826. &#39;plummer&#39; and &#39;hernquist&#39; are the supported values&#34;&#34;&#34;)
  827. r_vec = r[:,np.newaxis] * random_unit_vectors(size=particles)
  828. self.r_vec = np.concatenate([self.r_vec, r_vec], axis=0)
  829. # Orbit along random directions normal to the radial vector
  830. v_vec = np.cross(r_vec, random_unit_vectors(size=particles))
  831. self.v_vec = np.concatenate([self.v_vec, v_vec], axis=0)
  832. # Label the particles
  833. type_ = [[&#39;dark&#39;], 0]*particles
  834. self.type = np.concatenate([self.type, type_], axis=0)</code></pre>
  835. </details>
  836. </dd>
  837. <dt id="simulation.Galaxy.rotate"><code class="name flex">
  838. <span>def <span class="ident">rotate</span></span>(<span>self, theta, phi)</span>
  839. </code></dt>
  840. <dd>
  841. <section class="desc"><p>Rotates the galaxy so that its spin is along the (theta, phi)
  842. direction.</p>
  843. <h2 id="parameters">Parameters</h2>
  844. <dl>
  845. <dt><strong><code>theta</code></strong> :&ensp;<code>float</code></dt>
  846. <dd>polar angle.</dd>
  847. <dt><strong><code>phi</code></strong> :&ensp;<code>float</code></dt>
  848. <dd>azimuthal angle.</dd>
  849. </dl></section>
  850. <details class="source">
  851. <summary>Source code</summary>
  852. <pre><code class="python">def rotate(self, theta, phi):
  853. &#34;&#34;&#34;Rotates the galaxy so that its spin is along the (theta, phi)
  854. direction.
  855. Parameters:
  856. theta (float): polar angle.
  857. phi (float): azimuthal angle.
  858. &#34;&#34;&#34;
  859. M1 = np.array([[1, 0, 0],
  860. [0, np.cos(theta), np.sin(theta)],
  861. [0, -np.sin(theta), np.cos(theta)]])
  862. M2 = np.array([[np.cos(phi), np.sin(phi), 0],
  863. [-np.sin(phi), np.cos(phi), 0],
  864. [0, 0, 1]])
  865. M = np.matmul(M1, M2) # combine rotations
  866. self.r_vec = np.tensordot(self.r_vec, M, axes=[1, 0])
  867. self.v_vec = np.tensordot(self.v_vec, M, axes=[1, 0])</code></pre>
  868. </details>
  869. </dd>
  870. </dl>
  871. </dd>
  872. <dt id="simulation.Simulation"><code class="flex name class">
  873. <span>class <span class="ident">Simulation</span></span>
  874. </code></dt>
  875. <dd>
  876. <section class="desc"><p>"Main class for the gravitational simulation.</p>
  877. <h2 id="attributes">Attributes</h2>
  878. <dl>
  879. <dt><strong><code>r_vec</code></strong> :&ensp;<code>array</code></dt>
  880. <dd>position of the particles in the current timestep.
  881. Shape: (number of particles, 3)</dd>
  882. <dt><strong><code>rprev_vec</code></strong> :&ensp;<code>array</code></dt>
  883. <dd>position of the particles in the previous timestep.
  884. Shape: (number of particles, 3)</dd>
  885. <dt><strong><code>v_vec</code></strong> :&ensp;<code>array</code></dt>
  886. <dd>velocity in the current timestep.
  887. Shape: (number of particles, 3)</dd>
  888. <dt><strong><code>a_vec</code></strong> :&ensp;<code>array</code></dt>
  889. <dd>acceleration in the current timestep.
  890. Shape: (number of particles, 3)</dd>
  891. <dt><strong><code>mass</code></strong> :&ensp;<code>array</code></dt>
  892. <dd>mass of each particle in the simulation.
  893. Shape: (number of particles,)</dd>
  894. <dt><strong><code>type</code></strong> :&ensp;<code>array</code></dt>
  895. <dd>non-unique identifier for each particle.
  896. Shape: (number of particles,)</dd>
  897. <dt><strong><code>tracks</code></strong> :&ensp;<code>array</code></dt>
  898. <dd>list of positions through the simulation for central
  899. masses. Shape: (tracked particles, n+1, 3).</dd>
  900. <dt><strong><code>CONFIG</code></strong> :&ensp;<code>array</code></dt>
  901. <dd>configuration used to create the simulation.
  902. It will be saved along the state of the simulation.</dd>
  903. <dt><strong><code>dt</code></strong> :&ensp;<code>float</code></dt>
  904. <dd>timestep of the simulation</dd>
  905. <dt><strong><code>n</code></strong> :&ensp;<code>int</code></dt>
  906. <dd>current timestep. Initialized as n=0.</dd>
  907. <dt><strong><code>soft</code></strong> :&ensp;<code>float</code></dt>
  908. <dd>softening length used by the simulation.</dd>
  909. <dt><strong><code>verbose</code></strong> :&ensp;<code>boolean</code></dt>
  910. <dd>When True progress statements will be printed.</dd>
  911. </dl></section>
  912. <details class="source">
  913. <summary>Source code</summary>
  914. <pre><code class="python">class Simulation:
  915. &#34;&#34;&#34;&#34;Main class for the gravitational simulation.
  916. Attributes:
  917. r_vec (array): position of the particles in the current timestep.
  918. Shape: (number of particles, 3)
  919. rprev_vec (array): position of the particles in the previous timestep.
  920. Shape: (number of particles, 3)
  921. v_vec (array): velocity in the current timestep.
  922. Shape: (number of particles, 3)
  923. a_vec (array): acceleration in the current timestep.
  924. Shape: (number of particles, 3)
  925. mass (array): mass of each particle in the simulation.
  926. Shape: (number of particles,)
  927. type (array): non-unique identifier for each particle.
  928. Shape: (number of particles,)
  929. tracks (array): list of positions through the simulation for central
  930. masses. Shape: (tracked particles, n+1, 3).
  931. CONFIG (array): configuration used to create the simulation.
  932. It will be saved along the state of the simulation.
  933. dt (float): timestep of the simulation
  934. n (int): current timestep. Initialized as n=0.
  935. soft (float): softening length used by the simulation.
  936. verbose (boolean): When True progress statements will be printed.
  937. &#34;&#34;&#34;
  938. def __init__(self, dt, soft, verbose, CONFIG, method, **kwargs):
  939. &#34;&#34;&#34;Constructor for the Simulation class.
  940. Arguments:
  941. dt (float): timestep of the simulation
  942. n (int): current timestep. Initialized as n=0.
  943. soft (float): softening length used by the simulation.
  944. verbose (bool): When True progress statements will be printed.
  945. CONFIG (dict): configuration file used to create the simulation.
  946. method (string): Optional. Algorithm to use when computing the
  947. gravitational forces. One of &#39;bruteForce&#39;, &#39;bruteForce_numba&#39;,
  948. &#39;bruteForce_numbaopt&#39;, &#39;bruteForce_CPP&#39;, &#39;barnesHut_CPP&#39;.
  949. &#34;&#34;&#34;
  950. self.n = 0
  951. self.t = 0
  952. self.dt = dt
  953. self.soft = soft
  954. self.verbose = verbose
  955. self.CONFIG = CONFIG
  956. # Initialize empty arrays for all necessary properties
  957. self.r_vec = np.empty((0, 3))
  958. self.v_vec = np.empty((0, 3))
  959. self.a_vec = np.empty((0, 3))
  960. self.mass = np.empty((0,))
  961. self.type = np.empty((0, 2))
  962. algorithms = {
  963. &#39;bruteForce&#39;: acceleration.bruteForce,
  964. &#39;bruteForceNumba&#39;: acceleration.bruteForceNumba,
  965. &#39;bruteForceNumbaOptimized&#39;: acceleration.bruteForceNumbaOptimized,
  966. &#39;bruteForceCPP&#39;: acceleration.bruteForceCPP,
  967. &#39;barnesHutCPP&#39;: acceleration.barnesHutCPP
  968. }
  969. try:
  970. self.acceleration = algorithms[method]
  971. except: raise Exception(&#34;Method &#39;{}&#39; unknown&#34;.format(method))
  972. def add(self, body):
  973. &#34;&#34;&#34;Add a body to the simulation. It must expose the public attributes
  974. body.r_vec, body.v_vec, body.a_vec, body.type, body.mass.
  975. Arguments:
  976. body: Object to be added to the simulation (e.g. a Galaxy object)
  977. &#34;&#34;&#34;
  978. # Extend all relevant attributes by concatenating the body
  979. for name in [&#39;r_vec&#39;, &#39;v_vec&#39;, &#39;a_vec&#39;, &#39;type&#39;, &#39;mass&#39;]:
  980. simattr, bodyattr = getattr(self, name), getattr(body, name)
  981. setattr(self, name, np.concatenate([simattr, bodyattr], axis=0))
  982. # Order based on mass
  983. order = np.argsort(-self.mass)
  984. for name in [&#39;r_vec&#39;, &#39;v_vec&#39;, &#39;a_vec&#39;, &#39;type&#39;, &#39;mass&#39;]:
  985. setattr(self, name, getattr(self, name)[order])
  986. # Update the list of objects to keep track of
  987. self.tracks = np.empty((np.sum(self.type[:,0]==&#39;center&#39;), 0, 3))
  988. def step(self):
  989. &#34;&#34;&#34;Perform a single step of the simulation.
  990. Makes use of a 4th order Verlet integrator.
  991. &#34;&#34;&#34;
  992. # Calculate the acceleration
  993. self.a_vec = self.acceleration(self.r_vec, self.mass, soft=self.soft)
  994. # Update the state using the Verlet algorithm
  995. # (A custom algorithm is written mainly for learning purposes)
  996. self.r_vec, self.rprev_vec = (2*self.r_vec - self.rprev_vec
  997. + self.a_vec * self.dt**2, self.r_vec)
  998. self.n += 1
  999. # Update tracks
  1000. self.tracks = np.concatenate([self.tracks,
  1001. self.r_vec[self.type[:,0]==&#39;center&#39;][:,np.newaxis]], axis=1)
  1002. def run(self, tmax, saveEvery, outputFolder, **kwargs):
  1003. &#34;&#34;&#34;Run the galactic simulation.
  1004. Attributes:
  1005. tmax (float): Time to which the simulation will run to.
  1006. This is measured here since the start of the simulation,
  1007. not since pericenter.
  1008. saveEvery (int): The state is saved every saveEvery steps.
  1009. outputFolder (string): It will be saved to /data/outputFolder/
  1010. &#34;&#34;&#34;
  1011. # When the simulation starts, intialize self.rprev_vec
  1012. self.rprev_vec = self.r_vec - self.v_vec * self.dt
  1013. if self.verbose: print(&#39;Simulation starting. Bon voyage!&#39;)
  1014. while(self.t &lt; tmax):
  1015. self.step()
  1016. if(self.n % saveEvery == 0):
  1017. self.save(&#39;data/{}&#39;.format(outputFolder))
  1018. print(&#39;Simulation complete.&#39;)
  1019. def save(self, outputFolder):
  1020. &#34;&#34;&#34;Save the state of the simulation to the outputFolder.
  1021. Two files are saved:
  1022. sim{self.n}.pickle: serializing the state.
  1023. sim{self.n}.png: a simplified 2D plot of x, y.
  1024. &#34;&#34;&#34;
  1025. # Create the output folder if it doesn&#39;t exist
  1026. if not os.path.exists(outputFolder): os.makedirs(outputFolder)
  1027. # Compute some useful quantities
  1028. # v_vec is not required by the integrator, but useful
  1029. self.v_vec = (self.r_vec - self.rprev_vec)/self.dt
  1030. self.t = self.n * self.dt # prevents numerical rounding errors
  1031. # Serialize state
  1032. file = open(outputFolder+&#39;/data{}.pickle&#39;.format(self.n), &#34;wb&#34;)
  1033. pickle.dump({&#39;r_vec&#39;: self.r_vec, &#39;v_vec&#39;: self.v_vec,
  1034. &#39;type&#39;: self.type, &#39;mass&#39;: self.mass,
  1035. &#39;CONFIG&#39;: self.CONFIG, &#39;t&#39;: self.t,
  1036. &#39;tracks&#39;: self.tracks}, file)
  1037. # Save simplified plot of the current state.
  1038. # Its main use is to detect ill-behaved situations early on.
  1039. fig = plt.figure()
  1040. plt.xlim(-5, 5); plt.ylim(-5, 5); plt.axis(&#39;equal&#39;)
  1041. # Dark halo is plotted in red, disk in blue, bulge in green
  1042. PLTCON = [(&#39;dark&#39;, &#39;r&#39;, 0.3), (&#39;disk&#39;, &#39;b&#39;, 1.0), (&#39;bulge&#39;, &#39;g&#39;, 0.5)]
  1043. for type_, c, a in PLTCON:
  1044. plt.scatter(self.r_vec[self.type[:,0]==type_][:,0],
  1045. self.r_vec[self.type[:,0]==type_][:,1], s=0.1, c=c, alpha=a)
  1046. # Central mass as a magenta star
  1047. plt.scatter(self.r_vec[self.type[:,0]==&#39;center&#39;][:,0],
  1048. self.r_vec[self.type[:,0]==&#39;center&#39;][:,1], s=100, marker=&#34;*&#34;, c=&#39;m&#39;)
  1049. # Save to png file
  1050. fig.savefig(outputFolder+&#39;/sim{}.png&#39;.format(self.n), dpi=150)
  1051. plt.close(fig)
  1052. def project(self, theta, phi, view=0):
  1053. &#34;&#34;&#34;Projects the 3D simulation onto a plane as viewed from the
  1054. direction described by the (theta, phi, view). Angles in radians.
  1055. (This is used by the simulated annealing algorithm)
  1056. Parameters:
  1057. theta (float): polar angle.
  1058. phi (float): azimuthal angle.
  1059. view (float): rotation along line of sight.
  1060. &#34;&#34;&#34;
  1061. M1 = np.array([[np.cos(phi), np.sin(phi), 0],
  1062. [-np.sin(phi), np.cos(phi), 0],
  1063. [0, 0, 1]])
  1064. M2 = np.array([[1, 0, 0],
  1065. [0, np.cos(theta), np.sin(theta)],
  1066. [0, -np.sin(theta), np.cos(theta)]])
  1067. M3 = np.array([[np.cos(view), np.sin(view), 0],
  1068. [-np.sin(view), np.cos(view), 0],
  1069. [0, 0, 1]])
  1070. M = np.matmul(M1, np.matmul(M2, M3)) # combine rotations
  1071. r = np.tensordot(self.r_vec, M, axes=[1, 0])
  1072. return r[:,0:2]
  1073. def setOrbit(self, galaxy1, galaxy2, e, rmin, R0):
  1074. &#34;&#34;&#34;Sets the two galaxies galaxy1, galaxy2 in an orbit.
  1075. Parameters:
  1076. galaxy1 (Galaxy): 1st galaxy (main)
  1077. galaxy2 (Galaxy): 2nd galaxy (companion)
  1078. e: eccentricity of the orbit
  1079. rmin: distance of closest approach
  1080. R0: initial separation
  1081. &#34;&#34;&#34;
  1082. m1, m2 = np.sum(galaxy1.mass), np.sum(galaxy2.mass)
  1083. # Relevant formulae:
  1084. # $r_0 = r (1 + e) \cos(\phi)$, where $r_0$ ($\neq R_0$) is the semi-latus rectum
  1085. # $r_0 = r_\textup{min} (1 + e)$
  1086. # $v^2 = G M (2/r - 1/a)$, where a is the semimajor axis
  1087. # Solve the reduced two-body problem with reduced mass $\mu$ (mu)
  1088. M = m1 + m2
  1089. r0 = rmin * (1 + e)
  1090. try:
  1091. phi = np.arccos((r0/R0 - 1) / e) # inverting the orbit equation
  1092. phi = -np.abs(phi) # Choose negative (incoming) angle
  1093. ainv = (1 - e) / rmin # ainv = $1/a$, as a may be infinite
  1094. v = np.sqrt(M * (2/R0 - ainv))
  1095. # Finally, calculate the angle of motion. angle = tan(dy/dx)
  1096. # $dy/dx = ((dr/d\phi) sin(\phi) + r \cos(\phi))/((dr/d\phi) cos(\phi) - r \sin(\phi))$
  1097. dy = R0/r0 * e * np.sin(phi)**2 + np.cos(phi)
  1098. dx = R0/r0 * e * np.sin(phi) * np.cos(phi) - np.sin(phi)
  1099. vangle = np.arctan2(dy, dx)
  1100. except: raise Exception(&#39;&#39;&#39;The orbital parameters cannot be satisfied.
  1101. For elliptical orbits check that R0 is posible (&lt;rmax).&#39;&#39;&#39;)
  1102. # We now need the actual motion of each body
  1103. R_vec = np.array([[R0*np.cos(phi), R0*np.sin(phi), 0.]])
  1104. V_vec = np.array([[v*np.cos(vangle), v*np.sin(vangle), 0.]])
  1105. galaxy1.r_vec += m2/M * R_vec
  1106. galaxy1.v_vec += m2/M * V_vec
  1107. galaxy2.r_vec += -m1/M * R_vec
  1108. galaxy2.v_vec += -m1/M * V_vec
  1109. # Explicitely add the galaxies to the simulation
  1110. self.add(galaxy1)
  1111. self.add(galaxy2)
  1112. if self.verbose: print(&#39;Galaxies set in orbit.&#39;)</code></pre>
  1113. </details>
  1114. <h3>Methods</h3>
  1115. <dl>
  1116. <dt id="simulation.Simulation.__init__"><code class="name flex">
  1117. <span>def <span class="ident">__init__</span></span>(<span>self, dt, soft, verbose, CONFIG, method, **kwargs)</span>
  1118. </code></dt>
  1119. <dd>
  1120. <section class="desc"><p>Constructor for the Simulation class.</p>
  1121. <h2 id="arguments">Arguments</h2>
  1122. <dl>
  1123. <dt><strong><code>dt</code></strong> :&ensp;<code>float</code></dt>
  1124. <dd>timestep of the simulation</dd>
  1125. <dt><strong><code>n</code></strong> :&ensp;<code>int</code></dt>
  1126. <dd>current timestep. Initialized as n=0.</dd>
  1127. <dt><strong><code>soft</code></strong> :&ensp;<code>float</code></dt>
  1128. <dd>softening length used by the simulation.</dd>
  1129. <dt><strong><code>verbose</code></strong> :&ensp;<code>bool</code></dt>
  1130. <dd>When True progress statements will be printed.</dd>
  1131. <dt><strong><code>CONFIG</code></strong> :&ensp;<code>dict</code></dt>
  1132. <dd>configuration file used to create the simulation.</dd>
  1133. <dt><strong><code>method</code></strong> :&ensp;<code>string</code></dt>
  1134. <dd>Optional. Algorithm to use when computing the
  1135. gravitational forces. One of 'bruteForce', 'bruteForce_numba',
  1136. 'bruteForce_numbaopt', 'bruteForce_CPP', 'barnesHut_CPP'.</dd>
  1137. </dl></section>
  1138. <details class="source">
  1139. <summary>Source code</summary>
  1140. <pre><code class="python">def __init__(self, dt, soft, verbose, CONFIG, method, **kwargs):
  1141. &#34;&#34;&#34;Constructor for the Simulation class.
  1142. Arguments:
  1143. dt (float): timestep of the simulation
  1144. n (int): current timestep. Initialized as n=0.
  1145. soft (float): softening length used by the simulation.
  1146. verbose (bool): When True progress statements will be printed.
  1147. CONFIG (dict): configuration file used to create the simulation.
  1148. method (string): Optional. Algorithm to use when computing the
  1149. gravitational forces. One of &#39;bruteForce&#39;, &#39;bruteForce_numba&#39;,
  1150. &#39;bruteForce_numbaopt&#39;, &#39;bruteForce_CPP&#39;, &#39;barnesHut_CPP&#39;.
  1151. &#34;&#34;&#34;
  1152. self.n = 0
  1153. self.t = 0
  1154. self.dt = dt
  1155. self.soft = soft
  1156. self.verbose = verbose
  1157. self.CONFIG = CONFIG
  1158. # Initialize empty arrays for all necessary properties
  1159. self.r_vec = np.empty((0, 3))
  1160. self.v_vec = np.empty((0, 3))
  1161. self.a_vec = np.empty((0, 3))
  1162. self.mass = np.empty((0,))
  1163. self.type = np.empty((0, 2))
  1164. algorithms = {
  1165. &#39;bruteForce&#39;: acceleration.bruteForce,
  1166. &#39;bruteForceNumba&#39;: acceleration.bruteForceNumba,
  1167. &#39;bruteForceNumbaOptimized&#39;: acceleration.bruteForceNumbaOptimized,
  1168. &#39;bruteForceCPP&#39;: acceleration.bruteForceCPP,
  1169. &#39;barnesHutCPP&#39;: acceleration.barnesHutCPP
  1170. }
  1171. try:
  1172. self.acceleration = algorithms[method]
  1173. except: raise Exception(&#34;Method &#39;{}&#39; unknown&#34;.format(method))</code></pre>
  1174. </details>
  1175. </dd>
  1176. <dt id="simulation.Simulation.add"><code class="name flex">
  1177. <span>def <span class="ident">add</span></span>(<span>self, body)</span>
  1178. </code></dt>
  1179. <dd>
  1180. <section class="desc"><p>Add a body to the simulation. It must expose the public attributes
  1181. body.r_vec, body.v_vec, body.a_vec, body.type, body.mass.</p>
  1182. <h2 id="arguments">Arguments</h2>
  1183. <dl>
  1184. <dt><strong><code>body</code></strong></dt>
  1185. <dd>Object to be added to the simulation (e.g. a Galaxy object)</dd>
  1186. </dl></section>
  1187. <details class="source">
  1188. <summary>Source code</summary>
  1189. <pre><code class="python">def add(self, body):
  1190. &#34;&#34;&#34;Add a body to the simulation. It must expose the public attributes
  1191. body.r_vec, body.v_vec, body.a_vec, body.type, body.mass.
  1192. Arguments:
  1193. body: Object to be added to the simulation (e.g. a Galaxy object)
  1194. &#34;&#34;&#34;
  1195. # Extend all relevant attributes by concatenating the body
  1196. for name in [&#39;r_vec&#39;, &#39;v_vec&#39;, &#39;a_vec&#39;, &#39;type&#39;, &#39;mass&#39;]:
  1197. simattr, bodyattr = getattr(self, name), getattr(body, name)
  1198. setattr(self, name, np.concatenate([simattr, bodyattr], axis=0))
  1199. # Order based on mass
  1200. order = np.argsort(-self.mass)
  1201. for name in [&#39;r_vec&#39;, &#39;v_vec&#39;, &#39;a_vec&#39;, &#39;type&#39;, &#39;mass&#39;]:
  1202. setattr(self, name, getattr(self, name)[order])
  1203. # Update the list of objects to keep track of
  1204. self.tracks = np.empty((np.sum(self.type[:,0]==&#39;center&#39;), 0, 3))</code></pre>
  1205. </details>
  1206. </dd>
  1207. <dt id="simulation.Simulation.project"><code class="name flex">
  1208. <span>def <span class="ident">project</span></span>(<span>self, theta, phi, view=0)</span>
  1209. </code></dt>
  1210. <dd>
  1211. <section class="desc"><p>Projects the 3D simulation onto a plane as viewed from the
  1212. direction described by the (theta, phi, view). Angles in radians.
  1213. (This is used by the simulated annealing algorithm)</p>
  1214. <h2 id="parameters">Parameters</h2>
  1215. <dl>
  1216. <dt><strong><code>theta</code></strong> :&ensp;<code>float</code></dt>
  1217. <dd>polar angle.</dd>
  1218. <dt><strong><code>phi</code></strong> :&ensp;<code>float</code></dt>
  1219. <dd>azimuthal angle.</dd>
  1220. <dt><strong><code>view</code></strong> :&ensp;<code>float</code></dt>
  1221. <dd>rotation along line of sight.</dd>
  1222. </dl></section>
  1223. <details class="source">
  1224. <summary>Source code</summary>
  1225. <pre><code class="python">def project(self, theta, phi, view=0):
  1226. &#34;&#34;&#34;Projects the 3D simulation onto a plane as viewed from the
  1227. direction described by the (theta, phi, view). Angles in radians.
  1228. (This is used by the simulated annealing algorithm)
  1229. Parameters:
  1230. theta (float): polar angle.
  1231. phi (float): azimuthal angle.
  1232. view (float): rotation along line of sight.
  1233. &#34;&#34;&#34;
  1234. M1 = np.array([[np.cos(phi), np.sin(phi), 0],
  1235. [-np.sin(phi), np.cos(phi), 0],
  1236. [0, 0, 1]])
  1237. M2 = np.array([[1, 0, 0],
  1238. [0, np.cos(theta), np.sin(theta)],
  1239. [0, -np.sin(theta), np.cos(theta)]])
  1240. M3 = np.array([[np.cos(view), np.sin(view), 0],
  1241. [-np.sin(view), np.cos(view), 0],
  1242. [0, 0, 1]])
  1243. M = np.matmul(M1, np.matmul(M2, M3)) # combine rotations
  1244. r = np.tensordot(self.r_vec, M, axes=[1, 0])
  1245. return r[:,0:2]</code></pre>
  1246. </details>
  1247. </dd>
  1248. <dt id="simulation.Simulation.run"><code class="name flex">
  1249. <span>def <span class="ident">run</span></span>(<span>self, tmax, saveEvery, outputFolder, **kwargs)</span>
  1250. </code></dt>
  1251. <dd>
  1252. <section class="desc"><p>Run the galactic simulation.</p>
  1253. <h2 id="attributes">Attributes</h2>
  1254. <dl>
  1255. <dt><strong><code>tmax</code></strong> :&ensp;<code>float</code></dt>
  1256. <dd>Time to which the simulation will run to.
  1257. This is measured here since the start of the simulation,
  1258. not since pericenter.</dd>
  1259. <dt><strong><code>saveEvery</code></strong> :&ensp;<code>int</code></dt>
  1260. <dd>The state is saved every saveEvery steps.</dd>
  1261. <dt><strong><code>outputFolder</code></strong> :&ensp;<code>string</code></dt>
  1262. <dd>It will be saved to /data/outputFolder/</dd>
  1263. </dl></section>
  1264. <details class="source">
  1265. <summary>Source code</summary>
  1266. <pre><code class="python">def run(self, tmax, saveEvery, outputFolder, **kwargs):
  1267. &#34;&#34;&#34;Run the galactic simulation.
  1268. Attributes:
  1269. tmax (float): Time to which the simulation will run to.
  1270. This is measured here since the start of the simulation,
  1271. not since pericenter.
  1272. saveEvery (int): The state is saved every saveEvery steps.
  1273. outputFolder (string): It will be saved to /data/outputFolder/
  1274. &#34;&#34;&#34;
  1275. # When the simulation starts, intialize self.rprev_vec
  1276. self.rprev_vec = self.r_vec - self.v_vec * self.dt
  1277. if self.verbose: print(&#39;Simulation starting. Bon voyage!&#39;)
  1278. while(self.t &lt; tmax):
  1279. self.step()
  1280. if(self.n % saveEvery == 0):
  1281. self.save(&#39;data/{}&#39;.format(outputFolder))
  1282. print(&#39;Simulation complete.&#39;)</code></pre>
  1283. </details>
  1284. </dd>
  1285. <dt id="simulation.Simulation.save"><code class="name flex">
  1286. <span>def <span class="ident">save</span></span>(<span>self, outputFolder)</span>
  1287. </code></dt>
  1288. <dd>
  1289. <section class="desc"><p>Save the state of the simulation to the outputFolder.
  1290. Two files are saved:
  1291. sim{self.n}.pickle: serializing the state.
  1292. sim{self.n}.png: a simplified 2D plot of x, y.</p></section>
  1293. <details class="source">
  1294. <summary>Source code</summary>
  1295. <pre><code class="python">def save(self, outputFolder):
  1296. &#34;&#34;&#34;Save the state of the simulation to the outputFolder.
  1297. Two files are saved:
  1298. sim{self.n}.pickle: serializing the state.
  1299. sim{self.n}.png: a simplified 2D plot of x, y.
  1300. &#34;&#34;&#34;
  1301. # Create the output folder if it doesn&#39;t exist
  1302. if not os.path.exists(outputFolder): os.makedirs(outputFolder)
  1303. # Compute some useful quantities
  1304. # v_vec is not required by the integrator, but useful
  1305. self.v_vec = (self.r_vec - self.rprev_vec)/self.dt
  1306. self.t = self.n * self.dt # prevents numerical rounding errors
  1307. # Serialize state
  1308. file = open(outputFolder+&#39;/data{}.pickle&#39;.format(self.n), &#34;wb&#34;)
  1309. pickle.dump({&#39;r_vec&#39;: self.r_vec, &#39;v_vec&#39;: self.v_vec,
  1310. &#39;type&#39;: self.type, &#39;mass&#39;: self.mass,
  1311. &#39;CONFIG&#39;: self.CONFIG, &#39;t&#39;: self.t,
  1312. &#39;tracks&#39;: self.tracks}, file)
  1313. # Save simplified plot of the current state.
  1314. # Its main use is to detect ill-behaved situations early on.
  1315. fig = plt.figure()
  1316. plt.xlim(-5, 5); plt.ylim(-5, 5); plt.axis(&#39;equal&#39;)
  1317. # Dark halo is plotted in red, disk in blue, bulge in green
  1318. PLTCON = [(&#39;dark&#39;, &#39;r&#39;, 0.3), (&#39;disk&#39;, &#39;b&#39;, 1.0), (&#39;bulge&#39;, &#39;g&#39;, 0.5)]
  1319. for type_, c, a in PLTCON:
  1320. plt.scatter(self.r_vec[self.type[:,0]==type_][:,0],
  1321. self.r_vec[self.type[:,0]==type_][:,1], s=0.1, c=c, alpha=a)
  1322. # Central mass as a magenta star
  1323. plt.scatter(self.r_vec[self.type[:,0]==&#39;center&#39;][:,0],
  1324. self.r_vec[self.type[:,0]==&#39;center&#39;][:,1], s=100, marker=&#34;*&#34;, c=&#39;m&#39;)
  1325. # Save to png file
  1326. fig.savefig(outputFolder+&#39;/sim{}.png&#39;.format(self.n), dpi=150)
  1327. plt.close(fig)</code></pre>
  1328. </details>
  1329. </dd>
  1330. <dt id="simulation.Simulation.setOrbit"><code class="name flex">
  1331. <span>def <span class="ident">setOrbit</span></span>(<span>self, galaxy1, galaxy2, e, rmin, R0)</span>
  1332. </code></dt>
  1333. <dd>
  1334. <section class="desc"><p>Sets the two galaxies galaxy1, galaxy2 in an orbit.</p>
  1335. <h2 id="parameters">Parameters</h2>
  1336. <dl>
  1337. <dt><strong><code>galaxy1</code></strong> :&ensp;<a title="simulation.Galaxy" href="#simulation.Galaxy"><code>Galaxy</code></a></dt>
  1338. <dd>1st galaxy (main)</dd>
  1339. <dt><strong><code>galaxy2</code></strong> :&ensp;<a title="simulation.Galaxy" href="#simulation.Galaxy"><code>Galaxy</code></a></dt>
  1340. <dd>2nd galaxy (companion)</dd>
  1341. <dt><strong><code>e</code></strong></dt>
  1342. <dd>eccentricity of the orbit</dd>
  1343. <dt><strong><code>rmin</code></strong></dt>
  1344. <dd>distance of closest approach</dd>
  1345. <dt><strong><code>R0</code></strong></dt>
  1346. <dd>initial separation</dd>
  1347. </dl></section>
  1348. <details class="source">
  1349. <summary>Source code</summary>
  1350. <pre><code class="python">def setOrbit(self, galaxy1, galaxy2, e, rmin, R0):
  1351. &#34;&#34;&#34;Sets the two galaxies galaxy1, galaxy2 in an orbit.
  1352. Parameters:
  1353. galaxy1 (Galaxy): 1st galaxy (main)
  1354. galaxy2 (Galaxy): 2nd galaxy (companion)
  1355. e: eccentricity of the orbit
  1356. rmin: distance of closest approach
  1357. R0: initial separation
  1358. &#34;&#34;&#34;
  1359. m1, m2 = np.sum(galaxy1.mass), np.sum(galaxy2.mass)
  1360. # Relevant formulae:
  1361. # $r_0 = r (1 + e) \cos(\phi)$, where $r_0$ ($\neq R_0$) is the semi-latus rectum
  1362. # $r_0 = r_\textup{min} (1 + e)$
  1363. # $v^2 = G M (2/r - 1/a)$, where a is the semimajor axis
  1364. # Solve the reduced two-body problem with reduced mass $\mu$ (mu)
  1365. M = m1 + m2
  1366. r0 = rmin * (1 + e)
  1367. try:
  1368. phi = np.arccos((r0/R0 - 1) / e) # inverting the orbit equation
  1369. phi = -np.abs(phi) # Choose negative (incoming) angle
  1370. ainv = (1 - e) / rmin # ainv = $1/a$, as a may be infinite
  1371. v = np.sqrt(M * (2/R0 - ainv))
  1372. # Finally, calculate the angle of motion. angle = tan(dy/dx)
  1373. # $dy/dx = ((dr/d\phi) sin(\phi) + r \cos(\phi))/((dr/d\phi) cos(\phi) - r \sin(\phi))$
  1374. dy = R0/r0 * e * np.sin(phi)**2 + np.cos(phi)
  1375. dx = R0/r0 * e * np.sin(phi) * np.cos(phi) - np.sin(phi)
  1376. vangle = np.arctan2(dy, dx)
  1377. except: raise Exception(&#39;&#39;&#39;The orbital parameters cannot be satisfied.
  1378. For elliptical orbits check that R0 is posible (&lt;rmax).&#39;&#39;&#39;)
  1379. # We now need the actual motion of each body
  1380. R_vec = np.array([[R0*np.cos(phi), R0*np.sin(phi), 0.]])
  1381. V_vec = np.array([[v*np.cos(vangle), v*np.sin(vangle), 0.]])
  1382. galaxy1.r_vec += m2/M * R_vec
  1383. galaxy1.v_vec += m2/M * V_vec
  1384. galaxy2.r_vec += -m1/M * R_vec
  1385. galaxy2.v_vec += -m1/M * V_vec
  1386. # Explicitely add the galaxies to the simulation
  1387. self.add(galaxy1)
  1388. self.add(galaxy2)
  1389. if self.verbose: print(&#39;Galaxies set in orbit.&#39;)</code></pre>
  1390. </details>
  1391. </dd>
  1392. <dt id="simulation.Simulation.step"><code class="name flex">
  1393. <span>def <span class="ident">step</span></span>(<span>self)</span>
  1394. </code></dt>
  1395. <dd>
  1396. <section class="desc"><p>Perform a single step of the simulation.
  1397. Makes use of a 4th order Verlet integrator.</p></section>
  1398. <details class="source">
  1399. <summary>Source code</summary>
  1400. <pre><code class="python">def step(self):
  1401. &#34;&#34;&#34;Perform a single step of the simulation.
  1402. Makes use of a 4th order Verlet integrator.
  1403. &#34;&#34;&#34;
  1404. # Calculate the acceleration
  1405. self.a_vec = self.acceleration(self.r_vec, self.mass, soft=self.soft)
  1406. # Update the state using the Verlet algorithm
  1407. # (A custom algorithm is written mainly for learning purposes)
  1408. self.r_vec, self.rprev_vec = (2*self.r_vec - self.rprev_vec
  1409. + self.a_vec * self.dt**2, self.r_vec)
  1410. self.n += 1
  1411. # Update tracks
  1412. self.tracks = np.concatenate([self.tracks,
  1413. self.r_vec[self.type[:,0]==&#39;center&#39;][:,np.newaxis]], axis=1)</code></pre>
  1414. </details>
  1415. </dd>
  1416. </dl>
  1417. </dd>
  1418. </dl>
  1419. </section>
  1420. </article>
  1421. <nav id="sidebar">
  1422. <h1>Index</h1>
  1423. <div class="toc">
  1424. <ul></ul>
  1425. </div>
  1426. <ul id="index">
  1427. <li><h3><a href="#header-classes">Classes</a></h3>
  1428. <ul>
  1429. <li>
  1430. <h4><code><a title="simulation.Galaxy" href="#simulation.Galaxy">Galaxy</a></code></h4>
  1431. <ul class="">
  1432. <li><code><a title="simulation.Galaxy.__init__" href="#simulation.Galaxy.__init__">__init__</a></code></li>
  1433. <li><code><a title="simulation.Galaxy.addBulge" href="#simulation.Galaxy.addBulge">addBulge</a></code></li>
  1434. <li><code><a title="simulation.Galaxy.addDisk" href="#simulation.Galaxy.addDisk">addDisk</a></code></li>
  1435. <li><code><a title="simulation.Galaxy.addHalo" href="#simulation.Galaxy.addHalo">addHalo</a></code></li>
  1436. <li><code><a title="simulation.Galaxy.rotate" href="#simulation.Galaxy.rotate">rotate</a></code></li>
  1437. </ul>
  1438. </li>
  1439. <li>
  1440. <h4><code><a title="simulation.Simulation" href="#simulation.Simulation">Simulation</a></code></h4>
  1441. <ul class="two-column">
  1442. <li><code><a title="simulation.Simulation.__init__" href="#simulation.Simulation.__init__">__init__</a></code></li>
  1443. <li><code><a title="simulation.Simulation.add" href="#simulation.Simulation.add">add</a></code></li>
  1444. <li><code><a title="simulation.Simulation.project" href="#simulation.Simulation.project">project</a></code></li>
  1445. <li><code><a title="simulation.Simulation.run" href="#simulation.Simulation.run">run</a></code></li>
  1446. <li><code><a title="simulation.Simulation.save" href="#simulation.Simulation.save">save</a></code></li>
  1447. <li><code><a title="simulation.Simulation.setOrbit" href="#simulation.Simulation.setOrbit">setOrbit</a></code></li>
  1448. <li><code><a title="simulation.Simulation.step" href="#simulation.Simulation.step">step</a></code></li>
  1449. </ul>
  1450. </li>
  1451. </ul>
  1452. </li>
  1453. </ul>
  1454. </nav>
  1455. </main>
  1456. <footer id="footer">
  1457. <p>Generated by <a href="https://pdoc3.github.io/pdoc"><cite>pdoc</cite> 0.5.2</a>.</p>
  1458. </footer>
  1459. <script src="https://cdnjs.cloudflare.com/ajax/libs/highlight.js/9.12.0/highlight.min.js"></script>
  1460. <script>hljs.initHighlightingOnLoad()</script>
  1461. </body>
  1462. </html>