ForceBalance API  1.3
Automated optimization of force fields and empirical potentials
List of all members | Public Member Functions | Public Attributes
src.abinitio.AbInitio Class Reference

Subclass of Target for fitting force fields to ab initio data. More...

Inheritance diagram for src.abinitio.AbInitio:
[legend]
Collaboration diagram for src.abinitio.AbInitio:
[legend]

Public Member Functions

def __init__ (self, options, tgt_opts, forcefield)
 Initialization; define a few core concepts. More...
 
def build_invdist (self, mvals)
 
def compute_netforce_torque (self, xyz, force, QM=False)
 
def read_reference_data (self)
 Read the reference ab initio data from a file such as qdata.txt. More...
 
def indicate (self)
 
def energy_all (self)
 
def energy_force_all (self)
 
def energy_force_transform (self)
 
def energy_one (self, i)
 
def energy_force_one (self, i)
 
def energy_force_transform_one (self, i)
 
def get_energy_force (self, mvals, AGrad=False, AHess=False)
 LPW 7-13-2016. More...
 
def get_resp (self, mvals, AGrad=False, AHess=False)
 Electrostatic potential fitting. More...
 
def get (self, mvals, AGrad=False, AHess=False)
 

Public Attributes

 w_force
 Initialize the base class. More...
 
 savg
 Option for how much data to write to disk. More...
 
 energy_asymmetry
 
 asym
 
 loop_over_snapshots
 LPW 2018-02-11: This is set to True if the target calculates a single-point property over several existing snapshots. More...
 
 boltz_wts
 Boltzmann weights. More...
 
 eqm
 Reference (QM) energies. More...
 
 fqm
 Reference (QM) forces. More...
 
 espxyz
 ESP grid points. More...
 
 espval
 ESP values. More...
 
 qfnm
 The qdata.txt file that contains the QM energies and forces. More...
 
 e_err
 Qualitative Indicator: average energy error (in kJ/mol) More...
 
 e_err_pct
 
 f_err
 Qualitative Indicator: average force error (fractional) More...
 
 f_err_pct
 
 esp_err
 Qualitative Indicator: "relative RMS" for electrostatic potential. More...
 
 nf_err
 
 nf_err_pct
 
 tq_err_pct
 
 use_nft
 Whether to compute net forces and torques, or not. More...
 
 mol
 Read in the trajectory file. More...
 
 ns
 
 nparticles
 The number of (atoms + drude particles + virtual sites) More...
 
 engine
 Build keyword dictionaries to pass to engine. More...
 
 AtomLists
 Lists of atoms to do net force/torque fitting and excluding virtual sites. More...
 
 AtomMask
 
 buildx
 Read in the reference data. More...
 
 save_vmvals
 Save the mvals from the last time we updated the vsites. More...
 
 M_orig
 
 force_map
 
 nnf
 
 ntq
 
 smin
 
 qmatoms
 
 force
 
 fitatoms
 
 nesp
 
 nftqm
 
 fref
 
 w_energy
 
 w_netforce
 
 w_torque
 
 maxfatom
 
 maxfshot
 
 maxdf
 
 e_trm
 
 e_ctr
 
 e_ref
 
 f_trm
 
 f_ctr
 
 f_ref
 
 nf_trm
 
 nf_ctr
 
 nf_ref
 
 tq_trm
 
 tq_ctr
 
 tq_ref
 
 tq_err
 
 w_resp
 
 invdists
 
 esp_ref
 
 esp_err_pct
 
 respterm
 
 esp_trm
 
 esp_ctr
 
 objective
 

Detailed Description

Subclass of Target for fitting force fields to ab initio data.

Currently Gromacs-X2, Gromacs, Tinker, OpenMM and AMBER are supported.

We introduce the following concepts:

There are also these little details:

This subclass contains the 'get' method for building the objective function from any simulation software (a driver to run the program and read output is still required). The 'get' method can be overridden by subclasses like AbInitio_GMX.

Definition at line 84 of file abinitio.py.

Constructor & Destructor Documentation

