<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd"> <html xmlns="http://www.w3.org/1999/xhtml"> <head> <meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/> <meta http-equiv="X-UA-Compatible" content="IE=9"/> <meta name="generator" content="Doxygen 1.8.14"/> <meta name="viewport" content="width=device-width, initial-scale=1"/> <title>dcdfort: libdcdfort</title> <link href="tabs.css" rel="stylesheet" type="text/css"/> <script type="text/javascript" src="jquery.js"></script> <script type="text/javascript" src="dynsections.js"></script> <link href="search/search.css" rel="stylesheet" type="text/css"/> <script type="text/javascript" src="search/searchdata.js"></script> <script type="text/javascript" src="search/search.js"></script> <link href="doxygen.css" rel="stylesheet" type="text/css" /> </head> <body> <div id="top"><!-- do not remove this div, it is closed by doxygen! --> <div id="titlearea"> <table cellspacing="0" cellpadding="0"> <tbody> <tr style="height: 56px;"> <td id="projectalign" style="padding-left: 0.5em;"> <div id="projectname">dcdfort </div> <div id="projectbrief">Modern Fortran library for analyzing DCD trajectories</div> </td> </tr> </tbody> </table> </div> <!-- end header part --> <!-- Generated by Doxygen 1.8.14 --> <script type="text/javascript"> /* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */ var searchBox = new SearchBox("searchBox", "search",false,'Search'); /* @license-end */ </script> <script type="text/javascript" src="menudata.js"></script> <script type="text/javascript" src="menu.js"></script> <script type="text/javascript"> /* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */ $(function() { initMenu('',true,false,'search.php','Search'); $(document).ready(function() { init_search(); }); }); /* @license-end */</script> <div id="main-nav"></div> </div><!-- top --> <!-- window showing the filter options --> <div id="MSearchSelectWindow" onmouseover="return searchBox.OnSearchSelectShow()" onmouseout="return searchBox.OnSearchSelectHide()" onkeydown="return searchBox.OnSearchSelectKey(event)"> </div> <!-- iframe showing the search results (closed by default) --> <div id="MSearchResultsWindow"> <iframe src="javascript:void(0)" frameborder="0" name="MSearchResults" id="MSearchResults"> </iframe> </div> <div class="header"> <div class="headertitle"> <div class="title">libdcdfort </div> </div> </div><!--header--> <div class="contents"> <div class="textblock"><p>Copyright (C) 2017,2018 James W. Barnett</p> <p><a href="http://github.com/wesbarnett/dcdfort">http://github.com/wesbarnett/dcdfort</a></p> <h2>About</h2> <p>Fortran library for natively reading in DCD trajectory files generated from <a href="http://lammps.sandia.gov">LAMMPS</a> simulations for analysis. Uses an object-oriented style. For example, you can simply read in all simulation snapshots from a DCD file by adding the following lines:</p> <div class="fragment"><div class="line">using <a class="code" href="namespacedcdfort__trajectory.html">dcdfort_trajectory</a></div><div class="line"><span class="keywordtype">type</span>(Trajectory) :: trj</div><div class="line"></div><div class="line"><span class="keyword">call </span>trj%read(<span class="stringliteral">"mytrajectoryfile.dcd"</span>)</div></div><!-- fragment --><p>Now all information from the trajectory file (atom coordinates, box dimensions) is accessible via object getters. See the <a href="#api">API</a> below.</p> <p>This is similar to my other project <a href="https://github.com/wesbarnett/libgmxfort">libgmxfort</a>, except that this project does not require any plugins to read the binary trajectory files.</p> <p><b>Note:</b> DCD files generated from simulation packages other than LAMMPS will probably not work with this library. LAMMPS outputs DCD files as 32-bit CHARMM files with a unit cell and three dimensions, which is what this library can read in.</p> <p>Full API documentation is <a href="https://cdn.rawgit.com/wesbarnett/dcdfort/1.1/docs/html/index.html">here</a>.</p> <h2>Compilation</h2> <p>The GNU Fortran compiler version >= 7.1 is required, since non-constant error stop codes are used and we are using the <code>convert="swap"</code> GNU extension when reading in DCD files of different endiness.</p> <p>After cloning the repository, or extracting the release tarball, cd into the repository. Then:</p> <div class="fragment"><div class="line">meson build</div><div class="line">ninja -C build</div></div><!-- fragment --><p>If you want to change some of the compiler optimizations or other setttings, use <code>FFLAGS</code>:</p> <div class="fragment"><div class="line">FFLAGS="-O3" meson build</div><div class="line">ninja -C build</div></div><!-- fragment --><h2>Testing</h2> <p>To test your build, do:</p> <div class="fragment"><div class="line">ninja -C build test</div></div><!-- fragment --><p>If any tests do not pass, please file an issue.</p> <h2>Installation</h2> <p>The following will install the library to the location specified by the meson flag <code>--prefix</code>, which is <code>/usr/local</code> by default.</p> <div class="fragment"><div class="line">ninja -C build install</div></div><!-- fragment --><h2>Usage</h2> <p>Compile your Fortran trajectory analysis program with <code>-ldcdfort</code>. You may also need to use <code>-I</code> to point to where the modules files are even with all of the right environment variables set (by default at <code>/usr/local/include</code>).</p> <h3>pkg-config</h3> <p>A pkg-config file is included, so that it can be used in your program compilations. You may need to set the <code>PKG_CONFIG_PATH</code> environment variable to find the file (by default in the directory <code>/usr/local/lib/pkgconfig</code>). See <code>man 1 pkg-config</code> for more information.</p> <h3>API</h3> <p>Add <code>use <a class="el" href="namespacedcdfort__trajectory.html" title="Module that contains Trajectory type. ">dcdfort_trajectory</a></code> to your Fortran program in order to use the <code>Trajectory</code> class and <code>use <a class="el" href="namespacedcdfort__utils.html" title="Module that contains some useful utilities. ">dcdfort_utils</a></code> in order to use any of the other utilities. There is an example in the <code>example</code> folder on how to do this.</p> <p>Full API documentation is <a href="https://cdn.rawgit.com/wesbarnett/dcdfort/1.1/docs/html/index.html">here</a>.</p> <h4>Reading in trajectory and index files</h4> <p>Typically you will open a trajectory file (and optionally a corresponding GROMACS-style index file). Then you will read in the entire trajectory file at once, or you can read it in in chunks. Then you should close the trajectory file when done.</p> <p>The simplest way to use this library is to construct a <code>Trajectory</code> object and then use the <code>read()</code> method:</p> <div class="fragment"><div class="line"><span class="keywordtype">use </span><a class="code" href="namespacedcdfort__trajectory.html">dcdfort_trajectory</a></div><div class="line"><span class="keywordtype">implicit none</span></div><div class="line"><span class="keywordtype">type</span>(Trajectory) :: trj</div><div class="line"><span class="keyword">call </span>trj%read(<span class="stringliteral">"traj.dcd"</span>)</div></div><!-- fragment --><p>If you have a corresponding index file, you can add a second argument to open:</p> <div class="fragment"><div class="line"><span class="comment">call trj%read("traj.dcd", "index.ndx")</span></div></div><!-- fragment --><p>Now information regarding the index groups is stored in memory and can be used in some of the following methods.</p> <p>The <code>read()</code> method opens the dcd file, reads in all information, and then closes it. The <code>trj</code> object in this example now stores all of the coordinates and information from the .dcd file.</p> <p>To skip the first portion of a trajectory file with <code>read()</code> use the <code>skip</code> argument. The following skips the first <code>100</code> frames before reading them into memory.</p> <div class="fragment"><div class="line"><span class="comment">call trj%read("traj.dcd", "index.ndx", skip=100)</span></div></div><!-- fragment --><p>To only ready in every so many frames, use the <code>nevery</code> argument. The following reads in only every 10th snapshot into memory:</p> <div class="fragment"><div class="line"><span class="comment">call trj%read("traj.dcd", "index.ndx", nevery=10)</span></div></div><!-- fragment --><p>If you want to read in the trajectory file in frame-by-frame use <code>read_next()</code> instead of <code>read()</code>. To use this, you must additionally open and close the dcd file on your own. By default it reads in one frame:</p> <div class="fragment"><div class="line"><span class="keywordtype">integer</span> :: n</div><div class="line"><span class="keyword">call </span>trj%open(<span class="stringliteral">"traj.dcd"</span>, <span class="stringliteral">"index.ndx"</span>)</div><div class="line">n = trj%read_next()</div><div class="line"><span class="keyword">call </span>trj%close()</div></div><!-- fragment --><p>To read in more than one, specify an argument. The following reads in 10 frames:</p> <div class="fragment"><div class="line">n = trj%read_next(10)</div></div><!-- fragment --><p><code>read_next()</code> returns the number of frames actually read in. It is a function, and not a subroutine. This is useful for using it with a <code>do while</code> loop. For example:</p> <div class="fragment"><div class="line"><span class="keywordtype">use </span><a class="code" href="namespacedcdfort__trajectory.html">dcdfort_trajectory</a></div><div class="line"></div><div class="line"><span class="keywordtype">implicit none</span></div><div class="line"></div><div class="line"><span class="keywordtype">type</span>(Trajectory) :: trj</div><div class="line"><span class="keywordtype">integer</span> :: i, n</div><div class="line"></div><div class="line"><span class="keyword">call </span>trj%open(<span class="stringliteral">"traj.dcd"</span>, <span class="stringliteral">"index.ndx"</span>)</div><div class="line"></div><div class="line">n = trj%read_next(10)</div><div class="line"><span class="keywordflow">do</span> <span class="keywordflow">while</span> (n > 0)</div><div class="line"> <span class="keywordflow">do</span> i = 1, n</div><div class="line"> <span class="comment">! do some things with the frames read in</span></div><div class="line"><span class="keywordflow"> end do</span></div><div class="line"> n = trj%read_next(10)</div><div class="line"><span class="keywordflow">end do</span></div><div class="line"></div><div class="line"><span class="keyword">call </span>trj%close()</div></div><!-- fragment --><p>To skip a frame without reading it into memory use <code>skip_next()</code>. You can also pass an integer argument to indicate how many frames to skip. The function returns the actual number of frames skipped (you might be near the end of the file and not able to skip all you specified).</p> <h4>Getting simulation information</h4> <p>After calling <code>read()</code> or <code>read_next()</code> every atom's coordinates are accessible via the <code>x()</code> method. For example, to get the coordinates of the first atom in the first frame you would do the following. The frame is the first argument and the atom number is the second argument.</p> <div class="fragment"><div class="line"><span class="keywordtype">real</span> :: myatom(3)</div><div class="line"><span class="comment">! ...</span></div><div class="line">myatom = trj%x(1, 1)</div></div><!-- fragment --><p><b>Note</b>: Fortran uses one-based indexing, and that convention is retained here.</p> <p>If you read in an index file, you can get atom coordinates in relationship to that. The following gets the fifth atom in index group C in the 10th frame:</p> <div class="fragment"><div class="line">myatom = trj%x(10, 5, <span class="stringliteral">"C"</span>)</div></div><!-- fragment --><p>If the index group does not exist, then an error will be thrown, causing the program to stop.</p> <p><b>Note:</b> If you have more than one group in your index file with the same name, this will simply use the first group with that name. It's best not to repeat group names in your index file. The library will give you a warning if it finds that an index name is duplicated, but the program will continue.</p> <p>Note that when you use <code>x()</code> you will still have to give it the frame number as the first argument even if you only read in one frame with <code>read_next()</code>. You can always get the <em>total</em> number of frames in a trajectory file object with the <code>nframes</code> member:</p> <p>If you want direct access to the object storing a coordinate, do the following use <code>trjframeArray(i)xyz(j,k)</code> where <code>i</code> is the frame number, <code>j</code> are the x, y, and z coordinates (so <code>1</code>, <code>2</code>, and <code>3</code>), and <code>k</code> is the atom number. The <code>x()</code> method is just a convenient way to get this object.</p> <div class="fragment"><div class="line"><span class="keywordtype">integer</span> :: n</div><div class="line"><span class="comment">! ...</span></div><div class="line">n = trj%nframes</div></div><!-- fragment --><p>This is distinct from the number of frames read in using <code>read_next()</code>. The frame number passed to the <code>x()</code> method, and other methods here, is always in relationship to the number of frames read in, not the total number of frames in the file.</p> <p>To get the timestep corresponding with the first saved frame in the trajectory file do:</p> <div class="fragment"><div class="line"><span class="keywordtype">integer</span> :: istart</div><div class="line"><span class="comment">! ...</span></div><div class="line">istart = trj%istart</div></div><!-- fragment --><p>To get the timestep corresponding with the last saved frame in the trajectory file do:</p> <div class="fragment"><div class="line"><span class="keywordtype">integer</span> :: iend</div><div class="line"><span class="comment">! ...</span></div><div class="line">iend = trj%iend</div></div><!-- fragment --><p>To get how often frames were saved in your simulation to this trajectory file use the <code>nevery</code> object. This corresponds with the fifth column in a LAMMPS <code>dump dcd</code> line where you indicated to dump every this many timesteps. It is the column labeled <code>N</code> in the LAMMPS <a href="http://lammps.sandia.gov/doc/dump.html">dump manual page</a>.</p> <div class="fragment"><div class="line"><span class="keywordtype">real(8)</span> :: nevery</div><div class="line"><span class="comment">! ...</span></div><div class="line">nsavc = trj%nevery</div></div><!-- fragment --><p>To get the simulation timestep, use the <code>timestep</code> object. This corresponds to the <code>timestep</code> setting in LAMMPS.</p> <div class="fragment"><div class="line"><span class="keywordtype">real(8)</span> :: timestep</div><div class="line"><span class="comment">! ...</span></div><div class="line">delta = trj%timestep</div></div><!-- fragment --><p><b>Warning:</b> Some programs such as <em>catdcd</em> overwrite time step information. <em>dcdfort</em> outputs this information whenever it opens a file. If you intend on using this information in your analysis program, double check that it is correct. If you are only using LAMMPS output, you shouldn't have to worry about this.</p> <p>You can also get the number of atoms with the <code>natoms()</code> method:</p> <div class="fragment"><div class="line"><span class="keywordtype">integer</span> :: n</div><div class="line"><span class="comment">! ...</span></div><div class="line">n = trj%natoms()</div></div><!-- fragment --><p>If you want to know how many atoms are in an index group include the group name as an argument. In this example the group name is "C":</p> <div class="fragment"><div class="line">n = trj%natoms(<span class="stringliteral">"C"</span>)</div></div><!-- fragment --><p>If that index group does not exist, then the method will simply return 0.</p> <p>To get the box coordinates, use <code>box</code>. The following gets the box of the <code>2</code>nd frame:</p> <div class="fragment"><div class="line"><span class="keywordtype">real(8)</span> :: mybox(6)</div><div class="line"><span class="comment">! ...</span></div><div class="line">mybox = trj%box(2)</div></div><!-- fragment --><p>The first three elements of the array are the lengths of the unit cell. In an orthogonal simulation these are equivalent to the x, y, and z dimensions. With a triclinic box, these are the length of the unit cell vector along the x-axis (A), the length of the unit cell vector in the xy-plane (B), and the length of the unit cell vector in the yz-plane (C). The last three elements are the <b>cosine</b> of the box angles alpha, beta, and gamma. alpha is the angle between B and C, beta is the angle between A and C, and gamma is the angle between A and B.</p> <h4>Reading in specific groups only</h4> <p>As shown above, the most common use of this library is to use <code>read()</code> or <code>read_next()</code> to save all atom locations and then use getters like <code>x()</code> and <code>natoms()</code> to get information about them by specifying an index group as an argument.</p> <p>To save memory, you can save just a specific index group with read():</p> <div class="fragment"><div class="line">trj%read(xtcfile, ndxfile, <span class="stringliteral">"C"</span>)</div></div><!-- fragment --><p>If you do this, you only have access to the group above, and you should never pass an index group name to getters like x(), since only one group is available. If you do specify a group in a getter after already specifying it in <code>read()</code> or <code>read_next()</code>, you will get an error, and the program will stop.</p> <h4>Utilities</h4> <p>There are several functions and subroutines in the <code><a class="el" href="namespacedcdfort__utils.html" title="Module that contains some useful utilities. ">dcdfort_utils</a></code> module, including periodic boundary and distance calculations. Check out the source file for what is available.</p> <h2>License</h2> <p>This project is released under the following license.</p> <div class="fragment"><div class="line">libdcdfort</div><div class="line"></div><div class="line">Copyright (C) 2017,2018 James W. Barnett</div><div class="line"></div><div class="line">This program is free software; you can redistribute it and/or modify it under</div><div class="line">the terms of the GNU General Public License as published by the Free Software</div><div class="line">Foundation; either version 2 of the License, or (at your option) any later</div><div class="line">version.</div><div class="line"></div><div class="line">This program is distributed in the hope that it will be useful, but WITHOUT ANY</div><div class="line">WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A</div><div class="line">PARTICULAR PURPOSE. See the GNU General Public License for more details.</div><div class="line"></div><div class="line">You should have received a copy of the GNU General Public License along with</div><div class="line">this program; if not, write to the Free Software Foundation, Inc., 51 Franklin</div><div class="line">Street, Fifth Floor, Boston, MA 02110-1301 USA.</div></div><!-- fragment --><p>See the file <code>LICENSE</code> for full details. </p> </div></div><!-- contents --> <!-- start footer part --> <hr class="footer"/><address class="footer"><small> Generated by  <a href="http://www.doxygen.org/index.html"> <img class="footer" src="doxygen.png" alt="doxygen"/> </a> 1.8.14 </small></address> </body> </html>