Gromos Tutorial Pipeline - on a low level

NOTE This is the old way on how to do the Tutorial for Gromos! We want to remove this for the future, please checkout the new Tutorial!

[1]:
import os, sys
from pygromos.utils import bash
root_dir = os.getcwd()

#if package is not installed and path not set correct - this helps you out :)
sys.path.append(root_dir+"/..")


import pygromos
from pygromos.gromos.gromosPP import GromosPP
from pygromos.gromos.gromosXX import GromosXX
from pygromos.files.coord.cnf import Cnf
from pygromos.files.trajectory.trc import Trc

gromosPP_bin = None
gromosXX_bin = None
gromPP = GromosPP(gromosPP_bin)
gromXX = GromosXX(gromosXX_bin)

project_dir = os.path.abspath("../example_files/Tutorial_System")
input_dir = project_dir+"/input"

Build initial files

generate Topology

build single topologies

[2]:
from pygromos.data.ff import Gromos54A7

topo_dir = bash.make_folder(project_dir+'/a_topo')

## Make Cl-
sequence = "CL-"
solvent = "H2O"
top_cl = topo_dir+"/cl.top"

gromPP.make_top(in_building_block_lib_path=Gromos54A7.mtb,
                in_parameter_lib_path=Gromos54A7.ifp,
                in_sequence=sequence, in_solvent=solvent,out_top_path=top_cl)

## Make Peptide
sequence = "NH3+ VAL TYR ARG LYSH GLN COO-"
solvent = "H2O"
top_peptide = topo_dir+"/peptide.top"

gromPP.make_top(in_building_block_lib_path=Gromos54A7.mtb,
                in_parameter_lib_path=Gromos54A7.ifp,
                in_sequence=sequence,
                in_solvent=solvent,out_top_path=top_peptide)

[2]:
'/home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/a_topo/peptide.top'

combine topology

[3]:
top_system = topo_dir+"/vac_sys.top"
gromPP.com_top(in_topo_paths=[top_peptide, top_cl], topo_multiplier=[1,2], out_top_path=top_system)

[3]:
'/home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/a_topo/vac_sys.top'

generate coordinates

[4]:
coord_dir = bash.make_folder(project_dir+"/b_coord")
in_pdb = input_dir+"/peptide.pdb"
cnf_peptide = coord_dir+"/cnf_vacuum_peptide.cnf"

cnf_peptide = gromPP.pdb2gromos(in_pdb_path=in_pdb, in_top_path=top_peptide, out_cnf_path=cnf_peptide)
Cnf(cnf_peptide).view

add hydrogens

[5]:
cnf_hpeptide = coord_dir+"/vacuum_hpeptide.cnf"
cnf_hpeptide = gromPP.add_hydrogens(in_cnf_path=cnf_peptide, in_top_path=top_peptide, out_cnf_path=cnf_hpeptide)

Cnf(cnf_hpeptide).view

cnf to pdb

[6]:
out_pdb = coord_dir+"/vacuum_hpeptide.pdb"
out_pdb = gromPP.frameout(in_coord_path=cnf_hpeptide, in_top_path=top_peptide, out_file_path=out_pdb,
                               periodic_boundary_condition="v", out_file_format="pdb", time=0)

energy minimization - Vacuum

[7]:
from pygromos.data.simulation_parameters_templates import template_emin_vac
from pygromos.files.gromos_system import gromos_system

out_prefix = "vacuum_emin"
vacuum_emin_dir = bash.make_folder(project_dir+"/c_"+out_prefix)
os.chdir(vacuum_emin_dir)

grom_system = gromos_system.Gromos_System(work_folder=vacuum_emin_dir,
                                          system_name="in_"+out_prefix,
                                          in_top_path=top_peptide,
                                          in_cnf_path=cnf_hpeptide,
                                          in_imd_path=template_emin_vac,
                                         verbose=False)

grom_system.adapt_imd()
#del grom_system.imd.POSITIONRES
grom_system.imd.BOUNDCOND.NTB = 0
grom_system.write_files()


out_emin_vacuum = vacuum_emin_dir + "/" + out_prefix
gromXX.md_run(in_imd_path=grom_system.imd.path,
              in_topo_path=grom_system.top.path,
              in_coord_path=grom_system.cnf.path,
              out_prefix=out_emin_vacuum, verbose=True)
cnf_emin_vacuum = out_emin_vacuum+".cnf"
cnf_emin_vacuum

COMMAND:  export OMP_NUM_THREADS=1  &&  md @topo /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/a_topo/peptide.top @conf /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/b_coord/vacuum_hpeptide.cnf @input /home/mlehner/PyGromosTools/pygromos/data/simulation_parameters_templates/vac_emin.imd @fin /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/c_vacuum_emin/vacuum_emin.cnf > /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/c_vacuum_emin/vacuum_emin.omd

