Data Structures

The primary building blocks in an mBuild hierarchy inherit from the Compound class. Compounds maintain an ordered set of children which are other Compounds. In addition, an independent, ordered dictionary of labels is maintained through which users can reference any other Compound in the hierarchy via descriptive strings. Every Compound knows its parent Compound, one step up in the hierarchy, and knows which Compounds reference it in their labels. Ports are a special type of Compound which are used internally to connect different Compounds using the equivalence transformations described below.

Compounds at the bottom of an mBuild hierarchy, the leafs of the tree, are referred to as Particles and can be instantiated as foo = mb.Particle(name='bar'). Note however, that this merely serves to illustrate that this Compound is at the bottom of the hierarchy; Particle is simply an alias for Compound which can be used to clarify the intended role of an object you are creating. The method Compound.particles() traverses the hierarchy to the bottom and yields those Compounds. Compound.root() returns the compound at the top of the hierarchy.

Compound

class mbuild.compound.Compound(subcompounds=None, name=None, pos=None, charge=0.0, periodicity=None, port_particle=False)[source]

A building block in the mBuild hierarchy.

Compound is the superclass of all composite building blocks in the mBuild hierarchy. That is, all composite building blocks must inherit from compound, either directly or indirectly. The design of Compound follows the Composite design pattern (Gamma, Erich; Richard Helm; Ralph Johnson; John M. Vlissides (1995). Design Patterns: Elements of Reusable Object-Oriented Software. Addison-Wesley. p. 395. ISBN 0-201-63361-2.), with Compound being the composite, and Particle playing the role of the primitive (leaf) part, where Particle is in fact simply an alias to the Compound class.

Compound maintains a list of children (other Compounds contained within), and provides a means to tag the children with labels, so that the compounds can be easily looked up later. Labels may also point to objects outside the Compound’s containment hierarchy. Compound has built-in support for copying and deepcopying Compound hierarchies, enumerating particles or bonds in the hierarchy, proximity based searches, visualization, I/O operations, and a number of other convenience methods.

Parameters:
subcompounds : mb.Compound or list of mb.Compound, optional, default=None

One or more compounds to be added to self.

name : str, optional, default=self.__class__.__name__

The type of Compound.

pos : np.ndarray, shape=(3,), dtype=float, optional, default=[0, 0, 0]

The position of the Compound in Cartestian space

charge : float, optional, default=0.0

Currently not used. Likely removed in next release.

periodicity : np.ndarray, shape=(3,), dtype=float, optional, default=[0, 0, 0]

The periodic lengths of the Compound in the x, y and z directions. Defaults to zeros which is treated as non-periodic.

port_particle : bool, optional, default=False

Whether or not this Compound is part of a Port

Attributes:
bond_graph : mb.BondGraph

Graph-like object that stores bond information for this Compound

children : OrderedSet

Contains all children (other Compounds).

labels : OrderedDict

Labels to Compound/Atom mappings. These do not necessarily need not be in self.children.

parent : mb.Compound

The parent Compound that contains this part. Can be None if this compound is the root of the containment hierarchy.

referrers : set

Other compounds that reference this part with labels.

rigid_id : int, default=None

The ID of the rigid body that this Compound belongs to. Only Particles (the bottom of the containment hierarchy) can have integer values for rigid_id. Compounds containing rigid particles will always have rigid_id == None. See also contains_rigid.

boundingbox

Compute the bounding box of the compound.

center

The cartesian center of the Compound based on its Particles.

contains_rigid

Returns True if the Compound contains rigid bodies

max_rigid_id

Returns the maximum rigid body ID contained in the Compound.

n_particles

Return the number of Particles in the Compound.

n_bonds

Return the number of bonds in the Compound.

root

The Compound at the top of self’s hierarchy.

xyz

Return all particle coordinates in this compound.

xyz_with_ports

Return all particle coordinates in this compound including ports.

add(self, new_child, label=None, containment=True, replace=False, inherit_periodicity=True, reset_rigid_ids=True)[source]

Add a part to the Compound.

Note:
This does not necessarily add the part to self.children but may instead be used to add a reference to the part to self.labels. See ‘containment’ argument.
Parameters:
new_child : mb.Compound or list-like of mb.Compound

The object(s) to be added to this Compound.

