GenerateModel
GenerateModel is the primary user-facing module in QMzyme. GenerateModel
is used to load a starting structure, define QMzyme regions, and write calculation input.
To avoid unintended behavior, the initial structure must be pre-processed. I.e., ensure
hydrogens have been added, and the structure is representative of the system you hope
to study. If atomic charge information is not present in the input file(s), QMzyme
will guess atomic charges by referring to the residue names. Any residue name corresponding
to standard protein residue names, defined here,
are able to be parsed for their total charge. This library can also be found in configuration under
the dictionary protein_residues. If you have a non-protein residue QMzyme will assume its charge is 0. This
is important if you have a ligand with a non-zero charge that you will include in your calculation. After importing
QMzyme you can update the charge dictionary to add this information by adding to the residue_charges dictionary found
in configuration.
The starting structure is loaded in using MDAnalysis and converted to a Universe object.
There are a variety of ways to define QMzyme region(s), and once a region has been set it
can be further modified, i.e., through truncation schemes. Lastly, this module interfaces with
CalculateModel, Writer and
qprep to automate the creation of single- or multi-scale
calculation input files.
- class QMzyme.GenerateModel.GenerateModel(*args, name=None, universe=None, select_atoms='all', frame=0, pickle_file=None, **kwargs)
Bases:
QMzymeModelGenerateModel can be instantiated with an MDAnalysis Universe directly, or any combination of parameters that MDAnalysis.core.universe.Universe accepts to create a Universe. See https://userguide.mdanalysis.org/stable/universe.html for details.
- Parameters:
name (str, optional) -- Name to give to the QMzymeModel. This is used for default file naming purposes throughout the QMzyme package. If not provided, it will default to the base name of the universe filename attribute.
universe (MDAnalysis.Universe, default=None) -- MDAnalysis Universe object. If not specified, user will need to provide file(s) that MDAnalysis can use to create a Universe object.
select_atoms (str, default='all') -- MDAnalysis selection command to specify which atoms are stored in the universe.
frame (int, default=0) -- If trajectory was provided, specify a frame to extract coordinates from.
pickle_file (str, default=None) -- Provide name/path+file of previously pickled QMzymeModel object to initialize
- Usage:
To instantiate a model from a PDB file called "filename.pdb":
model = QMzyme.GenerateModel("filename.pdb")
If "filename.pdb" contains any components you know you do not want included in your model, you can initialize the GenerateModel instance from a selection of atoms by using the select_atoms argument:
model = QMzyme.GenerateModel("filename.pdb", select_atoms="not resname WAT")
You can also initialize the model from a topology and trajectory file(s) and specify what frame to take coordinates from:
model = QMzyme.GenerateModel("filename.pdb", "filename.dcd", frame=100)
- set_catalytic_center(selection)
This method calls the set_region method and forces the region name to be 'catalytic_center'. See
set_region().
- set_region(selection, name=None, **kwargs)
Method to create a QMzymeRegion and add to the QMzymeModel regions list.
- Parameters:
selection (See options below, required) --
Determines what atoms are included in the region. A variety of input options are accepted:
str that can be interpreted by MDAnalysis selection commands
an MDAnalysis AtomGroup
any concrete class of
SelectionScheme, i.e.,DistanceCutoff. Options can be found inSelectionSchemes.
name (str, optional) -- Name of the resulting region.
kwargs -- Keyword arguments that might be needed if a
SelectionSchemeis used. For example, the parameter cutoff is required to use theDistanceCutoffscheme.
- Usage:
model.set_region(selection="resid 10 or resid 15", name="two_residues")
from QMzyme.SelectionSchemes import DistanceCutoff model.set_region(selection=DistanceCutoff, cutoff=5)
Note
When using a
SelectionSchemethe scheme class must be imported.
- truncate(scheme=<class 'QMzyme.TruncationSchemes.TerminalAlphaCarbon'>, name=None)
Method to truncate QMzymeModel. All QMzymeModel regions with assigned methods will be combined and truncated according to the specified scheme. The resulting region will be saved as the QMzymeModel truncated attribute.
- Parameters:
scheme (
TruncationSchemeconcrete class, default=:class:~QMzyme.TruncationSchemes.TerminalAlphaCarbon) -- Specifies the truncation scheme to use. Options can be found inTruncationSchemes.name (str, optional) -- Name to give the truncated model. If None, the original region name will be the combination of calculation methods and the suffix '_combined_region_truncated'.
- write_input(filename=None, memory='24GB', nprocs=12, reset_calculation=False)
Method to write calculation file input. The code will automatically detect what type of calculation file to prepare based on the calculation methods that have been assigned to the model region(s). Once this is called the QMzymeModel object will automatically be serialized using the pickle library and saved under the filename {self.name+'.pkl'} in the current working directory.
- Parameters:
filename (str, optional) -- Name to use for resulting file. If not specified, the file will be named according to the region(s) name. The file format ending does not need to be included.
memory (str, optional) -- Amount of memory to specify in the input file.
nprocs (int, optional) -- Number of processors to specify in the input file.
Note
A
QM_Methodmust first be assigned to a region.