Modeller¶
- class openmm.app.modeller.Modeller(topology, positions)¶
- Modeller provides tools for editing molecular models, such as adding water or missing hydrogens. - To use it, create a Modeller object, specifying the initial Topology and atom positions. You can then call various methods to change the model in different ways. Each time you do, a new Topology and list of coordinates is created to represent the changed model. Finally, call getTopology() and getPositions() to get the results. - __init__(topology, positions)¶
- Create a new Modeller object - Parameters
- topology (Topology) – the initial Topology of the model 
- positions (list) – the initial atomic positions 
 
 
 - Methods - __init__(topology, positions)- Create a new Modeller object - add(addTopology, addPositions)- Add chains, residues, atoms, and bonds to the model. - addExtraParticles(forcefield[, …])- Add missing extra particles to the model that are required by a force field. - addHydrogens([forcefield, pH, variants, …])- Add missing hydrogens to the model. - addMembrane(forcefield[, lipidType, …])- Add a lipid membrane to the model. - addSolvent(forcefield[, model, boxSize, …])- Add solvent (both water and ions) to the model to fill a periodic box. - convertWater([model])- Convert all water molecules to a different water model. - delete(toDelete)- Delete chains, residues, atoms, and bonds from the model. - Delete all water molecules from the model. - Get the atomic positions. - Get the Topology of the model. - loadHydrogenDefinitions(file)- Load an XML file containing definitions of hydrogens that should be added by addHydrogens(). - getTopology()¶
- Get the Topology of the model. 
 - getPositions()¶
- Get the atomic positions. 
 - add(addTopology, addPositions)¶
- Add chains, residues, atoms, and bonds to the model. - Specify what to add by providing a new Topology object and the corresponding atomic positions. All chains, residues, atoms, and bonds contained in the Topology are added to the model. - Parameters
- addTopology (Topology) – a Topology whose contents should be added to the model 
- addPositions (list) – the positions of the atoms to add 
 
 
 - delete(toDelete)¶
- Delete chains, residues, atoms, and bonds from the model. - You can specify objects to delete at any granularity: atoms, residues, or chains. Passing in an Atom object causes that Atom to be deleted. Passing in a Residue object causes that Residue and all Atoms it contains to be deleted. Passing in a Chain object causes that Chain and all Residues and Atoms it contains to be deleted. - In all cases, when an Atom is deleted, any bonds it participates in are also deleted. You also can specify a bond (as a tuple of Atom objects) to delete just that bond without deleting the Atoms it connects. - Parameters
- toDelete (list) – a list of Atoms, Residues, Chains, and bonds (specified as tuples of Atoms) to delete 
 
 - deleteWater()¶
- Delete all water molecules from the model. 
 - convertWater(model='tip3p')¶
- Convert all water molecules to a different water model. - Deprecated - Use addExtraParticles() instead. It performs the same function but in a more general way. - Parameters
- model (string='tip3p') – the water model to convert to. Supported values are ‘tip3p’, ‘spce’, ‘tip4pew’, and ‘tip5p’. 
 
 - addSolvent(forcefield, model='tip3p', boxSize=None, boxVectors=None, padding=None, numAdded=None, boxShape='cube', positiveIon='Na+', negativeIon='Cl-', ionicStrength=0 M, neutralize=True, residueTemplates={})¶
- Add solvent (both water and ions) to the model to fill a periodic box. - The algorithm works as follows: - Water molecules are added to fill the box. 
- Water molecules are removed if their distance to any solute atom is less than the sum of their van der Waals radii. 
- If the solute is charged and neutralize=True, enough positive or negative ions are added to neutralize it. Each ion is added by randomly selecting a water molecule and replacing it with the ion. 
- Ion pairs are added to give the requested total ionic strength. Note that only monovalent ions are currently supported. 
 - The box size can be specified in any of several ways: - You can explicitly give the vectors defining the periodic box to use. 