label : str, optional

A descriptive string for the part.

containment : bool, optional, default=True

Add the part to self.children.

replace : bool, optional, default=True

Replace the label if it already exists.

inherit_periodicity : bool, optional, default=True

Replace the periodicity of self with the periodicity of the Compound being added

reset_rigid_ids : bool, optional, default=True

If the Compound to be added contains rigid bodies, reset the rigid_ids such that values remain distinct from rigid_ids already present in self. Can be set to False if attempting to add Compounds to an existing rigid body.

add_bond(self, particle_pair)[source]

Add a bond between two Particles.

Parameters:
particle_pair : indexable object, length=2, dtype=mb.Compound

The pair of Particles to add a bond between

all_ports(self)[source]

Return all Ports referenced by this Compound and its successors

Returns:
list of mb.Compound

A list of all Ports referenced by this Compound and its successors

ancestors(self)[source]

Generate all ancestors of the Compound recursively.

Yields:
mb.Compound

The next Compound above self in the hierarchy

available_ports(self)[source]

Return all unoccupied Ports referenced by this Compound.

Returns:
list of mb.Compound

A list of all unoccupied ports referenced by the Compound

bonds(self)[source]

Return all bonds in the Compound and sub-Compounds.

Yields:
tuple of mb.Compound

The next bond in the Compound

See also

bond_graph.edges_iter
Iterates over all edges in a BondGraph
boundingbox

Compute the bounding box of the compound.

Returns:
mb.Box

The bounding box for this Compound

center

The cartesian center of the Compound based on its Particles.

Returns:
np.ndarray, shape=(3,), dtype=float

The cartesian center of the Compound based on its Particles

contains_rigid

Returns True if the Compound contains rigid bodies

If the Compound contains any particle with a rigid_id != None then contains_rigid will return True. If the Compound has no children (i.e. the Compound resides at the bottom of the containment hierarchy) then contains_rigid will return False.

Returns:
bool

True if the Compound contains any particle with a rigid_id != None

Notes

The private variable ‘_check_if_contains_rigid_bodies’ is used to help cache the status of ‘contains_rigid’. If ‘_check_if_contains_rigid_bodies’ is False, then the rigid body containment of the Compound has not changed, and the particle tree is not traversed, boosting performance.

energy_minimize(self, forcefield='UFF', steps=1000, **kwargs)[source]

Perform an energy minimization on a Compound