◆ __init__()

def src.abinitio.AbInitio.__init__ (   self,
  options,
  tgt_opts,
  forcefield 
)

Initialization; define a few core concepts.

Todo:
Obtain the number of true atoms (or the particle -> atom mapping) from the force field.

Definition at line 94 of file abinitio.py.

Here is the call graph for this function:

Member Function Documentation

◆ build_invdist()

def src.abinitio.AbInitio.build_invdist (   self,
  mvals 
)

Definition at line 216 of file abinitio.py.

◆ compute_netforce_torque()

def src.abinitio.AbInitio.compute_netforce_torque (   self,
  xyz,
  force,
  QM = False 
)

Definition at line 245 of file abinitio.py.

◆ energy_all()

def src.abinitio.AbInitio.energy_all (   self)

Definition at line 505 of file abinitio.py.

◆ energy_force_all()

def src.abinitio.AbInitio.energy_force_all (   self)

Definition at line 512 of file abinitio.py.

◆ energy_force_one()

def src.abinitio.AbInitio.energy_force_one (   self,
  i 
)

Definition at line 544 of file abinitio.py.

◆ energy_force_transform()

def src.abinitio.AbInitio.energy_force_transform (   self)

Definition at line 519 of file abinitio.py.

Here is the call graph for this function:

◆ energy_force_transform_one()

def src.abinitio.AbInitio.energy_force_transform_one (   self,
  i 
)

Definition at line 551 of file abinitio.py.

Here is the call graph for this function:

◆ energy_one()

def src.abinitio.AbInitio.energy_one (   self,
  i 
)

Definition at line 537 of file abinitio.py.

◆ get()

def src.abinitio.AbInitio.get (   self,
  mvals,
  AGrad = False,
  AHess = False 
)

Definition at line 1090 of file abinitio.py.

Here is the call graph for this function:

◆ get_energy_force()

def src.abinitio.AbInitio.get_energy_force (   self,
  mvals,
  AGrad = False,
  AHess = False 
)

LPW 7-13-2016.

This code computes the least squares objective function for the energy and force. The most recent revision simplified the code to make it easier to maintain, and to remove the covariance matrix and dual-weight features.

The equations have also been simplified. Previously I normalizing each force component (or triples of components belonging to an atom) separately. I was also computing the objective function separately from the "indicators", which led to confusion regarding why they resulted in different values. The code was structured around computing the objective function by multiplying an Applequist polytensor containing both energy and force with an inverse covariance matrix, but it became apparent over time and experience that the more complicated approach was not worth it, given the opacity that it introduced into how things were computed.

In the new code, the objective function is computed in a simple way. For the energies we compute a weighted sum of squared differences between E_MM and E_QM, minus (optionally) the mean energy gap, for the numerator. The weighted variance of the QM energies <E_QM^2>-<E_QM>^2 is the denominator. The user-supplied w_energy option is a prefactor that multiplies this term. If w_normalize is set to True (no longer the default), the prefactor is further divided by the sum of all of the weights. The indicators are set to the square roots of the numerator and denominator above.

For the forces we compute the same weighted sum, where each term in the sum is the norm-squared of F_MM - F_QM, with no option to subtract the mean. The denominator is computed in an analogous way using the norm-squared F_QM, and the prefactor is w_force. The same approach is applied if the user asks for net forces and/or torques. The indicators are computed from the square roots of the numerator and denominator, divided by the number of atoms (molecules) for forces (net forces / torques).

In equation form, the objective function is given by:

[ = {W_E}[ {{{( {{i {N_s}} {{w_i}{{( {E_i^{MM} - E_i^{QM}} )}^2}} } ) - {{( {{{ E}^{MM}} - {{ E}^{QM}}} )}^2}}} {{{i {N_s}} {{w_i}{{( {E_i^{QM} - {{ E}^{QM}}} )}^2}} }}} ] + {W_F}[ {{{{i {N_s}} {{w_i}{a {N_a}} {{{| {{{F}}_{i,a}^{MM} - {{F}}_{i,a}^{QM}} |}^2}} } }} {{{i {N_s}} {{w_i}{a {N_a}} {{{| {{{F}}_{i,a}^{QM}} |}^2}} } }}} ]]