- Alternatively, for a rectangular box you can simply give the dimensions of the unit cell. 
- You can give a padding distance. A bounding sphere containing the solute is determined, and the box size is set to (sphere diameter)+(padding). This guarantees no atom in the solute will come closer than the padding distance to any atom of another periodic copy. If the sphere diameter is less than the padding distance, the box size is set to 2*(padding) to ensure no atom is closer than the padding distance to two periodic copies of any other atom. 
- You can specify the total number of molecules (both waters and ions) to add. A box is then created whose size is just large enough hold the specified amount of solvent. 
- Finally, if none of the above options is specified, the existing Topology’s box vectors are used. 
 - When specifying either a padding distance or a number of molecules, you can specify a shape for the periodic box: cubic, rhombic dodecahedron, or truncated octahedron. Using a non-rectangular box allows the same distance between periodic copies to be achieved with a smaller box. The most compact option is a rhombic dodecahedron, for which the box volume is 70.7% the volume of a cubic box with the same amount of padding. - There exist many different water models, many of which are very similar to each other. This method creates preequilibrated water boxes for a limited set of water models. In most cases, they work equally well for other models that involve the same number of particles. For example, to simulate a box of TIP3P-FB water, use this method to create a box of TIP3P water, construct a System using TIP3P-FB parameters, and perform a local energy minimization to correct for the small differences between the models. Likewise, a box of TIP4P-Ew water can be used for most four site water models. - Parameters
- forcefield (ForceField) – the ForceField to use for determining van der Waals radii and atomic charges 
- model (str='tip3p') – the water model to use. Supported values are ‘tip3p’, ‘spce’, ‘tip4pew’, ‘tip5p’, and ‘swm4ndp’ (polarizable). 
- boxSize (Vec3=None) – the size of the box to fill with water 
- boxVectors (tuple of Vec3=None) – the vectors defining the periodic box to fill with water 
- padding (distance=None) – the padding distance to use 
- numAdded (int=None) – the total number of molecules (waters and ions) to add 
- boxShape (str='cube') – the box shape to use. Allowed values are ‘cube’, ‘dodecahedron’, and ‘octahedron’. If padding and numAdded are both None, this is ignored. 
- positiveIon (string='Na+') – the type of positive ion to add. Allowed values are ‘Cs+’, ‘K+’, ‘Li+’, ‘Na+’, and ‘Rb+’ 
- negativeIon (string='Cl-') – the type of negative ion to add. Allowed values are ‘Cl-’, ‘Br-’, ‘F-’, and ‘I-’. Be aware that not all force fields support all ion types. 
- ionicStrength (concentration=0*molar) – the total concentration of ions (both positive and negative) to add. This does not include ions that are added to neutralize the system. Note that only monovalent ions are currently supported. 
- neutralize (bool=True) – whether to add ions to neutralize the system 
- residueTemplates (dict=dict()) – specifies which template the ForceField should use for particular residues. The keys should be Residue objects from the Topology, and the values should be the names of the templates to use for them. This is useful when a ForceField contains multiple templates that can match the same residue (e.g Fe2+ and Fe3+ templates in the ForceField for a monoatomic iron ion in the Topology). 
 
 
 - static loadHydrogenDefinitions(file)¶
- Load an XML file containing definitions of hydrogens that should be added by addHydrogens(). - The built in hydrogens.xml file containing definitions for standard amino acids and nucleotides is loaded automatically. This method can be used to load additional definitions for other residue types. They will then be used in subsequent calls to addHydrogens(). - Parameters
- file (string or file) – An XML file containing hydrogen definitions. It may be either an absolute file path, a path relative to the current working directory, a path relative to this module’s data subdirectory (for built in sets of definitions), or an open file-like object with a read() method from which the data can be loaded. 
 
 - addHydrogens(forcefield=None, pH=7.0, variants=None, platform=None, residueTemplates={})¶