Default behavior utilizes Open Babel (http://openbabel.org/docs/dev/) to perform an energy minimization/geometry optimization on a Compound by applying a generic force field

Can also utilize OpenMM (http://openmm.org/) to energy minimize after atomtyping a Compound using Foyer (https://github.com/mosdef-hub/foyer) to apply a forcefield XML file that contains valid SMARTS strings.

This function is primarily intended to be used on smaller components, with sizes on the order of 10’s to 100’s of particles, as the energy minimization scales poorly with the number of particles.

Parameters:
steps : int, optional, default=1000

The number of optimization iterations

forcefield : str, optional, default=’UFF’

The generic force field to apply to the Compound for minimization. Valid options are ‘MMFF94’, ‘MMFF94s’, ‘’UFF’, ‘GAFF’, and ‘Ghemical’. Please refer to the Open Babel documentation (http://open-babel. readthedocs.io/en/latest/Forcefields/Overview.html) when considering your choice of force field. Utilizing OpenMM for energy minimization requires a forcefield XML file with valid SMARTS strings. Please refer to (http://docs. openmm.org/7.0.0/userguide/application.html#creating-force-fields) for more information.

Keyword Arguments
————
algorithm : str, optional, default=’cg’

The energy minimization algorithm. Valid options are ‘steep’, ‘cg’, and ‘md’, corresponding to steepest descent, conjugate gradient, and equilibrium molecular dynamics respectively. For _energy_minimize_openbabel

scale_bonds : float, optional, default=1

Scales the bond force constant (1 is completely on). For _energy_minimize_openmm

scale_angles : float, optional, default=1

Scales the angle force constant (1 is completely on) For _energy_minimize_openmm

scale_torsions : float, optional, default=1

Scales the torsional force constants (1 is completely on) For _energy_minimize_openmm Note: Only Ryckaert-Bellemans style torsions are currently supported

scale_nonbonded : float, optional, default=1

Scales epsilon (1 is completely on) For _energy_minimize_openmm

References

If using _energy_minimize_openmm(), please cite: .. [R92550878d20b-1] P. Eastman, M. S. Friedrichs, J. D. Chodera, R. J. Radmer,

C. M. Bruns, J. P. Ku, K. A. Beauchamp, T. J. Lane, L.-P. Wang, D. Shukla, T. Tye, M. Houston, T. Stich, C. Klein, M. R. Shirts, and V. S. Pande. “OpenMM 4: A Reusable, Extensible, Hardware Independent Library for High Performance Molecular Simulation.” J. Chem. Theor. Comput. 9(1): 461-469. (2013).

If using _energy_minimize_openbabel(), please cite: .. [R92550878d20b-1] O’Boyle, N.M.; Banck, M.; James, C.A.; Morley, C.;

Vandermeersch, T.; Hutchison, G.R. “Open Babel: An open chemical toolbox.” (2011) J. Cheminf. 3, 33
[2]Open Babel, version X.X.X http://openbabel.org, (installed Month Year)

If using the ‘MMFF94’ force field please also cite the following: .. [R92550878d20b-3] T.A. Halgren, “Merck molecular force field. I. Basis, form,

scope, parameterization, and performance of MMFF94.” (1996) J. Comput. Chem. 17, 490-519
[4]T.A. Halgren, “Merck molecular force field. II. MMFF94 van der Waals and electrostatic parameters for intermolecular interactions.” (1996) J. Comput. Chem. 17, 520-552
[5]T.A. Halgren, “Merck molecular force field. III. Molecular geometries and vibrational frequencies for MMFF94.” (1996) J. Comput. Chem. 17, 553-586
[6]T.A. Halgren and R.B. Nachbar, “Merck molecular force field. IV. Conformational energies and geometries for MMFF94.” (1996) J. Comput. Chem. 17, 587-615
[7]T.A. Halgren, “Merck molecular force field. V. Extension of MMFF94 using experimental data, additional computational data, and empirical rules.” (1996) J. Comput. Chem. 17, 616-641

If using the ‘MMFF94s’ force field please cite the above along with: .. [R92550878d20b-8] T.A. Halgren, “MMFF VI. MMFF94s option for energy minimization

studies.” (1999) J. Comput. Chem. 20, 720-729

If using the ‘UFF’ force field please cite the following: .. [R92550878d20b-3] Rappe, A.K., Casewit, C.J., Colwell, K.S., Goddard, W.A. III,

Skiff, W.M. “UFF, a full periodic table force field for molecular mechanics and molecular dynamics simulations.” (1992) J. Am. Chem. Soc. 114, 10024-10039

If using the ‘GAFF’ force field please cite the following: .. [R92550878d20b-3] Wang, J., Wolf, R.M., Caldwell, J.W., Kollman, P.A., Case, D.A.

“Development and testing of a general AMBER force field” (2004) J. Comput. Chem. 25, 1157-1174

If using the ‘Ghemical’ force field please cite the following: .. [R92550878d20b-3] T. Hassinen and M. Perakyla, “New energy terms for reduced

protein models implemented in an off-lattice force field” (2001) J. Comput. Chem. 22, 1229-1242
from_parmed(self, structure, coords_only=False, infer_hierarchy=True)[source]

Extract atoms and bonds from a pmd.Structure.

Will create sub-compounds for every chain if there is more than one and sub-sub-compounds for every residue.

Parameters:
structure : pmd.Structure

The structure to load.

coords_only : bool

Set preexisting atoms in compound to coordinates given by structure.

infer_hierarchy : bool, optional, default=True

If true, infer compound hierarchy from chains and residues

from_pybel(self, pybel_mol, use_element=True, coords_only=False, infer_hierarchy=True)[source]

Create a Compound from a Pybel.Molecule

pybel_mol: pybel.Molecule use_element : bool, default True

If True, construct mb Particles based on the pybel Atom’s element. If False, construcs mb Particles based on the pybel Atom’s type
coords_only : bool, default False
Set preexisting atoms in compound to coordinates given by structure. Note: Not yet implemented, included only for parity with other conversion functions
infer_hierarchy : bool, optional, default=True
If True, infer hierarchy from residues
from_trajectory(self, traj, frame=-1, coords_only=False, infer_hierarchy=True)[source]

Extract atoms and bonds from a md.Trajectory.

Will create sub-compounds for every chain if there is more than one and sub-sub-compounds for every residue.

Parameters:
traj : mdtraj.Trajectory

The trajectory to load.

frame : int, optional, default=-1 (last)

The frame to take coordinates from.

coords_only : bool, optional, default=False

Only read coordinate information

infer_hierarchy : bool, optional, default=True

If True, infer compound hierarchy from chains and residues

generate_bonds(self, name_a, name_b, dmin, dmax)[source]

Add Bonds between all pairs of types a/b within [dmin, dmax].

Parameters:
name_a : str

The name of one of the Particles to be in each bond

name_b : str

The name of the other Particle to be in each bond

dmin : float

The minimum distance between Particles for considering a bond

dmax : float

The maximum distance between Particles for considering a bond

get_smiles(self)[source]

Get SMILES string for compound

Bond order is guessed with pybel and may lead to incorrect SMILES strings.

Returns:
smiles_string: str
label_rigid_bodies(self, discrete_bodies=None, rigid_particles=None)[source]

Designate which Compounds should be treated as rigid bodies

If no arguments are provided, this function will treat the compound as a single rigid body by providing all particles in self with the same rigid_id. If discrete_bodies is not None, each instance of a Compound with a name found in discrete_bodies will be treated as a unique rigid body. If rigid_particles is not None, only Particles (Compounds at the bottom of the containment hierarchy) matching this name will be considered part of the rigid body.

Parameters:
discrete_bodies : str or list of str, optional, default=None

Name(s) of Compound instances to be treated as unique rigid bodies. Compound instances matching this (these) name(s) will be provided with unique rigid_ids

rigid_particles : str or list of str, optional, default=None

Name(s) of Compound instances at the bottom of the containment hierarchy (Particles) to be included in rigid bodies. Only Particles matching this (these) name(s) will have their rigid_ids altered to match the rigid body number.

Examples

Creating a rigid benzene

>>> import mbuild as mb
>>> from mbuild.utils.io import get_fn
>>> benzene = mb.load(get_fn('benzene.mol2'))
>>> benzene.label_rigid_bodies()

Creating a semi-rigid benzene, where only the carbons are treated as a rigid body

>>> import mbuild as mb
>>> from mbuild.utils.io import get_fn
>>> benzene = mb.load(get_fn('benzene.mol2'))
>>> benzene.label_rigid_bodies(rigid_particles='C')

Create a box of rigid benzenes, where each benzene has a unique rigid body ID.

>>> import mbuild as mb
>>> from mbuild.utils.io import get_fn
>>> benzene = mb.load(get_fn('benzene.mol2'))
>>> benzene.name = 'Benzene'
>>> filled = mb.fill_box(benzene,
...                      n_compounds=10,
...                      box=[0, 0, 0, 4, 4, 4])
>>> filled.label_rigid_bodies(distinct_bodies='Benzene')

Create a box of semi-rigid benzenes, where each benzene has a unique rigid body ID and only the carbon portion is treated as rigid.

>>> import mbuild as mb
>>> from mbuild.utils.io import get_fn
>>> benzene = mb.load(get_fn('benzene.mol2'))
>>> benzene.name = 'Benzene'
>>> filled = mb.fill_box(benzene,
...                      n_compounds=10,
...                      box=[0, 0, 0, 4, 4, 4])
>>> filled.label_rigid_bodies(distinct_bodies='Benzene',
...                           rigid_particles='C')
max_rigid_id

Returns the maximum rigid body ID contained in the Compound.

This is usually used by compound.root to determine the maximum rigid_id in the containment hierarchy.

Returns:
int or None

The maximum rigid body ID contained in the Compound. If no rigid body IDs are found, None is returned

min_periodic_distance(self, xyz0, xyz1)[source]

Vectorized distance calculation considering minimum image.

Parameters:
xyz0 : np.ndarray, shape=(3,), dtype=float

Coordinates of first point

xyz1 : np.ndarray, shape=(3,), dtype=float

Coordinates of second point

Returns:
float

Vectorized distance between the two points following minimum image convention

n_bonds

Return the number of bonds in the Compound.

Returns:
int

The number of bonds in the Compound

n_particles

Return the number of Particles in the Compound.

Returns:
int

The number of Particles in the Compound

particles(self, include_ports=False)[source]

Return all Particles of the Compound.

Parameters:
include_ports : bool, optional, default=False

Include port particles

Yields:
mb.Compound

The next Particle in the Compound

particles_by_name(self, name)[source]

Return all Particles of the Compound with a specific name

Parameters:
name : str

Only particles with this name are returned

Yields:
mb.Compound

The next Particle in the Compound with the user-specified name

particles_in_range(self, compound, dmax, max_particles=20, particle_kdtree=None, particle_array=None)[source]

Find particles within a specified range of another particle.

Parameters:
compound : mb.Compound

Reference particle to find other particles in range of

dmax : float

Maximum distance from ‘compound’ to look for Particles

max_particles : int, optional, default=20

Maximum number of Particles to return

particle_kdtree : mb.PeriodicCKDTree, optional

KD-tree for looking up nearest neighbors. If not provided, a KD- tree will be generated from all Particles in self

particle_array : np.ndarray, shape=(n,), dtype=mb.Compound, optional

Array of possible particles to consider for return. If not provided, this defaults to all Particles in self

Returns:
np.ndarray, shape=(n,), dtype=mb.Compound

Particles in range of compound according to user-defined limits

See also

periodic_kdtree.PerioidicCKDTree
mBuild implementation of kd-trees
scipy.spatial.ckdtree
Further details on kd-trees
referenced_ports(self)[source]

Return all Ports referenced by this Compound.

Returns:
list of mb.Compound

A list of all ports referenced by the Compound

remove(self, objs_to_remove)[source]

Cleanly remove children from the Compound.

Parameters:
objs_to_remove : mb.Compound or list of mb.Compound

The Compound(s) to be removed from self

remove_bond(self, particle_pair)[source]

Deletes a bond between a pair of Particles

Parameters:
particle_pair : indexable object, length=2, dtype=mb.Compound

The pair of Particles to remove the bond between

rigid_particles(self, rigid_id=None)[source]

Generate all particles in rigid bodies.

If a rigid_id is specified, then this function will only yield particles with a matching rigid_id.

Parameters:
rigid_id : int, optional

Include only particles with this rigid body ID

Yields:
mb.Compound

The next particle with a rigid_id that is not None, or the next particle with a matching rigid_id if specified

root

The Compound at the top of self’s hierarchy.

Returns:
mb.Compound

The Compound at the top of self’s hierarchy

rotate(self, theta, around)[source]

Rotate Compound around an arbitrary vector.

Parameters:
theta : float

The angle by which to rotate the Compound, in radians.

around : np.ndarray, shape=(3,), dtype=float

The vector about which to rotate the Compound.

save(self, filename, show_ports=False, forcefield_name=None, forcefield_files=None, forcefield_debug=False, box=None, overwrite=False, residues=None, combining_rule='lorentz', foyer_kwargs=None, **kwargs)[source]

Save the Compound to a file.

Parameters:
filename : str

Filesystem path in which to save the trajectory. The extension or prefix will be parsed and control the format. Supported extensions are: ‘hoomdxml’, ‘gsd’, ‘gro’, ‘top’, ‘lammps’, ‘lmp’, ‘mcf’

show_ports : bool, optional, default=False

Save ports contained within the compound.

forcefield_files : str, optional, default=None

Apply a forcefield to the output file using a forcefield provided by the foyer package.

forcefield_name : str, optional, default=None

Apply a named forcefield to the output file using the foyer package, e.g. ‘oplsaa’. Forcefields listed here: https://github.com/mosdef-hub/foyer/tree/master/foyer/forcefields

forcefield_debug : bool, optional, default=False

Choose level of verbosity when applying a forcefield through foyer. Specifically, when missing atom types in the forcefield xml file, determine if the warning is condensed or verbose.

box : mb.Box, optional, default=self.boundingbox (with buffer)

Box information to be written to the output file. If ‘None’, a bounding box is used with 0.25nm buffers at each face to avoid overlapping atoms.

overwrite : bool, optional, default=False

Overwrite if the filename already exists

residues : str of list of str

Labels of residues in the Compound. Residues are assigned by checking against Compound.name.

combining_rule : str, optional, default=’lorentz’

Specify the combining rule for nonbonded interactions. Only relevant when the foyer package is used to apply a forcefield. Valid options are ‘lorentz’ and ‘geometric’, specifying Lorentz-Berthelot and geometric combining rules respectively.

foyer_kwargs : dict, optional, default=None

Keyword arguments to provide to foyer.Forcefield.apply.

**kwargs

Depending on the file extension these will be passed to either write_gsd, write_hoomdxml, write_lammpsdata, write_mcf, or parmed.Structure.save. See https://parmed.github.io/ParmEd/html/structobj/parmed.structure.Structure.html#parmed.structure.Structure.save

Other Parameters:
 
ref_distance : float, optional, default=1.0

Normalization factor used when saving to .gsd and .hoomdxml formats for converting distance values to reduced units.

ref_energy : float, optional, default=1.0

Normalization factor used when saving to .gsd and .hoomdxml formats for converting energy values to reduced units.

ref_mass : float, optional, default=1.0

Normalization factor used when saving to .gsd and .hoomdxml formats for converting mass values to reduced units.

atom_style: str, default=’full’

Defines the style of atoms to be saved in a LAMMPS data file. The following atom styles are currently supported: ‘full’, ‘atomic’, ‘charge’, ‘molecular’ see http://lammps.sandia.gov/doc/atom_style.html for more information on atom styles.

unit_style: str, default=’real’

Defines to unit style to be save in a LAMMPS data file. Defaults to ‘real’ units. Current styles are supported: ‘real’, ‘lj’ see https://lammps.sandia.gov/doc/99/units.html for more information on unit styles

See also

formats.gsdwrite.write_gsd
Write to GSD format
formats.hoomdxml.write_hoomdxml
Write to Hoomd XML format
formats.xyzwriter.write_xyz
Write to XYZ format
formats.lammpsdata.write_lammpsdata
Write to LAMMPS data format
formats.cassandramcf.write_mcf
Write to Cassandra MCF format
formats.json_formats.compound_to_json
Write to a json file

Notes

When saving the compound as a json, only the following arguments are used:
  • filename
  • show_ports
spin(self, theta, around)[source]

Rotate Compound in place around an arbitrary vector.

Parameters:
theta : float

The angle by which to rotate the Compound, in radians.

around : np.ndarray, shape=(3,), dtype=float

The axis about which to spin the Compound.

successors(self)[source]

Yield Compounds below self in the hierarchy.

Yields:
mb.Compound

The next Particle below self in the hierarchy

to_intermol(self, molecule_types=None)[source]

Create an InterMol system from a Compound.

Parameters:
molecule_types : list or tuple of subclasses of Compound
Returns:
intermol_system : intermol.system.System
to_networkx(self, names_only=False)[source]

Create a NetworkX graph representing the hierarchy of a Compound.

Parameters:
names_only : bool, optional, default=False
Store only the names of the

compounds in the graph, appended with their IDs, for distinction even if they have the same name. When set to False, the default behavior, the nodes are the compounds themselves.

Returns:
G : networkx.DiGraph

See also

mbuild.bond_graph

Notes

This digraph is not the bondgraph of the compound.

to_parmed(self, box=None, title='', residues=None, show_ports=False, infer_residues=False)[source]

Create a ParmEd Structure from a Compound.

Parameters:
box : mb.Box, optional, default=self.boundingbox (with buffer)

Box information to be used when converting to a Structure. If ‘None’, a bounding box is used with 0.25nm buffers at each face to avoid overlapping atoms, unless self.periodicity is not None, in which case those values are used for the box lengths.

title : str, optional, default=self.name

Title/name of the ParmEd Structure

residues : str of list of str

Labels of residues in the Compound. Residues are assigned by checking against Compound.name.

show_ports : boolean, optional, default=False

Include all port atoms when converting to a Structure.

infer_residues : bool, optional, default=False

Attempt to assign residues based on names of children.

Returns:
parmed.structure.Structure

ParmEd Structure object converted from self

See also

parmed.structure.Structure
Details on the ParmEd Structure object
to_pybel(self, box=None, title='', residues=None, show_ports=False, infer_residues=False)[source]

Create a pybel.Molecule from a Compound

box : mb.Box, def None title : str, optional, default=self.name

Title/name of the ParmEd Structure
residues : str of list of str
Labels of residues in the Compound. Residues are assigned by checking against Compound.name.
show_ports : boolean, optional, default=False
Include all port atoms when converting to a Structure.
infer_residues : bool, optional, default=False
Attempt to assign residues based on names of children

pybel.Molecule

Notes

Most of the mb.Compound is first converted to openbabel.OBMol And then pybel creates a pybel.Molecule from the OBMol Bond orders are assumed to be 1 OBMol atom indexing starts at 1, with spatial dimension Angstrom

to_trajectory(self, show_ports=False, chains=None, residues=None, box=None)[source]

Convert to an md.Trajectory and flatten the compound.

Parameters:
show_ports : bool, optional, default=False

Include all port atoms when converting to trajectory.

chains : mb.Compound or list of mb.Compound

Chain types to add to the topology

residues : str of list of str

Labels of residues in the Compound. Residues are assigned by checking against Compound.name.

box : mb.Box, optional, default=self.boundingbox (with buffer)

Box information to be used when converting to a Trajectory. If ‘None’, a bounding box is used with a 0.5nm buffer in each dimension. to avoid overlapping atoms, unless self.periodicity is not None, in which case those values are used for the box lengths.

Returns:
trajectory : md.Trajectory

See also

_to_topology
translate(self, by)[source]

Translate the Compound by a vector

Parameters:
by : np.ndarray, shape=(3,), dtype=float
translate_to(self, pos)[source]

Translate the Compound to a specific position

Parameters:
pos : np.ndarray, shape=3(,), dtype=float
unlabel_rigid_bodies(self)[source]

Remove all rigid body labels from the Compound

update_coordinates(self, filename, update_port_locations=True)[source]

Update the coordinates of this Compound from a file.

Parameters:
filename : str

Name of file from which to load coordinates. Supported file types are the same as those supported by load()

update_port_locations : bool, optional, default=True

Update the locations of Ports so that they are shifted along with their anchor particles. Note: This conserves the location of Ports with respect to the anchor Particle, but does not conserve the orientation of Ports with respect to the molecule as a whole.

See also

load
Load coordinates from a file
visualize(self, show_ports=False, backend='py3dmol', color_scheme={})[source]

Visualize the Compound using py3dmol (default) or nglview.

Allows for visualization of a Compound within a Jupyter Notebook.

Parameters:
show_ports : bool, optional, default=False

Visualize Ports in addition to Particles

backend : str, optional, default=’py3dmol’

Specify the backend package to visualize compounds Currently supported: py3dmol, nglview

color_scheme : dict, optional

Specify coloring for non-elemental particles keys are strings of the particle names values are strings of the colors i.e. {‘_CGBEAD’: ‘blue’}

xyz

Return all particle coordinates in this compound.

Returns:
pos : np.ndarray, shape=(n, 3), dtype=float

Array with the positions of all particles.

xyz_with_ports

Return all particle coordinates in this compound including ports.

Returns:
pos : np.ndarray, shape=(n, 3), dtype=float

Array with the positions of all particles and ports.

Port

class mbuild.port.Port(anchor=None, orientation=None, separation=0)[source]

A set of four ghost Particles used to connect parts.

Parameters:
anchor : mb.Particle, optional, default=None

A Particle associated with the port. Used to form bonds.

orientation : array-like, shape=(3,), optional, default=[0, 1, 0]

Vector along which to orient the port

separation : float, optional, default=0

Distance to shift port along the orientation vector from the anchor particle position. If no anchor is provided, the port will be shifted from the origin.

Attributes:
anchor : mb.Particle, optional, default=None

A Particle associated with the port. Used to form bonds.

up : mb.Compound

Collection of 4 ghost particles used to perform equivalence transforms. Faces the opposite direction as self[‘down’].

down : mb.Compound

Collection of 4 ghost particles used to perform equivalence transforms. Faces the opposite direction as self[‘up’].

used : bool

Status of whether a port has been occupied following an equivalence transform.

access_labels

List of labels used to access the Port

Returns:
list of str

Strings that can be used to access this Port relative to self.root

center

The cartesian center of the Port

direction

The unit vector pointing in the ‘direction’ of the Port