In the previous code (ForTune, 2011 and previous) this subroutine used analytic first derivatives of the energy and force to build the derivatives of the objective function. Here I will take a simplified approach, because building the analytic derivatives require substantial modifications of the engine code, which is unsustainable. We use a finite different calculation of the first derivatives of the energies and forces to get the exact first derivative and approximate second derivative of the objective function..

Parameters
[in]mvalsMathematical parameter values
[in]AGradSwitch to turn on analytic gradient
[in]AHessSwitch to turn on analytic Hessian
Returns
Answer Contribution to the objective function

Definition at line 626 of file abinitio.py.

◆ get_resp()

def src.abinitio.AbInitio.get_resp (   self,
  mvals,
  AGrad = False,
  AHess = False 
)

Electrostatic potential fitting.

Implements the RESP objective function. (In Python so obviously not optimized.) Return the charges acting on the system.

Definition at line 975 of file abinitio.py.

◆ indicate()

def src.abinitio.AbInitio.indicate (   self)

Definition at line 464 of file abinitio.py.

Here is the call graph for this function:

◆ read_reference_data()

def src.abinitio.AbInitio.read_reference_data (   self)

Read the reference ab initio data from a file such as qdata.txt.

Todo:
Add an option for picking any slice out of qdata.txt, helpful for cross-validation
     After reading in the information from qdata.txt, it is converted
     into the GROMACS energy units (kind of an arbitrary choice);
     forces (kind of a misnomer in qdata.txt) are multipled by -1
     to convert gradients to forces.

     We typically subtract out the mean energies of all energy arrays
     because energy/force matching does not account for zero-point
     energy differences between MM and QM (i.e. energy of electrons
     in core orbitals).

     A 'hybrid' ensemble is possible where we use 50% MM and 50% QM
     weights.  Please read more in LPW and Troy Van Voorhis, JCP
     Vol. 133, Pg. 231101 (2010), doi:10.1063/1.3519043.  In the
     updated version of the code (July 13 2016), this feature is
     currently not implemented due to disuse, but it is easy to
     re-implement if desired.

     Finally, note that using non-Boltzmann weights degrades the
     statistical information content of the snapshots.  This
     problem will generally become worse if the ensemble to which
     we're reweighting is dramatically different from the one we're
     sampling from.  We end up with a set of Boltzmann weights like
     [1e-9, 1e-9, 1.0, 1e-9, 1e-9 ... ] and this is essentially just
     one snapshot.

     Here, we have a measure for the information content of our snapshots,
     which comes easily from the definition of information entropy:

     S = -1*Sum_i(P_i*log(P_i))
     InfoContent = exp(-S)

     With uniform weights, InfoContent is equal to the number of snapshots;
     with horrible weights, InfoContent is closer to one.

Definition at line 360 of file abinitio.py.

Member Data Documentation

◆ asym

src.abinitio.AbInitio.asym

Definition at line 148 of file abinitio.py.

◆ AtomLists

src.abinitio.AbInitio.AtomLists

Lists of atoms to do net force/torque fitting and excluding virtual sites.

Definition at line 204 of file abinitio.py.

◆ AtomMask

src.abinitio.AbInitio.AtomMask

Definition at line 205 of file abinitio.py.

◆ boltz_wts

src.abinitio.AbInitio.boltz_wts

Boltzmann weights.

Definition at line 162 of file abinitio.py.

◆ buildx

src.abinitio.AbInitio.buildx

Read in the reference data.

The below two options are related to whether we want to rebuild virtual site positions. Rebuild the distance matrix if virtual site positions have changed

Definition at line 210 of file abinitio.py.

◆ e_ctr

src.abinitio.AbInitio.e_ctr

Definition at line 942 of file abinitio.py.

◆ e_err

src.abinitio.AbInitio.e_err

