gromacs - parser for Gromacs files

This Python module allows one to use GROMACS data files as the input to an ESPResSo++ simulation, set interactions for given particle types and convert GROMACS potential tables into ESPResSo++ tables. It containts functions: read(), setInteractions(), convertTable()

Some tips for using the gromacs parser:

Tip 1.

topol.top includes solvent via #include statements

If the included .itp file ONLY contains the solvent molecule you’re using (e.g. spc/e water using spce.itp) then this is okay.

But if the .itp file contains info about many molecules (e.g. you want to use one ion from ions.itp), then gromacs.py will just take the first one listed. You must edit your topol.top file to explicitly include the solvent molecule you’re using.

e.g. replace:

; Include topology for ions
#include "amber03.ff/ions.itp"

by:

; Include topology for ions
[ moleculetype ]
; molname       nrexcl
CL              1

[ atoms ]
; id    at type         res nr  residu name     at name  cg nr  charge
1       Cl              1       CL              CL       1      -1.00000

Tip 2. impropers

impropers in the topol.top file (function type 4) need to be labelled ‘[ impropers ]’, not ‘[ dihedrals ]’ as in standard gromacs format”

Also, the dihedrals should be listed before the impropers (this is usuall the case by default in gromacs-format files).

Tip 3.

For rigid SPC/E water using Settle, spce.itp file should look like this:

[ moleculetype ]
; molname       nrexcl
SOL             2

[ atoms ]
; id  at type     res nr  res name  at name  cg nr  charge    mass
  1   OW_spc      1       SOL       OW       1      -0.8476   15.99940
  2   HW_spc      1       SOL       HW1      1       0.4238    1.00800
  3   HW_spc      1       SOL       HW2      1       0.4238    1.00800

[ bonds ]
; i     j       funct   length  force.c.
1       2       1       0.1     345000  0.1     345000
1       3       1       0.1     345000  0.1     345000

[ angles ]
; i     j       k       funct   angle   force.c.
2       1       3       1       109.47  383     109.47  383

The bonds section is used to generate exclusions, but bond and angle parameters are not relevant if the Settle extension is used. The geometry is that specified in the python script when adding the Settle extension

Include modified spce file in topol.top, e.g. replace

#include “amber03.ff/spce.itp”

by

#include “amber03.ff/spce-for-espressopp.itp”

Tip 4.

Use absolute paths for any include files which are not in the standard gromacs topology directory ($GMXLIB)

e.g. replace

#include “mynewresidue.itp”

by

#include “path/to/mynewres/file/mynewresidue.itp”

Tip 5.

The parser won’t work if the particles ids in the include files conflict with the particle ids in the topol.top file itself, and the bonded interaction parameters in the itp file need to be looked up via particle type in the standard gromacs topology directory ($GMXLIB)

i.e. Okay for an itp file like spce.itp above, where the bonds and angles parameters are given in the itp file, as in:

[ bonds ]
; i     j       funct   length  force.c.
1       2       1       0.1     345000  0.1     345000

Not okay for an itp file containing lines like:

[ bonds ]
; i     j       funct   length  force.c.
1       2       1
espressopp.tools.gromacs.convertTable(gro_in_file, esp_out_file, sigma=1.0, epsilon=1.0, c6=1.0, c12=1.0)

Convert GROMACS tabulated file into ESPResSo++ tabulated file (new file is created). First column of input file can be either distance or angle. For non-bonded files, c6 and c12 can be provided. Default value for sigma, epsilon, c6 and c12 is 1.0. Electrostatics are not taken into account (f and fd columns).

Keyword arguments: gro_in_file – the GROMACS tabulated file name (bonded, nonbonded, angle or dihedral). esp_out_file – filename of the ESPResSo++ tabulated file to be written. sigma – optional, depending on whether you want to convert units or not. epsilon – optional, depending on whether you want to convert units or not. c6 – optional c12 – optional

espressopp.tools.gromacs.read(gro_file, top_file='', doRegularExcl=True)

Read GROMACS data files.

Arguments: :param gro_file: – contains coordinates of all particles, the number of particles, velocities and box size. :type gro_file: string :param top_file: – contains topology information. Included topology files (.itp) are also read :type gro_file: string :param doRegularExcl: – if True, exclusions are generated automatically based on the nregxcl parameter (see gromacs manual) :type doRegularExcl: bool

espressopp.tools.gromacs.setLennardJones14Interactions(system, defaults, atomtypeparams, onefourlist, cutoff)

Set lennard jones interactions which were read from gromacs based on the atomypes

espressopp.tools.gromacs.setLennardJonesInteractions(system, defaults, atomtypeparams, verletlist, cutoff, hadress=False, adress=False, ftpl=None)

Set lennard jones interactions which were read from gromacs based on the atomypes

espressopp.tools.gromacs.setLennardJonesInteractionsTI(system, defaults, atomtypeparams, verletlist, cutoff, epsilonB, sigmaSC, alphaSC, powerSC, lambdaTI, pidlist, annihilate=True, hadress=False, adress=False, ftpl=None)

Set lennard jones interactions which were read from gromacs based on the atomypes

espressopp.tools.gromacs.setTabulatedInteractions(potentials, particleTypes, system, interaction)

Set interactions for all given particle types. Return value is a system with all interactions added.

Keyword arguments: potentials – is a dictionary where key is a string composed of two particle types and value is a potential. example: {“A_A”:potAA, “A_B”:potAB, “B_B”:potBB} particleTypes – is a dictionary where key is the particle type, and value is a list of particles of that type. example: {“A”:[“A1m”, “A2m”],”B”:[“B1u”,”B2u”]} system – is the system to which the interaction will be added interaction – is the interaction to which to add the potentials