CGModel class

class cg_openmm.cg_model.cgmodel.CGModel(particle_type_list={}, charges={}, monomer_types=None, sequence=None, positions=None, bond_lengths={}, bond_force_constants={}, bond_angle_force_constants={}, torsion_force_constants={}, torsion_periodicities={}, equil_bond_angles={}, torsion_phase_angles={}, binary_interaction_parameters={}, go_model=False, go_repulsive_epsilon=None, constrain_bonds=False, include_nonbonded_forces=True, include_bond_forces=True, include_bond_angle_forces=True, include_torsion_forces=True, angle_style='harmonic', nonbond_repulsive_exp=12, nonbond_attractive_exp=6, exclusions={}, hbonds={}, rosetta_functional_form=False, check_energy_conservation=True, use_structure_library=False, random_positions=False, system=None, topology=None, simulation=None)[source]

Coarse-grained model class object containing:

  • particle and residue definitions

  • monomer sequence

  • bonded force field parameters

  • nonbonded force field parameters

  • initial particle positions

Example

from openmm import unit
from cg_openmm.cg_model.cgmodel import CGModel

# Specify backbone (bb) and sidechain (sc) particle parameters:
sigma = 0.3 * unit.nanometer
epsilon = 2 * unit.kilojoule_per_mole
mass = 100 * unit.amu

bb = {"particle_type_name": "bb", "sigma": sigma, "epsilon": epsilon, "mass": mass}
sc = {"particle_type_name": "sc", "sigma": sigma, "epsilon": epsilon, "mass": mass}

# Define monomer (residue):
A = {
    "monomer_name": "A",
    "particle_sequence": [bb, sc],
    "bond_list": [[0, 1]],
    "start": 0,
    "end": 0}

# Specify bonded parameters:
bond_lengths = {
    "default_bond_length": 0.35 * unit.nanometer,
    "bb_bb_bb_bond_length": 0.40 * unit.nanometer}

bond_force_constants = {
    "default_bond_force_constant": 1000 * unit.kilojoule_per_mole / unit.nanometer**2}

equil_bond_angles = {
    "default_equil_bond_angle": 120.0 * unit.degrees,
    "bb_bb_bb_equil_bond_angle": 150.0 * unit.degrees}

bond_angle_force_constants = {
    "default_bond_angle_force_constant": 100.0 * unit.kilojoule_per_mole / unit.radian**2}

torsion_phase_angles = {
    "default_torsion_phase_angle": 150 * unit.degrees}

torsion_force_constants = {
    "default_torsion_force_constant": 2.0 * unit.kilojoule_per_mole,
    "bb_bb_bb_bb_torsion_force_constant": 5.0 * unit.kilojoule_per_mole}

torsion_periodicities = {
    "default_torsion_periodicity": 1}

# Define oligomer sequence:
sequence = 12 * [A]

# Initial particle positions determined from random builder

cgmodel = CGModel(
    particle_type_list=[bb, sc],
    bond_lengths=bond_lengths,
    bond_force_constants=bond_force_constants,
    bond_angle_force_constants=bond_angle_force_constants,
    torsion_force_constants=torsion_force_constants,
    equil_bond_angles=equil_bond_angles,
    torsion_phase_angles=torsion_phase_angles,
    torsion_periodicities=torsion_periodicities,
    include_nonbonded_forces=True,
    include_bond_forces=True,
    include_bond_angle_forces=True,
    include_torsion_forces=True,
    constrain_bonds=False,
    sequence=sequence,
    monomer_types=[A],
)
__init__(particle_type_list={}, charges={}, monomer_types=None, sequence=None, positions=None, bond_lengths={}, bond_force_constants={}, bond_angle_force_constants={}, torsion_force_constants={}, torsion_periodicities={}, equil_bond_angles={}, torsion_phase_angles={}, binary_interaction_parameters={}, go_model=False, go_repulsive_epsilon=None, constrain_bonds=False, include_nonbonded_forces=True, include_bond_forces=True, include_bond_angle_forces=True, include_torsion_forces=True, angle_style='harmonic', nonbond_repulsive_exp=12, nonbond_attractive_exp=6, exclusions={}, hbonds={}, rosetta_functional_form=False, check_energy_conservation=True, use_structure_library=False, random_positions=False, system=None, topology=None, simulation=None)[source]

Initialize definitions for all of the attributes of a coarse-grained model

Note

Default definitions are applied to force field parameters not explicitly defined.

Parameters
  • positions (Quantity() ( np.array( [cgmodel.num_beads,3] ), simtk.unit )) – Positions for the particles in the coarse-grained model (default = None)

  • particle_type_list (list) – list of all particle types (default = None)

  • monomer_types – A list of dictionary objects containing the properties for unique monomer types used to construct a heteropolymeric coarse-grained model (default = None)

  • sequence – The sequence from which to build a heteropolymer. Defined using a list of ‘monomer_types’, each of which contains the properties for that monomer (default = None (Homopolymer))

  • bond_lengths – Equilibrium bond lengths for all bond types (default = None)

  • bond_angle_force_constants – Harmonic bond angle-bending force constants for all bond types (default = None)

  • bond_force_constants – Harmonic bond-stretching force constants for all bond types (default = None)

  • equil_bond_angles – Equilibrium bond bending angle for all angle types (default = None)

  • torsion_force_constants – Torsion force constants for all torsion types (default = None)

  • torsion_phase_angles – Periodic torsion phase angle for all unique torsion angle definitions (default = 0)

  • binary_interaction_parameters (dict( 'type_name1_type_name2_binary_interaction': float )) – Binary interaction parameters used to scale nonbonded interactions between unlike particles (default=None)

  • go_model (bool) – If True, the binary interaction parameters will be applied only to the attractive component of the nonbonded potential. Otherwise, binary_interaction_parameters will be applied to the total potential (default=False)

  • go_repulsive_epsilon (Quantity ( float*unit.kilojoule_per_mole )) – If not None and go_model=True, use a fixed value for repulsive interactions, applied only to pairs that also have a binary interaction parameters < 1. (default=None)

  • constrain_bonds (Bool) – Option to use rigid bond constaints during a simulation of the energy for the system (default = False)

  • include_bond_forces (Bool) – Include bond stretching potentials when calculating the potential energy (default = True)

  • include_nonbonded_forces (Bool) – Include nonbonded interactions when calculating the potential energy (default = True)

  • include_bond_angle_forces (Bool) – Include contributions from bond angle forces when calculating the potential energy (default = True)

  • include_torsion_forces (Bool) – Include contributions from torsions when calculating the potential energy (default = True)

  • angle_style (str) – Functional form to use for bond-bending angle potential (‘harmonic’, ‘restricted’, or ‘cosine’) (default=’harmonic’)

  • nonbond_repulsive_exp (float) – Repulsive exponent for custom nonbonded Mie potential. For now this same exponent is applied to all pair types. (default=12)

  • nonbond_attractive_exp (float) – Attractive exponent for custom nonbonded Mie potential. For now this same exponent is applied to all pair types. (default=6)

  • exclusions (dict( list( int ) )) – Nonbonded weights for [1-2, 1-3, 1-4] interactions (default = [0,0,1])

  • hbonds (dict( 'donors': list(); 'acceptors': list(); 'epsilon_hb': Quantity * unit.kilojoule_per_mole; 'sigma_hb': Quantity * unit.angstrom; 'theta_d': Quantity * unit.degrees; 'theta_a': Quantity * unit.degrees )) – Dictionary containing directional CustomHbondedForce potential information.

  • rosetta_functional_form (Bool) – Option to use nonbonded exclusions consistent with Rosetta

  • check_energy_conservation (Bool) – Flag designating whether or not to perform a test OpenMM simulation with this coarse-grained model (default = True).

  • use_structure_library (Bool) – Flag designating whether or not to use a structure from the foldamers ensemble as the initial positions for the particles in the coarse-grained model (default = False)

  • random_positions (Bool) – Flag designating whether or not to generate a set of random coordinates for the coarse-grained model (default = None)

  • system (System()) – OpenMM System() object, which stores the forces for the coarse grained model (default = None)

  • topology (Topology()) – OpenMM Topology() object, which stores bonds, angles, and torsions coarse-grained model (default = None)

  • simulation (Simulation()) – OpenMM Simulation() object (default = None)

build_polymer(sequence)[source]

Used to build a polymer, or reset the properties for a polymer after parameters such as the polymer_length or sequence have been modified.

Parameters

sequence – The sequence from which to build a heteropolymer. Defined using a list of ‘monomer_types’, each of which contains the properties for that monomer (default = None (Homopolymer))

Returns

  • self.nonbonded_interaction_list ( List( List( int, int ) ) ) - A list of the nonbonded interactions (which don’t violate exclusion rules) in the coarse-grained model

create_particle_list()[source]

Get a list of particles, where the indices correspond to those in the system/topology.

Returns

  • particle_list ( List( str ) ) - A list of unique particles in the coarse-grained model

export(filename)[source]

Export a cgmodel to a pickle file.

Parameters

filename (str) – filename for exported cgmodel

get_bond_angle_force_constant(angle)[source]

Determines the correct bond angle force constant for a bond angle between three particles, given their indices within the coarse-grained model

Parameters

angle (List ( int )) – A list of the indices for the particles in an angle

Returns

  • bond_angle_force_constant ( Quantity() ) - The assigned bond angle force constant for the provided particles

get_bond_angle_list()[source]

Construct a list of bond angles, which can be used to build bond angle potentials for the coarse-grained model

Parameters

CGModel (class) – CGModel() class object

Returns

  • bond_angles ( List( List( int, int, int ) ) ) - A list of indices for all of the bond angles in the coarse-grained model

get_bond_force_constant(bond)[source]

Determines the correct bond force constant for two particles, given their indices

Parameters

bond (List ( int )) – A list of the indices for the particles in a bond

Returns

  • bond_force_constant ( Quantity() ) - The assigned bond force constant for the provided particles

get_bond_length(bond)[source]

Determines the correct bond length for two particles, given their indices.

Parameters

bond (List ( int )) – A list of the indices for the particles in a bond

Returns

  • bond_length ( Quantity() ) - The assigned bond length for the provided particles

get_bond_list()[source]

Construct a bond list for the coarse-grained model

Returns

  • bond_list ( List( List( int, int ) ) ) - A list of the bonds in the coarse-grained model.

get_equil_bond_angle(angle)[source]

Determines the correct equilibrium bond angle between three particles, given their indices within the coarse-grained model

Parameters

angle (List ( int )) – A list of the indices for the particles in an angle

Returns

  • equil_bond_angle (float) - The assigned equilibrium bond angle for the provided particles

get_exclusions(pair)[source]

Gets the exclusion rules applied to this specific pair of particles.

Parameters

pair (List ( int )) – A list of 2 indices defining a pair

Returns

  • exclusion_rules (list(int)) - List of exclusions for [1-2, 1-3, 1-4] nonbonded interactions (0 = no interaction)

get_nonbonded_exclusion_list(rosetta_functional_form=False)[source]

Get a list of the nonbonded interaction exclusions, which are assigned if two particles are separated by less than three bonds

Parameters

rosetta_functional_form (Bool) – Option to use nonbonded exclusions consistent with Rosetta

Returns

  • exclusion_list ( List( List( int, int ) ) ) - A list of the nonbonded particle interaction exclusions for the coarse-grained model

get_nonbonded_interaction_list()[source]

Construct a nonbonded interaction list for the coarse-grained model

Returns

  • interaction_list ( List( List( int, int ) ) ) - A list of the nonbonded interactions (which don’t violate exclusion rules) in the coarse-grained model

get_num_beads()[source]

Calculate the number of beads in a coarse-grained model class object

Returns

  • num_beads (int) - The total number of beads in the coarse-grained model

get_particle_attribute(particle, attribute)[source]

Get various attributes of a particle, given either the index or the particle dictionary

Parameters
  • particle (int or dict()) – Index of the particle of interest OR particle dictionary

  • attribute (str) – options are “monomer”, “monomer_type”, “name”, “index”, “type”, “mass”, “charge”, “epsilon”, “sigma”, “particle_type_name”

