Jekyll2018-08-13T07:26:14+00:00https://ayushsuhane.github.io/Ayush SuhaneGraduate StudentFinal Benchmarks (Part II)2018-08-12T00:00:00+00:002018-08-12T00:00:00+00:00https://ayushsuhane.github.io/FinalBenchmarks-PartII<p>In the final set of blog posts, lets have a look at how all the changes in MDAnalysis resulted in superior performance of distance based evaluations. In particular, we will look into three different applications (1) Radial Distribution Function (RDF) (2) Identifying Bonds (3) Distance based Selections. Actual Performance comparisons with all the previous commits can be cheked <a href="https://www.mdanalysis.org/benchmarks/#/">here</a>. Below mentioned benchmarks are local performance comparisons which were made during the incorporation of the code, and some might require a little bit of extra efforts to replicate as the original code will not be in the repository anymore (particularly distance based selections). In any case, I would appreciate to check the commits <a href="https://github.com/MDAnalysis/mdanalysis/pull/2035">here</a> to use the previously implemented methods.</p>
<h2 id="rdf">RDF</h2>
<p>Radial Distribution Function describes the variation of density as a function of distance between two sets of particles. More details on RDF can be found <a href="https://en.wikipedia.org/wiki/Radial_distribution_function">here</a>. RDF is implemented in <code class="highlighter-rouge">analysis.rdf</code> module of MDAnalysis, where it builds a sparse matrix of distances and then bins the distances to generate a RDF plot. As discussed before, sparse matrices are fast for smaller datasets but are inefficient in terms of memory and time especially when the size/ number of particles increase beyond 10K in a machine with 4 GB RAM. Therefore, the current implementation fails to calculate RDF for larger dataset. This functionality can be extended by the use of <code class="highlighter-rouge">capped_function</code> along with the promise of superior performance than the current implementation. The main gist of rdf calculation lies in <code class="highlighter-rouge">_single_frame</code> method in <code class="highlighter-rouge">analysis.rdf.InterRDF</code> function:</p>
<div class="language-python highlighter-rouge"><div class="highlight"><pre class="highlight"><code>
<span class="n">distances</span><span class="o">.</span><span class="n">distance_array</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">g1</span><span class="o">.</span><span class="n">positions</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">g2</span><span class="o">.</span><span class="n">positions</span><span class="p">,</span>
<span class="n">box</span><span class="o">=</span><span class="bp">self</span><span class="o">.</span><span class="n">u</span><span class="o">.</span><span class="n">dimensions</span><span class="p">,</span> <span class="n">result</span><span class="o">=</span><span class="bp">self</span><span class="o">.</span><span class="n">_result</span><span class="p">)</span>
<span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_exclusion_mask</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_exclusion_mask</span><span class="p">[:]</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">_maxrange</span>
<span class="n">count</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">histogram</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_result</span><span class="p">,</span> <span class="o">**</span><span class="bp">self</span><span class="o">.</span><span class="n">rdf_settings</span><span class="p">)[</span><span class="mi">0</span><span class="p">]</span>
</code></pre></div></div>
<p>which can be replaced by:</p>
<div class="language-python highlighter-rouge"><div class="highlight"><pre class="highlight"><code>
<span class="n">pairs</span><span class="p">,</span> <span class="n">dist</span> <span class="o">=</span> <span class="n">capped_distance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">g1</span><span class="o">.</span><span class="n">positions</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">g2</span><span class="o">.</span><span class="n">positions</span><span class="p">,</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_maxrange</span><span class="p">,</span>
<span class="n">box</span><span class="o">=</span><span class="bp">self</span><span class="o">.</span><span class="n">u</span><span class="o">.</span><span class="n">dimensions</span><span class="p">)</span>
<span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_exclusion_block</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span>
<span class="n">A</span><span class="p">,</span> <span class="n">B</span> <span class="o">=</span> <span class="n">pairs</span><span class="p">[:,</span> <span class="mi">0</span><span class="p">]</span><span class="o">//</span><span class="bp">self</span><span class="o">.</span><span class="n">_exclusion_block</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">pairs</span><span class="p">[:,</span> <span class="mi">1</span><span class="p">]</span><span class="o">//</span><span class="bp">self</span><span class="o">.</span><span class="n">_exclusion_block</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span>
<span class="n">C</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">where</span><span class="p">(</span><span class="n">A</span> <span class="o">!=</span> <span class="n">B</span><span class="p">)[</span><span class="mi">0</span><span class="p">]</span>
<span class="n">dist</span> <span class="o">=</span> <span class="n">dist</span><span class="p">[</span><span class="n">C</span><span class="p">]</span>
<span class="n">count</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">histogram</span><span class="p">(</span><span class="n">dist</span><span class="p">,</span> <span class="o">**</span><span class="bp">self</span><span class="o">.</span><span class="n">rdf_settings</span><span class="p">)[</span><span class="mi">0</span><span class="p">]</span>
</code></pre></div></div>
<p>where we first calculate the pairs and distance for each pairs which can be used to construct a histogram. To deal with exclusion blocks, a simple mask using strides in <code class="highlighter-rouge">g1, g2</code> can be used in a way such that if the indices of atom in <code class="highlighter-rouge">g1</code> and <code class="highlighter-rouge">g2</code> are in the multiple of exclusion block indices, they should be ignored. More information about the details and performance can be found in this <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/RDF_Comparison.ipynb">notebook</a>. The performance is tested against the datafiles (TPR, XTC) which are already present in <code class="highlighter-rouge">MDAnalysisTests</code> . Below is the performance of older RDF function and the current implemenation of RDF function which uses <code class="highlighter-rouge">capped_distance</code> function.</p>
<p><img src="/images/120818_rdf.PNG" alt="RDF Comparison Benchmarks" /></p>
<p>As can be seen, the performance of RDF has been improved significantly especially for larger set of particles. The x-axis in the graph represents the number of atoms in <code class="highlighter-rouge">g1, g2</code> each. The improvements in time are of the factor of ~ 2-3 for single core computations. However, there might be more improvements possible due to parallellization and SIMD instructions. More discussions on the possible ways to improve RDF calculations can also be found <a href="https://github.com/MDAnalysis/mdanalysis/pull/2013">here</a>.</p>
<h2 id="bonds">Bonds</h2>
<p>One of the important requirement in Analysis of Molecular DYnamics simulations is to be able to distinguish bonds between all the atoms. From an algorithmic perspective, it is a problem of fixed neighbor search where two particles are said to be mutually bonded if the distance between them is smaller than the sum of their individual radius. However, this requires finding neighbors of every atoms, which is time consuming specially when tackled with bruteforce approach. <code class="highlighter-rouge">self_capped_distance</code> is an ideal candidate for such applications. The main approach is to identify all the pairs which are within the distance equivalent to the maximum diameter among all the atoms followed by checking the individual distances between them.</p>
<div class="language-python highlighter-rouge"><div class="highlight"><pre class="highlight"><code>
<span class="n">pairs</span><span class="p">,</span> <span class="n">dist</span> <span class="o">=</span> <span class="n">distances</span><span class="o">.</span><span class="n">self_capped_distance</span><span class="p">(</span><span class="n">coords</span><span class="p">,</span>
<span class="n">max_cutoff</span><span class="o">=</span><span class="mf">2.0</span><span class="o">*</span><span class="n">max_vdw</span><span class="p">,</span>
<span class="n">min_cutoff</span><span class="o">=</span><span class="n">lower_bound</span><span class="p">,</span>
<span class="n">box</span><span class="o">=</span><span class="n">box</span><span class="p">)</span>
<span class="k">for</span> <span class="n">idx</span><span class="p">,</span> <span class="p">(</span><span class="n">i</span><span class="p">,</span> <span class="n">j</span><span class="p">)</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">pairs</span><span class="p">):</span>
<span class="n">d</span> <span class="o">=</span> <span class="p">(</span><span class="n">vdwradii</span><span class="p">[</span><span class="n">atomtypes</span><span class="p">[</span><span class="n">i</span><span class="p">]]</span> <span class="o">+</span> <span class="n">vdwradii</span><span class="p">[</span><span class="n">atomtypes</span><span class="p">[</span><span class="n">j</span><span class="p">]])</span><span class="o">*</span><span class="n">fudge_factor</span>
<span class="k">if</span> <span class="p">(</span><span class="n">dist</span><span class="p">[</span><span class="n">idx</span><span class="p">]</span> <span class="o"><</span> <span class="n">d</span><span class="p">):</span>
<span class="n">bonds</span><span class="o">.</span><span class="n">append</span><span class="p">((</span><span class="n">atoms</span><span class="p">[</span><span class="n">i</span><span class="p">]</span><span class="o">.</span><span class="n">index</span><span class="p">,</span> <span class="n">atoms</span><span class="p">[</span><span class="n">j</span><span class="p">]</span><span class="o">.</span><span class="n">index</span><span class="p">))</span>
</code></pre></div></div>
<p>The notebook with actual benchmarks can be found <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/guessbonds_benchmark.ipynb">here</a>, which shows a performance improvement by an order of magnitude for the case under consideration with 8284 bonds. The file used for benchmark is <code class="highlighter-rouge">smal.gro</code> and can be found <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/small.gro">here</a>. Apparently, the factor of time improvement also improves with increase in number of datasets primarily due to the scaling of time complexity of brute force and <code class="highlighter-rouge">periodic KDTree</code>/<code class="highlighter-rouge">nsgrid</code>.</p>
<p>In this case the time taken by old implementation of <code class="highlighter-rouge">topology.guessers.guess_bonds</code> took 3.26 secs, while using <code class="highlighter-rouge">capped_function</code> improved the time to 394 ms.</p>
<h2 id="distance-based-selections">Distance Based Selections</h2>
<p>Before <code class="highlighter-rouge">capped_distance</code>, distance based selections utilized two different methods <code class="highlighter-rouge">_apply_kdtree</code> and <code class="highlighter-rouge">_apply_distmat</code> which relied on user’s judgement and required user to have apriori knowledge on both the methods. Distance based selectiion was also memory intensive for PBC calculations as <code class="highlighter-rouge">_apply_distmat</code> was enforced for periodic selections. <code class="highlighter-rouge">capped_distance</code> deals with these problems in an elegant way by hiding the technical details of the method from the user to focus on the analysis rather than spending time understanding the different algorithms of distance calculations. ‘around’ selection was benchmarked with GRO datafile for the current use case. For the benchmark case of <code class="highlighter-rouge">around 5.0 resid 1</code>, previous and current implementation was checked for around selections. The <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/AroundBenchmark.ipynb">notebook</a> contains the benchmark informations. It can be seen that it took 44 ms with the implemented bruteforce approach for non-PBC calculations. The time can be improved by using KDtree to 37 ms. While the improvement is not that significant with kdtree, <code class="highlighter-rouge">nsgrid</code> improves the performance by an order of magnitude. <code class="highlighter-rouge">NSGrid</code> took 5 ms for the similar task which is a significant improvement in distance based selections (around a factor of 10).</p>
<p>Overall, three cases which depended on distance evaluations were selected and significant improvement in performance is demonstrated. The next task is to implement capped distance in other distance based analysis such as contact maps, cumulative RDF etc.</p>
<p>Till then, Adios Amigos!!!</p>Improvements in RDF, Guessbonds, Distance Based SelectionsFinal Benchmarks (Part I)2018-08-11T00:00:00+00:002018-08-11T00:00:00+00:00https://ayushsuhane.github.io/FinalBenchmarks<p>Lets begin by reiterating the objectives of the project followed by a bird’s eye view to the approach and finally the quantification of performance improvements in distance evaluations.</p>
<p>The objective of the project is to improve the performance of distance based calculations and their corresponding applications in MDAnalysis. Brute force calculations are highly space inefficient and have a time complexity of O(N^2), which becomes a bottleneck in Molecular Dynamics simulations especially when the size of datapoints is in millions. As discussed before, multiple data structures and algorithms exist which can be more efficient depending on the size and distribution of the datapoints. Binary tree, cell-lists are few among many data structures, which can be used for faster computations to identify nearest neighbours/fixed radius neighbors.</p>
<p>One of the primary aim was to introduce a framework which can automatically select the fastest method based on the size of data, cutoff distance and domain size. This would keep the methods hidden from the user and choose the fastest available method internally, and user doesn’t have to worry about the technical details of the method. However, an extra handle would be provided to the user to override this automatic selection, if desired, and to choose a method of choice. Another advantage of this framework would be an easy extensibility to other faster algorithms/data structures. Addition of another fast algorithm would require just a separate definition and its correspoding rules in one single file and all the applications could reflect the superior performance without much effort.</p>
<p>To introduce the framework, a <code class="highlighter-rouge">capped_distance</code> is defined as discussed in previous posts, which takes in query coordinates, search coordinates, cutoff distance as mandatory arguments and <code class="highlighter-rouge">min_cutoff</code>, <code class="highlighter-rouge">box</code>, and <code class="highlighter-rouge">method</code> as optional arguments. Currently <code class="highlighter-rouge">method</code> supports <code class="highlighter-rouge">bruteforce</code>, <code class="highlighter-rouge">pkdtree</code> and <code class="highlighter-rouge">nsgrid</code>, which direct to naive brute force implementation, periodic KDTree( a wrapper around <code class="highlighter-rouge">scipy.spatial</code> ) and <code class="highlighter-rouge">cell-list</code> data structure motivated by <code class="highlighter-rouge">GROMACS</code> (more precisely <code class="highlighter-rouge">nsgrid.c <https://github.com/gromacs/gromacs/commits/master/src/mdlib/nsgrid.c>_</code>) respectively. To automatically select the method, an internal method <code class="highlighter-rouge">_determine_method</code> is used, which includes rules to choose the method. At present, the rules are based on the optimized cutoff distance and size of data points, while in principle, different parameters such as distribution can also be used here. Every method is defined by a separate function, for instance a brute force method inside a <code class="highlighter-rouge">capped_distance</code> has same signature as that of <code class="highlighter-rouge">capped_distance</code> except the <code class="highlighter-rouge">method</code> argument and is defined as <code class="highlighter-rouge">_bruteforce_capped</code> as an internal function. Similarly, a self_capped_distance with identical signature is also defined to deal with all the distance evaluations related to all the atoms in the system.</p>
<p>The next step is to optimize individual functions. While one can evaluate the distance matrix, which is fast but limited by sytem memory, another way is to loop over every query and check all the neighbors by searching for neighbors of one query at a time. Even though looping over numpy array is a slower, but this approach is not limited by system memory atleast for the size dealt in MD simulations. However, since we have multiple methods implemented in <code class="highlighter-rouge">capped_distance</code>, it would be intuitive to use brute-force algorithm to the limit of system memory and then switch over to other fast methods as a first iteration. Similarly, for KDTree, the current implementation involved using <code class="highlighter-rouge">Biopython.KDTree</code> which can be replaced with more stable <code class="highlighter-rouge">scipy.spatial</code>. The current implementation of periodic boundary conditions involves looping over every coordinate and evaluating the distance with relevant periodic images to find the minimum distance. As explained in previous post, this operation of looping can be replaced by a more faster cython version, where all the relevant periodic images of particles are concatenated to the original particles. This is followed by a non-periodic search over the new combined dataset which should yield the fixed radius neighbors following Periodic Boudary conditions. This is included in <code class="highlighter-rouge">augment_coordinates</code> in <code class="highlighter-rouge">lib._augment.pyx</code> in MDAanalysis. Augmenting coordinates also has another advantage which allows easy extensibility to any near neighbor search algorithm to handle Periodic Boundary condition with a little overhead cost. Another method <code class="highlighter-rouge">nsgrid</code> which is also discussed in previous blogposts is also added to the <code class="highlighter-rouge">capped_distance</code> function.</p>
<p>Now moving on to the performance comparisons of these methods, the identified parameters are cutoff distance, datasize, and pbc/no pbc conditions. Almost all of the Benchmarks referred in this blog series can be found <a href="https://github.com/ayushsuhane/Benchmarks_Distance/tree/master/Notebooks">here</a>, but the benchmarks for <code class="highlighter-rouge">capped_distance</code> can be found in this <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/BM_CappedNS.ipynb">notebook</a>. Below are few of the important snippets to compare the performances for <code class="highlighter-rouge">capped_function</code> with different methods. To begin with, lets compare a single query search for all the implemented methods, its variation with cutoff radius and its sensivity with PBC/no-PBC condition.</p>
<table>
<thead>
<tr>
<th style="text-align: center">Cutoff= 4 units</th>
<th style="text-align: center">DataSize = 100K</th>
<th style="text-align: center">NoPBC Datasize = 100K</th>
</tr>
</thead>
<tbody>
<tr>
<td style="text-align: center"><img src="/images/110718_sqnum.PNG" width="200" /></td>
<td style="text-align: center"><img src="/images/110718_sqcutpbc.PNG" width="200" /></td>
<td style="text-align: center"><img src="/images/110718_sqcutnopbc.PNG" width="200" /></td>
</tr>
</tbody>
</table>
<p>As anticipated, bruteforce should be the fastest for this case, as the cost of construction of data structure is larger than querying all the combinations of distances. However, for larger number of queries, the other two methods are expected to be more faster than the naive <code class="highlighter-rouge">bruteforce</code> implementation. As a next comparison, lets change the number of query points keeping the size of the data constant i.e. searching varying number of query points against a fixed data set and assess the performance.</p>
<table>
<thead>
<tr>
<th style="text-align: center">Cutoff= 4 units</th>
<th style="text-align: center">DataSize = 100K</th>
<th style="text-align: center">NoPBC Datasize = 100K</th>
</tr>
</thead>
<tbody>
<tr>
<td style="text-align: center"><img src="/images/110718_mqnum.PNG" width="200" /></td>
<td style="text-align: center"><img src="/images/110718_mqcut.PNG" width="200" /></td>
<td style="text-align: center"><img src="/images/110718_mqcutnopbc.PNG" width="200" /></td>
</tr>
</tbody>
</table>
<p>As it can be seen, <code class="highlighter-rouge">NSGrid</code> is the fastest for multiple query searches where size is greater than 10 particles. Furthermore, one can contemplate from the trends that <code class="highlighter-rouge">NSGrid</code> slows down if the search radius is above a critical value, i.e. greater than 20-30% of the box size in this case.</p>
<p>Another aspect is searching all the particles within a certain cutoff distance, which is defined in <code class="highlighter-rouge">lib.distances.self_capped_distance</code>. The performance benchmarks below are for PBC and there is a similar trend for noPBC calculations as well.</p>
<table>
<thead>
<tr>
<th style="text-align: center">Cutoff= 4 units</th>
<th style="text-align: center">DataSize = 100K</th>
<th style="text-align: center">Datasize = 100K</th>
</tr>
</thead>
<tbody>
<tr>
<td style="text-align: center"><img src="/images/110718_selfnum.PNG" width="200" /></td>
<td style="text-align: center"><img src="/images/110718_selfcut.PNG" width="200" /></td>
<td style="text-align: center"><img src="/images/110718_selfcut10k.PNG" width="200" /></td>
</tr>
</tbody>
</table>
<p>As it can be seen here, PeriodicKDTree shows superior performance, if the desired cutoff distance is below a critical value (typically below 3% of the domain size) and <code class="highlighter-rouge">nsgrid</code> becomes superior for larger distances.</p>
<p>The above benchmarks can be used in <code class="highlighter-rouge">_determine_method</code> and <code class="highlighter-rouge">_determine_method_self</code> to chose the method automatically without user’s intervention. These rules can then be applied for multiple distance based applications in MDAnalysis to gain performance over the previous methods.</p>
<p>The next blogpost will cover the actual improvements for few of the applications in MDAnalysis such as Radial Distribution Function, Distance based selections and identification of bonds in a given dataset. See you there!!</p>Benchmarks for automatic selectionCell Lists and advantages2018-08-08T00:00:00+00:002018-08-08T00:00:00+00:00https://ayushsuhane.github.io/Cell-List-Algorithm<p>We had discussed earlier that cell-lists have algorithmic advantage over naive brute force method as well as tree based data structures for uniform distributed datasets. To reiterate, KDTree require o(NlogN) time for construction of data structure and o(logN) for the search. Cell-list algorithm, on the contrary, require O(N) time for the construction of the data structure, and O(1) time for search for a best case. The best case for cell-list is achievable when the particles are well-dispered in the domain, and the minimum required distance between neighbors is such there exist only one particle in a cell. However, such a case ceases to exist almost always in practical applications.</p>
<p>At its core, cell-list data structure includes splitting of the domain into smaller subdomains called cells. This is followed by sorting the atoms into these cells. Once sorted, for any given query, its cell-location can be identified (<code class="highlighter-rouge">cellindex</code>) and naive brute force can be used to calculate the distances between query and all the particles residing in nearest neighbor cells. As can be contemplated, a competition between cell-list and KDtree should exist since increasing the number of particles per cell would eventually lead to increase in evaluations (<code class="highlighter-rouge">O(N^2)</code>) as compared to KDTree (<code class="highlighter-rouge">O(logN)</code>).</p>
<p>From an implementation perspective, the basic implementation of cell-list is quite straightforward and can also be found in Appendix F, Page 552 of
<code class="highlighter-rouge">Understanding Molecular Dynamics: From Algorithm to Applications</code> by Frenkel and Smit. Provided the cutoff radius, cell-size is typically chosen as the cutoff radius. All the particles are then distributed into their respective cells based on their position. Any query for a nearest neighbour / fixed radius search, identifies the cell of the query and searches all the neighboring cells only i.e. 27 cells. This reduces the computational cost tremendously as compared to naive brute force calculations.</p>
<p>While its implementation is also relatively straightforward for orthogonal boxes, another challenge remains is to make cell-list work with triclinic boxes. The challenge of different systems arise due to the requirement to handle Periodic boundary conditions in Molecular dynamics. To explain in brief, a periodic boundary condition in terms of distance calculations means a minimum distance between two particles and their periodic images should be considered as the distance between the particles. Inorder to satisfy the minimum image convention, the particles must be translated in all three directions equivalent to the box vectors, and corresponding distances should be evaluated to get the minimum distance. Inorder to calculate the distances in triclinic systems, an equivalent brick shaped box is constructed and all the particles are wrapped around the brick shaped box. The dimension of the brick-shaped box are defined as the component of box vectors in respective directions i.e. <code class="highlighter-rouge">[[ax, ay, az], [bx, by, bz], [cx, cy, cz]]</code> as the box vectors, then the boundary of the brick shaped box is defined by <code class="highlighter-rouge">[-ax/2, ax/2], [-bx/2, by/2], [-cz/2, cz/2]]</code>. To evaluate the minimum distances, every displacement vector is modified by the addition/substraction of box vectors such that the distance becomes less than half the box width in every direction. In python, it can be written as follows:</p>
<p>.. code-block ::python:</p>
<div class="highlighter-rouge"><div class="highlight"><pre class="highlight"><code>for i in range(natoms):
for m in range(DIM - 1, -1, -1):
while bbox_coords[i, m] < 0:
for d in range(m+1):
bbox_coords[i, d] += self.c_pbcbox.box[m][d]
while bbox_coords[i, m] >= self.c_pbcbox.box[m][m]:
for d in range(m+1):
bbox_coords[i, d] -= self.c_pbcbox.box[m][d]
</code></pre></div></div>
<p>Another practical problem is that the cell-list data structure can become huge for small cutoff radius. This is prevented by defining a <code class="highlighter-rouge">max_gridsize</code> which limits the number of cells within a domain. A class <code class="highlighter-rouge">FastNS</code> is defined which contains all the serach functions, while other helper Classes like <code class="highlighter-rouge">PBCBox</code> handles the brick-shaped box and evaluation of minimum distance for periodic/non-periodic calculations, <code class="highlighter-rouge">NSResults</code> is a class which initializes and maintains containers to access the results of search query, and <code class="highlighter-rouge">NSGrid</code> contains all the grid related functions to distribute atoms in cells, define cell index and more. For reference, the actual file <code class="highlighter-rouge">nsgrid</code> can be found <a href="https://github.com/MDAnalysis/mdanalysis/blob/develop/package/MDAnalysis/lib/nsgrid.pyx">here</a>. Few of the relevant snippets are shown below:</p>
<p>In <code class="highlighter-rouge">NSGrid</code> object, a linear mapping of ever coordinate to cellindex in defined as :</p>
<p>.. code-block ::python:</p>
<div class="highlighter-rouge"><div class="highlight"><pre class="highlight"><code>cdef int coord2cellid(self, vec coord) nogil:
return <int> (coord[ZZ] / cellsize[ZZ]) * (cell_offsets[ZZ]) +\
<int> (coord[YY] / cellsize[YY]) * cell_offsets[YY] + \
<int> (coord[XX] / cellsize[XX])
</code></pre></div></div>
<p>However, this is valid only under the assumption that all the <code class="highlighter-rouge">coord</code> are within the brick shaped box spanned by cells. Inorder to distribute the coordinates to their respective cell</p>
<p>.. code-block ::python:</p>
<div class="highlighter-rouge"><div class="highlight"><pre class="highlight"><code>for i in range(ncoords):
# Add bead to grid cell
cellindex = self.cellids[i]
self.beadids[cellindex * self.nbeads_per_cell + beadcounts[cellindex]] = i
beadcounts[cellindex] += 1
</code></pre></div></div>
<p>This snippet, finds the <code class="highlighter-rouge">cellid</code> of every coordinate, and traverses the <code class="highlighter-rouge">beadids</code> array which is a list of all the atoms present in a particluar cell. A fixed size of <code class="highlighter-rouge">nbeads_per_cell</code> is used to make a uniform array, where <code class="highlighter-rouge">nbeads_per_cell</code> is the maximum number of atoms in any cell. Once the particles are distributed in the grid, search can be performed using search functions defined in <code class="highlighter-rouge">FastNS</code> class.</p>
<p>Primarily, two methods are introduced in the <code class="highlighter-rouge">FastNS</code> class with different objectives. The <code class="highlighter-rouge">search</code> method searches the already populated grid with a list of queries to return all the pairs between list of query coordinates and populated coordinates within a fixed radius. <code class="highlighter-rouge">self_search</code> searches the coordinates within the grid to find all the neighbours. While in first look, <code class="highlighter-rouge">self_search</code> might look as a subset of <code class="highlighter-rouge">search</code> method, but time can be halved directly for this special case as distance between every pair needs to be calculated just once. This small trick saves a significant amount of time from the <code class="highlighter-rouge">search</code> method.</p>
<p>The <code class="highlighter-rouge">search</code> method finds the pairs and distances by creating another grid with the <code class="highlighter-rouge">NSGrid</code> object followed by searching querines from one grid against initially populated grid. In contrast, the <code class="highlighter-rouge">self_search</code> doesn’t involve creating another grid, but searches within itself. This is the layout of the Grid-Search method which is supposed to increase the performance of distance evaluations. We shall compare the performance results in the next blog post.</p>
<p>Stay tuned to see the performance improvements!!</p>Cell-list implementation in MDAnalysisAugment Coordinates for Implementation of Periodic Boundary Conditions2018-07-06T00:00:00+00:002018-07-06T00:00:00+00:00https://ayushsuhane.github.io/Augment-Coordinates-for-Periodic-Boundary-Conditions<p>Periodic Boundary conditions are one of the essential requirements of Molecular Dynamics Simulations. Majority of the MD simulations are executed with PBC in all three directions. Therefore, all the distance computations in MDAnalysis should also support PBC in triclinic simulation boxes. In simpler terms, a particle exiting from one of the face of a box will also enter simulataneously from the corresponding opposite face. In three dimensions, this creates 26 exact replicas of any particle in the neighbourhood of the box under consideration. Inorder to calculate the distance between two particles under the influence of periodic boundary conditions, a minimum image convention is followed. Minimum image distance is the smallest possible distance between two particles in the box as well as between the images in the surrounding duplicate boxes. This naive implementation of PBC for distance calculation increases the computation cost by 27 times.</p>
<p>Present implementation to handle PBC using KDTree in MDAnalysis is a more advanced approach. Rather than evaluating all the images, only the relevant images are generated and checked based on the relative position of both the particles. For a given box dimensions, its reciprocal vectors are evaluated. If <code class="highlighter-rouge">[[a1, a2, a3], [b1, b2, b3], [c1, c2, c3]]</code> are the dimensions of the box such that <code class="highlighter-rouge">[a, b, c]</code> are the three vectors of the box, the reciprocal vectors <code class="highlighter-rouge">[ra, rb, rc]</code> are evaluated as:</p>
<pre><code class="language-math">ra = b X c
rb = c X a
rc = a X b
</code></pre>
<p>The reciprocal vectors represent the normal vectors to the plane formed by the its member vectors i.e. <code class="highlighter-rouge">ra</code> is the plane normal vector of a plane formed by <code class="highlighter-rouge">b</code> and <code class="highlighter-rouge">c</code> vectors. Since the problem is inherently three dimensional, both the particles in consideration have a location defined by their respective three dimensional vectors. A dot product with the reciprocal vectors yield the distance of the coordinate from the corresponding face along the direction of normal (reciprocal) vector. This concept can be exploited to limit the generation of images which are within a certain distance from the box faces. The cutoff distance of <code class="highlighter-rouge">[a/2, b/2, c/2]</code> would be the condition of equivalency of this approach with the naive PBC implementation described above.</p>
<p>However, typical problems in MDAnalysis doesnot require to go to the extreme of choosing half the box distance as the cutoff distance. A conditional evaluation of image distance only of the coordinate lie inside the cutoff distance saves a significant amount of time. The current algorithm in MDAnalysis is pure python implementation and can be described as follows:</p>
<ul>
<li>Evaluate the reciprocal vectors.</li>
<li>Find the distance of the coordinates from each of the 6 faces.</li>
<li>If any of the distance is less than the specified distance, generate images and check the distance based on the minimum image convention.</li>
<li>Return the minimum distance.</li>
</ul>
<p>This has two drawbacks though :</p>
<ul>
<li>Pure python implementation is slower as the evaluation of minimum distance requires multiple loops.</li>
<li>Other algorithms for near neighbour searching which donot have specific implementation of PBC cannot be handled direclty.</li>
</ul>
<p>To make the algorithm more flexible, an approach of augment coordinates is adopted. In this approach, all the particles within the cutoff distance are duplicated as per the Periodic Boundary Conditions. Another mapping of these particles with the original particles is also constructed to undo the augmentation. These duplicated particles along with the original particles can now be operated with any algorithm to evaluate the distances without any special consideratio to PBC. All the indices/coordinates can be mapped back to the original set of particles using the mapping described above. A cython implementation is implemented for triclinic boxes for fast looping based on this <a href="http://www.richardjgowers.com/2018/06/28/make_halos.html">notebook</a>. Furthermore, memoryviews are used here based on this <a href="https://jakevdp.github.io/blog/2012/08/08/memoryview-benchmarks/">observation</a>.</p>
<p>Similar algorithm described above is used here for cython implementation. To evaluate the distance of particle from each of the faces, dot product operation of two vectors i.e. <code class="highlighter-rouge">[x, y, z]</code> (particle coordinate) and <code class="highlighter-rouge">[x, y, z] - ([a1, a2, a3] + [b1, b2, b3] + [c1, c2, c3])</code> i.e distance of coordinate from the opposite diagonal box is performed. Furthermore, duplicates in the diagonal images of box are generated based on the condition if the particle distance evaluated from more than one faces in less than the specified distance. For example, if the particle lies close to the xy plane but far from the yz plane, then corresponding duplicate will be generated in xy plane near the opposite plane. The function to augment the coordinates and undo the arguments can be written as</p>
<div class="language-python highlighter-rouge"><div class="highlight"><pre class="highlight"><code> <span class="o">%%</span><span class="n">cython</span> <span class="o">--</span><span class="n">annotate</span>
<span class="n">cimport</span> <span class="n">cython</span>
<span class="kn">import</span> <span class="nn">cython</span>
<span class="n">cimport</span> <span class="n">numpy</span> <span class="k">as</span> <span class="n">np</span>
<span class="kn">import</span> <span class="nn">numpy</span> <span class="k">as</span> <span class="n">np</span>
<span class="k">from</span> <span class="n">libc</span><span class="o">.</span><span class="n">math</span> <span class="n">cimport</span> <span class="n">sqrt</span>
<span class="nd">@cython.boundscheck</span><span class="p">(</span><span class="bp">False</span><span class="p">)</span>
<span class="nd">@cython.wraparound</span><span class="p">(</span><span class="bp">False</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">make_halo</span><span class="p">(</span><span class="nb">float</span><span class="p">[:,</span> <span class="p">::</span><span class="mi">1</span><span class="p">]</span> <span class="n">coordinates</span><span class="p">,</span> <span class="nb">float</span><span class="p">[:,::</span><span class="mi">1</span><span class="p">]</span> <span class="n">box</span><span class="p">,</span> <span class="nb">float</span><span class="p">[:,::</span><span class="mi">1</span><span class="p">]</span> <span class="n">reciprocal</span><span class="p">,</span> <span class="nb">float</span> <span class="n">r</span><span class="p">):</span>
<span class="s">"""Calculate augmented coordinate set
Parameters
----------
coordinates : np.ndarray
coordinates to augment
box : np.ndarray
size of box
rm : np.ndarray
reciprocal matrix
r : float
thickness of halo region to buffer by
Returns
-------
augmented : np.ndarray
coordinates of the new augmented coordinates
indices : np.ndarray
original indices of the augmented coordinates
"""</span>
<span class="n">cdef</span> <span class="n">bint</span> <span class="n">lo_x</span><span class="p">,</span> <span class="n">hi_x</span><span class="p">,</span> <span class="n">lo_y</span><span class="p">,</span> <span class="n">hi_y</span><span class="p">,</span> <span class="n">lo_z</span><span class="p">,</span> <span class="n">hi_z</span>
<span class="n">cdef</span> <span class="nb">int</span> <span class="n">i</span><span class="p">,</span> <span class="n">j</span><span class="p">,</span> <span class="n">p</span><span class="p">,</span> <span class="n">N</span>
<span class="n">cdef</span> <span class="nb">float</span> <span class="n">shiftX</span><span class="p">[</span><span class="mi">3</span><span class="p">]</span>
<span class="n">cdef</span> <span class="nb">float</span> <span class="n">shiftY</span><span class="p">[</span><span class="mi">3</span><span class="p">]</span>
<span class="n">cdef</span> <span class="nb">float</span> <span class="n">shiftZ</span><span class="p">[</span><span class="mi">3</span><span class="p">]</span>
<span class="n">cdef</span> <span class="nb">float</span> <span class="n">coord</span><span class="p">[</span><span class="mi">3</span><span class="p">]</span>
<span class="n">cdef</span> <span class="nb">float</span> <span class="n">end</span><span class="p">[</span><span class="mi">3</span><span class="p">]</span>
<span class="n">cdef</span> <span class="nb">float</span> <span class="n">other</span><span class="p">[</span><span class="mi">3</span><span class="p">]</span>
<span class="n">cdef</span> <span class="nb">int</span> <span class="n">dim</span>
<span class="n">dim</span> <span class="o">=</span> <span class="n">coordinates</span><span class="o">.</span><span class="n">shape</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span>
<span class="c"># room for adding triclinic support by using (3,) vectors</span>
<span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="n">dim</span><span class="p">):</span>
<span class="n">shiftX</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">box</span><span class="p">[</span><span class="mi">0</span><span class="p">,</span> <span class="n">i</span><span class="p">]</span>
<span class="n">shiftY</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">box</span><span class="p">[</span><span class="mi">1</span><span class="p">,</span> <span class="n">i</span><span class="p">]</span>
<span class="n">shiftZ</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">box</span><span class="p">[</span><span class="mi">2</span><span class="p">,</span> <span class="n">i</span><span class="p">]</span>
<span class="n">end</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">box</span><span class="p">[</span><span class="mi">0</span><span class="p">,</span> <span class="n">i</span><span class="p">]</span> <span class="o">+</span> <span class="n">box</span><span class="p">[</span><span class="mi">1</span><span class="p">,</span> <span class="n">i</span><span class="p">]</span> <span class="o">+</span> <span class="n">box</span><span class="p">[</span><span class="mi">2</span><span class="p">,</span><span class="n">i</span><span class="p">]</span>
<span class="n">N</span> <span class="o">=</span> <span class="n">coordinates</span><span class="o">.</span><span class="n">shape</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span>
<span class="n">p</span> <span class="o">=</span> <span class="mi">0</span> <span class="c"># output counter</span>
<span class="c"># allocate output arrays</span>
<span class="c"># could be more conservative with this</span>
<span class="c"># or use C++ vectors + push etc</span>
<span class="n">cdef</span> <span class="nb">float</span><span class="p">[:,</span> <span class="p">:]</span> <span class="n">output</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">zeros</span><span class="p">((</span><span class="n">N</span><span class="p">,</span> <span class="mi">3</span><span class="p">),</span> <span class="n">dtype</span><span class="o">=</span><span class="n">np</span><span class="o">.</span><span class="n">float32</span><span class="p">)</span>
<span class="n">cdef</span> <span class="nb">int</span><span class="p">[:]</span> <span class="n">indices</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">zeros</span><span class="p">(</span><span class="n">N</span><span class="p">,</span> <span class="n">dtype</span><span class="o">=</span><span class="n">np</span><span class="o">.</span><span class="n">int32</span><span class="p">)</span>
<span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">0</span><span class="p">,</span><span class="n">N</span><span class="p">):</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coordinates</span><span class="p">[</span><span class="n">i</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span>
<span class="n">other</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">end</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">coordinates</span><span class="p">[</span><span class="n">i</span><span class="p">,</span><span class="n">j</span><span class="p">]</span>
<span class="c"># identify the condition </span>
<span class="n">lo_x</span> <span class="o">=</span> <span class="n">_dot</span><span class="p">(</span><span class="o">&</span><span class="n">coord</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="o">&</span><span class="n">reciprocal</span><span class="p">[</span><span class="mi">0</span><span class="p">,</span><span class="mi">0</span><span class="p">])</span> <span class="o"><=</span> <span class="n">r</span>
<span class="n">hi_x</span> <span class="o">=</span> <span class="n">_dot</span><span class="p">(</span><span class="o">&</span><span class="n">other</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="o">&</span><span class="n">reciprocal</span><span class="p">[</span><span class="mi">0</span><span class="p">,</span><span class="mi">0</span><span class="p">])</span> <span class="o"><=</span> <span class="n">r</span>
<span class="n">lo_y</span> <span class="o">=</span> <span class="n">_dot</span><span class="p">(</span><span class="o">&</span><span class="n">coord</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="o">&</span><span class="n">reciprocal</span><span class="p">[</span><span class="mi">1</span><span class="p">,</span><span class="mi">0</span><span class="p">])</span> <span class="o"><=</span> <span class="n">r</span>
<span class="n">hi_y</span> <span class="o">=</span> <span class="n">_dot</span><span class="p">(</span><span class="o">&</span><span class="n">other</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="o">&</span><span class="n">reciprocal</span><span class="p">[</span><span class="mi">1</span><span class="p">,</span><span class="mi">0</span><span class="p">])</span> <span class="o"><=</span> <span class="n">r</span>
<span class="n">lo_z</span> <span class="o">=</span> <span class="n">_dot</span><span class="p">(</span><span class="o">&</span><span class="n">coord</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="o">&</span><span class="n">reciprocal</span><span class="p">[</span><span class="mi">2</span><span class="p">,</span><span class="mi">0</span><span class="p">])</span> <span class="o"><=</span> <span class="n">r</span>
<span class="n">hi_z</span> <span class="o">=</span> <span class="n">_dot</span><span class="p">(</span><span class="o">&</span><span class="n">other</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="o">&</span><span class="n">reciprocal</span><span class="p">[</span><span class="mi">2</span><span class="p">,</span><span class="mi">0</span><span class="p">])</span> <span class="o"><=</span> <span class="n">r</span>
<span class="k">if</span> <span class="n">lo_x</span><span class="p">:</span>
<span class="c"># if X, face piece</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="c"># add to output</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="c"># keep record of which index this augmented position was created from</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_y</span><span class="p">:</span>
<span class="c"># if X&Y, edge piece</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="c"># if X&Y&Z, corner piece</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_y</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_x</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_y</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_y</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftX</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_y</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_y</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftY</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">if</span> <span class="n">lo_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">+</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">elif</span> <span class="n">hi_z</span><span class="p">:</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">output</span><span class="p">[</span><span class="n">p</span><span class="p">,</span> <span class="n">j</span><span class="p">]</span> <span class="o">=</span> <span class="n">coord</span><span class="p">[</span><span class="n">j</span><span class="p">]</span> <span class="o">-</span> <span class="n">shiftZ</span><span class="p">[</span><span class="n">j</span><span class="p">]</span>
<span class="n">indices</span><span class="p">[</span><span class="n">p</span><span class="p">]</span> <span class="o">=</span> <span class="n">i</span>
<span class="n">p</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="k">return</span> <span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">(</span><span class="n">output</span><span class="p">[:</span><span class="n">p</span><span class="p">]),</span> <span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">(</span><span class="n">indices</span><span class="p">[:</span><span class="n">p</span><span class="p">])</span>
<span class="nd">@cython.boundscheck</span><span class="p">(</span><span class="bp">False</span><span class="p">)</span>
<span class="nd">@cython.wraparound</span><span class="p">(</span><span class="bp">False</span><span class="p">)</span>
<span class="n">cdef</span> <span class="nb">float</span> <span class="n">_dot</span><span class="p">(</span><span class="nb">float</span> <span class="o">*</span><span class="n">a</span><span class="p">,</span> <span class="nb">float</span> <span class="o">*</span><span class="n">b</span><span class="p">):</span>
<span class="s">"""Return dot product of two sequences in range."""</span>
<span class="n">cdef</span> <span class="n">ssize_t</span> <span class="n">n</span>
<span class="n">cdef</span> <span class="nb">float</span> <span class="n">sum1</span>
<span class="n">cdef</span> <span class="n">ssize_t</span> <span class="n">dim</span>
<span class="n">dim</span><span class="o">=</span><span class="mi">3</span>
<span class="n">sum1</span> <span class="o">=</span> <span class="mf">0.0</span>
<span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="n">dim</span><span class="p">):</span>
<span class="n">sum1</span> <span class="o">+=</span> <span class="n">a</span><span class="p">[</span><span class="n">n</span><span class="p">]</span> <span class="o">*</span> <span class="n">b</span><span class="p">[</span><span class="n">n</span><span class="p">]</span>
<span class="k">return</span> <span class="n">sum1</span>
<span class="nd">@cython.boundscheck</span><span class="p">(</span><span class="bp">False</span><span class="p">)</span>
<span class="nd">@cython.wraparound</span><span class="p">(</span><span class="bp">False</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">undo_augment</span><span class="p">(</span><span class="nb">int</span><span class="p">[:]</span> <span class="n">results</span><span class="p">,</span> <span class="nb">int</span><span class="p">[:]</span> <span class="n">translation</span><span class="p">,</span> <span class="nb">int</span> <span class="n">nreal</span><span class="p">):</span>
<span class="s">"""Translate augmented indices back to originals
Note: modifies results in place!
Parameters
----------
results : ndarray of ints
indices of coordinates, including "augmented" indices
translation : ndarray of ints
original indices of augmented coordinates
nreal : int
number of real coordinates, ie values in results equal or larger than this
need to be translated to their real counterpart
Returns
-------
results : ndarray of ints
"""</span>
<span class="n">cdef</span> <span class="nb">int</span> <span class="n">N</span>
<span class="n">N</span> <span class="o">=</span> <span class="n">results</span><span class="o">.</span><span class="n">shape</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span>
<span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="n">N</span><span class="p">):</span>
<span class="k">if</span> <span class="n">results</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">>=</span> <span class="n">nreal</span><span class="p">:</span>
<span class="n">results</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">translation</span><span class="p">[</span><span class="n">results</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">-</span> <span class="n">nreal</span><span class="p">]</span>
<span class="k">return</span> <span class="n">results</span>
</code></pre></div></div>
<p>Followed by the definition of this function, the API for distance evaluations which returns all the pairs within a specified distance using KDTree can be written as:</p>
<div class="language-python highlighter-rouge"><div class="highlight"><pre class="highlight"><code><span class="k">def</span> <span class="nf">augment_kdtree</span><span class="p">(</span><span class="n">coords</span><span class="p">,</span> <span class="n">cutoff</span><span class="p">,</span> <span class="n">box</span><span class="p">):</span>
<span class="s">"""
box of type [A,B,C,a1,a2,a3]
"""</span>
<span class="n">box</span> <span class="o">=</span> <span class="n">triclinic_vectors</span><span class="p">(</span><span class="n">box</span><span class="p">)</span>
<span class="c"># Define reciprocal matrix</span>
<span class="n">rm</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">zeros</span><span class="p">(</span><span class="mi">9</span><span class="p">,</span> <span class="n">dtype</span><span class="o">=</span><span class="n">np</span><span class="o">.</span><span class="n">float32</span><span class="p">)</span><span class="o">.</span><span class="n">reshape</span><span class="p">(</span><span class="mi">3</span><span class="p">,</span> <span class="mi">3</span><span class="p">)</span>
<span class="n">rm</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">cross</span><span class="p">(</span><span class="n">box</span><span class="p">[</span><span class="mi">1</span><span class="p">],</span> <span class="n">box</span><span class="p">[</span><span class="mi">2</span><span class="p">])</span>
<span class="n">rm</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">cross</span><span class="p">(</span><span class="n">box</span><span class="p">[</span><span class="mi">2</span><span class="p">],</span> <span class="n">box</span><span class="p">[</span><span class="mi">0</span><span class="p">])</span>
<span class="n">rm</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">cross</span><span class="p">(</span><span class="n">box</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">box</span><span class="p">[</span><span class="mi">1</span><span class="p">])</span>
<span class="k">for</span> <span class="n">i</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">3</span><span class="p">):</span>
<span class="n">rm</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">/=</span> <span class="n">norm</span><span class="p">(</span><span class="n">rm</span><span class="p">[</span><span class="n">i</span><span class="p">])</span> <span class="c"># normalize</span>
<span class="n">aug</span><span class="p">,</span> <span class="n">idx</span> <span class="o">=</span> <span class="n">make_halo</span><span class="p">(</span><span class="n">coords</span><span class="p">,</span> <span class="n">box</span><span class="p">,</span> <span class="n">rm</span><span class="p">,</span> <span class="n">cutoff</span><span class="p">)</span>
<span class="n">aug_coord</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">concatenate</span><span class="p">([</span><span class="n">coords</span><span class="p">,</span> <span class="n">aug</span><span class="p">])</span>
<span class="n">kdtree</span> <span class="o">=</span> <span class="n">spatial</span><span class="o">.</span><span class="n">cKDTree</span><span class="p">(</span><span class="n">aug_coord</span><span class="p">)</span>
<span class="n">pairs</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">(</span><span class="nb">list</span><span class="p">(</span><span class="n">kdtree</span><span class="o">.</span><span class="n">query_pairs</span><span class="p">(</span><span class="n">cutoff</span><span class="p">)),</span> <span class="n">dtype</span><span class="o">=</span><span class="n">np</span><span class="o">.</span><span class="n">int32</span><span class="p">)</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">pairs</span><span class="p">)</span> <span class="o">></span> <span class="mi">1</span><span class="p">:</span>
<span class="n">undo_augment</span><span class="p">(</span><span class="n">pairs</span><span class="p">[:,</span> <span class="mi">0</span><span class="p">],</span> <span class="n">idx</span><span class="p">,</span> <span class="nb">len</span><span class="p">(</span><span class="n">coords</span><span class="p">))</span>
<span class="n">undo_augment</span><span class="p">(</span><span class="n">pairs</span><span class="p">[:,</span> <span class="mi">1</span><span class="p">],</span> <span class="n">idx</span><span class="p">,</span> <span class="nb">len</span><span class="p">(</span><span class="n">coords</span><span class="p">))</span>
<span class="n">pairs</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">unique</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">sort</span><span class="p">(</span><span class="n">pairs</span><span class="p">),</span> <span class="n">axis</span><span class="o">=</span><span class="mi">0</span><span class="p">)</span>
<span class="k">return</span> <span class="n">pairs</span>
</code></pre></div></div>
<p>Similarly, other distance evaluations such as atom selections can be performed using similar methods. To benchmark the cython implementation of KDtree and already implemented PKDtree in MDAnalysis, a use case to identify all the bonds i.e. distance between two particles less than the sum of radius of each particle is chosen. Two different real case are also chosen with ~12k and ~1000k particles to check the extremeties of performance. While more detailed analysis is provided <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/Augment_Functionality.ipynb">here</a> and <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/Augment_Functionality.ipynb">here</a>, major highlights of the benchmark is provided below.</p>
<p>For the two cases described above the guess_bonds function in <code class="highlighter-rouge">MDAnalysis.topology.guessers.py</code> took <code class="highlighter-rouge">2.86 sec</code> for smaller system, while the <code class="highlighter-rouge">PKDTree</code> in <code class="highlighter-rouge">MDAnalysis.lib.pkdtree</code> took <code class="highlighter-rouge">2.69 sec</code>, whereas the function defined here took <code class="highlighter-rouge">920 ms</code> i.e increase in performance by 3 times. For the larger system, the computation time for ` guess_bonds<code class="highlighter-rouge"> was more than </code>30 min<code class="highlighter-rouge"> whereas the current implementation of PKDTree in MDAnalysis took around </code>7 mins<code class="highlighter-rouge">. For this case, the method described above took </code>2 mins` i.e. the performance boost of 3 times than already implemented PKDtree (but not in guess bonds), whereas the performance boost is significantly high compared to brute force approach currently implemented in MDAnalysis for guessing the bonds.</p>
<p>The next step is to include this functionality in MDAnalysis and replacement of Bio.KDtree, which is obsolete in BioPython, with scipy.spatial.cKDTree.</p>
<p>See you next time!!!</p>Handling Periodic boundary conditions with augmenting particle coordinatesOptimization of distance evaluations2018-06-23T00:00:00+00:002018-06-23T00:00:00+00:00https://ayushsuhane.github.io/Optimization%20of%20distance%20evaluations<p>In the last blog post, we discussed the advantages of gridsearch algorithm for three dimensional datasets, but it can be memory expensive based on the required cutoff radius. As KDTree and brute force are implemented and deeply rooted in MDAnalysis, an optimization of selection of algorithm i.e. choosing the algorithm based on certain empirical rules hidden from the user, would be helpful in fast distance evaluations. For this reason, the primary objective is to maintain a function which can handle different algorithms and switch between them based on the data provided to it.
The skeleton of the function included the definition of the function which takes two sets of coordinates i.e. atom positions accessed by <code class="highlighter-rouge">Universe.atoms.positions</code>, maximum and minimum cutoff radius where minimum cutoff radius is kept as optional, box which can be accessed from <code class="highlighter-rouge">universe.dimensions</code> but is kept as optional along with an optional handle provided to the user to override the automatic selection of the method. The implememtation of periodic boundary conditions are automatically considered when a box is supplied. Currently, KDTree and Bruteforce method are supported in this function, but the flexible structure allows adding other methods i.e. grid search method once it becomes available.</p>
<p>The main functionality is handled by <code class="highlighter-rouge">cappedfunction(reference, configuration, max_cutoff, min_cutoff=None, box=None, method=None)</code> which resides at <code class="highlighter-rouge">MDAnalysis.lib.distances.py</code>. The function returns array of all the indices of the pairs as well as their corresponding distances which are within the distance specified by <code class="highlighter-rouge">max_cutoff</code> and <code class="highlighter-rouge">min_cutoff</code>. Unless the method is specified by the user, a function (<code class="highlighter-rouge">_determine_method</code>) which automatically choses the method based on the number of particles in either set of coordinates and/or search radius is called by the <code class="highlighter-rouge">cappeddistance</code> function. This function (<code class="highlighter-rouge">determine_method</code>) returns the method object which points towards an efficient algorithm. Each algorithm has a different method defined i.e. KDTree has <code class="highlighter-rouge">_pkdtree_capped</code> and bruteforce has <code class="highlighter-rouge">_bruteforce_capped</code> method defined in the function as:</p>
<div class="language-python highlighter-rouge"><div class="highlight"><pre class="highlight"><code><span class="k">def</span> <span class="nf">_bruteforce_capped</span><span class="p">(</span><span class="n">reference</span><span class="p">,</span> <span class="n">configuration</span><span class="p">,</span> <span class="n">max_cutoff</span><span class="p">,</span> <span class="n">min_cutoff</span><span class="o">=</span><span class="bp">None</span><span class="p">,</span> <span class="n">box</span><span class="o">=</span><span class="bp">None</span><span class="p">):</span>
<span class="n">pairs</span><span class="p">,</span> <span class="n">distance</span> <span class="o">=</span> <span class="p">[],</span> <span class="p">[]</span>
<span class="n">reference</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">asarray</span><span class="p">(</span><span class="n">reference</span><span class="p">,</span> <span class="n">dtype</span><span class="o">=</span><span class="n">np</span><span class="o">.</span><span class="n">float32</span><span class="p">)</span>
<span class="n">configuration</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">asarray</span><span class="p">(</span><span class="n">configuration</span><span class="p">,</span> <span class="n">dtype</span><span class="o">=</span><span class="n">np</span><span class="o">.</span><span class="n">float32</span><span class="p">)</span>
<span class="k">if</span> <span class="n">reference</span><span class="o">.</span><span class="n">shape</span> <span class="o">==</span> <span class="p">(</span><span class="mi">3</span><span class="p">,</span> <span class="p">):</span>
<span class="n">reference</span> <span class="o">=</span> <span class="n">reference</span><span class="p">[</span><span class="bp">None</span><span class="p">,</span> <span class="p">:]</span>
<span class="k">if</span> <span class="n">configuration</span><span class="o">.</span><span class="n">shape</span> <span class="o">==</span> <span class="p">(</span><span class="mi">3</span><span class="p">,</span> <span class="p">):</span>
<span class="n">configuration</span> <span class="o">=</span> <span class="n">configuration</span><span class="p">[</span><span class="bp">None</span><span class="p">,</span> <span class="p">:]</span>
<span class="n">_check_array</span><span class="p">(</span><span class="n">reference</span><span class="p">,</span> <span class="s">'reference'</span><span class="p">)</span>
<span class="n">_check_array</span><span class="p">(</span><span class="n">configuration</span><span class="p">,</span> <span class="s">'configuration'</span><span class="p">)</span>
<span class="k">for</span> <span class="n">i</span><span class="p">,</span> <span class="n">coords</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="n">reference</span><span class="p">):</span>
<span class="n">dist</span> <span class="o">=</span> <span class="n">distance_array</span><span class="p">(</span><span class="n">coords</span><span class="p">[</span><span class="bp">None</span><span class="p">,</span> <span class="p">:],</span> <span class="n">configuration</span><span class="p">,</span> <span class="n">box</span><span class="o">=</span><span class="n">box</span><span class="p">)[</span><span class="mi">0</span><span class="p">]</span>
<span class="k">if</span> <span class="n">min_cutoff</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span>
<span class="n">idx</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">where</span><span class="p">((</span><span class="n">dist</span> <span class="o"><=</span> <span class="n">max_cutoff</span><span class="p">)</span> <span class="o">&</span> <span class="p">(</span><span class="n">dist</span> <span class="o">></span> <span class="n">min_cutoff</span><span class="p">))[</span><span class="mi">0</span><span class="p">]</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">idx</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">where</span><span class="p">((</span><span class="n">dist</span> <span class="o"><=</span> <span class="n">max_cutoff</span><span class="p">))[</span><span class="mi">0</span><span class="p">]</span>
<span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="n">idx</span><span class="p">:</span>
<span class="n">pairs</span><span class="o">.</span><span class="n">append</span><span class="p">((</span><span class="n">i</span><span class="p">,</span> <span class="n">j</span><span class="p">))</span>
<span class="n">distance</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="n">dist</span><span class="p">[</span><span class="n">j</span><span class="p">])</span>
<span class="k">return</span> <span class="n">pairs</span><span class="p">,</span> <span class="n">distance</span>
</code></pre></div></div>
<p>and similarly for KDTree. To remove the circular imports within <code class="highlighter-rouge">pkdtree.py</code> and <code class="highlighter-rouge">distances.py</code>, local import to the modules are performed in the <code class="highlighter-rouge">_pkdtree_capped</code> function. The rules are selected such as if either of the particle set consisted of less than 5000 particles or if the number of particles is greater than 100k and the maximum cutoff distance is greater than 10% of the box size naive brute force is selected, whereas kdtree is selected in all the other cases. A notebook with the benchmarks to evaluate these rules can be found <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/Capped_Distance.ipynb">here</a>.</p>
<p>For any new method, it needs to be registered into the <code class="highlighter-rouge">_determine_method</code> function dictionary i.e</p>
<div class="language-python highlighter-rouge"><div class="highlight"><pre class="highlight"><code><span class="n">methods</span> <span class="o">=</span> <span class="p">{</span><span class="s">'bruteforce'</span><span class="p">:</span> <span class="n">_bruteforce_capped</span><span class="p">,</span>
<span class="s">'pkdtree'</span><span class="p">:</span> <span class="n">_pkdtree_capped</span><span class="p">}</span>
</code></pre></div></div>
<p>and a corresponding function needs to be defined as well. The corresponding rules can be inserted in the <code class="highlighter-rouge">_determine_method</code> without disturbing the main algorithm for distance evaluations.</p>
<p>This is the brief discussion of the method which has various advantages i.e. it can be used to replace a lot of analysis functions such as guess_bonds, radial distribution function analysis, particle selections etc in an optimized way and totally hidden from the user (along with the flexibility to give user the control to override), and resulting in a more maintainable and clear code.</p>
<p>This is all there is to it now. We shall discuss it later as well, after the implementation of gridsearch method where optimized rules will be added for gridsearch method in this function as well.</p>Inclusion of capped functionInitial Benchmarks (Continued ..)2018-06-09T00:00:00+00:002018-06-09T00:00:00+00:00https://ayushsuhane.github.io/Initial-Benchmarks-Continued<p>This is the continuation of previous post, and includes more benchmarks for different use cases. Last post concluded few initial benchmarks comparing KDtree, Cellgrid and Brute force method for two use case (1)Contact searches (2) Near neighbour selection. An obvious conclusion that cell list and KDtree are better alternatives as compared to brute force for large number of particles was quantified in terms of execution time. However, brute force is a better option for very less number of particles. Other data structures were also compared for single point queries and it was concluded that neighbour search as implemented in FATSLiM and Octree are viable candidates for distance calculations.</p>
<p>Since it was established that <a href="https://github.com/PointCloudLibrary">Octree</a> and neighbour search in <a href="http://fatslim.github.io/">FATSLiM</a> could be good replacements for fast distance evaluations for single point queries (fixed radius selection and near neighbour selections), a more detailed analysis of pair contact searches is covered in this post. For contact searches, grid based methods become more and more useful as usage of sparse matrix can become handy when the density is high as opposed to traversing the tree to identify each particle’s node in tree based data structures. This usage of sparse matrix reduces the query time, which becomes a huge time saving factor for large number of particles. As a next step to previous benchmarks, the Cellgrid package is optimized and used here to search for the contact pairs(see <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/CellGrid_Optimization.ipynb">here</a>). As a modification, a optional argument was added to the <code class="highlighter-rouge">capped_self_distance_array</code> which calculates the cellsize corresponding to 30 particles per cell assuming a uniform distribution of particles. The cellsize has a lower bound of cutoff distance, since decreasing the cellsize below the cutoff distance would require to extend the search to cells beyond the immediate neihbourhood of the query point.</p>
<p>CellGrid package, KDTree and brute force methods are first compared for the use case to find all the pairs within a fixed distance (see <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/BM_PairContact_MOD.ipynb">here</a>). For brute force calculations, <code class="highlighter-rouge">distance_array</code> is used in place of <code class="highlighter-rouge">self_distance_array</code> due to its high memory consumption. To avoid the repetatitve calculations due to each pair, care is taken to evaluate distance only once for each pair. Similarly, individual queries is pruned in KDtree and a list with all the pairs is maintained. For the case of Cellgrid, a distance matrix containing the particles within its neighbourhood is constructed. This matrix is masked based on the distances. To study the scaling behaviour of all three methods, a cutoff distance of 10 distance units is fixed and execution time is registered for different particle densities in a fixed box.</p>
<p><img src="/images/090518_paircon_rad10.PNG" alt="alt text" /></p>
<p>Time to evaluate particle within a distance of 10 distance units</p>
<p>Since, cutoff radius is an important factor for performance, similar studies were performed to characterize the effect of cutoff radius for different particle densities. Although the trend with cutoff distance is continuously increasing for complete range of particle density i.e. for number of particles ranging from 100 to 100k. However, a significant relative shift can be seen from low to high number of particles. While brute force takes notably less time as compared to other data structures, the transition to cellgrid becomes evident at higher particle density.</p>
<p><img src="/images/090518_paircon_n100.PNG" alt="Variation of execution time for different cutoff distances for 100 particles" /></p>
<p>Variation of execution time for different cutoff distances for 100 particles</p>
<p><img src="/images/090518_paircon_n17k.PNG" alt="Variation of execution time for different cutoff distances for 17k particles" /></p>
<p>Variation of execution time for different cutoff distances for 17k particles</p>
<p><img src="/images/090518_paircon_n100.PNG" alt="Variation of execution time for different cutoff distances for 100k particles" /></p>
<p>Variation of execution time for different cutoff distances for 100k particles</p>
<p>It can be seen that Cell-lists become advantageous as the number of particles increase. Since FATSLiM have a more efficient implementation of grid search, The next step is to check the timings of PBC aware Neighbour search module for pair contact searches. For this benchmark, we chose a particular implementation of bonds_guess, which is implemented in MDAnalysis at <code class="highlighter-rouge">MDAnalysis.topology.guessers</code>. The goal of this function is to identify the bonds between atoms by identifying the neighbouring atoms and checking the distance between them relative to the sum of their radius. Current implementation evaluates all the pair contacts for each particle using naive distance algorithm. It is anticipated that this algorithm is very costly and can be replaced with other data structures to improve the performance. It can be seen that while tree and grid structures both are advantageous at very large particle densities, Neighbour search (linear grid search) is more advantageous at intermediate particle densities as well. As expected a transition from brute force to KDTree/ Cell-List is achieved at lower particle densities i.e. around 1k for cell-list and around 6k for Periodic KDTree (see <a href="http://localhost:8888/notebooks/GuessBonds.ipynb">here</a>).</p>
<p><img src="/images/090518_bondsguess.PNG" alt="Benchmarks for Guessing bonds between atoms in a static dataset of atoms" /></p>
<p>Benchmarks for Guessing bonds between atoms in a static dataset of atoms</p>
<p>Apart from the increase in performance, another trick which only calculates distances from half the neighbours (as implemented in Cellgrid) is not implemented in Neighbour search routine. It is anticipated that performance of this method can be further increased by adopting this procedure.</p>
<p>In conclusion, for fixed radius neighbor selection, cython optimized neighbour search routine is a robust algorithm which can be used in tandem with brute force method and /KDTree in few cases depending on the application. Another recommendation would be to switch from Biopython KDtree to scipy cKDTree, which is found to be faster than the former(see <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/OctreeComparison.ipynb">here</a>). The order of magnitude depends on the query point, but a performance improvement of the order of 5 times is observed as a best case i.e. for uniform distribution. For contact searches, it is already shown in this blog that, a cutoff of around 1000 particles can be used to switch between brute force and neighbour search routine from FATSLiM. Modified script for the neihbour search routing can be found <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/GuessBonds/core_ns.pyx">here</a>.</p>
<p>Thanks to Sébastian to introduce the cython optimized Neighbour search routine along with the benchmarks (see <a href="https://github.com/seb-buch/Benchmarks_Distance/blob/master/Notebooks/CythonNS.ipynb">here</a>).</p>
<p>Till Next Time!!!</p>Initial comparison with cell-listInitial Benchmarks of Distance Calculations2018-05-23T00:00:00+00:002018-05-23T00:00:00+00:00https://ayushsuhane.github.io/First%20Post-Initial-Benchmarks<p><em>Using Jupyter Notebooks for the first time, and they are awesome!!!</em></p>
<p>This is the first post after the start of GSoC. A reflection on the past few weeks and summary of the things I have tried, failed and learned should be a good start to this blog.</p>
<h1 id="journey-so-far">Journey so far</h1>
<hr />
<p>Community bonding period started as soon as the good news arrived, but it ended quickly before I could notice it. Although exams and few other chores also had their fair share in consuming the time. With the remaining time at hand, I focussed on the first part of the project.
To give a brief introduction, the project focusses on speeding up the distance calculations, which is one of the fundamental requirement in the MDAnalysis package and Molecular Dynamics Analysis packages in general. More information on the project can be found <a href="https://summerofcode.withgoogle.com/projects/#5050592943144960">here</a>.</p>
<p>Over the course of summer, the goal is to identify the suitable candidate for efficient distance calculations. Typically, the number of particles in a Molecular dynamics ranges from few hundred atoms to couple of millions of atoms. The solution to efficient distance based calculations vary significantly based on different factors such as number of particles, their distribution and size of particles. While various data structures exist to handle and manipulate large number of datasets, an efficient solution depends on the problem at hand due to different advantages and disadvantages of these data structures. For instance, for the case of few hundred particles, naive calculation of distances are efficient as opposed to different data strucutures which consume a huge amount of preprocessing time. However, for a relatively large static dataset, it is wise to spend the time in constructing a data-structures to reduce the query times as naive calculations scale with O(N^2) whereas other data structures have query times far lower than the brute force algorithm.</p>
<p>Having discussed the main gist of the project, the remaining blog includes the dissection of the projects into sub-modules which I will be focussing during the course of next few months. Since, at its core, the main problem is the identification of the correct data structure which is robust enough to handle different particle distribution and sizes, the first part of the project is to compare and benchmark different data structures against certain use cases. For the sake of brevity, uniform distribution of point sized particles are selected to be used, unless specified, for the benchmark studies. Two main use cases are identified for distance calculations, which can be listed as follows :</p>
<ul>
<li>Typically selections of atoms form a major class of problems required in the analysis of MD trajectories. An around distance selection statement like <code class="highlighter-rouge">around protein 3.5 </code> selects all the atoms not beloning to the protein within a distance of 3.5 units from the protein particles. As a use case, we consider selection of all the particles within a certain distance from a query point.</li>
<li>Another major class which depends heavily on distance calculations is the evaluation of all possible pairs of particles within a definite distance from each other. While this forms the crux of the problem, similar analysis like contact maps, bond formation also form the subset of this problem.</li>
</ul>
<p>Based on the above description, one might think selection use case as a subset of pair calculations such that similar query operation is executed for every particle. While it is correct, but due to large number of computations involved in pair distance evaluations, a small modification in accessing the particles might lead to a significant change in overall time. For instance, in a hierarchical tree based structure, every query operation needs to traverse the tree to find its position within the tree which itself is time consuming. If the time to find the position can be made independent of the number of particles, this seemingly small modification leads to huge reductions in total time.</p>
<p>Since the solution revolves around different data structures, it is best to list down different data strcutures and their current state. While doing so, it should also be kept in mind that a typical boundary condition frequently used by Molecular dynamics community, Periodic boundary conditions, should be supported by the data structure.</p>
<ul>
<li>Naive Distance - This is not a data structure, per se, but a method to evaluate the distances and is already implemented in MDAnalysis. The algorithmic running time for this is O(N^2) and more so the space complexity also scales with O(N^2) if distance between all the particles are evaluated. Furthermore, the algortihm allows easy implementation of periodic boundary condiions.</li>
<li>Grid Structures - These are linear data structures, where the domain is divided into linear grids. Preprocessing time for this data structure scales as O(N), as every particle is alloted a cell in the grid. For distance calculations, only the particles in the neighbourhood are used for distance evaluations. This specific step reduces a huge amount of time. The query time for such data structure is also dependent on the distribution of particles as opposed to independent naive distance calculations. The query time for nearest neighbour for this class of data structure is independent of the number of particles, provided the particles are arranged in a sorted array. For the case useful for MDAnalysis, the query time for an around selection depends on the number of particles in every cell i.e. depends on the density of particles irrespective of the number. Another factor which affects the performance is the cellsize. Smaller the cell size, higher is the cell checking with very low number of cells actually having a particle whereas large size of cell tends to perform poorly and close to naive distance algorithm. Typically, a cellsize equivalent to cutoff radius is used, which is classified as cell-list algorithm. Various algorithms exist for grid structures. Some of them are Cell-list, Combined cell-list and verlet, Full-list, twin grid, excelfile etc. An excellent review is provided <a href="https://nanopdf.com/download/project-proposal-andrew-noske-homesite_pdf">here</a>. A particular advantage of this method is easy parallellization and independent query time along with linear scaling for construction time. Due to these advantages, cell-list combined with verlet is generally used in Molecular dynamics codes. However, the advantage of verlet lists is pronounced in dynamic data sets, when the atom position donot change much in consecutive data sets. For MDAnalysis, these advantages become less pronounced due to static datasets and therefore this algorithm is not pursued further. Due to the similarity to naive distance calculations, this data structures also allows easy implementation of Periodic boundary conditions.</li>
<li>Tree Structures - Some of the often used tree structures are KDTree, Octree and R-tree. KDtree is a bindary tree which is often used for high dimensional searches. KDtree is a space partitioning data structure for organizing points in a k dimensional space. Typically, these are used for range searches and near neighbour searches in k dimensional space. For the requirements of MDAnalysis, KDtree has been implememented to generate tree structure in 3 dimensions. As opposed to a typical tree, where every leaf node corresponds to a single atom, leaf node in a current implementation can hold a list of particles. This step significantly reduces the time of distance calculations from a query for both fixed radius and nearest neighbour searches. More information on KDtree can be found <a href="https://en.wikipedia.org/wiki/K-d_tree">here</a>. Tree structure usually require huge amount of construction time which scales as O(NlogN) for a bucket size of 1. The time changes with the bucket size as well and typically a bucket size of 10-30 particles is preferred. Query time for KDtree scales with O(logN), where most of the time is spent in traversing the tree to locate the node of query point. Implementation of periodic boundary conditions is little tricky in tree-structures. However, it is already implemented in MDAnalysis. Octree is another structure which is based on recursive division of a cube into its eight octants. This structure shares the same advantages and disadvantages of KDtree. While KDtree is a good choice for manipulation of high dimensional datasets, octree is specifically suited for three dimensional dataset. Due to an organized structure, accessibility of different nodes in faster in this structure. However, since every node has eight childs with each child having equal probability of containing the query point, the query time is dependent on the order of tree traversal. Inspite of the clear advantages of Octree over KDtree, I didn’t come across any implementation of periodic boundary conditions in OCtree. It is more likely that the scaling of Octree will follow a similar trend to the periodic implemetation of KDtree. The construction time for OCtree scales as O(N) and query time scales as O(klogN), where k depends on the number of particles in each cell and total number of bins. More information on Octree can be found <a href="https://en.wikipedia.org/wiki/Octree">here</a>.</li>
</ul>
<p>Regarding the current status of different data structures, the brute force/ naive distance algorithm currelty resides in <code class="highlighter-rouge">MDAnalysis.lib.distances</code>, which is parallellization enabled as well as pbc aware. In terms of grid structure, a pure python implementation is used in GNM module of MDAnalysis, a CellGrid package also offers a python implementation of cell-list structure. Another package of fast distance calculations based on the cell-list data structure resides in FATSLiM package. More information can be found <a href="http://fatslim.github.io/">here</a>. All the mentioned cell-list structures are PBC aware. KDtree is also implemented in MDAanlysis, which is PBC aware and sits in <code class="highlighter-rouge">MDAnalysis.lib.pkdtree</code>. Octree is not implemented to handle periodic boundary conditions (PBC), but a fast implemetation for non-PBC is implemented in Point Cloud Library (<a href="https://github.com/PointCloudLibrary">PCL</a>). Other than these, non PBC KDtree data structures are also implemented in <a href="https://biopython.org/">Biopython</a> and <a href="https://docs.scipy.org/doc/scipy/reference/generated/scipy.spatial.cKDTree.html">Scipy</a>.</p>
<p>The task now is to analyze the methods and their suitable applications. The plan is to provide benchmarks for different use cases and establish the best and robust method for distance evaluations. Although there might not be a single method to serve this purpose, another goal is to assess the limiting conditions when one method becomes superior over the other in terms of performance. The tasks can be dissected into following mini-tasks:</p>
<ol>
<li>Benchmarks for build and query time for single queries.</li>
<li>Benchmarks for pair contacts for bond analysis.</li>
<li>Benchmarks for KDtree to compare Scipy, Biopython and current implementation of Periodic KDtree.</li>
<li>Relevant benchmark for Octree and cython implementation of Nearest neighbour search implemented in FATSliM.</li>
<li>The cellgrid module in MDAnalysis is not optimized in terms of grid size and cutoff distances. The next task is to optimize the Cellgrid API and Benchmark it against the previous cases.</li>
</ol>
<h1 id="benchmark-studies">Benchmark Studies</h1>
<hr />
<p>Here, I present the initial benchmarks, a condensed version of conclusions that can be derived from initial benchmarks is also presented. For all the studies presented below, I am using Intel Core i5-5200 processor clocked at 2.20 GHz and RAM of 8GB. The related ipython notebooks are available <a href="https://github.com/ayushsuhane/Benchmarks_Distance">here</a>. To collect the timings, an inbuilt magic method <code class="highlighter-rouge">timeit</code> is used. For any function timeit can be used as <code class="highlighter-rouge">%timeit function(arguments)</code> in ipython notebook. Similarly, for memory details, <code class="highlighter-rouge">memit</code> magic method is used here.</p>
<ul>
<li>For Build and Single query timings related to different data structures (particularly, <a href="https://github.com/MDAnalysis/cellgrid">Cellgrid</a>, <a href="https://github.com/MDAnalysis/mdanalysis/blob/develop/package/MDAnalysis/lib/pkdtree.py">PeriodicKDTree</a>, <a href="https://docs.scipy.org/doc/scipy/reference/generated/scipy.spatial.cKDTree.html">KDTree Scipy</a>, and pure python implementation of Cell-list) are used <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/SingleQuery.ipynb">here</a> for benchmarks. For single point queries, it can be concluded that both the query and build time are larger for the case of Cell-list case shown here. However, it is anticipated that an efficient implementation of cell-ist should consume less time as opposed to KDtree for single point query.</li>
</ul>
<p><img src="/images/230518_single_bm.PNG" alt="alt text" title="Build and Query time for single queries using Cell-Lists and Periodic KDtree" /></p>
<ul>
<li>For pair contacts, naive distances, KDtree and Cellgrid packages are used for the benchmarks. For the sake of parametric study, a cutoff distance of 10 units is chosen i.e the motive is to time different functions to find all particle pair within a distance of 10 distance units. For naive distances, <code class="highlighter-rouge">self_distance_array</code>, as implemented in MDAnalysis is used which constructs an array and stores the distance of each pari of particles, one at a time. This causes sudden increase in memory usage which is also documented <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/BM_PairContact.ipynb">here</a>. Timings for differetn cutoff radius are also noted to increase with increase in cutoff distances for a fixed number of particles.
The goal, here, was to characterize the CellGrid package against already implemented techniques. It also involves the variation of time due to different cutoff radius. It is observed that smaller cutoff distance consumes a huge amount of time. Since, cellsize is other parameter which controls the time. A more interesting question was whether a significant reduction in timing can be achieved by increasing the cellsize for a given cutoff distance. As a result, density was identified as a parameter and particles per unit cell was optimized for lower calculation time. The difference in time was found to be significant for small number of particles. A particle density of 30 particles per cell is identified as an optimized parameter (see <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/CellGrid_Optimization.ipynb">here</a>). The time for optimized cutoff was also compared against tree structures, but KDtree took less time than the optimized Cellgrid method.</li>
</ul>
<p><img src="/images/230518_cellgrid_opt.PNG" alt="alt text" title="CellGrid optimization with particles per cell as the parameter for cutoff distance of 5 units" />
<img src="/images/230518_cellgrid_opt2.PNG" alt="alt text" title="CellGrid optimization with particles per cell as the parameter for cutoff distance of 8 units" /></p>
<ul>
<li>A use case of selections with different particle type and distribution also needs to be checked. As most of the MDAnalysis community requires selection of one group of particles around another group of particles. Two different geometries i.e. Slab distribution and spherical protein with solvent distribution are chosen for benchmark studies. For spherical protein particle surrounded by a solvent, KDtree aced among the cellgrid and brute force method for large number of particles. (see <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/BM_Selection_SphProtein.ipynb">here</a>.</li>
</ul>
<p><img src="/images/230518_spherical_bm.PNG" alt="alt text" title="Comparison of spherical selection for different cutoff radius" /></p>
<p><img src="/images/230518_spherical_bm2.PNG" alt="alt text" title="Benchmark for selection arounf a spherical protein particle with number of particles" /></p>
<ul>
<li>Single queries are benchmarked against non-PBC and PBC aware methods by choosing the query point deep inside the box such that periodic images are never required for computations in any of the methods. Here, Octree and cython based neighbour search method is also used for benchmarking. The library implemented in FATSLiM showed the fastest build and query time for large number of particles. However, for small number of particles, non-periodic octree evaluated the selections in shortest time.(see <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/OctreeComparison.ipynb">here</a>). Similarly, spherical selection was also evaluated for all the data structures. As anticipated, Octree and Neighbour search in FATSLiM shows the best performance(see <a href="https://github.com/ayushsuhane/Benchmarks_Distance/blob/master/Notebooks/SphericalParticleSelection.ipynb">here</a>).</li>
</ul>
<p><img src="/images/230518_octree_bm.PNG" alt="alt text" title="Single point query benchmarks for PBC and non PBC structures and scaling behavious with number of particles" /></p>
<p><img src="/images/230518_octree_bm2.PNG" alt="alt text" title="Single point benchmarks with different cutoff distance" /></p>
<p>Based on these benchmarks, it can be concluded that Neighbour search (FATSLiM) and Octree are definitely good alternatives to reduce the time of distance calculations. Although it was not possible to evaluate the pair contacts of all the particles using octree, but the scaling can be anticipated to be similar to the periodic KDTree.</p>
<p>All this for now. See you in the next blog post.</p>Using Jupyter Notebooks for the first time, and they are awesome!!!Up and running!2018-04-25T00:00:00+00:002018-04-25T00:00:00+00:00https://ayushsuhane.github.io/Hello-World<p>Only if the gas molecule knew, where the Brownian motion is taking it.</p>
<p>Just completed a setup of the blogpost by forking the amazing and simple <a href="http://github.com/barryclark/jekyll-now/">Jekyll Now</a>. Anyone who wants to setup a minimal blog post quickly, I would highly recommend it.</p>Only if the gas molecule knew, where the Brownian motion is taking it.