gromacs_py.tools package

Submodules

gromacs_py.tools.monitor module

Collection of function to monitor a simulation in real time.

gromacs_py.tools.monitor.extract_log_dict(func_input_dict, tail_line_num=20)

Read the last lines of a gromacs .log file and return a dictionnary containing time, step and all energetic terms.

gromacs_py.tools.monitor.isnotebook()

Return if the command is launch from a notebook or not Taken from: https://stackoverflow.com/questions/15411967/ how-can-i-check-if-code-is-executed-in-the-ipython-notebook

Example

>>> isnotebook()
False
gromacs_py.tools.monitor.print_log_file(proc, func_input_dict, tail_line_num=20)

Monitor .log file information. The func_input_dict should contains several keys:

  • terms: list of energetic terms to extract, eg. [‘Potential’,

    ‘Temperature’]

  • log: path of the log file (Defined in os_command.run_background())

  • refresh_time: time interval to refresh log extract (default=1.0 s)

Parameters
  • proc (subprocess object) – running subprocess

  • func_input_dict (dict) – dictionnary containing parameters for log extract

  • tail_line_num (int, default=20) – number of line to read at the end of .log file

Example:

>>> TEST_OUT = str(getfixture('tmpdir'))
>>> import sys
>>> sys.path.insert(0, os.path.abspath(os.path.join(MONITOR_LIB_DIR, '../..')))
>>> from gromacs_py import gmx 

...
>>> prot = gmx.GmxSys(name='1y0m', coor_file=TEST_PATH+'/1y0m.pdb')
>>> ###################################
>>> ####   Create the topologie:   ###
>>> ###################################
>>> prot.prepare_top(out_folder=os.path.join(TEST_OUT, 'top_SH3'), vsite='hydrogens') 
pdb2pqr30... --ff CHARMM --ffout CHARMM --keep-chain --titration-state-method=propka --with-ph=7.00 tmp_pdb2pqr.pdb 00_1y0m.pqr
gmx pdb2gmx -f 01_1y0m_good_his.pdb -o 1y0m_pdb2gmx.pdb -p 1y0m_pdb2gmx.top -i 1y0m_posre.itp -water tip3p -ff charmm36-jul2017 -ignh -vsite hydrogens
>>> ######################################
>>> ### Monitor an energy minimisation ###
>>> ######################################
>>> monitor = {'function': print_log_file,           'terms':['Potential'],           'file_check_ext':'log'}
>>> prot.em(out_folder=os.path.join(TEST_OUT, 'em_SH3'), nsteps=50,    constraints='none', create_box_flag=True, monitor=monitor, nstlog=1)    
gmx editconf -f .../top_SH3/1y0m_pdb2gmx.pdb -o .../top_SH3/1y0m_pdb2gmx_box.pdb -bt dodecahedron -d 1.0
gmx grompp -f 1y0m.mdp -c ../top_SH3/1y0m_pdb2gmx_box.pdb -r ../top_SH3/1y0m_pdb2gmx_box.pdb -p ../top_SH3/1y0m_pdb2gmx.top -po out_1y0m.mdp -o 1y0m.tpr -maxwarn 1
gmx mdrun -s 1y0m.tpr -deffnm 1y0m -nt 0 -ntmpi 0 -nsteps -2 -nocopyright...
gromacs_py.tools.monitor.progress_bar(proc, func_input_dict, tail_line_num=20)

Monitor .log file timestep. The func_input_dict should contains several keys:

  • nsteps: Total number of steps during the simulation

  • log: path of the log file (Defined in os_command.run_background())

  • refresh_time: time interval to refresh log extract (default=1.0 s)

Parameters
  • proc (subprocess object) – running subprocess

  • func_input_dict (dict) – dictionnary containing parameters for log extract

  • tail_line_num (int, default=20) – number of line to read at the end of .log file

Example:

>>> TEST_OUT = str(getfixture('tmpdir'))
>>> import sys
>>> sys.path.insert(0, os.path.abspath(os.path.join(MONITOR_LIB_DIR, '../..')))
>>> from gromacs_py import gmx 