[7]:
'/home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/c_vacuum_emin/vacuum_emin.cnf'
[8]:
from pygromos.simulations.modules.preset_simulation_modules import emin

emin_gromos_system = emin(in_gromos_system=grom_system)
emin_gromos_system

################################################################################

emin
################################################################################

============================================================
FOUND RESULT: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/c_vacuum_emin/emin/analysis/data/emin.cnf
 GOING TO SKIPT THIS SUBMISSION!
[8]:

GROMOS SYSTEM: emin
################################################################################
WORKDIR: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/c_vacuum_emin
LAST CHECKPOINT: None

GromosXX_bin: None
GromosPP_bin: None
FILES:
        imd: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/c_vacuum_emin/emin/input/emin.imd
        top: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/c_vacuum_emin/emin/input/emin.top
        cnf: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/c_vacuum_emin/emin/analysis/data/emin.cnf
FUTURE PROMISE: False
SYSTEM:
        PROTEIN:        protein  nresidues: 5 natoms: 71
        LIGANDS:        []  resID: []  natoms: 0
        Non-LIGANDS:    []  nmolecules: 0  natoms: 0
        SOLVENT:        []  nmolecules: 0  natoms: 0


[9]:
emin_gromos_system.cnf.view

Solvatistation and Solvent Energy Minimization

build box system

[10]:
from pygromos.data.solvent_coordinates import spc
out_prefix = "box"
box_dir = bash.make_folder(project_dir+"/d_"+out_prefix)

cnf_box = gromPP.sim_box(in_top_path=top_peptide, in_cnf_path=cnf_emin_vacuum,in_solvent_cnf_file_path=spc,
                          out_cnf_path=box_dir+"/"+out_prefix+".cnf",
                         periodic_boundary_condition="r", minwall=0.8, threshold=0.23, rotate=True)

Cnf(cnf_box).view

to pdb

[11]:
out_pdb = box_dir+"/"+out_prefix+".pdb"

out_pdb = gromPP.frameout(in_coord_path=cnf_box, in_top_path=top_peptide, out_file_path=out_pdb,
                               periodic_boundary_condition="r", out_file_format="pdb", include="ALL", time=0)

Add Ions

[12]:
out_prefix = "ion"
cnf_ion = gromPP.ion(in_cnf_path=cnf_box,
                     in_top_path=top_peptide,
                     out_cnf_path=box_dir+"/"+out_prefix+".cnf",
                     negative=[2, "CL-"],verbose=True           )

Cnf(cnf_ion).view
ion @topo /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/a_topo/peptide.top @pos /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/d_box/box.cnf @pbc v @potential 0.8 @mindist 0.8 @negative 2 CL-

Energy Minimization BOX

[13]:
from pygromos.data.simulation_parameters_templates import template_emin
from pygromos.files.gromos_system import gromos_system

out_prefix = "box_emin"
box_emin_dir = bash.make_folder(project_dir+"/e_"+out_prefix)
os.chdir(box_emin_dir)

grom_system = gromos_system.Gromos_System(work_folder=box_emin_dir,
                                          system_name="in_"+out_prefix,
                                          in_top_path=top_system,
                                          in_cnf_path=cnf_ion,
                                          in_imd_path=template_emin)

grom_system.adapt_imd()
grom_system.imd.STEP.NSTLIM = 3000
grom_system.imd.PRINTOUT.NTPR = 300
grom_system.write_files()

cnf_reference_position = grom_system.cnf.write_refpos(box_emin_dir+"/"+out_prefix+"_refpos.rpf")
cnf_position_restraint = grom_system.cnf.write_possrespec(box_emin_dir+"/"+out_prefix+"_posres.pos", residues=list(filter(lambda x: x != "SOLV", grom_system.cnf.get_residues())))

out_emin_box = box_emin_dir + "/" + out_prefix
gromXX.md_run(in_imd_path=grom_system.imd.path,
              in_topo_path=grom_system.top.path,
              in_coord_path=grom_system.cnf.path,
              in_refpos_path=cnf_reference_position,
              in_posresspec_path=cnf_position_restraint,
              out_prefix=out_emin_box, verbose=True)

cnf_emin_box =out_emin_box+".cnf"
cnf_emin_box = gromPP.frameout(in_coord_path=cnf_emin_box, in_top_path=top_system, out_file_path=cnf_emin_box,
                               periodic_boundary_condition="r cog", out_file_format="cnf", include="ALL", time=0)

POSRES /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/e_box_emin/box_emin_posres.pos
COMMAND:  export OMP_NUM_THREADS=1  &&  md @topo /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/a_topo/vac_sys.top @conf /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/d_box/ion.cnf @input /home/mlehner/PyGromosTools/pygromos/data/simulation_parameters_templates/emin.imd @posresspec /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/e_box_emin/box_emin_posres.pos @refpos /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/e_box_emin/box_emin_refpos.rpf @fin /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/e_box_emin/box_emin.cnf > /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/e_box_emin/box_emin.omd