- Add missing hydrogens to the model. - Some residues can exist in multiple forms depending on the pH and properties of the local environment. These variants differ in the presence or absence of particular hydrogens. In particular, the following variants are supported: - Aspartic acid:
- ASH: Neutral form with a hydrogen on one of the delta oxygens ASP: Negatively charged form without a hydrogen on either delta oxygen 
- Cysteine:
- CYS: Neutral form with a hydrogen on the sulfur CYX: No hydrogen on the sulfur (either negatively charged, or part of a disulfide bond) 
- Glutamic acid:
- GLH: Neutral form with a hydrogen on one of the epsilon oxygens GLU: Negatively charged form without a hydrogen on either epsilon oxygen 
- Histidine:
- HID: Neutral form with a hydrogen on the ND1 atom HIE: Neutral form with a hydrogen on the NE2 atom HIP: Positively charged form with hydrogens on both ND1 and NE2 HIN: Negatively charged form without a hydrogen on either ND1 or NE2 
- Lysine:
- LYN: Neutral form with two hydrogens on the zeta nitrogen LYS: Positively charged form with three hydrogens on the zeta nitrogen 
 - The variant to use for each residue is determined by the following rules: - The most common variant at the specified pH is selected. 
- Any Cysteine that participates in a disulfide bond uses the CYX variant regardless of pH. 
- For a neutral Histidine residue, the HID or HIE variant is selected based on which one forms a better hydrogen bond. 
 - You can override these rules by explicitly specifying a variant for any residue. To do that, provide a list for the ‘variants’ parameter, and set the corresponding element to the name of the variant to use. - A special case is when the model already contains a hydrogen that should not be present in the desired variant. If you explicitly specify a variant using the ‘variants’ parameter, the residue will be modified to match the desired variant, removing hydrogens if necessary. On the other hand, for residues whose variant is selected automatically, this function will only add hydrogens. It will never remove ones that are already present in the model, regardless of the specified pH. - In all cases, the positions of existing atoms (including existing hydrogens) are not modified. - Definitions for standard amino acids and nucleotides are built in. You can call loadHydrogenDefinitions() to load additional definitions for other residue types. - Parameters
- forcefield (ForceField=None) – the ForceField to use for determining the positions of hydrogens. If this is None, positions will be picked which are generally reasonable but not optimized for any particular ForceField. 
- pH (float=7.0) – the pH based on which to select variants 
- variants (list=None) – an optional list of variants to use. If this is specified, its length must equal the number of residues in the model. variants[i] is the name of the variant to use for residue i (indexed starting at 0). If an element is None, the standard rules will be followed to select a variant for that residue. Alternatively, an element may specify exactly which hydrogens to add. In that case, variants[i] should be a list of tuples [(name1, parent1), (name2, parent2), …]. Each tuple specifies the name of a hydrogen and the name of the parent atom it should be bonded to. 
- platform (Platform=None) – the Platform to use when computing the hydrogen atom positions. If this is None, the default Platform will be used. 
- residueTemplates (dict=dict()) – specifies which template the ForceField should use for particular residues. The keys should be Residue objects from the Topology, and the values should be the names of the templates to use for them. This is useful when a ForceField contains multiple templates that can match the same residue (e.g Fe2+ and Fe3+ templates in the ForceField for a monoatomic iron ion in the Topology). 
 
- Returns
- a list of what variant was actually selected for each residue, in the same format as the variants parameter 
- Return type
- list 
 
 - addExtraParticles(forcefield, ignoreExternalBonds=False, residueTemplates={})¶
- Add missing extra particles to the model that are required by a force field. - Some force fields use “extra particles” that do not represent actual atoms, but still need to be included in the System. Examples include lone pairs, Drude particles, and the virtual sites used in some water models to adjust the charge distribution. Extra particles can be recognized by the fact that their element is None. - This method is primarily used to add extra particles, but it can also remove them. It tries to match every residue in the Topology to a template in the force field. If there is no match, it will both add and remove extra particles as necessary to make it match. - Parameters
- forcefield (ForceField) – the ForceField defining what extra particles should be present 
- ignoreExternalBonds (boolean=False) – If true, ignore external bonds when matching residues to templates. This is useful when the Topology represents one piece of a larger molecule, so chains are not terminated properly. 
- residueTemplates (dict=dict()) – specifies which template the ForceField should use for particular residues. The keys should be Residue objects from the Topology, and the values should be the names of the templates to use for them. This is useful when a ForceField contains multiple templates that can match the same residue (e.g Fe2+ and Fe3+ templates in the ForceField for a monoatomic iron ion in the Topology). 
 
 
 - addMembrane(forcefield, lipidType='POPC', membraneCenterZ=0 nm, minimumPadding=1 nm, positiveIon='Na+', negativeIon='Cl-', ionicStrength=0 M, neutralize=True, residueTemplates={}, platform=None)¶