...
>>> prot = gmx.GmxSys(name='1y0m', coor_file=TEST_PATH+'/1y0m.pdb')
>>> ###################################
>>> ####   Create the topologie:   ###
>>> ###################################
>>> prot.prepare_top(out_folder=os.path.join(TEST_OUT, 'top_SH3'), vsite='hydrogens') 
pdb2pqr30... --ff CHARMM --ffout CHARMM --keep-chain --titration-state-method=propka --with-ph=7.00 tmp_pdb2pqr.pdb 00_1y0m.pqr
gmx pdb2gmx -f 01_1y0m_good_his.pdb -o 1y0m_pdb2gmx.pdb -p 1y0m_pdb2gmx.top -i 1y0m_posre.itp -water tip3p -ff charmm36-jul2017 -ignh -vsite hydrogens
>>> ######################################
>>> ### Monitor an energy minimisation ###
>>> ######################################
>>> monitor = PROGRESS_BAR
>>> prot.em(out_folder=os.path.join(TEST_OUT, 'em_SH3'), nsteps=50,    constraints='none', create_box_flag=True, monitor=monitor, nstlog=1)    
gmx editconf -f .../top_SH3/1y0m_pdb2gmx.pdb -o .../top_SH3/1y0m_pdb2gmx_box.pdb -bt dodecahedron -d 1.0
gmx grompp -f 1y0m.mdp -c ../top_SH3/1y0m_pdb2gmx_box.pdb -r ../top_SH3/1y0m_pdb2gmx_box.pdb -p ../top_SH3/1y0m_pdb2gmx.top -po out_1y0m.mdp -o 1y0m.tpr -maxwarn 1
gmx mdrun -s 1y0m.tpr -deffnm 1y0m -nt 0 -ntmpi 0 -nsteps -2 -nocopyright
gromacs_py.tools.monitor.read_xvg(xvg_file)

Read a .xvg file and return a pandas dataframe.

Parameters
  • xvg_file (string (Default: 'time')) – path of the xvg file

  • x_axis – name of first column

Example

>>> xvg_file = os.path.join(TEST_PATH, 'volume.xvg')
>>> vol_df = read_xvg(xvg_file)
>>> vol_df.head()
   Time (ps)      Volume
0        0.0  171.237213
1        5.0  135.081039
2       10.0   94.224792
3       15.0   59.942383
4       20.0   58.125397
gromacs_py.tools.monitor.simulation_plot(proc, func_input_dict, refresh_time=1.0)

This function is used for monitoring a simulation in real time. Function can be excecuted by the gromacs.tools.os_command.run_background() function. The function monitors a trajectory file, and launch the analysis if the file has been modified. It can plot as function of time an analysis of a simulation. Analysis is passed as input function.

Warning

Need to add the following lines to be run in jupyter notebook:

  • %matplotlib notebook

Example:

>>> TEST_OUT = str(getfixture('tmpdir'))
>>> import sys
>>> # print(os.path.abspath(os.path.join(MONITOR_LIB_DIR, '../..')))
>>> sys.path.insert(0, os.path.abspath(os.path.join(MONITOR_LIB_DIR, '../..')))
>>> from gromacs_py import gmx 
>>> prot = gmx.GmxSys(name='1y0m', coor_file=TEST_PATH+'/1y0m.pdb')
>>> ###################################
>>> ####   Create the topologie:   ###
>>> ###################################
>>> prot.prepare_top(out_folder=os.path.join(TEST_OUT, 'top_SH3'), vsite='hydrogens') 
pdb2pqr30... --ff CHARMM --ffout CHARMM --keep-chain --titration-state-method=propka --with-ph=7.00 tmp_pdb2pqr.pdb 00_1y0m.pqr
gmx pdb2gmx -f 01_1y0m_good_his.pdb -o 1y0m_pdb2gmx.pdb -p 1y0m_pdb2gmx.top -i 1y0m_posre.itp -water tip3p -ff charmm36-jul2017 -ignh -vsite hydrogens
>>> ######################################
>>> ### Monitor an energy minimisation ###
>>> ######################################
>>> monitor = {'function': simulation_plot,           'extract_func': [{'func': extract_log_dict,                             'term': 'Potential'},                           {'func': extract_log_dict,                             'term': 'Temperature'}],           'file_check_ext':'log'}
>>> prot.em(out_folder=os.path.join(TEST_OUT, 'em_SH3'), nsteps=50,    constraints='none', create_box_flag=True, monitor=monitor, nstlog=10)    
gmx editconf -f ...top_SH3/1y0m_pdb2gmx.pdb -o ...top_SH3/1y0m_pdb2gmx_box.pdb -bt dodecahedron -d 1.0
gmx grompp -f 1y0m.mdp -c .../top_SH3/1y0m_pdb2gmx_box.pdb -r .../top_SH3/1y0m_pdb2gmx_box.pdb -p .../top_SH3/1y0m_pdb2gmx.top -po out_1y0m.mdp -o 1y0m.tpr -maxwarn 1
gmx mdrun -s 1y0m.tpr -deffnm 1y0m -nt 0 -ntmpi 0 -nsteps -2 -nocopyright

Module contents