[14]:
Cnf(cnf_emin_box).view

Simulation

Equilibration NVP

[15]:
from pygromos.data.simulation_parameters_templates import template_md_tut as template_md
from pygromos.files.gromos_system import gromos_system

out_prefix = "eq_NVP"
eq_NVP_dir = bash.make_folder(project_dir+"/f_"+out_prefix)
os.chdir(eq_NVP_dir)

grom_system = gromos_system.Gromos_System(work_folder=eq_NVP_dir,
                                          system_name="in_"+out_prefix,
                                          in_top_path=top_system,
                                          in_cnf_path=cnf_emin_box,
                                          in_imd_path=template_md)

grom_system.adapt_imd(not_ligand_residues="CL-")
grom_system.imd.STEP.NSTLIM = 1000
grom_system.imd.WRITETRAJ.NTWX = 10
grom_system.imd.WRITETRAJ.NTWE = 10
grom_system.imd.INITIALISE.NTIVEL = 1
grom_system.imd.INITIALISE.NTISHK = 1
grom_system.imd.INITIALISE.NTISHI = 1
grom_system.imd.INITIALISE.NTIRTC = 1

grom_system.imd.randomize_seed()
grom_system.rebase_files()
grom_system.write_files()

out_eq_NVP = eq_NVP_dir + "/" + out_prefix
gromXX.md_run(in_imd_path=grom_system.imd.path,
              in_topo_path=grom_system.top.path,
              in_coord_path=grom_system.cnf.path,
              out_tre=True, out_trc=True,
              out_prefix=out_eq_NVP)

cnf_eq_NVP = out_eq_NVP+".cnf"
cnf_eq_NVP

[15]:
'/home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/f_eq_NVP/eq_NVP.cnf'
[16]:
cnf = Cnf(cnf_eq_NVP)
trc = Trc(traj_path=out_eq_NVP+".trc", in_cnf=cnf)
trc.view

MD NVP

[17]:
grom_system
[17]:

GROMOS SYSTEM: in_eq_NVP
################################################################################
WORKDIR: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/f_eq_NVP
LAST CHECKPOINT: None

GromosXX_bin: None
GromosPP_bin: None
FILES:
        imd: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/f_eq_NVP/in_eq_NVP.imd
        top: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/f_eq_NVP/in_eq_NVP.top
        cnf: /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/f_eq_NVP/in_eq_NVP.cnf
FUTURE PROMISE: False
SYSTEM:
        PROTEIN:        protein  nresidues: 5 natoms: 71
        LIGANDS:        []  resID: []  natoms: 0
        Non-LIGANDS:    ['CL-']  nmolecules: 0  natoms: 2
        SOLVENT:        SOLV  nmolecules: 930  natoms: 2790


[18]:
from pygromos.data.simulation_parameters_templates import template_md
from pygromos.files.gromos_system import gromos_system

out_prefix = "md"
md_dir = bash.make_folder(project_dir+"/g_"+out_prefix)
os.chdir(md_dir)

grom_system = gromos_system.Gromos_System(work_folder=md_dir,
                                          system_name="in_"+out_prefix,
                                          in_top_path=top_system,
                                          in_cnf_path=cnf_eq_NVP,
                                          in_imd_path=template_md)

grom_system.adapt_imd(not_ligand_residues="CL-")
grom_system.imd.STEP.NSTLIM = 1000
grom_system.imd.WRITETRAJ.NTWX = 10
grom_system.imd.WRITETRAJ.NTWE = 10
grom_system.imd.INITIALISE.NTIVEL = 0

grom_system.rebase_files()
grom_system.write_files()

out_md = md_dir + "/" + out_prefix
gromXX.md_run(in_imd_path=grom_system.imd.path,
              in_topo_path=grom_system.top.path,
              in_coord_path=grom_system.cnf.path,
              out_tre=True, out_trc=True,
              out_prefix=out_md, verbose=True)

cnf_md = out_md+".cnf"
cnf_md
COMMAND:  export OMP_NUM_THREADS=1  &&  md @topo /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/g_md/in_md.top @conf /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/g_md/in_md.cnf @input /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/g_md/in_md.imd @fin /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/g_md/md.cnf @trc /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/g_md/md.trc @tre /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/g_md/md.tre > /home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/g_md/md.omd

[18]:
'/home/mlehner/PyGromosTools/docs/sphinx_project/Examples/example_files/Tutorial_System/g_md/md.cnf'

Analysis

[19]:
cnf = Cnf(cnf_md)
cnf.view

#shows pbc :)
[20]:
trc = Trc(traj_path=out_md+".trc", in_cnf=cnf)
trc.TITLE = "\n".join(trc.TITLE)
trc.image_molecules()
trc.view
[21]:
out_prefix = "ana"
md_dir = bash.make_folder(project_dir+"/h_"+out_prefix)
os.chdir(md_dir)
[ ]: