build.cg_build submodule
- cg_openmm.build.cg_build.add_force(cgmodel, force_type=None, rosetta_functional_form=False)[source]
Given a ‘cgmodel’ and ‘force_type’ as input, this function adds the OpenMM force corresponding to ‘force_type’ to ‘cgmodel.system’.
- Parameters
cgmodel – CGModel() class object.
type – class
force_type (str) – Designates the kind of ‘force’ provided. (Valid options include: “Bond”, “Nonbonded”, “Angle”, and “Torsion”)
- Returns
cgmodel (class) - ‘foldamers’ CGModel() class object
force (class) - An OpenMM Force() object.
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> force_type = "Bond" >>> cgmodel,force = add_force(cgmodel,force_type=force_type)
- cg_openmm.build.cg_build.add_new_elements(cgmodel)[source]
Add coarse grained particle types to OpenMM.
- Parameters
cgmodel (class) – CGModel object (contains all attributes for a coarse grained model).
- Returns
new_particles (list) - a list of the particle names that were added to OpenMM’s ‘Element’ List.
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> particle_types = add_new_elements(cgmodel)
Warning
If the particle names were user defined, and any of the names conflict with existing element names in OpenMM, OpenMM will issue an error exit.
- cg_openmm.build.cg_build.add_rosetta_exception_parameters(cgmodel, nonbonded_force, particle_index_1, particle_index_2)[source]
- cg_openmm.build.cg_build.build_system(cgmodel, rosetta_functional_form=False, verify=True)[source]
Builds an OpenMM System() object, given a CGModel() as input.
- Parameters
cgmodel (class) – CGModel() class object
- Returns
system ( System() ) - OpenMM System() object
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> system = build_system(cgmodel) >>> cgmodel.system = system
- cg_openmm.build.cg_build.build_topology(cgmodel, use_pdbfile=False, pdbfile=None)[source]
Construct an OpenMM Topology() class object for our coarse grained model,
- Parameters
cgmodel (class) – CGModel() class object
use_pdbfile (Logical) – Determines whether or not to use a PDB file in order to generate the Topology().
pdbfile (str) – Name of a PDB file to use when building the topology.
- Returns
topology (Topology() ) - OpenMM Topology() object
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> from foldamers.util.iotools import write_pdbfile_without_topology >>> input_pdb = "top.pdb" >>> cgmodel = CGModel() >>> write_pdbfile_without_topology(cgmodel,input_pdb) >>> topology = build_topology(cgmodel,use_pdbfile=True,pdbfile=input_pdb) >>> cgmodel.topology = topology
Warning
When ‘use_pdbfile’=True, this function will use the PDBFile() class object from OpenMM to build the Topology(). In order for this approach to function correctly, the particle names in the PDB file must match the particle names in the coarse grained model.
- cg_openmm.build.cg_build.check_force(cgmodel, force, force_type=None)[source]
Given an OpenMM Force(), this function determines if there are any problems with its configuration.
- Parameters
cgmodel (class) – CGModel() class object.
force – An OpenMM Force() object.
force_type (str) – Designates the kind of ‘force’ provided. (Valid options include: “Nonbonded”)
- Returns
‘success’ (Logical) - a variable indicating if the force test passed.
- Example
>>> from simtk.openmm.openmm import NonbondedForce >>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> force = NonbondedForce() >>> force_type = "Nonbonded" >>> test_result = check_force(cgmodel,force,force_type="Nonbonded")
- cg_openmm.build.cg_build.check_forces(cgmodel)[source]
Given a cgmodel that contains positions and an an OpenMM System() object, this function tests the forces for cgmodel.system.
More specifically, this function confirms that the model does not have any “NaN” or unphysically large forces.
- Parameters
cgmodel – CGModel() class object.
type – class
- Returns
success (Logical) - Indicates if this cgmodel has unphysical forces.
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> pass_forces_test = check_forces(cgmodel)
- cg_openmm.build.cg_build.get_num_forces(cgmodel)[source]
Given a CGModel() class object, this function determines how many forces we are including when evaluating the energy.
- Parameters
cgmodel (class) – CGModel() class object
- Returns
total_forces (int) - Number of forces in the coarse grained model
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> total_number_forces = get_num_forces(cgmodel)
- cg_openmm.build.cg_build.verify_system(cgmodel)[source]
Given a CGModel() class object, this function confirms that its OpenMM System() object is configured correctly.
- Parameters
cgmodel (class) – CGModel() class object
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> verify_system(cgmodel)
Warning
The function will force an error exit if the system is invalid, and will proceed as normal if the system is valid.
- cg_openmm.build.cg_build.verify_topology(cgmodel)[source]
Given a coarse grained model that contains a Topology() (cgmodel.topology), this function verifies the validity of the topology.
- Parameters
cgmodel (class) – CGModel() class object.
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> verify_topology(cgmodel)
Warning
The function will force an error exit if the topology is invalid, and will proceed as normal if the topology is valid.
- cg_openmm.build.cg_build.write_xml_file(cgmodel, xml_file_name)[source]
Write an XML-formatted forcefield file for a coarse grained model.
- Parameters
cgmodel (class) – CGModel() class object.
xml_file_name (str) – Path to XML output file.
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> cgmodel = CGModel() >>> xml_file_name = "openmm_cgmodel.xml" >>> write_xml_file(cgmodel,xml_file_name)