- Add a lipid membrane to the model. - This method actually adds both a membrane and a water box. It is best to build them together, both to avoid adding waters inside the membrane and to ensure that lipid head groups are properly solvated. For that reason, this method includes many of the same arguments as addSolvent(). - The membrane is added in the XY plane, and the existing protein is assumed to already be oriented and positioned correctly. When possible, it is recommended to start with a model from the Orientations of Proteins in Membranes (OPM) database at http://opm.phar.umich.edu. Otherwise, it is up to you to select the protein position yourself. - The algorithm is based on the one described in Wolf et al., J. Comp. Chem. 31, pp. 2169-2174 (2010). It begins by tiling copies of a pre-equilibrated membrane patch to create a membrane of the desired size. Next it scales down the protein by 50% along the X and Y axes. Any lipid within a cutoff distance of the scaled protein is removed. It also ensures that equal numbers of lipids are removed from each leaf of the membrane. Finally, it runs molecular dynamics to let the membrane relax while gradually scaling the protein back up to its original size. - The size of the membrane and water box are determined by the minimumPadding argument. All pre-existing atoms are guaranteed to be at least this far from any edge of the periodic box. It is also possible for the periodic box to have more padding than requested. In particular, it only adds whole copies of the pre-equilibrated membrane patch, so the box dimensions will always be integer multiples of the patch size. That may lead to a larger membrane than what you requested. - This method has built in support for POPC, POPE, DLPC, DLPE, DMPC, DOPC and DPPC lipids. You can also build other types of membranes by providing a pre-equilibrated, solvated membrane patch that can be tiled in the XY plane to form the membrane. - Parameters
- forcefield (ForceField) – the ForceField to use for determining atomic charges and for relaxing the membrane 
- lipidType (string or object) – the type of lipid to use. Supported string values are ‘POPC’, ‘POPE’, ‘DLPC’, ‘DLPE’, ‘DMPC’, ‘DOPC’, and ‘DPPC’. For other types of lipids, provide a PDBFile or PDBxFile object (or any other object with “topology” and “positions” fields) containing a membrane patch. 
- membraneCenterZ (distance=0*nanometer) – the position along the Z axis of the center of the membrane 
- minimumPadding (distance=1*nanometer) – the padding distance to use 
- positiveIon (string='Na+') – the type of positive ion to add. Allowed values are ‘Cs+’, ‘K+’, ‘Li+’, ‘Na+’, and ‘Rb+’ 
- negativeIon (string='Cl-') – the type of negative ion to add. Allowed values are ‘Cl-’, ‘Br-’, ‘F-’, and ‘I-’. Be aware that not all force fields support all ion types. 
- ionicStrength (concentration=0*molar) – the total concentration of ions (both positive and negative) to add. This does not include ions that are added to neutralize the system. Note that only monovalent ions are currently supported. 
- neutralize (bool=True) – whether to add ions to neutralize the system 
- residueTemplates (dict=dict()) – specifies which template the ForceField should use for particular residues. The keys should be Residue objects from the Topology, and the values should be the names of the templates to use for them. This is useful when a ForceField contains multiple templates that can match the same residue (e.g Fe2+ and Fe3+ templates in the ForceField for a monoatomic iron ion in the Topology). 
- platform (Platform=None) – the Platform to use when computing the hydrogen atom positions. If this is None, the default Platform will be used. 
 
 
 