Qualitative Indicator: average energy error (in kJ/mol)

Definition at line 174 of file abinitio.py.

◆ e_err_pct

src.abinitio.AbInitio.e_err_pct

Definition at line 175 of file abinitio.py.

◆ e_ref

src.abinitio.AbInitio.e_ref

Definition at line 944 of file abinitio.py.

◆ e_trm

src.abinitio.AbInitio.e_trm

Definition at line 941 of file abinitio.py.

◆ energy_asymmetry

src.abinitio.AbInitio.energy_asymmetry

Definition at line 147 of file abinitio.py.

◆ engine

src.abinitio.AbInitio.engine

Build keyword dictionaries to pass to engine.

Create engine object.

Definition at line 202 of file abinitio.py.

◆ eqm

src.abinitio.AbInitio.eqm

Reference (QM) energies.

Definition at line 164 of file abinitio.py.

◆ esp_ctr

src.abinitio.AbInitio.esp_ctr

Definition at line 1081 of file abinitio.py.

◆ esp_err

src.abinitio.AbInitio.esp_err

Qualitative Indicator: "relative RMS" for electrostatic potential.

Definition at line 180 of file abinitio.py.

◆ esp_err_pct

src.abinitio.AbInitio.esp_err_pct

Definition at line 1061 of file abinitio.py.

◆ esp_ref

src.abinitio.AbInitio.esp_ref

Definition at line 1060 of file abinitio.py.

◆ esp_trm

src.abinitio.AbInitio.esp_trm

Definition at line 1080 of file abinitio.py.

◆ espval

src.abinitio.AbInitio.espval

ESP values.

Definition at line 170 of file abinitio.py.

◆ espxyz

src.abinitio.AbInitio.espxyz

ESP grid points.

Definition at line 168 of file abinitio.py.

◆ f_ctr

src.abinitio.AbInitio.f_ctr

Definition at line 949 of file abinitio.py.

◆ f_err

src.abinitio.AbInitio.f_err

Qualitative Indicator: average force error (fractional)

Definition at line 177 of file abinitio.py.

◆ f_err_pct

src.abinitio.AbInitio.f_err_pct

Definition at line 178 of file abinitio.py.

◆ f_ref

src.abinitio.AbInitio.f_ref

Definition at line 951 of file abinitio.py.

◆ f_trm

src.abinitio.AbInitio.f_trm

Definition at line 948 of file abinitio.py.

◆ fitatoms

src.abinitio.AbInitio.fitatoms

Definition at line 407 of file abinitio.py.

◆ force

src.abinitio.AbInitio.force

Definition at line 399 of file abinitio.py.

◆ force_map

src.abinitio.AbInitio.force_map

Definition at line 254 of file abinitio.py.

◆ fqm

src.abinitio.AbInitio.fqm

Reference (QM) forces.

Definition at line 166 of file abinitio.py.

◆ fref

src.abinitio.AbInitio.fref

Definition at line 455 of file abinitio.py.

◆ invdists

src.abinitio.AbInitio.invdists

Definition at line 983 of file abinitio.py.

◆ loop_over_snapshots

src.abinitio.AbInitio.loop_over_snapshots

LPW 2018-02-11: This is set to True if the target calculates a single-point property over several existing snapshots.

Definition at line 160 of file abinitio.py.

◆ M_orig

src.abinitio.AbInitio.M_orig

Definition at line 214 of file abinitio.py.

◆ maxdf

src.abinitio.AbInitio.maxdf

Definition at line 711 of file abinitio.py.

◆ maxfatom

src.abinitio.AbInitio.maxfatom

Definition at line 709 of file abinitio.py.

◆ maxfshot

src.abinitio.AbInitio.maxfshot

Definition at line 710 of file abinitio.py.

◆ mol

src.abinitio.AbInitio.mol

Read in the trajectory file.

Set the number of snapshots.

Definition at line 188 of file abinitio.py.

◆ nesp

src.abinitio.AbInitio.nesp

Definition at line 429 of file abinitio.py.

◆ nf_ctr