Returns

  • attribute of interest

get_particle_charge(particle)[source]

Returns the charge for a particle, given its index within the coarse-grained model, or the dict

Parameters

particle_index (int) – Index of the particle for which we would like to determine the type

Returns

  • particle_charge ( Quantity() ) - The charge for the provided particle index

get_particle_epsilon(particle)[source]

Returns the Lennard-Jones epsilon value for a particle, given its index within the coarse-grained model.

Parameters

particle_index (int) – Index of the particle for which we would like to determine the type

Returns

  • epsilon ( Quantity() ) - The assigned Lennard-Jones epsilon value for the provided particle index

get_particle_index(particle)[source]

Returns the index of a particle, given its index within the model or the particle dictionary. Obviously, kind of redundant if using the index instead of the particle dictionary

Parameters

particle_index (int) – Index of the particle for which we would like to determine the type

Returns

  • particle_name ( str ) - The name of the particle

get_particle_mass(particle)[source]

Returns the mass of a particle, given its index within the coarse-grained model or the particle dictionary

Parameters

particle (int or dict()) – Index of the particle for which we would like to determine the type, or dict()

Returns

  • epsilon ( Quantity() ) - The assigned Lennard-Jones epsilon value for the provided particle index

get_particle_monomer(particle)[source]

Indicates which monomer index a particle belongs to

Parameters

particle_index (int) – Index of the particle for which we would like to determine the monomer type

Returns

  • monomer_type (dict) : monomer type

get_particle_monomer_type(particle)[source]

Indicates which type of monomer a particle belongs to

Parameters

particle_index (int) – Index of the particle for which we would like to determine the monomer type

Returns

  • monomer_type (dict) : monomer type

get_particle_name(particle)[source]

Returns the name of a particle, given its index within the model

Parameters

particle_index (int) – Index of the particle for which we would like to determine the type

Returns

  • particle_name ( str ) - The name of the particle

get_particle_sigma(particle)[source]

Returns the Lennard-Jones potential sigma value for a particle, given the particle index

Parameters

particle_index (int) – Index of the particle for which we would like to determine the type

Returns

  • sigma ( Quantity() ) - The assigned Lennard-Jones sigma value for the provided particle index

get_particle_type(particle)[source]

Gives the type of a particle (a dictionary)

Parameters

particle (int or dict()) – Index of the particle for which we would like to determine the type OR particle dictionary

Returns

  • particle_type (str):

get_particle_type_name(particle)[source]

Gives the type name of a particle.

Parameters

particle_index (int or dict()) – Index of the particle for which we would like to determine the type name OR particle dictionary

Returns

  • particle_type_name (str):

get_torsion_force_constant(torsion)[source]

Determines the correct torsion force constant for a torsion (bond angle involving four particles), given their indices within the coarse-grained model. For sums of periodic torsions, this returns a list.

Parameters

torsion (List( int )) – A list of the indices for the particles in a torsion

Returns

  • torsion_force_constant ( Quantity() ) - The assigned torsion force constant for the provided particles

get_torsion_list()[source]

Construct a list of particle indices from which to define torsions for the coarse-grained model

Returns

  • torsions ( List( List( int, int, int, int ) ) ) - A list of the particle indices for the torsions in the coarse-grained model

get_torsion_periodicity(torsion)[source]

Determines the periodicity for a torsion, given a quadruplet of particle indices. For sums of periodic torsions, this returns a list.

Parameters

torsion (List( int )) – A list of the indices for the particles in a torsion

Returns

  • torsion_periodicity ( int ) - The periodicity for the input torsion

get_torsion_phase_angle(torsion)[source]

Determines the phase_angle for a torsion, given indices of the 4 particles within the coarse-grained model For sums of periodic torsions, this returns a list.

Parameters

torsion (List( int )) – A list of the indices for the particles in a torsion

Returns

  • torsion_phase_angle (float) - The assigned periodic torsion phase angle for the provided particles

process_monomer_types()[source]

Clean up a list of ‘monomer_types’ for all unique monomers.