ForceBalance API
1.3
Automated optimization of force fields and empirical potentials
|
Force field class. More...
Public Member Functions | |
def | __init__ (self, options, verbose=True, printopt=True) |
Instantiation of force field class. More... | |
def | fromfile (cls, fnm) |
def | __getstate__ (self) |
def | __setstate__ (self, state) |
def | addff (self, ffname, xmlScript=False) |
Parse a force field file and add it to the class. More... | |
def | check_dupes (self, pid, ffname, ln, pfld) |
Check to see whether a given parameter ID already exists, and provide an alternate if needed. More... | |
def | addff_txt (self, ffname, fftype, xmlScript) |
Parse a text force field and create several important instance variables. More... | |
def | addff_xml (self, ffname) |
Parse an XML force field file and create important instance variables. More... | |
def | make (self, vals=None, use_pvals=False, printdir=None, precision=12) |
Create a new force field using provided parameter values. More... | |
def | make_redirect (self, mvals) |
def | find_spacings (self) |
def | create_pvals (self, mvals) |
Converts mathematical to physical parameters. More... | |
def | create_mvals (self, pvals) |
Converts physical to mathematical parameters. More... | |
def | rsmake (self, printfacs=True) |
Create the rescaling factors for the coordinate transformation in parameter space. More... | |
def | make_rescale (self, scales, mvals=None, G=None, H=None, multiply=True, verbose=False) |
Obtain rescaled versions of the inputs according to dictionary values in "scales" (i.e. More... | |
def | mktransmat (self) |
Create the transformation matrix to rescale and rotate the mathematical parameters. More... | |
def | list_map (self) |
Create the plist, which is like a reversed version of the parameter map. More... | |
def | print_map (self, vals=None, precision=4) |
Prints out the (physical or mathematical) parameter indices, IDs and values in a visually appealing way. More... | |
def | sprint_map (self, vals=None, precision=4) |
Prints out the (physical or mathematical) parameter indices, IDs and values to a string. More... | |
def | assign_p0 (self, idx, val) |
Assign physical parameter values to the 'pvals0' array. More... | |
def | assign_field (self, idx, pid, fnm, ln, pfld, mult, cmd=None) |
Record the locations of a parameter in a txt file; [[file name, line number, field number, and multiplier]]. More... | |
def | plot_parameter_tree (fnm) |
def | get_mathid (self, pid) |
For a given unique parameter identifier, return all of the mvals that it is connected to. More... | |
def | __eq__ (self, other) |
Public Attributes | |
ffdata | |
As these options proliferate, the force field class becomes less standalone. More... | |
ffdata_isxml | |
map | |
The mapping of interaction type -> parameter number. More... | |
plist | |
The listing of parameter number -> interaction types. More... | |
patoms | |
A listing of parameter number -> atoms involved. More... | |
pfields | |
A list where pfields[i] = [pid,'file',line,field,mult,cmd], basically a new way to modify force field files; when we modify the force field file, we go to the specific line/field in a given file and change the number. More... | |
pTree | |
Improved representation of pfields as a networkx graph. More... | |
offxml_unit_strs | |
rs | |
List of rescaling factors. More... | |
tm | |
The transformation matrix for mathematical -> physical parameters. More... | |
tmI | |
The transpose of the transformation matrix. More... | |
excision | |
Indices to exclude from optimization / Hessian inversion. More... | |
np | |
The total number of parameters. More... | |
pvals0 | |
Initial value of physical parameters. More... | |
Readers | |
A dictionary of force field reader classes. More... | |
atomnames | |
A list of atom names (this is new, for ESP fitting) More... | |
FFAtomTypes | |
WORK IN PROGRESS ## This is a dictionary of {'AtomType':{'Mass' : float, 'Charge' : float, 'ParticleType' : string ('A', 'S', or 'D'), 'AtomicNumber' : int}}. More... | |
FFMolecules | |
redirect | |
Creates plist from map. More... | |
linedestroy_save | |
Destruction dictionary (experimental). More... | |
prmdestroy_save | |
linedestroy_this | |
prmdestroy_this | |
tinkerprm | |
openmmxml | |
offxml | |
openff_forcefield | |
amber_mol2 | |
amber_frcmod | |
rs_ord | |
Takes the dictionary 'BONDS':{3:'B', 4:'K'}, 'VDW':{4:'S', 5:'T'}, and turns it into a list of term types ['BONDSB','BONDSK','VDWS','VDWT']. More... | |
rs_type | |
qmap | |
qid | |
qid2 | |
Force field class.
This class contains all methods for force field manipulation. To create an instance of this class, an input file is required containing the list of force field file names. Everything else inside this class pertaining to force field generation is self-contained.
For details on force field parsing, see the detailed documentation for addff.
Definition at line 208 of file forcefield.py.
def src.forcefield.FF.__init__ | ( | self, | |
options, | |||
verbose = True , |
|||
printopt = True |
|||
) |
Instantiation of force field class.
Many variables here are initialized to zero, but they are filled out by methods like addff, rsmake, and mktransmat.
Definition at line 216 of file forcefield.py.
def src.forcefield.FF.__eq__ | ( | self, | |
other | |||
) |
Definition at line 1557 of file forcefield.py.
def src.forcefield.FF.__getstate__ | ( | self | ) |
Definition at line 352 of file forcefield.py.
def src.forcefield.FF.__setstate__ | ( | self, | |
state | |||
) |
Definition at line 361 of file forcefield.py.
def src.forcefield.FF.addff | ( | self, | |
ffname, | |||
xmlScript = False |
|||
) |
Parse a force field file and add it to the class.
First, figure out the type of force field file. This is done either by explicitly specifying the type using for example, ffname force_field.xml:openmm
or we figure it out by looking at the file extension.
Next, parse the file. Currently we support two classes of files - text and XML. The two types are treated very differently; for XML we use the parsers in libxml (via the python lxml module), and for text files we have our own in-house parsing class. Within text files, there is also a specialized GROMACS and TINKER parser as well as a generic text parser.
The job of the parser is to determine the following things: 1) Read the user-specified selection of parameters being fitted 2) Build a mapping (dictionary) of parameter identifier -> index in parameter vector
3) Build a list of physical parameter values 4) Figure out where to replace the parameter values in the force field file when the values are changed 5) Figure out which parameters need to be repeated or sign-flipped
Generally speaking, each parameter value in the force field file has a unique parameter identifier
. The identifier consists of three parts - the interaction type, the parameter subtype (within that interaction type), and the atoms involved.
— If XML: —
The force field file is read in using the lxml Python module. Specify which parameter you want to fit using by adding a 'parameterize' element to the end of the force field XML file, like so.
In this example, the parameter identifier would look like
Vdw/74/epsilon
.
— If GROMACS (.itp) or TINKER (.prm) : —
Follow the rules in the ITP_Reader or Tinker_Reader derived class. Read the documentation in the class documentation or the 'feed' method to learn more. In all cases the parameter is tagged using
# PRM 3
(where # denotes a comment, the word PRM stays the same, and 3 is the field number starting from zero.)
— If normal text : —
The parameter identifier is simply built using the file name, line number, and field. Thus, the identifier is unique but completely noninformative (which is not ideal for our purposes, but it works.)
— Endif —
[in] | ffname | Name of the force field file |
Definition at line 440 of file forcefield.py.
def src.forcefield.FF.addff_txt | ( | self, | |
ffname, | |||
fftype, | |||
xmlScript | |||
) |
Parse a text force field and create several important instance variables.
Each line is processed using the 'feed' method as implemented in the reader class. This essentially allows us to create the correct parameter identifier (pid), because the pid comes from more than the current line, it also depends on the section that we're in.
When 'PRM' or 'RPT' is encountered, we do several things:
Additionally, when 'PRM' is encountered:
When 'RPT' is encountered we don't expand the parameter vector because this parameter is a copy of an existing one. If the parameter identifier is preceded by MINUS_, we chop off the prefix but remember that the sign needs to be flipped.
Definition at line 564 of file forcefield.py.
def src.forcefield.FF.addff_xml | ( | self, | |
ffname | |||
) |
Parse an XML force field file and create important instance variables.
This was modeled after addff_txt, but XML and text files are fundamentally different, necessitating two different methods.
We begin with an _ElementTree object. We search through the tree for the 'parameterize' and 'parameter_repeat' keywords. Each time the keyword is encountered, we do the same four actions that I describe in addff_txt.
It's hard to specify precisely the location in an XML file to change a force field parameter. I can create a list of tree elements (essentially pointers to elements within a tree), but this method breaks down when I copy the tree because I have no way to refer to the copied tree elements. Fortunately, lxml gives me a way to represent a tree using a flat list, and my XML file 'locations' are represented using the positions in the list.
Definition at line 684 of file forcefield.py.
def src.forcefield.FF.assign_field | ( | self, | |
idx, | |||
pid, | |||
fnm, | |||
ln, | |||
pfld, | |||
mult, | |||
cmd = None |
|||
) |
Record the locations of a parameter in a txt file; [[file name, line number, field number, and multiplier]].
Note that parameters can have multiple locations because of the repetition functionality.
LPW 2019-07-08: Build a graph representation of the parameters where each node is a specific parameter in the force field file, and contains attributes such as the mathematical ID, file name, etc. This is an improved representation of the data contained in self.pfields.
[in] | idx | The (not necessarily unique) index of the parameter. |
[in] | pid | The unique parameter name. |
[in] | fnm | The file name of the parameter field. |
[in] | ln | The line number within the file (or the node index in the flattened xml) |
[in] | pfld | The field within the line (or the name of the attribute in the xml) |
[in] | mult | The multiplier (this is usually 1.0) |
Definition at line 1508 of file forcefield.py.
def src.forcefield.FF.assign_p0 | ( | self, | |
idx, | |||
val | |||
) |
Assign physical parameter values to the 'pvals0' array.
[in] | idx | The index to which we assign the parameter value. |
[in] | val | The parameter value to be inserted. |
Definition at line 1485 of file forcefield.py.
def src.forcefield.FF.check_dupes | ( | self, | |
pid, | |||
ffname, | |||
ln, | |||
pfld | |||
) |
Check to see whether a given parameter ID already exists, and provide an alternate if needed.
Definition at line 511 of file forcefield.py.
def src.forcefield.FF.create_mvals | ( | self, | |
pvals | |||
) |
Converts physical to mathematical parameters.
We create the inverse transformation matrix using SVD.
[in] | pvals | The physical parameters |
Definition at line 1052 of file forcefield.py.
def src.forcefield.FF.create_pvals | ( | self, | |
mvals | |||
) |
Converts mathematical to physical parameters.
First, mathematical parameters are rescaled and rotated by multiplying by the transformation matrix, followed by adding the original physical parameters.
[in] | mvals | The mathematical parameters |
Definition at line 1013 of file forcefield.py.
def src.forcefield.FF.find_spacings | ( | self | ) |
def src.forcefield.FF.fromfile | ( | cls, | |
fnm | |||
) |
Definition at line 346 of file forcefield.py.
def src.forcefield.FF.get_mathid | ( | self, | |
pid | |||
) |
For a given unique parameter identifier, return all of the mvals that it is connected to.
Definition at line 1547 of file forcefield.py.
def src.forcefield.FF.list_map | ( | self | ) |
Create the plist, which is like a reversed version of the parameter map.
More convenient for printing.
Definition at line 1453 of file forcefield.py.
def src.forcefield.FF.make | ( | self, | |
vals = None , |
|||
use_pvals = False , |
|||
printdir = None , |
|||
precision = 12 |
|||
) |
Create a new force field using provided parameter values.
This big kahuna does a number of things: 1) Creates the physical parameters from the mathematical parameters 2) Creates force fields with physical parameters substituted in 3) Prints the force fields to the specified file.
It does NOT store the mathematical parameters in the class state (since we can only hold one set of parameters).
[in] | printdir | The directory that the force fields are printed to; as usual this is relative to the project root directory. |
[in] | vals | Input parameters. I previously had an option where it uses stored values in the class state, but I don't think that's a good idea anymore. |
[in] | use_pvals | Switch for whether to bypass the coordinate transformation and use physical parameters directly. |
[in] | precision | Number of decimal points to print out |
Definition at line 780 of file forcefield.py.
def src.forcefield.FF.make_redirect | ( | self, | |
mvals | |||
) |
def src.forcefield.FF.make_rescale | ( | self, | |
scales, | |||
mvals = None , |
|||
G = None , |
|||
H = None , |
|||
multiply = True , |
|||
verbose = False |
|||
) |
Obtain rescaled versions of the inputs according to dictionary values in "scales" (i.e.
a replacement or multiplicative factor on self.rs_ord). Note that self.rs and self.rs_ord are not updated in this function. You need to do that outside.
The purpose of this function is to simulate the effect of changing the parameter scale factors in the force field. If the scale factor is N, then a unit change in the mathematical parameters produces a change of N in the physical parameters. Thus, for a given point in the physical parameter space, the mathematical parameters are proportional to 1/N, the gradient is proportional to N, and the Hessian is proportional to N^2.
mvals : numpy.ndarray Parameters to be transformed, if desired. Must be same length as number of parameters. G : numpy.ndarray Gradient to be transformed, if desired. Must be same length as number of parameters. H : numpy.ndarray Hessian to be transformed, if desired. Must be square matrix with side-length = number of parameters. scales : OrderedDict A dictionary with the same keys as self.rs_ord and floating point values. These represent the values with which to multiply the existing scale factors (if multiply == True) or replace them (if multiply == False). Pro Tip: Create this variable from a copy of self.rs_ord multiply : bool When set to True, the new scale factors are the existing scale factors verbose : bool Loud and noisy
answer : OrderedDict Output dictionary containing : 'rs' : New parameter scale factors (multiplied by scales if multiply=True, or replaced if multiply=False) 'rs_ord' : New parameter scale factor dictionary 'mvals' : New parameter values (if mvals is provided) 'G' : New gradient (if G is provided) 'H' : New Hessian (if H is provided)
Definition at line 1195 of file forcefield.py.
def src.forcefield.FF.mktransmat | ( | self | ) |
Create the transformation matrix to rescale and rotate the mathematical parameters.
For point charge parameters, project out perturbations that change the total charge.
First build these:
'qmap' : Just a list of parameter indices that point to charges.
'qid' : For each parameter in the qmap, a list of the affected atoms :) A potential target for the molecule-specific thang.
Then make this:
'qtrans2' : A transformation matrix that rotates the charge parameters. The first row is all zeros (because it corresponds to increasing the charge on all atoms) The other rows correspond to changing one of the parameters and decreasing all of the others equally such that the overall charge is preserved.
'qmat2' : An identity matrix with 'qtrans2' pasted into the right place
'transmat': 'qmat2' with rows and columns scaled using self.rs
'excision': Parameter indices that need to be 'cut out' because they are irrelevant and mess with the matrix diagonalization
Only project out changes in total charge of a molecule, and perhaps generalize to fragments of molecules or other types of parameters.
The AMOEBA selection of charge depends not only on the atom type, but what that atom is bonded to.
Build the matrix that ensures the net charge does not change.
Definition at line 1272 of file forcefield.py.
def src.forcefield.FF.plot_parameter_tree | ( | fnm | ) |
Definition at line 1533 of file forcefield.py.
def src.forcefield.FF.print_map | ( | self, | |
vals = None , |
|||
precision = 4 |
|||
) |
Prints out the (physical or mathematical) parameter indices, IDs and values in a visually appealing way.
Definition at line 1465 of file forcefield.py.
def src.forcefield.FF.rsmake | ( | self, | |
printfacs = True |
|||
) |
Create the rescaling factors for the coordinate transformation in parameter space.
The proper choice of rescaling factors (read: prior widths in maximum likelihood analysis) is still a black art. This is a topic of current research.
@param[in] printfacs List for printing out the resecaling factors
Definition at line 1072 of file forcefield.py.
def src.forcefield.FF.sprint_map | ( | self, | |
vals = None , |
|||
precision = 4 |
|||
) |
Prints out the (physical or mathematical) parameter indices, IDs and values to a string.
Definition at line 1473 of file forcefield.py.
src.forcefield.FF.amber_frcmod |
Definition at line 472 of file forcefield.py.
src.forcefield.FF.amber_mol2 |
Definition at line 468 of file forcefield.py.
src.forcefield.FF.atomnames |
A list of atom names (this is new, for ESP fitting)
Definition at line 287 of file forcefield.py.
src.forcefield.FF.excision |
Indices to exclude from optimization / Hessian inversion.
Some customized constraints here.
Quadrupoles must be traceless
Definition at line 279 of file forcefield.py.
src.forcefield.FF.FFAtomTypes |
WORK IN PROGRESS ## This is a dictionary of {'AtomType':{'Mass' : float, 'Charge' : float, 'ParticleType' : string ('A', 'S', or 'D'), 'AtomicNumber' : int}}.
Definition at line 297 of file forcefield.py.
src.forcefield.FF.ffdata |
As these options proliferate, the force field class becomes less standalone.
I need to think of a good solution here... The root directory of the project File names of force fields Directory containing force fields, relative to project directory Priors given by the user :) Whether to constrain the charges. Whether to constrain the charges. Switch for AMOEBA direct or mutual. AMOEBA mutual dipole convergence tolerance. Switch for rigid water molecules Switch for constraining bonds involving hydrogen Bypass the transformation and use physical parameters directly Allow duplicate parameter names (internally construct unique names) The content of all force field files are stored in memory
Definition at line 255 of file forcefield.py.
src.forcefield.FF.ffdata_isxml |
Definition at line 256 of file forcefield.py.
src.forcefield.FF.FFMolecules |
Definition at line 310 of file forcefield.py.
src.forcefield.FF.linedestroy_save |
Destruction dictionary (experimental).
Definition at line 338 of file forcefield.py.
src.forcefield.FF.linedestroy_this |
Definition at line 340 of file forcefield.py.
src.forcefield.FF.map |
The mapping of interaction type -> parameter number.
Definition at line 258 of file forcefield.py.
src.forcefield.FF.np |
The total number of parameters.
Definition at line 281 of file forcefield.py.
src.forcefield.FF.offxml |
Definition at line 464 of file forcefield.py.
src.forcefield.FF.offxml_unit_strs |
Definition at line 271 of file forcefield.py.
src.forcefield.FF.openff_forcefield |
Definition at line 465 of file forcefield.py.
src.forcefield.FF.openmmxml |
Definition at line 456 of file forcefield.py.
src.forcefield.FF.patoms |
A listing of parameter number -> atoms involved.
Definition at line 262 of file forcefield.py.
src.forcefield.FF.pfields |
A list where pfields[i] = [pid,'file',line,field,mult,cmd], basically a new way to modify force field files; when we modify the force field file, we go to the specific line/field in a given file and change the number.
Definition at line 267 of file forcefield.py.
src.forcefield.FF.plist |
The listing of parameter number -> interaction types.
Definition at line 260 of file forcefield.py.
src.forcefield.FF.prmdestroy_save |
Definition at line 339 of file forcefield.py.
src.forcefield.FF.prmdestroy_this |
Definition at line 341 of file forcefield.py.
src.forcefield.FF.pTree |
Improved representation of pfields as a networkx graph.
Definition at line 269 of file forcefield.py.
src.forcefield.FF.pvals0 |
Initial value of physical parameters.
Definition at line 283 of file forcefield.py.
src.forcefield.FF.qid |
Definition at line 1274 of file forcefield.py.
src.forcefield.FF.qid2 |
Definition at line 1275 of file forcefield.py.
src.forcefield.FF.qmap |
Definition at line 1273 of file forcefield.py.
src.forcefield.FF.Readers |
A dictionary of force field reader classes.
Definition at line 285 of file forcefield.py.
src.forcefield.FF.redirect |
Creates plist from map.
Prints the plist to screen. Make the rescaling factors. Make the transformation matrix. Redirection dictionary (experimental).
Definition at line 336 of file forcefield.py.
src.forcefield.FF.rs |
List of rescaling factors.
The array of rescaling factors.
Definition at line 273 of file forcefield.py.
src.forcefield.FF.rs_ord |
Takes the dictionary 'BONDS':{3:'B', 4:'K'}, 'VDW':{4:'S', 5:'T'}, and turns it into a list of term types ['BONDSB','BONDSK','VDWS','VDWT'].
Definition at line 1130 of file forcefield.py.
src.forcefield.FF.rs_type |
Definition at line 1133 of file forcefield.py.
src.forcefield.FF.tinkerprm |
Definition at line 449 of file forcefield.py.
src.forcefield.FF.tm |
The transformation matrix for mathematical -> physical parameters.
Definition at line 275 of file forcefield.py.
src.forcefield.FF.tmI |
The transpose of the transformation matrix.
Definition at line 277 of file forcefield.py.