src.abinitio.AbInitio.nf_ctr

Definition at line 956 of file abinitio.py.

◆ nf_err

src.abinitio.AbInitio.nf_err

Definition at line 181 of file abinitio.py.

◆ nf_err_pct

src.abinitio.AbInitio.nf_err_pct

Definition at line 182 of file abinitio.py.

◆ nf_ref

src.abinitio.AbInitio.nf_ref

Definition at line 957 of file abinitio.py.

◆ nf_trm

src.abinitio.AbInitio.nf_trm

Definition at line 955 of file abinitio.py.

◆ nftqm

src.abinitio.AbInitio.nftqm

Definition at line 451 of file abinitio.py.

◆ nnf

src.abinitio.AbInitio.nnf

Definition at line 315 of file abinitio.py.

◆ nparticles

src.abinitio.AbInitio.nparticles

The number of (atoms + drude particles + virtual sites)

Definition at line 197 of file abinitio.py.

◆ ns

src.abinitio.AbInitio.ns

Definition at line 195 of file abinitio.py.

◆ ntq

src.abinitio.AbInitio.ntq

Definition at line 316 of file abinitio.py.

◆ objective

src.abinitio.AbInitio.objective

Definition at line 1111 of file abinitio.py.

◆ qfnm

src.abinitio.AbInitio.qfnm

The qdata.txt file that contains the QM energies and forces.

Definition at line 172 of file abinitio.py.

◆ qmatoms

src.abinitio.AbInitio.qmatoms

Definition at line 396 of file abinitio.py.

◆ respterm

src.abinitio.AbInitio.respterm

Definition at line 1072 of file abinitio.py.

◆ save_vmvals

src.abinitio.AbInitio.save_vmvals

Save the mvals from the last time we updated the vsites.

Definition at line 212 of file abinitio.py.

◆ savg

src.abinitio.AbInitio.savg

Option for how much data to write to disk.

Whether to do energy and force calculations for the whole trajectory, or to do one calculation per snapshot. OpenMM-only option - whether to run the energies and forces internally. Whether we have virtual sites (set at the global option level) Attenuate the weights as a function of energy What is the energy denominator? (Valid for 'attenuate') Set upper cutoff energy Assign a greater weight to MM snapshots that underestimate the QM energy (surfaces referenced to QM absolute minimum)

Definition at line 147 of file abinitio.py.

◆ smin

src.abinitio.AbInitio.smin

Definition at line 386 of file abinitio.py.

◆ tq_ctr

src.abinitio.AbInitio.tq_ctr

Definition at line 961 of file abinitio.py.

◆ tq_err

src.abinitio.AbInitio.tq_err

Definition at line 963 of file abinitio.py.

◆ tq_err_pct

src.abinitio.AbInitio.tq_err_pct

Definition at line 183 of file abinitio.py.

◆ tq_ref

src.abinitio.AbInitio.tq_ref

Definition at line 962 of file abinitio.py.

◆ tq_trm

src.abinitio.AbInitio.tq_trm

Definition at line 960 of file abinitio.py.

◆ use_nft

src.abinitio.AbInitio.use_nft

Whether to compute net forces and torques, or not.

Definition at line 185 of file abinitio.py.

◆ w_energy

src.abinitio.AbInitio.w_energy

Definition at line 651 of file abinitio.py.

◆ w_force

src.abinitio.AbInitio.w_force

Initialize the base class.

Number of snapshots Whether to match Absolute Energies (make sure you know what you're doing) Number of atoms that we are fitting Whether to fit Energies. Whether to fit Forces. Whether to fit Electrostatic Potential. Weights for the three components.

Definition at line 121 of file abinitio.py.

◆ w_netforce

src.abinitio.AbInitio.w_netforce

Definition at line 651 of file abinitio.py.

◆ w_resp

src.abinitio.AbInitio.w_resp

Definition at line 976 of file abinitio.py.

◆ w_torque

src.abinitio.AbInitio.w_torque

Definition at line 651 of file abinitio.py.


The documentation for this class was generated from the following file: