ForceBalance API  1.3
Automated optimization of force fields and empirical potentials
tinkerio.py
Go to the documentation of this file.
1 """ @package forcebalance.tinkerio TINKER input/output.
2 
3 This serves as a good template for writing future force matching I/O
4 modules for other programs because it's so simple.
5 
6 @author Lee-Ping Wang
7 @date 01/2012
8 """
9 from __future__ import division
10 from __future__ import print_function
11 
12 from builtins import str
13 from builtins import zip
14 from builtins import range
15 import os, shutil
16 from re import match, sub
17 from forcebalance.nifty import *
18 from forcebalance.nifty import _exec
19 import time
20 import numpy as np
21 import networkx as nx
22 from copy import deepcopy
23 from forcebalance import BaseReader
24 from subprocess import Popen, PIPE
25 from forcebalance.engine import Engine
26 from forcebalance.abinitio import AbInitio
27 from forcebalance.vibration import Vibration
28 from forcebalance.moments import Moments
29 from forcebalance.liquid import Liquid
30 from forcebalance.molecule import Molecule, BuildLatticeFromLengthsAngles
31 from forcebalance.binding import BindingEnergy
32 from forcebalance.interaction import Interaction
33 from forcebalance.finite_difference import in_fd
34 from collections import OrderedDict
35 
36 # All TINKER force field parameter types, which should eventually go into pdict
37 # at some point (for full compatibility).
38 allp = ['atom', 'vdw', 'vdw14', 'vdwpr', 'hbond', 'bond', 'bond5', 'bond4',
39  'bond3', 'electneg', 'angle', 'angle5', 'angle4', 'angle3', 'anglef',
40  'strbnd', 'ureybrad', 'angang', 'opbend', 'opdist', 'improper', 'imptors',
41  'torsion', 'torsion5', 'torsion4', 'pitors', 'strtors', 'tortors', 'charge',
42  'dipole', 'dipole5', 'dipole4', 'dipole3', 'multipole', 'polarize', 'piatom',
43  'pibond', 'pibond5', 'pibond4', 'metal', 'biotype', 'mmffvdw', 'mmffbond',
44  'mmffbonder', 'mmffangle', 'mmffstrbnd', 'mmffopbend', 'mmfftorsion', 'mmffbci',
45  'mmffpbci', 'mmffequiv', 'mmffdefstbn', 'mmffcovrad', 'mmffprop', 'mmffarom']
46 
47 # All possible output from analyze's energy component breakdown.
48 eckeys = ['Angle-Angle', 'Angle Bending', 'Atomic Multipoles', 'Bond Stretching', 'Charge-Charge',
49  'Charge-Dipole', 'Dipole-Dipole', 'Extra Energy Terms', 'Geometric Restraints', 'Implicit Solvation',
50  'Improper Dihedral', 'Improper Torsion', 'Metal Ligand Field', 'Out-of-Plane Bend', 'Out-of-Plane Distance',
51  'Pi-Orbital Torsion', 'Polarization', 'Reaction Field', 'Stretch-Bend', 'Stretch-Torsion',
52  'Torsional Angle', 'Torsion-Torsion', 'Urey-Bradley', 'Van der Waals']
53 
54 from forcebalance.output import getLogger
55 logger = getLogger(__name__)
56 
57 pdict = {'VDW' : {'Atom':[1], 2:'S',3:'T',4:'D'}, # Van der Waals distance, well depth, distance from bonded neighbor?
58  'BOND' : {'Atom':[1,2], 3:'K',4:'B'}, # Bond force constant and equilibrium distance (Angstrom)
59  'ANGLE' : {'Atom':[1,2,3], 4:'K',5:'B'}, # Angle force constant and equilibrium angle
60  'STRBND' : {'Atom':[1,2,3], 4:'K1',5:'K2'}, # Two stretch-bend force constants (usually same)
61  'OPBEND' : {'Atom':[1,2,3,4], 5:'K'}, # Out-of-plane bending force constant
62  'UREYBRAD' : {'Atom':[1,2,3], 4:'K',5:'B'}, # Urey-Bradley force constant and equilibrium distance (Angstrom)
63  'TORSION' : ({'Atom':[1,2,3,4], 5:'1K', 6:'1B',
64  8:'2K', 9:'2B', 11:'3K', 12:'3B'}), # Torsional force constants and equilibrium phi-angles
65  'PITORS' : {'Atom':[1,2], 3:'K'}, # Pi-torsion force constants (usually 6.85 ..)
66  'CHARGE' : {'Atom':[1], 2:''}, # Atomic partial charge (OPLS style)
67  # Note torsion-torsion (CMAP) not implemented at this time.
68  'MCHARGE' : {'Atom':[1,2,3], 4:''}, # Atomic charge
69  'DIPOLE' : {0:'X',1:'Y',2:'Z'}, # Dipole moment in local frame
70  'QUADX' : {0:'X'}, # Quadrupole moment, X component
71  'QUADY' : {0:'X',1:'Y'}, # Quadrupole moment, Y component
72  'QUADZ' : {0:'X',1:'Y',2:'Z'}, # Quadrupole moment, Y component
73  'POLARIZE' : {'Atom':[1], 2:'A',3:'T'}, # Atomic dipole polarizability
74  'BOND-CUBIC' : {'Atom':[], 0:''}, # Below are global parameters.
75  'BOND-QUARTIC' : {'Atom':[], 0:''},
76  'ANGLE-CUBIC' : {'Atom':[], 0:''},
77  'ANGLE-QUARTIC': {'Atom':[], 0:''},
78  'ANGLE-PENTIC' : {'Atom':[], 0:''},
79  'ANGLE-SEXTIC' : {'Atom':[], 0:''},
80  'DIELECTRIC' : {'Atom':[], 0:''},
81  'POLAR-SOR' : {'Atom':[], 0:''}
82  # Ignored for now: stretch/bend coupling, out-of-plane bending,
83  # torsional parameters, pi-torsion, torsion-torsion
84  }
85 
86 class Tinker_Reader(BaseReader):
87  """Finite state machine for parsing TINKER force field files.
88 
89  This class is instantiated when we begin to read in a file. The
90  feed(line) method updates the state of the machine, informing it
91  of the current interaction type. Using this information we can
92  look up the interaction type and parameter type for building the
93  parameter ID.
94 
95  """
96 
97  def __init__(self,fnm):
98  super(Tinker_Reader,self).__init__(fnm)
99 
100  self.pdict = pdict
101 
102  self.atom = []
103 
104  def feed(self, line):
105  """ Given a line, determine the interaction type and the atoms involved (the suffix).
106 
107  TINKER generally has stuff like this:
108 
109  @verbatim
110  bond-cubic -2.55
111  bond-quartic 3.793125
112 
113  vdw 1 3.4050 0.1100
114  vdw 2 2.6550 0.0135 0.910 # PRM 4
115 
116  multipole 2 1 2 0.25983
117  -0.03859 0.00000 -0.05818
118  -0.03673
119  0.00000 -0.10739
120  -0.00203 0.00000 0.14412
121  @endverbatim
122 
123  The '#PRM 4' has no effect on TINKER but it indicates that we
124  are tuning the fourth field on the line (the 0.910 value).
125 
126  @todo Put the rescaling factors for TINKER parameters in here.
127  Currently we're using the initial value to determine the
128  rescaling factor which is not very good.
129 
130  Every parameter line is prefaced by the interaction type
131  except for 'multipole' which is on multiple lines. Because
132  the lines that come after 'multipole' are predictable, we just
133  determine the current line using the previous line.
134 
135  Random note: Unit of force is kcal / mole / angstrom squared.
136 
137  """
138  s = line.split()
139  self.ln += 1
140  # No sense in doing anything for an empty line or a comment line.
141  if len(s) == 0 or match('^#',line): return None, None
142  # From the line, figure out the interaction type. If this line doesn't correspond
143  # to an interaction, then do nothing.
144  if s[0].upper() in pdict:
145  self.itype = s[0].upper()
146  # This is kind of a silly hack that allows us to take care of the 'multipole' keyword,
147  # because of the syntax that the TINKER .prm file uses.
148  elif s[0].upper() == 'MULTIPOLE':
149  self.itype = 'MCHARGE'
150  elif self.itype == 'MCHARGE':
151  self.itype = 'DIPOLE'
152  elif self.itype == 'DIPOLE':
153  self.itype = 'QUADX'
154  elif self.itype == 'QUADX':
155  self.itype = 'QUADY'
156  elif self.itype == 'QUADY':
157  self.itype = 'QUADZ'
158  else:
159  self.itype = None
160 
161  if self.itype in pdict:
162  if 'Atom' in pdict[self.itype]:
163  # List the atoms in the interaction.
164  self.atom = [s[i] for i in pdict[self.itype]['Atom']]
165  # The suffix of the parameter ID is built from the atom #
166  # types/classes involved in the interaction.
167  self.suffix = '/'+'.'.join(self.atom)
168 
169 def write_key(fout, options, fin=None, defaults={}, verbose=False, prmfnm=None, chk=[]):
170  """
171  Create or edit a TINKER .key file.
172  @param[in] fout Output file name, can be the same as input file name.
173  @param[in] options Dictionary containing .key options. Existing options are replaced, new options are added at the end.
174  Passing None causes options to be deleted. To pass an option without an argument, use ''.
175  @param[in] fin Input file name.
176  @param[in] defaults Default options to add to the mdp only if they don't already exist.
177  @param[in] verbose Print out all modifications to the file.
178  @param[in] prmfnm TINKER parameter file name.
179  @param[in] chk Crash if the key file does NOT have these options by the end.
180  """
181  # Make sure that the keys are lowercase, and the values are all strings.
182  options = OrderedDict([(key.lower(), str(val) if val is not None else None) for key, val in options.items()])
183  if 'parameters' in options and prmfnm is not None:
184  logger.error("Please pass prmfnm or 'parameters':'filename.prm' in options but not both.\n")
185  raise RuntimeError
186  elif 'parameters' in options:
187  prmfnm = options['parameters']
188 
189  if prmfnm is not None:
190  # Account for both cases where the file name may end with .prm
191  prms = [prmfnm]
192  if prms[0].endswith('.prm'):
193  prms.append(prmfnm[:-4])
194  else:
195  prms = []
196 
197  # Options that cause the program to crash if they are overwritten
198  clashes = []
199  # List of lines in the output file.
200  out = []
201  # List of options in the output file.
202  haveopts = []
203  skip = 0
204  prmflag = 0
205  if fin is not None and os.path.isfile(fin):
206  for line0 in open(fin).readlines():
207  line1 = line0.replace('\n','').expandtabs()
208  line = line0.strip().expandtabs()
209  s = line.split()
210  if skip:
211  out.append(line1)
212  skip -= 1
213  continue
214  # Skip over these three cases:
215  # 1) Empty lines get appended to the output and skipped.
216  if len(line) == 0 or set(line).issubset([' ']):
217  out.append('')
218  continue
219  # 2) Lines that start with comments are skipped as well.
220  if line.startswith("#"):
221  out.append(line1)
222  continue
223  # 3) Lines that correspond to force field parameters are skipped
224  if s[0].lower() in allp:
225  out.append(line1)
226  # 3a) For AMOEBA multipole parameters, skip four additional lines
227  if s[0].lower() == "multipole":
228  skip = 4
229  continue
230  # Now split by the comment character.
231  s = line.split('#',1)
232  data = s[0]
233  comms = s[1] if len(s) > 1 else None
234  # Now split off the key and value fields at the space.
235  ds = data.split(' ',1)
236  keyf = ds[0]
237  valf = ds[1] if len(ds) > 1 else ''
238  key = keyf.strip().lower()
239  haveopts.append(key)
240  if key == 'parameters':
241  val0 = valf.strip()
242  if val0 == '':
243  warn_press_key("Expected a parameter file name but got none")
244  # This is the case where "parameters" correctly corresponds to optimize.in
245  prmflag = 1
246  if prmfnm is None or val0 in prms or val0[:-4] in prms:
247  out.append(line1)
248  continue
249  else:
250  logger.info(line + '\n')
251  warn_press_key("The above line was found in %s, but we expected something like 'parameters %s'; replacing." % (line,prmfnm))
252  options['parameters'] = prmfnm
253  if key in options:
254  # This line replaces the line in the .key file with the value provided in the dictionary.
255  val = options[key]
256  val0 = valf.strip()
257  if key in clashes and val != val0:
258  logger.error("write_key tried to set %s = %s but its original value was %s = %s\n" % (key, val, key, val0))
259  raise RuntimeError
260  # Passing None as the value causes the option to be deleted
261  if val is None:
262  continue
263  if len(val) < len(valf):
264  valf = ' ' + val + ' '*(len(valf) - len(val)-1)
265  else:
266  valf = ' ' + val + ' '
267  lout = [keyf, ' ', valf]
268  if comms is not None:
269  lout += ['#',comms]
270  out.append(''.join(lout))
271  else:
272  out.append(line1)
273  # Options that don't already exist are written at the bottom.
274  for key, val in options.items():
275  key = key.lower()
276  if val is None: continue
277  if key not in haveopts:
278  haveopts.append(key)
279  out.append("%-20s %s" % (key, val))
280  # Fill in default options.
281  for key, val in defaults.items():
282  key = key.lower()
283  options[key] = val
284  if key not in haveopts:
285  haveopts.append(key)
286  out.append("%-20s %s" % (key, val))
287  # If parameters are not specified, they are printed at the top.
288  if not prmflag and prmfnm is not None:
289  out.insert(0,"parameters %s" % prmfnm)
290  options["parameters"] = prmfnm
291  elif not prmflag:
292  if not os.path.exists('%s.prm' % os.path.splitext(fout)[0]):
293  logger.error('No parameter file detected, this will cause TINKER to crash\n')
294  raise RuntimeError
295  for i in chk:
296  if i not in haveopts:
297  logger.error('%s is expected to be in the .key file, but not found\n' % i)
298  raise RuntimeError
299  # Finally write the key file.
300  file_out = wopen(fout)
301  for line in out:
302  print(line, file=file_out)
303  if verbose:
304  printcool_dictionary(options, title="%s -> %s with options:" % (fin, fout))
305  file_out.close()
306 
307 class TINKER(Engine):
308 
309  """ Engine for carrying out general purpose TINKER calculations. """
310 
311  def __init__(self, name="tinker", **kwargs):
312 
313  self.valkwd = ['tinker_key', 'tinkerpath', 'tinker_prm']
314  self.warn_vn = False
315  super(TINKER,self).__init__(name=name, **kwargs)
317  def setopts(self, **kwargs):
318 
319  """ Called by __init__ ; Set TINKER-specific options. """
320 
321 
322  if 'tinkerpath' in kwargs:
323  self.tinkerpath = kwargs['tinkerpath']
324  if not os.path.exists(os.path.join(self.tinkerpath,"dynamic")):
325  warn_press_key("The 'dynamic' executable indicated by %s doesn't exist! (Check tinkerpath)" \
326  % os.path.join(self.tinkerpath,"dynamic"))
327  else:
328  warn_once("The 'tinkerpath' option was not specified; using default.")
329  if which('mdrun') == '':
330  warn_press_key("Please add TINKER executables to the PATH or specify tinkerpath.")
331  self.tinkerpath = which('dynamic')
332 
333  def readsrc(self, **kwargs):
334 
335  """ Called by __init__ ; read files from the source directory. """
336 
337  self.key = onefile(kwargs.get('tinker_key'), 'key')
338  self.prm = onefile(kwargs.get('tinker_prm'), 'prm')
339  if 'mol' in kwargs:
340  self.mol = kwargs['mol']
341  else:
342  crdfile = onefile(kwargs.get('coords'), 'arc', err=True)
343  self.mol = Molecule(crdfile)
345  def calltinker(self, command, stdin=None, print_to_screen=False, print_command=False, **kwargs):
346 
347  """ Call TINKER; prepend the tinkerpath to calling the TINKER program. """
348 
349  csplit = command.split()
350  # Sometimes the engine changes dirs and the key goes missing, so we link it.
351  if "%s.key" % self.name in csplit and not os.path.exists("%s.key" % self.name):
352  LinkFile(self.abskey, "%s.key" % self.name)
353  prog = os.path.join(self.tinkerpath, csplit[0])
354  csplit[0] = prog
355  o = _exec(' '.join(csplit), stdin=stdin, print_to_screen=print_to_screen, print_command=print_command, rbytes=1024, **kwargs)
356  # Determine the TINKER version number.
357  for line in o[:10]:
358  if "Version" in line:
359  vw = line.split()[2]
360  if len(vw.split('.')) <= 2:
361  vn = float(vw)
362  else:
363  vn = float(vw.split('.')[:2])
364  vn_need = 6.3
365  try:
366  if vn < vn_need:
367  if self.warn_vn:
368  warn_press_key("ForceBalance requires TINKER %.1f - unexpected behavior with older versions!" % vn_need)
369  self.warn_vn = True
370  except:
371  logger.error("Unable to determine TINKER version number!\n")
372  raise RuntimeError
373  for line in o[-10:]:
374  # Catch exceptions since TINKER does not have exit status.
375  if "TINKER is Unable to Continue" in line:
376  for l in o:
377  logger.error("%s\n" % l)
378  time.sleep(1)
379  logger.error("TINKER may have crashed! (See above output)\nThe command was: %s\nThe directory was: %s\n" % (' '.join(csplit), os.getcwd()))
380  raise RuntimeError
381  break
382  for line in o:
383  if 'D+' in line:
384  logger.info(line+'\n')
385  warn_press_key("TINKER returned a very large floating point number! (See above line; will give error on parse)")
386  return o
387 
388  def prepare(self, pbc=False, **kwargs):
389 
390  """ Called by __init__ ; prepare the temp directory and figure out the topology. """
391 
392  # Call TINKER but do nothing to figure out the version number.
393  o = self.calltinker("dynamic", persist=1, print_error=False)
394 
395  self.rigid = False
396 
397 
398  tk_chk = []
399  tk_opts = OrderedDict([("digits", "10"), ("archive", "")])
400  tk_defs = OrderedDict()
401 
402  prmtmp = False
403 
404  if hasattr(self,'FF'):
405  if not os.path.exists(self.FF.tinkerprm):
406  # If the parameter files don't already exist, create them for the purpose of
407  # preparing the engine, but then delete them afterward.
408  prmtmp = True
409  self.FF.make(np.zeros(self.FF.np))
410  if self.FF.rigid_water:
411  tk_opts["rattle"] = "water"
412  self.rigid = True
413  if self.FF.amoeba_pol == 'mutual':
414  tk_opts['polarization'] = 'mutual'
415  if self.FF.amoeba_eps is not None:
416  tk_opts['polar-eps'] = str(self.FF.amoeba_eps)
417  else:
418  tk_defs['polar-eps'] = '1e-6'
419  elif self.FF.amoeba_pol == 'direct':
420  tk_opts['polarization'] = 'direct'
421  else:
422  warn_press_key("Using TINKER without explicitly specifying AMOEBA settings. Are you sure?")
423  self.prm = self.FF.tinkerprm
424  prmfnm = self.FF.tinkerprm
425  elif self.prm:
426  prmfnm = self.prm
427  else:
428  prmfnm = None
429 
430  # Periodic boundary conditions may come from the TINKER .key file.
431  keypbc = False
432  minbox = 1e10
433  if self.key:
434  for line in open(os.path.join(self.srcdir, self.key)).readlines():
435  s = line.split()
436  if len(s) > 0 and s[0].lower() == 'a-axis':
437  keypbc = True
438  minbox = float(s[1])
439  if len(s) > 0 and s[0].lower() == 'b-axis' and float(s[1]) < minbox:
440  minbox = float(s[1])
441  if len(s) > 0 and s[0].lower() == 'c-axis' and float(s[1]) < minbox:
442  minbox = float(s[1])
443  if keypbc and (not pbc):
444  warn_once("Deleting PBC options from the .key file.")
445  tk_opts['a-axis'] = None
446  tk_opts['b-axis'] = None
447  tk_opts['c-axis'] = None
448  tk_opts['alpha'] = None
449  tk_opts['beta'] = None
450  tk_opts['gamma'] = None
451  if pbc:
452  if 'boxes' in self.mol.Data:
453  minbox = min([self.mol.boxes[0].a, self.mol.boxes[0].b, self.mol.boxes[0].c])
454  if (not keypbc) and 'boxes' not in self.mol.Data:
455  logger.error("Periodic boundary conditions require either (1) a-axis to be in the .key file or (b) boxes to be in the coordinate file.\n")
456  raise RuntimeError
457  self.pbc = pbc
458  if pbc:
459  tk_opts['ewald'] = ''
460  if minbox <= 10:
461  warn_press_key("Periodic box is set to less than 10 Angstroms across")
462  # TINKER likes to use up to 7.0 Angstrom for PME cutoffs
463  rpme = 0.5*(float(int(minbox - 1))) if minbox <= 15 else 7.0
464  if 'nonbonded_cutoff' in kwargs:
465  rpme = kwargs['nonbonded_cutoff']
466  if rpme > 0.5*(float(int(minbox - 1))):
467  warn_press_key("nonbonded_cutoff = %.1f should be smaller than half the box size = %.1f Angstrom" % (rpme, minbox))
468  tk_defs['ewald-cutoff'] = "%f" % rpme
469  # TINKER likes to use up to 9.0 Angstrom for vdW cutoffs
470  rvdw = 0.5*(float(int(minbox - 1))) if minbox <= 19 else 9.0
471  if 'nonbonded_cutoff' in kwargs and 'vdw_cutoff' not in kwargs:
472  warn_press_key('AMOEBA detected and nonbonded_cutoff is set, but not vdw_cutoff (so it will be set equal to nonbonded_cutoff)')
473  rvdw = kwargs['nonbonded_cutoff']
474  if 'vdw_cutoff' in kwargs:
475  rvdw = kwargs['vdw_cutoff']
476  if rvdw > 0.5*(float(int(minbox - 1))):
477  warn_press_key("vdw_cutoff = %.1f should be smaller than half the box size = %.1f Angstrom" % (rvdw, minbox))
478  tk_defs['vdw-cutoff'] = "%f" % rvdw
479  if (minbox*0.5 - rpme) > 2.5 and (minbox*0.5 - rvdw) > 2.5:
480  tk_defs['neighbor-list'] = ''
481  elif (minbox*0.5 - rpme) > 2.5:
482  tk_defs['mpole-list'] = ''
483  else:
484  if 'nonbonded_cutoff' in kwargs or 'vdw_cutoff' in kwargs:
485  warn_press_key('No periodic boundary conditions, your provided nonbonded_cutoff and vdw_cutoff will not be used')
486  tk_opts['ewald'] = None
487  tk_opts['ewald-cutoff'] = None
488  tk_opts['vdw-cutoff'] = None
489  # This seems to have no effect on the kinetic energy.
490  # tk_opts['remove-inertia'] = '0'
491 
492  write_key("%s.key" % self.name, tk_opts, os.path.join(self.srcdir, self.key) if self.key else None, tk_defs, verbose=False, prmfnm=prmfnm)
493  self.abskey = os.path.abspath("%s.key")
494 
495  self.mol[0].write(os.path.join("%s.xyz" % self.name), ftype="tinker")
496 
497 
498  self.mol.require('tinkersuf')
499 
500 
501  o = self.calltinker("analyze %s.xyz P,C" % (self.name), stdin="ALL")
503 
504  mode = 0
505  self.AtomMask = []
506  self.AtomLists = defaultdict(list)
507  ptype_dict = {'atom': 'A', 'vsite': 'D'}
508  G = nx.Graph()
509  for line in o:
510  s = line.split()
511  if len(s) == 0: continue
512  # TINKER 8.2.1 -> 8.7.1 changed printout to "Atom Definition Parameters"
513  if "Atom Type Definition Parameters" in line or "Atom Definition Parameters" in line:
514  mode = 1
515  if mode == 1:
516  if isint(s[0]): mode = 2
517  if mode == 2:
518  if isint(s[0]):
519  G.add_node(int(s[0]))
520  mass = float(s[5])
521  self.AtomLists['Mass'].append(mass)
522  if mass < 1.0:
523  # Particles with mass less than one count as virtual sites.
524  self.AtomLists['ParticleType'].append('D')
525  else:
526  self.AtomLists['ParticleType'].append('A')
527  self.AtomMask.append(mass >= 1.0)
528  else:
529  mode = 0
530  if "List of 1-2 Connected Atomic Interactions" in line:
531  mode = 3
532  if mode == 3:
533  if isint(s[0]): mode = 4
534  if mode == 4:
535  if isint(s[0]):
536  a = int(s[0])
537  b = int(s[1])
538  G.add_edge(a, b)
539  else: mode = 0
540  # Use networkx to figure out a list of molecule numbers.
541  if len(list(G.nodes())) > 0:
542  # The following code only works in TINKER 6.2
543  gs = [G.subgraph(c).copy() for c in nx.connected_components(G)]
544  # Deprecated in networkx 2.2
545  # gs = list(nx.connected_component_subgraphs(G))
546  tmols = [gs[i] for i in np.argsort(np.array([min(list(g.nodes())) for g in gs]))]
547  mnodes = [list(m.nodes()) for m in tmols]
548  self.AtomLists['MoleculeNumber'] = [[i+1 in m for m in mnodes].index(1) for i in range(self.mol.na)]
549  else:
550  grouped = [i.L() for i in self.mol.molecules]
551  self.AtomLists['MoleculeNumber'] = [[i in g for g in grouped].index(1) for i in range(self.mol.na)]
552  if prmtmp:
553  for f in self.FF.fnms:
554  os.unlink(f)
555 
556  def optimize(self, shot, method="newton", align=True, crit=1e-4):
557 
558  """ Optimize the geometry and align the optimized geometry to the starting geometry. """
559 
560  if os.path.exists('%s.xyz_2' % self.name):
561  os.unlink('%s.xyz_2' % self.name)
562 
563  self.mol[shot].write('%s.xyz' % self.name, ftype="tinker")
564 
565  if method == "newton":
566  if self.rigid: optprog = "optrigid"
567  else: optprog = "optimize"
568  elif method == "bfgs":
569  if self.rigid: optprog = "minrigid"
570  else: optprog = "minimize"
571  else:
572  raise RuntimeError("Incorrect choice of method for TINKER.optimize()")
573 
574  # Actually run the minimize program
575  o = self.calltinker("%s %s.xyz %f" % (optprog, self.name, crit))
576 
577 
597 
598  # Silently align the optimized geometry.
599  M12 = Molecule("%s.xyz" % self.name, ftype="tinker") + Molecule("%s.xyz_2" % self.name, ftype="tinker")
600  if not self.pbc:
601  M12.align(center=False)
602  M12[1].write("%s.xyz_2" % self.name, ftype="tinker")
603  # print("LPW copying %s.xyz2" % self.name)
604  # shutil.copy2("%s.xyz_2" % self.name, "bak.xyz")
605  rmsd = M12.ref_rmsd(0)[1]
606  cnvgd = 0
607  mode = 0
608  for line in o:
609  s = line.split()
610  if len(s) == 0: continue
611  if "Optimally Conditioned Variable Metric Optimization" in line: mode = 1
612  if "Limited Memory BFGS Quasi-Newton Optimization" in line: mode = 1
613  if mode == 1 and isint(s[0]): mode = 2
614  if mode == 2:
615  if isint(s[0]): E = float(s[1])
616  else: mode = 0
617  if "Normal Termination" in line:
618  cnvgd = 1
619  if not cnvgd:
620  for line in o:
621  logger.info(str(line) + '\n')
622  logger.info("The minimization did not converge in the geometry optimization - printout is above.\n")
623  return E, rmsd, M12[1]
624 
625  def evaluate_(self, xyzin, force=False, dipole=False):
626 
627  """
628  Utility function for computing energy, and (optionally) forces and dipoles using TINKER.
629 
630  Inputs:
631  xyzin: TINKER .xyz file name.
632  force: Switch for calculating the force.
633  dipole: Switch for calculating the dipole.
634 
635  Outputs:
636  Result: Dictionary containing energies, forces and/or dipoles.
637  """
638 
639  Result = OrderedDict()
640  # If we want the dipoles (or just energies), analyze is the way to go.
641  if dipole or (not force):
642  oanl = self.calltinker("analyze %s -k %s" % (xyzin, self.name), stdin="G,E,M", print_to_screen=False)
643  # Read potential energy and dipole from file.
644  eanl = []
645  dip = []
646  for line in oanl:
647  s = line.split()
648  if 'Total Potential Energy : ' in line:
649  eanl.append(float(s[4]) * 4.184)
650  if dipole:
651  if 'Dipole X,Y,Z-Components :' in line:
652  dip.append([float(s[i]) for i in range(-3,0)])
653  Result["Energy"] = np.array(eanl)
654  Result["Dipole"] = np.array(dip)
655  # If we want forces, then we need to call testgrad.
656  if force:
657  E = []
658  F = []
659  Fi = []
660  o = self.calltinker("testgrad %s -k %s y n n" % (xyzin, self.name))
661  i = 0
662  ReadFrc = 0
663  for line in o:
664  s = line.split()
665  if "Total Potential Energy" in line:
666  E.append(float(s[4]) * 4.184)
667  if "Cartesian Gradient Breakdown over Individual Atoms" in line:
668  ReadFrc = 1
669  if ReadFrc and len(s) == 6 and all([s[0] == 'Anlyt',isint(s[1]),isfloat(s[2]),isfloat(s[3]),isfloat(s[4]),isfloat(s[5])]):
670  ReadFrc = 2
671  if self.AtomMask[i]:
672  Fi += [-1 * float(j) * 4.184 * 10 for j in s[2:5]]
673  i += 1
674  if ReadFrc == 2 and len(s) < 6:
675  ReadFrc = 0
676  F.append(Fi)
677  Fi = []
678  i = 0
679  Result["Energy"] = np.array(E)
680  Result["Force"] = np.array(F)
681  return Result
682 
683  def get_charges(self):
684  logger.error('TINKER engine does not have get_charges (should be easy to implement however.)')
685  raise NotImplementedError
686 
687  def energy_force_one(self, shot):
688 
689  """ Computes the energy and force using TINKER for one snapshot. """
690 
691  self.mol[shot].write("%s.xyz" % self.name, ftype="tinker")
692  Result = self.evaluate_("%s.xyz" % self.name, force=True)
693  return np.hstack((Result["Energy"].reshape(-1,1), Result["Force"]))
695  def energy(self):
696 
697  """ Computes the energy using TINKER over a trajectory. """
698 
699  if hasattr(self, 'md_trajectory') :
700  x = self.md_trajectory
701  else:
702  x = "%s.xyz" % self.name
703  self.mol.write(x, ftype="tinker")
704  return self.evaluate_(x)["Energy"]
705 
706  def energy_force(self):
707 
708  """ Computes the energy and force using TINKER over a trajectory. """
710  if hasattr(self, 'md_trajectory') :
711  x = self.md_trajectory
712  else:
713  x = "%s.xyz" % self.name
714  self.mol.write(x, ftype="tinker")
715  Result = self.evaluate_(x, force=True)
716  return np.hstack((Result["Energy"].reshape(-1,1), Result["Force"]))
717 
718  def energy_dipole(self):
719 
720  """ Computes the energy and dipole using TINKER over a trajectory. """
722  if hasattr(self, 'md_trajectory') :
723  x = self.md_trajectory
724  else:
725  x = "%s.xyz" % self.name
726  self.mol.write(x, ftype="tinker")
727  Result = self.evaluate_(x, dipole=True)
728  return np.hstack((Result["Energy"].reshape(-1,1), Result["Dipole"]))
729 
730  def normal_modes(self, shot=0, optimize=True):
731  # This line actually runs TINKER
732  if optimize:
733  self.optimize(shot, crit=1e-6)
734  o = self.calltinker("vibrate %s.xyz_2 a" % (self.name))
735  else:
736  warn_once("Asking for normal modes without geometry optimization?")
737  self.mol[shot].write('%s.xyz' % self.name, ftype="tinker")
738  o = self.calltinker("vibrate %s.xyz a" % (self.name))
739  # Read the TINKER output. The vibrational frequencies are ordered.
740  # The six modes with frequencies closest to zero are ignored
741  readev = False
742  calc_eigvals = []
743  calc_eigvecs = []
744  for line in o:
745  s = line.split()
746  if "Vibrational Normal Mode" in line:
747  freq = float(s[-2])
748  readev = False
749  calc_eigvals.append(freq)
750  calc_eigvecs.append([])
751  elif "Atom" in line and "Delta X" in line:
752  readev = True
753  elif readev and len(s) == 4 and all([isint(s[0]), isfloat(s[1]), isfloat(s[2]), isfloat(s[3])]):
754  calc_eigvecs[-1].append([float(i) for i in s[1:]])
755  calc_eigvals = np.array(calc_eigvals)
756  calc_eigvecs = np.array(calc_eigvecs)
757  # Sort by frequency absolute value and discard the six that are closest to zero
758  calc_eigvecs = calc_eigvecs[np.argsort(np.abs(calc_eigvals))][6:]
759  calc_eigvals = calc_eigvals[np.argsort(np.abs(calc_eigvals))][6:]
760  # Sort again by frequency
761  calc_eigvecs = calc_eigvecs[np.argsort(calc_eigvals)]
762  calc_eigvals = calc_eigvals[np.argsort(calc_eigvals)]
763  os.system("rm -rf *.xyz_* *.[0-9][0-9][0-9]")
764  return calc_eigvals, calc_eigvecs
765 
766  def multipole_moments(self, shot=0, optimize=True, polarizability=False):
767 
768  """ Return the multipole moments of the 1st snapshot in Debye and Buckingham units. """
769 
770  # This line actually runs TINKER
771  if optimize:
772  self.optimize(shot, crit=1e-6)
773  o = self.calltinker("analyze %s.xyz_2 M" % (self.name))
774  else:
775  self.mol[shot].write('%s.xyz' % self.name, ftype="tinker")
776  o = self.calltinker("analyze %s.xyz M" % (self.name))
777  # Read the TINKER output.
778  qn = -1
779  ln = 0
780  for line in o:
781  s = line.split()
782  if "Dipole X,Y,Z-Components" in line:
783  dipole_dict = OrderedDict(zip(['x','y','z'], [float(i) for i in s[-3:]]))
784  elif "Quadrupole Moment Tensor" in line:
785  qn = ln
786  quadrupole_dict = OrderedDict([('xx',float(s[-3]))])
787  elif qn > 0 and ln == qn + 1:
788  quadrupole_dict['xy'] = float(s[-3])
789  quadrupole_dict['yy'] = float(s[-2])
790  elif qn > 0 and ln == qn + 2:
791  quadrupole_dict['xz'] = float(s[-3])
792  quadrupole_dict['yz'] = float(s[-2])
793  quadrupole_dict['zz'] = float(s[-1])
794  ln += 1
795 
796  calc_moments = OrderedDict([('dipole', dipole_dict), ('quadrupole', quadrupole_dict)])
797 
798  if polarizability:
799  if optimize:
800  o = self.calltinker("polarize %s.xyz_2" % (self.name))
801  else:
802  o = self.calltinker("polarize %s.xyz" % (self.name))
803  # Read the TINKER output.
804  pn = -1
805  ln = 0
806  polarizability_dict = OrderedDict()
807  for line in o:
808  s = line.split()
809  if "Molecular Polarizability Tensor" in line:
810  pn = ln
811  elif pn > 0 and ln == pn + 2:
812  polarizability_dict['xx'] = float(s[-3])
813  polarizability_dict['yx'] = float(s[-2])
814  polarizability_dict['zx'] = float(s[-1])
815  elif pn > 0 and ln == pn + 3:
816  polarizability_dict['xy'] = float(s[-3])
817  polarizability_dict['yy'] = float(s[-2])
818  polarizability_dict['zy'] = float(s[-1])
819  elif pn > 0 and ln == pn + 4:
820  polarizability_dict['xz'] = float(s[-3])
821  polarizability_dict['yz'] = float(s[-2])
822  polarizability_dict['zz'] = float(s[-1])
823  ln += 1
824  calc_moments['polarizability'] = polarizability_dict
825  os.system("rm -rf *.xyz_* *.[0-9][0-9][0-9]")
826  return calc_moments
827 
828  def energy_rmsd(self, shot=0, optimize=True):
829 
830  """ Calculate energy of the selected structure (optionally minimize and return the minimized energy and RMSD). In kcal/mol. """
831 
832  rmsd = 0.0
833  # This line actually runs TINKER
834  # xyzfnm = sysname+".xyz"
835  if optimize:
836  E_, rmsd, _ = self.optimize(shot)
837  o = self.calltinker("analyze %s.xyz_2 E" % self.name)
838  #----
839  # Two equivalent ways to get the RMSD, here for reference.
840  #----
841  # M1 = Molecule("%s.xyz" % self.name, ftype="tinker")
842  # M2 = Molecule("%s.xyz_2" % self.name, ftype="tinker")
843  # M1 += M2
844  # rmsd = M1.ref_rmsd(0)[1]
845  #----
846  # oo = self.calltinker("superpose %s.xyz %s.xyz_2 1 y u n 0" % (self.name, self.name))
847  # for line in oo:
848  # if "Root Mean Square Distance" in line:
849  # rmsd = float(line.split()[-1])
850  #----
851  os.system("rm %s.xyz_2" % self.name)
852  else:
853  o = self.calltinker("analyze %s.xyz E" % self.name)
854  # Read the TINKER output.
855  E = None
856  for line in o:
857  if "Total Potential Energy" in line:
858  E = float(line.split()[-2].replace('D','e'))
859  if E is None:
860  logger.error("Total potential energy wasn't encountered when calling analyze!\n")
861  raise RuntimeError
862  if optimize and abs(E-E_) > 0.1:
863  warn_press_key("Energy from optimize and analyze aren't the same (%.3f vs. %.3f)" % (E, E_))
864  return E, rmsd
865 
866  def interaction_energy(self, fraga, fragb):
867 
868  """ Calculate the interaction energy for two fragments. """
869 
870  self.A = TINKER(name="A", mol=self.mol.atom_select(fraga), tinker_key="%s.key" % self.name, tinkerpath=self.tinkerpath)
871  self.B = TINKER(name="B", mol=self.mol.atom_select(fragb), tinker_key="%s.key" % self.name, tinkerpath=self.tinkerpath)
872 
873  # Interaction energy needs to be in kcal/mol.
874  return (self.energy() - self.A.energy() - self.B.energy()) / 4.184
875 
876  def molecular_dynamics(self, nsteps, timestep, temperature=None, pressure=None, nequil=0, nsave=1000, minimize=True, anisotropic=False, threads=1, verbose=False, **kwargs):
877 
878  """
879  Method for running a molecular dynamics simulation.
880 
881  Required arguments:
882  nsteps = (int) Number of total time steps
883  timestep = (float) Time step in FEMTOSECONDS
884  temperature = (float) Temperature control (Kelvin)
885  pressure = (float) Pressure control (atmospheres)
886  nequil = (int) Number of additional time steps at the beginning for equilibration
887  nsave = (int) Step interval for saving and printing data
888  minimize = (bool) Perform an energy minimization prior to dynamics
889  threads = (int) Specify how many OpenMP threads to use
890 
891  Returns simulation data:
892  Rhos = (array) Density in kilogram m^-3
893  Potentials = (array) Potential energies
894  Kinetics = (array) Kinetic energies
895  Volumes = (array) Box volumes
896  Dips = (3xN array) Dipole moments
897  EComps = (dict) Energy components
898  """
899 
900  md_defs = OrderedDict()
901  md_opts = OrderedDict()
902  # Print out averages only at the end.
903  md_opts["printout"] = nsave
904  md_opts["openmp-threads"] = threads
905  # Langevin dynamics for temperature control.
906  if temperature is not None:
907  md_defs["integrator"] = "stochastic"
908  else:
909  md_defs["integrator"] = "beeman"
910  md_opts["thermostat"] = None
911  # Periodic boundary conditions.
912  if self.pbc:
913  md_opts["vdw-correction"] = ''
914  if temperature is not None and pressure is not None:
915  md_defs["integrator"] = "beeman"
916  md_defs["thermostat"] = "bussi"
917  md_defs["barostat"] = "montecarlo"
918  if anisotropic:
919  md_opts["aniso-pressure"] = ''
920  elif pressure is not None:
921  warn_once("Pressure is ignored because temperature is turned off.")
922  else:
923  if pressure is not None:
924  warn_once("Pressure is ignored because pbc is set to False.")
925  # Use stochastic dynamics for the gas phase molecule.
926  # If we use the regular integrators it may miss
927  # six degrees of freedom in calculating the kinetic energy.
928  md_opts["barostat"] = None
929 
930  eq_opts = deepcopy(md_opts)
931  if self.pbc and temperature is not None and pressure is not None:
932  eq_opts["integrator"] = "beeman"
933  eq_opts["thermostat"] = "bussi"
934  eq_opts["barostat"] = "berendsen"
935 
936  if minimize:
937  if verbose: logger.info("Minimizing the energy...")
938  self.optimize(0, method="bfgs", crit=1)
939  os.system("mv %s.xyz_2 %s.xyz" % (self.name, self.name))
940  if verbose: logger.info("Done\n")
941 
942  # Run equilibration.
943  if nequil > 0:
944  write_key("%s-eq.key" % self.name, eq_opts, "%s.key" % self.name, md_defs)
945  if verbose: printcool("Running equilibration dynamics", color=0)
946  if self.pbc and pressure is not None:
947  self.calltinker("dynamic %s -k %s-eq %i %f %f 4 %f %f" % (self.name, self.name, nequil, timestep, float(nsave*timestep)/1000,
948  temperature, pressure), print_to_screen=verbose)
949  else:
950  self.calltinker("dynamic %s -k %s-eq %i %f %f 2 %f" % (self.name, self.name, nequil, timestep, float(nsave*timestep)/1000,
951  temperature), print_to_screen=verbose)
952  os.system("rm -f %s.arc" % (self.name))
953 
954  # Run production.
955  if verbose: printcool("Running production dynamics", color=0)
956  write_key("%s-md.key" % self.name, md_opts, "%s.key" % self.name, md_defs)
957  if self.pbc and pressure is not None:
958  odyn = self.calltinker("dynamic %s -k %s-md %i %f %f 4 %f %f" % (self.name, self.name, nsteps, timestep, float(nsave*timestep/1000),
959  temperature, pressure), print_to_screen=verbose)
960  else:
961  odyn = self.calltinker("dynamic %s -k %s-md %i %f %f 2 %f" % (self.name, self.name, nsteps, timestep, float(nsave*timestep/1000),
962  temperature), print_to_screen=verbose)
963 
964  # Gather information.
965  os.system("mv %s.arc %s-md.arc" % (self.name, self.name))
966  self.md_trajectory = "%s-md.arc" % self.name
967  edyn = []
968  kdyn = []
969  temps = []
970  for line in odyn:
971  s = line.split()
972  if 'Current Potential' in line:
973  edyn.append(float(s[2]))
974  if 'Current Kinetic' in line:
975  kdyn.append(float(s[2]))
976  if len(s) > 0 and s[0] == 'Temperature' and s[2] == 'Kelvin':
977  temps.append(float(s[1]))
978 
979  # Potential and kinetic energies converted to kJ/mol.
980  edyn = np.array(edyn) * 4.184
981  kdyn = np.array(kdyn) * 4.184
982  temps = np.array(temps)
983 
984  if verbose: logger.info("Post-processing to get the dipole moments\n")
985  oanl = self.calltinker("analyze %s-md.arc" % self.name, stdin="G,E,M", print_to_screen=False)
986 
987  # Read potential energy and dipole from file.
988  eanl = []
989  dip = []
990  mass = 0.0
991  ecomp = OrderedDict()
992  havekeys = set()
993  first_shot = True
994  for ln, line in enumerate(oanl):
995  strip = line.strip()
996  s = line.split()
997  if 'Total System Mass' in line:
998  mass = float(s[-1])
999  if 'Total Potential Energy : ' in line:
1000  eanl.append(float(s[4]))
1001  if 'Dipole X,Y,Z-Components :' in line:
1002  dip.append([float(s[i]) for i in range(-3,0)])
1003  if first_shot:
1004  for key in eckeys:
1005  if strip.startswith(key):
1006  if key in ecomp:
1007  ecomp[key].append(float(s[-2])*4.184)
1008  else:
1009  ecomp[key] = [float(s[-2])*4.184]
1010  if key in havekeys:
1011  first_shot = False
1012  havekeys.add(key)
1013  else:
1014  for key in havekeys:
1015  if strip.startswith(key):
1016  if key in ecomp:
1017  ecomp[key].append(float(s[-2])*4.184)
1018  else:
1019  ecomp[key] = [float(s[-2])*4.184]
1020  for key in ecomp:
1021  ecomp[key] = np.array(ecomp[key])
1022  ecomp["Potential Energy"] = edyn
1023  ecomp["Kinetic Energy"] = kdyn
1024  ecomp["Temperature"] = temps
1025  ecomp["Total Energy"] = edyn+kdyn
1026 
1027  # Energies in kilojoules per mole
1028  eanl = np.array(eanl) * 4.184
1029  # Dipole moments in debye
1030  dip = np.array(dip)
1031  # Volume of simulation boxes in cubic nanometers
1032  # Conversion factor derived from the following:
1033  # In [22]: 1.0 * gram / mole / (1.0 * nanometer)**3 / AVOGADRO_CONSTANT_NA / (kilogram/meter**3)
1034  # Out[22]: 1.6605387831627252
1035  conv = 1.6605387831627252
1036  if self.pbc:
1037  vol = np.array([BuildLatticeFromLengthsAngles(*[float(j) for j in line.split()]).V \
1038  for line in open("%s-md.arc" % self.name).readlines() \
1039  if (len(line.split()) == 6 and isfloat(line.split()[1]) \
1040  and all([isfloat(i) for i in line.split()[:6]]))]) / 1000
1041  rho = conv * mass / vol
1042  else:
1043  vol = None
1044  rho = None
1045  prop_return = OrderedDict()
1046  prop_return.update({'Rhos': rho, 'Potentials': edyn, 'Kinetics': kdyn, 'Volumes': vol, 'Dips': dip, 'Ecomps': ecomp})
1047  return prop_return
1048 
1049 class Liquid_TINKER(Liquid):
1050  """ Condensed phase property matching using TINKER. """
1051  def __init__(self,options,tgt_opts,forcefield):
1052  # Number of threads in running "dynamic".
1053  self.set_option(tgt_opts,'md_threads')
1054  # Number of threads in running "dynamic".
1055  self.set_option(options,'tinkerpath')
1056  # Name of the liquid coordinate file.
1057  self.set_option(tgt_opts,'liquid_coords',default='liquid.xyz',forceprint=True)
1058  # Name of the gas coordinate file.
1059  self.set_option(tgt_opts,'gas_coords',default='gas.xyz',forceprint=True)
1060  # Class for creating engine object.
1061  self.engine_ = TINKER
1062  # Name of the engine to pass to npt.py.
1063  self.engname = "tinker"
1064  # Command prefix.
1065  self.nptpfx = ""
1066  # Extra files to be linked into the temp-directory.
1067  self.nptfiles = ['%s.key' % os.path.splitext(f)[0] for f in [self.liquid_coords, self.gas_coords]]
1068  # Set some options for the polarization correction calculation.
1069  self.gas_engine_args = {'tinker_key' : '%s.key' % os.path.splitext(self.gas_coords)[0]}
1070  # Scripts to be copied from the ForceBalance installation directory.
1071  self.scripts = []
1072  # Initialize the base class.
1073  super(Liquid_TINKER,self).__init__(options,tgt_opts,forcefield)
1074  # Error checking.
1075  for i in self.nptfiles:
1076  if not os.path.exists(os.path.join(self.root, self.tgtdir, i)):
1077  logger.error('Please provide %s; it is needed to proceed.\n' % i)
1078  raise RuntimeError
1079  # Send back the trajectory file.
1080  self.extra_output = ['liquid.dyn']
1081  if self.save_traj > 0:
1082  self.extra_output += ['liquid-md.arc']
1083  # Dictionary of .dyn files used to restart simulations.
1084  self.DynDict = OrderedDict()
1085  self.DynDict_New = OrderedDict()
1086  # These functions need to be called after self.nptfiles is populated
1087  self.post_init(options)
1088 
1089  def npt_simulation(self, temperature, pressure, simnum):
1090  """ Submit a NPT simulation to the Work Queue. """
1091  if self.goodstep and (temperature, pressure) in self.DynDict_New:
1092  self.DynDict[(temperature, pressure)] = self.DynDict_New[(temperature, pressure)]
1093  if (temperature, pressure) in self.DynDict:
1094  dynsrc = self.DynDict[(temperature, pressure)]
1095  dyndest = os.path.join(os.getcwd(), 'liquid.dyn')
1096  logger.info("Copying .dyn file: %s to %s\n" % (dynsrc, dyndest))
1097  shutil.copy2(dynsrc,dyndest)
1098  self.nptfiles.append(dyndest)
1099  self.DynDict_New[(temperature, pressure)] = os.path.join(os.getcwd(),'liquid.dyn')
1100  super(Liquid_TINKER, self).npt_simulation(temperature, pressure, simnum)
1101  self.last_traj = [i for i in self.last_traj if '.dyn' not in i]
1102 
1103 class AbInitio_TINKER(AbInitio):
1104  """ Subclass of Target for force and energy matching using TINKER. """
1105  def __init__(self,options,tgt_opts,forcefield):
1106 
1107  self.set_option(tgt_opts,'coords',default="all.arc")
1108  self.set_option(tgt_opts,'tinker_key',default="shot.key")
1109  self.engine_ = TINKER
1110 
1111  super(AbInitio_TINKER,self).__init__(options,tgt_opts,forcefield)
1112 
1113 class BindingEnergy_TINKER(BindingEnergy):
1114  """ Binding energy matching using TINKER. """
1115  def __init__(self,options,tgt_opts,forcefield):
1116  self.engine_ = TINKER
1117 
1118  super(BindingEnergy_TINKER,self).__init__(options,tgt_opts,forcefield)
1119 
1120 class Interaction_TINKER(Interaction):
1121  """ Subclass of Target for interaction matching using TINKER. """
1122  def __init__(self,options,tgt_opts,forcefield):
1123 
1124  self.set_option(tgt_opts,'coords',default="all.arc")
1125  self.set_option(tgt_opts,'tinker_key',default="shot.key")
1126  self.engine_ = TINKER
1128  super(Interaction_TINKER,self).__init__(options,tgt_opts,forcefield)
1129 
1130 class Moments_TINKER(Moments):
1131  """ Subclass of Target for multipole moment matching using TINKER. """
1132  def __init__(self,options,tgt_opts,forcefield):
1133 
1134  self.set_option(tgt_opts,'coords',default="input.xyz")
1135  self.set_option(tgt_opts,'tinker_key',default="input.key")
1136  self.engine_ = TINKER
1138  super(Moments_TINKER,self).__init__(options,tgt_opts,forcefield)
1140 class Vibration_TINKER(Vibration):
1141  """ Vibrational frequency matching using TINKER. """
1142  def __init__(self,options,tgt_opts,forcefield):
1143 
1144  self.set_option(tgt_opts,'coords',default="input.xyz")
1145  self.set_option(tgt_opts,'tinker_key',default="input.key")
1146  self.engine_ = TINKER
1147 
1148  super(Vibration_TINKER,self).__init__(options,tgt_opts,forcefield)
1149 
def energy_rmsd(self, shot=0, optimize=True)
Calculate energy of the selected structure (optionally minimize and return the minimized energy and R...
Definition: tinkerio.py:846
Engine for carrying out general purpose TINKER calculations.
Definition: tinkerio.py:313
def readsrc(self, kwargs)
Called by init ; read files from the source directory.
Definition: tinkerio.py:341
def energy_dipole(self)
Computes the energy and dipole using TINKER over a trajectory.
Definition: tinkerio.py:734
Vibrational mode fitting module.
Nifty functions, intended to be imported by any module within ForceBalance.
def BuildLatticeFromLengthsAngles(a, b, c, alpha, beta, gamma)
This function takes in three lattice lengths and three lattice angles, and tries to return a complete...
Definition: molecule.py:437
Ab-initio fitting module (energies, forces, resp).
def energy_force(self)
Computes the energy and force using TINKER over a trajectory.
Definition: tinkerio.py:721
def write_key(fout, options, fin=None, defaults={}, verbose=False, prmfnm=None, chk=[])
Create or edit a TINKER .key file.
Definition: tinkerio.py:184
def LinkFile(src, dest, nosrcok=False)
Definition: nifty.py:1313
valkwd
Keyword args that aren&#39;t in this list are filtered out.
Definition: tinkerio.py:318
def molecular_dynamics(self, nsteps, timestep, temperature=None, pressure=None, nequil=0, nsave=1000, minimize=True, anisotropic=False, threads=1, verbose=False, kwargs)
Method for running a molecular dynamics simulation.
Definition: tinkerio.py:916
def isint(word)
ONLY matches integers! If you have a decimal point? None shall pass!
Definition: molecule.py:401
def optimize(self, shot, method="newton", align=True, crit=1e-4)
Optimize the geometry and align the optimized geometry to the starting geometry.
Definition: tinkerio.py:567
def __init__(self, options, tgt_opts, forcefield)
Definition: tinkerio.py:1138
def which(fnm)
Definition: nifty.py:1362
def npt_simulation(self, temperature, pressure, simnum)
Submit a NPT simulation to the Work Queue.
Definition: tinkerio.py:1111
def energy_force_one(self, shot)
Computes the energy and force using TINKER for one snapshot.
Definition: tinkerio.py:700
def multipole_moments(self, shot=0, optimize=True, polarizability=False)
Return the multipole moments of the 1st snapshot in Debye and Buckingham units.
Definition: tinkerio.py:783
Binding energy fitting module.
Subclass of Target for force and energy matching using TINKER.
Definition: tinkerio.py:1126
Subclass of Target for interaction matching using TINKER.
Definition: tinkerio.py:1145
tinkerpath
The directory containing TINKER executables (e.g.
Definition: tinkerio.py:329
def calltinker(self, command, stdin=None, print_to_screen=False, print_command=False, kwargs)
Call TINKER; prepend the tinkerpath to calling the TINKER program.
Definition: tinkerio.py:354
AtomMask
If the coordinates do not come with TINKER suffixes then throw an error.
Definition: tinkerio.py:514
atom
The atom numbers in the interaction (stored in the TINKER parser)
Definition: tinkerio.py:104
def energy(self)
Computes the energy using TINKER over a trajectory.
Definition: tinkerio.py:709
Vibrational frequency matching using TINKER.
Definition: tinkerio.py:1167
engine_
Default file names for coordinates and key file.
Definition: tinkerio.py:1131
def __init__(self, options, tgt_opts, forcefield)
Definition: tinkerio.py:1127
def warn_press_key(warning, timeout=10)
Definition: nifty.py:1599
def printcool_dictionary(Dict, title="Dictionary Keys : Values", bold=False, color=2, keywidth=25, topwidth=50, center=True, leftpad=0)
See documentation for printcool; this is a nice way to print out keys/values in a dictionary...
Definition: nifty.py:366
def __init__(self, options, tgt_opts, forcefield)
Definition: tinkerio.py:1146
Subclass of Target for multipole moment matching using TINKER.
Definition: tinkerio.py:1156
def printcool(text, sym="#", bold=False, color=2, ansi=None, bottom='-', minwidth=50, center=True, sym2="=")
Cool-looking printout for slick formatting of output.
Definition: nifty.py:321
def feed(self, line)
Given a line, determine the interaction type and the atoms involved (the suffix). ...
Definition: tinkerio.py:140
def get_charges(self)
Definition: tinkerio.py:694
Interaction energy fitting module.
def onefile(fnm=None, ext=None, err=False)
Definition: nifty.py:1119
rigid
Attempt to set some TINKER options.
Definition: tinkerio.py:404
Finite state machine for parsing TINKER force field files.
Definition: tinkerio.py:97
Binding energy matching using TINKER.
Definition: tinkerio.py:1137
def prepare(self, pbc=False, kwargs)
Called by init ; prepare the temp directory and figure out the topology.
Definition: tinkerio.py:398
def warn_once(warning, warnhash=None)
Prints a warning but will only do so once in a given run.
Definition: nifty.py:1611
pdict
The parameter dictionary (defined in this file)
Definition: tinkerio.py:102
def interaction_energy(self, fraga, fragb)
Calculate the interaction energy for two fragments.
Definition: tinkerio.py:885
engine_
Default file names for coordinates and key file.
Definition: tinkerio.py:1150
def normal_modes(self, shot=0, optimize=True)
Definition: tinkerio.py:745
Multipole moment fitting module.
def evaluate_(self, xyzin, force=False, dipole=False)
Utility function for computing energy, and (optionally) forces and dipoles using TINKER.
Definition: tinkerio.py:647
def isfloat(word)
Matches ANY number; it can be a decimal, scientific notation, integer, or what have you...
Definition: molecule.py:406
def setopts(self, kwargs)
Called by init ; Set TINKER-specific options.
Definition: tinkerio.py:324
def __init__(self, name="tinker", kwargs)
Definition: tinkerio.py:316
def wopen(dest, binary=False)
If trying to write to a symbolic link, remove it first.
Definition: nifty.py:1304
Condensed phase property matching using TINKER.
Definition: tinkerio.py:1070
Matching of liquid bulk properties.
def __init__(self, options, tgt_opts, forcefield)
Definition: tinkerio.py:1071
def __init__(self, fnm)
Definition: tinkerio.py:99