OpenSMOG

The OpenSMOG classes perform molecular dynamics using Structure-Based Models (SBM) for biomolecular simulations. OpenSMOG uses force fields generated by SMOG 2, and it allows the simulations of a wide variety of potential forms, including commonly employed C-alpha and all-atom variants. Details about the default models in SMOG 2 can be found in the following resources:

  • SMOG server: https://smog-server.org/smog2/

  • C-alpha: Clementi, C., Nymeyer, H. and Onuchic, J.N., 2000. Topological and energetic factors: what determines the structural details of the transition state ensemble and “en-route” intermediates for protein folding? An investigation for small globular proteins. Journal of molecular biology, 298(5), pp.937-953.

  • All-Atom: Whitford, P.C., Noel, J.K., Gosavi, S., Schug, A., Sanbonmatsu, K.Y. and Onuchic, J.N., 2009. An all‐atom structure‐based potential for proteins: bridging minimal models with all‐atom empirical forcefields. Proteins: Structure, Function, and Bioinformatics, 75(2), pp.430-441.

class OpenSMOG.SBM(time_step, collision_rate, r_cutoff, temperature, name='OpenSMOG')[source]

Bases: object

The SBM class performs Molecular dynamics simulations using structure-based models to investigate a broad range of biomolecular dynamics, including domain rearrangements in proteins, folding and ligand binding in RNA, and large-scale rearrangements in ribonucleoprotein assemblies. In its simplest form, a structure-based model defines a particular structure (usually obtained from X-ray, or NMR, methods) as the energetic global minimum.

The SBM sets the environment to start the molecular dynamics simulations.

Parameters
  • time_step (float, required) – Simulation time step in units of \(\tau\).

  • collision_rate (float, required) – Friction/Damping constant in units of reciprocal time (\(1/\tau\)).

  • r_cutoff (float, required) – Cutoff distance to consider non-bonded interactions in units of nanometers.

  • temperature (float, required) – Temperature in reduced units.

  • name (str) – Name used in the output files. (Default value: OpenSMOG).

createReporters(trajectory=True, trajectoryName=None, energies=True, energiesName=None, energy_components=False, energy_componentsName=None, interval=1000)[source]

Creates the reporters to provide the output data.

Parameters
  • trajectory (bool, optional) – Whether to save the trajectory .dcd file containing the position of the atoms as a function of time. (Default value: True).

  • energies (bool, optional) – Whether to save the energies in a .txt file containing five columns, comma-delimited. The header of the files shows the information of each collum: #”Step”,”Potential Energy (kJ/mole)”,”Kinetic Energy (kJ/mole)”,”Total Energy (kJ/mole)”,”Temperature (K)”. (Default value: True).

  • forces (bool, optional) – Whether to save the potential energy for each applied force in a .txt file containing several columns, comma-delimited. The header of the files shows the information of each column. An example of the header is: #”Step”,”eletrostastic”,”Non-Contacts”,”Bonds”,”Angles”,”Dihedrals”,”contact_1-10-12”. (Default value: False).

  • interval (int, required) – Frequency to write the data to the output files. (Default value: 10**3)

createSimulation()[source]

Creates the simulation context and loads into the OpenMM platform.

loadGro(Grofile)[source]

Loads the .gro file format in the OpenMM system platform. The input files are generated using SMOG2 software with the flag -OpenSMOG. Details on how to create the files can be found in the SMOG2 User Manual. A tutorial on how to generate the inputs files for the default all-atom and C-alpha models can be found here.

Parameters

Grofile (file, required) – Initial structure for the MD simulations in .gro file format generated by SMOG2 software with the flag -OpenSMOG. (Default value: None).

loadSystem(Grofile, Topfile, Xmlfile)[source]

Loads the input files in the OpenMM system platform. The input files are generated using SMOG2 software with the flag -OpenSMOG. Details on how to create the files can be found in the SMOG2 User Manual. A tutorial on how to generate the inputs files for default all-atom and C-alpha models can be found here.

Parameters
  • Grofile (file, required) – Initial structure for the MD simulations in .gro file format generated by SMOG2 software with the flag -OpenSMOG. (Default value: None).

  • Topfile (file, required) – Topology .top file format generated by SMOG2 software with the flag -OpenSMOG. The topology file lists the interactions between the system atoms except for the “Native Contacts” potential that is provided to OpenSMOG in a .xml file. (Default value: None).

  • Xmlfile (file, required) – The .xml file that contains the all information that defines the “Contact” potential. The .xml file is generated by SMOG2 software with the flag -OpenSMOG, which support custom potential functions. (Default value: None).

loadTop(Topfile)[source]

Loads the .top file format in the OpenMM system platform. The input files are generated using SMOG2 software with the flag -OpenSMOG. Details on how to create the files can be found in the SMOG2 User Manual. A tutorial on how to generate the inputs files for the default all-atom and C-alpha models can be found here.

Parameters

Topfile (file, required) – Topology .top file format generated by SMOG2 software with the flag -OpenSMOG. The topology file defines the interactions between atoms, except for the “Native Contacts” potential that is provided to OpenSMOG in the form of a .xml file. (Default value: None).

loadXml(Xmlfile)[source]

Loads the .xml file format in the OpenMM system platform. The input files are generated using SMOG2 software with the flag -OpenSMOG. Details on how to create the files can be found in the SMOG2 User Manual. A tutorial on how to generate the inputs files for default all-atom and C-alpha models can be found here.

Parameters

Xmlfile (file, required) – The .xml file that contains all information that defines the “Contact” potentials. The .xml file is generated by SMOG2 software with the flag -OpenSMOG, which support custom potential energy functions. (Default value: None).

printHeader()[source]
run(nsteps, report=True, interval=10000)[source]

Run the molecular dynamics simulation.

Parameters
  • nsteps (int, required) – Number of steps to be performed in the simulation. (Default value: 10**7)

  • report (bool, optional) – Whether to print the simulation progress. (Default value: True).

  • interval (int, required) – Frequency to print the simulation progress. (Default value: 10**4)

saveFolder(folder)[source]

Sets the folder path to save data.

Parameters

folder (str, optional) – Folder path to save the simulation data. If the folder path does not exist, the function will create the directory.

setup_openmm(platform='opencl', precision='single', GPUindex='default', integrator='langevin')[source]

Sets up the parameters of the simulation OpenMM platform.

Parameters
  • platform (str, optional) – Platform to use in the simulations. Opitions are CUDA, OpenCL, HIP, CPU, Reference. (Default value: OpenCL).

  • precision (str, optional) – Numerical precision type of the simulation. Options are single, mixed, double. (Default value: single). For details check the OpenMM Documentation.

  • GPUIndex (str, optional) – Set of Platform device index IDs. Ex: 0,1,2 for the system to use the devices 0, 1 and 2. (Use only when GPU != default).

  • integrator (str) – Integrator to use in the simulations. Options are langevin, variableLangevin, verlet, variableVerlet and, brownian. (Default value: langevin).