Previous topic

Gromacs parameter MDP file format

Next topic

Gromacs topology file (ITP) parser

This Page

Gromacs NDX index file format

The .ndx file contains lists of atom indices that are grouped in sections by group names. The classes NDX and uniqueNDX can parse such ndx files and provide convenient access to the individual groups.

class gromacs.fileformats.ndx.NDX(filename=None, **kwargs)

Gromacs index file.

Represented as a ordered dict where the keys are index group names and values are numpy arrays of atom numbers.

Use the NDX.read() and NDX.write() methods for I/O. Access groups by name via the NDX.get() and NDX.set() methods.

Alternatively, simply treat the NDX instance as a dictionary. Setting a key automatically transforms the new value into a integer 1D numpy array (not a set, as would be the make_ndx behaviour).

Note

The index entries themselves are ordered and can contain duplicates so that output from NDX can be easily used for g_dih and friends. If you need set-like behaviour you will have do use gromacs.formats.uniqueNDX or gromacs.cbook.IndexBuilder (which uses make_ndx throughout).

Example

Read index file, make new group and write to disk:

ndx = NDX()
ndx.read('system.ndx')
print ndx['Protein']
ndx['my_group'] = [2, 4, 1, 5]   # add new group
ndx.write('new.ndx')

Or quicker (replacing the input file system.ndx):

ndx = NDX('system')          # suffix .ndx is automatically added
ndx['chi1'] = [2, 7, 8, 10]
ndx.write()
format = '%6d'

standard ndx file format: ‘%6d’

get(name)

Return index array for index group name.

groups

Return a list of all groups.

ncol = 15

standard ndx file format: 15 columns

ndxlist

Return a list of groups in the same format as gromacs.cbook.get_ndx_groups().

Format:
[ {‘name’: group_name, ‘natoms’: number_atoms, ‘nr’: # group_number}, ....]
read(filename=None)

Read and parse index file filename.

set(name, value)

Set or add group name as a 1D numpy array.

size(name)

Return number of entries for group name.

sizes

Return a dict with group names and number of entries,

write(filename=None, ncol=15, format='%6d')

Write index file to filename (or overwrite the file that the index was read from)

class gromacs.fileformats.ndx.uniqueNDX(filename=None, **kwargs)

Index that behaves like make_ndx, i.e. entries behaves as sets, not lists.

The index lists behave like sets: - adding sets with ‘+’ is equivalent to a logical OR: x + y == “x | y” - subtraction ‘-‘ is AND: x - y == “x & y” - see join() for ORing multiple groups (x+y+z+...)

Example

I = uniqueNDX('system.ndx')
I['SOLVENT'] = I['SOL'] + I['NA+'] + I['CL-']
join(*groupnames)

Return an index group that contains atoms from all groupnames.

The method will silently ignore any groups that are not in the index.

Example

Always make a solvent group from water and ions, even if not all ions are present in all simulations:

I['SOLVENT'] = I.join('SOL', 'NA+', 'K+', 'CL-')
class gromacs.fileformats.ndx.IndexSet

set which defines ‘+’ as union (OR) and ‘-‘ as intersection (AND).