Reference Guide

Modules

An overview is provived of all modules included in QuickFF. The information below is based on the documentation provided in the source code.

quickff.cost – Least-Square cost function for hessian matching

class quickff.cost.HessianFCCost(system, ai, valence, fit_indices, ffrefs=[], do_mass_weighting=True)[source]

A class to implement the least-square cost function to fit the force field hessian to the ab initio hessian.

Arguments

system
a Yaff system object
ai
an instance of the Reference representing the ab initio input
valence
A ValenceFF object containing all valence terms.
fit_indices
a list of indices indicating the terms for which the force constants should be determined.

Optional Arguments

ffrefs
a list of Reference instances representing possible a priori determined contributions to the force field (such as eg. electrostatics and van der Waals)

do_mass_weighting

estimate(init=None, lower=None, upper=None, do_svd=False, svd_rcond=0.0)[source]

Estimate the force constants by minimizing the cost function

quickff.io – Input/Ouput methods

class quickff.io.VASPRun(filename, field_labels=[])[source]

Load information from a vasprun.xml file

Arguments

filename
Filename of vasprun.xml

Optional Arguments

field_labels
List of things we want to read. If an empty list is provided, only numbers, masses, initial positions and initial cell vectors are read.
quickff.io.read_abinitio(fn, do_hess=True)[source]

Wrapper to read all information from an ab initio calculation that QuickFF needs. Currently Gaussian .fchk and VASP .xml files are supported.

Optional Arguments

do_hess
Extract the hessian from the ab initio output. For qff-input-ei.py, it is interesting to be able to switch this off.
quickff.io.dump_charmm22_prm(valence, fn)[source]

Dump supported parameters in a CHARMM22 parameter file.

Arguments

valence
Instance of ValenceFF, which defines the force field.
fn
The filename to write to.
quickff.io.dump_charmm22_psf(system, valence, fn)[source]

Dump supported internal coordinates in a CHARMM psf file.

Arguments

system
Instance of yaff.System class
valence
Instance of ValenceFF, which defines the force field.
fn
The filename to write to.

quickff.settings – QuickFF configuration settings

class quickff.settings.Settings(fn=None, **kwargs)[source]

Class to control the behaviour of a Quickff run

Keyword Arguments

fn

file name of a config file from which settings can be read. Each line contains a single setting (except the lines starting with a #) and should have the following syntax

key: value
kwargs
each setting can also be parsed directly through the use of keyword arguments in the __init__ constructor with the syntax key=value. The settings parsed through the use of kwargs overwrite those in the config file.

quickff.perturbation – Perturbation trajectories for internal coordinates

class quickff.perturbation.Trajectory(term, start, end, numbers, nsteps=11)[source]

A class to store a perturbation trajectory

Arguments

term
an instance of the Term class for which the perturbation trajectory will be computed
numbers
list of atom numbers required for dumping xyz coords
start
a float defining the lower limit of the perturbation value of the given ic. If not given, a standard value is choosen according to the kind of internal coordinate.
end
a float defining the upper limit of the perturbations value of the given ic. If not given, a standard value is choosen according to the kind of internal coordinate.

Optional Arguments

nsteps
an integer defining the number of steps in the perturbation trajectory. The default value is 11 steps.
plot(ai, ffrefs=[], valence=None, fn='default', eunit='kjmol', suffix='')[source]

Method to plot the energy contributions along a perturbation trajectory associated to a given ic. This method assumes that the coords, fc and rv attributes have been computed and assigned.

Arguments

ai
an instance of the Reference representing the ab initio input

Optional Arguments

ffrefs
a list of Reference instances representing possible a priori determined contributions to the force field (such as eg. electrostatics and van der Waals)
valence
an instance of ValenceFF which will be used to plot the covalent contribution. If not given, only the contribution of the IC corresponding to the trajectory will be plotted using the values of fc and rv
fn
a string defining the name of the figure
eunit
a string describing the conversion of the unit of energy. More info regarding possible strings can be found in the MolMod documentation.
suffix
a string to be added to the filename at the end. Is overwritten when fn is specified.
to_xyz(fn=None)[source]

Method to write the given trajectory to an XYZ file. This method assumes that the coords attribute has been assigned.

Optional Arguments

fn
a string defining the name of the output file

quickff.program – Program classes implementing a sequence of fitting steps

class quickff.program.BaseProgram(system, ai, settings, ffrefs=[])[source]

Base program which implements all possible steps of a force field fitting program. The actual sequence of the steps are defined in the deriving classes.

Arguments

system
a Yaff System instance defining the system
ai
a Reference instance corresponding to the ab initio input data
settings
a Settings instance defining all QuickFF settings

Optional Arguments

ffrefs
a list of Reference instances defining the a-priori force field contributions.
average_pars()[source]

Average force field parameters over master and slaves.

do_bendclin(thresshold=0.08726646259971647)[source]

No Harmonic bend can have a rest value equal to are large than 180 deg - thresshold. If a master (or its slaves) has such a rest value, convert master and all slaves to BendCLin (which corresponds to a chebychev1 potential with sign=+1):

0.5*K*[1+cos(theta)]
do_cross_init()[source]

Add cross terms to the valence list and initialize parameters.

do_eq_setrv(tasks, logger_level=3)[source]

Set the rest values to their respective AI equilibrium values.

do_hc_estimatefc(tasks, logger_level=3, do_svd=False, svd_rcond=0.0, do_mass_weighting=True)[source]

Refine force constants using Hessian Cost function.

Arguments

tasks
A list of strings identifying which terms should have their force constant estimated from the hessian cost function. Using such a flag, one can distinguish between for example force constant refinement (flag=HC_FC_DIAG) of the diagonal terms and force constant estimation of the cross terms (flag=HC_FC_CROSS). If the string ‘all’ is present in tasks, all fc’s will be estimated.

Optional Arguments

logger_level
print level at which the resulting parameters should be dumped to the logger. By default, the parameters will only be dumped at the highest log level.
do_svd
whether or not to do an SVD decomposition before solving the set of equations and explicitly throw out the degrees of freedom that correspond to the lowest singular values.
do_mass_weighting
whether or not to apply mass weighing to the ab initio hessian and the force field contributions before doing the fitting.
do_pt_estimate(do_valence=False, energy_noise=None, logger_level=3)[source]

Estimate force constants and rest values from the perturbation trajectories

Optional Arguments

do_valence
if set to True, the current valence force field will be used to estimate the contribution of all other valence terms.
do_pt_generate()[source]

Generate perturbation trajectories.

do_pt_postprocess()[source]

Do some first post processing of the ff parameters estimated from the perturbation trajectories including:

  • detecting bend patterns with rest values of 90 and 180 deg
  • detecting bend patterns with rest values only close to 180 deg
  • transforming SqOopDist with rv=0.0 to OopDist
  • averaging parameters
do_sqoopdist_to_oopdist(thresshold=0.0001889726133921252)[source]

Transform a SqOopdist term with a rest value that has been set to zero, to a term Oopdist (harmonic in Oopdist instead of square of Oopdist) with a rest value of 0.0 A.

do_squarebend(thresshold=0.3490658503988659)[source]

Identify bend patterns in which 4 atoms of type A surround a central atom of type B with A-B-A angles of 90/180 degrees. A simple harmonic pattern will not be adequate since a rest value of 90 and 180 degrees is possible for the same A-B-A term. Therefore, a cosine term with multiplicity of 4 is used (which corresponds to a chebychev4 potential with sign=-1):

V = K/2*[1-cos(4*theta)]

To identify the patterns, it is assumed that the rest values have already been estimated from the perturbation trajectories. For each master and slave of a BENDAHARM term, its rest value is computed and checked if it lies either the interval [90-thresshold,90+thresshold] or [180-thresshold,180]. If this is the case, the new cosine term is used.

Optional arguments

thresshold
the (half) the width of the interval around 180 deg (90 degrees) to check if a square BA4
make_output()[source]

Dump Yaff parameters, Yaff system, plot energy contributions along perturbation trajectories and dump perturbation trajectories to XYZ files.

plot_trajectories(do_valence=False, suffix='')[source]

Plot energy contributions along perturbation trajectories and

print_system()[source]

dump overview of atoms (and associated parameters) in the system

reset_system()[source]

routine to reset the system coords to the ai equilbrium

run()[source]

Sequence of instructions, should be implemented in the inheriting classes. The various inheriting classes distinguish themselves by means of the instructions implemented in this routine.

update_cross_pars()[source]

Set the rest values of cross terms to the rest values of the corresponding diagonal terms. Set the force constants to zero.

update_trajectory_terms()[source]

Routine to make self.valence.terms and the term attribute of each trajectory in self.trajectories consistent again. This is usefull if the trajectory were read from a file and the valenceFF instance was modified.

write_trajectories()[source]

Write perturbation trajectories to XYZ files.

class quickff.program.MakeTrajectories(system, ai, settings, ffrefs=[])[source]

Construct the perturbation trajectories and store them. This program does not derive the force field.

Arguments

system
a Yaff System instance defining the system
ai
a Reference instance corresponding to the ab initio input data
settings
a Settings instance defining all QuickFF settings

Optional Arguments

ffrefs
a list of Reference instances defining the a-priori force field contributions.
class quickff.program.PlotTrajectories(system, ai, settings, ffrefs=[])[source]

Read the perturbation trajectories, dump to XYZ files and plot the energy contributions. This program does not derive the force field.

Arguments

system
a Yaff System instance defining the system
ai
a Reference instance corresponding to the ab initio input data
settings
a Settings instance defining all QuickFF settings

Optional Arguments

ffrefs
a list of Reference instances defining the a-priori force field contributions.
class quickff.program.DeriveFF(system, ai, settings, ffrefs=[])[source]

Derive a force field for the given system. After the hessian fit of the force constants, the rest values are refined by revisiting the perturbation trajectories with an extra a priori term representing the current valence contribution. Finally, the force constants are refined by means of a final hessian fit.

Arguments

system
a Yaff System instance defining the system
ai
a Reference instance corresponding to the ab initio input data
settings
a Settings instance defining all QuickFF settings

Optional Arguments

ffrefs
a list of Reference instances defining the a-priori force field contributions.

quickff.reference – Representation of the ab initio reference data

class quickff.reference.SecondOrderTaylor(name, coords=None, energy=0.0, grad=None, hess=None, pbc=[0, 0, 0])[source]

Second-order Taylor expansion model, can be used for the ab initio input.

energy(coords)[source]

Compute the energy for the given positions

gradient(coords)[source]

Compute the gradient for the given positions

hessian(coords)[source]

Compute the hessian for the given positions

update(coords=None, grad=None, hess=None, pbc=None)[source]

Method to update one or more of the attributes. The meaning of the optional arguments is the same as with the constructor __init__

class quickff.reference.YaffForceField(name, ff)[source]

A model object based on a YAFF force field. Such a model can be used to account for non-covalent interactions or residual covalent interactions.

quickff.reference.get_ei_ff(name, system, charges, scales, radii=None, average=True, pbc=[0, 0, 0])[source]

A routine to construct a Yaff force field for the electrostatics

Arguments

name
A string for the name of the force field. This name will show in possible plots visualizing contributions along perturbation trajectories
system
A Yaff System instance representing the system
charges
A numpy array containing the charge of each separate atom
scales
A list of 4 floats, representing scale1, scale2, scale3, scale4, i.e. the electrostatic scaling factors

Optional Arguments

radii
A numpy array containing the gaussian radii of each separate atom. If this argument is omitted, point charges are used.
average
If set to True, the charges and radii will first be averaged over atom types. This is True by default.

quickff.tools – Various convenient tools

quickff.tools.global_translation(coords)[source]

A function to generate vectors that represent global translations of a system.

Arguments

coords
a (N,3) numpy array describing the system that has to be translated
quickff.tools.global_rotation(coords)[source]

A function to generate vectors that represent global translations of a system. Rx is a matrix of rotatino around the x-axis minus the identity matrix. VRx is a vector of rotation around x-axis.

Arguments

coords
a (N,3) numpy array describing the system that has to be translated
quickff.tools.fitpar(xs, ys, rcond=-1)[source]

Fit a parabola to the samples (xs, ys):

ys[:] = a*xs[:]^2 + b*xs[:] + c

Returns the parabola parameters a, b and c.

Arguments

xs
a (N) numpy array containing the x values of the samples
ys
a (N) numpy array containing the x values of the samples
quickff.tools.boxqp(A, B, bndl, bndu, x0, threshold=1e-09, status=False)[source]

Minimize the function

1/2*xT.A.x - B.x

subject to

bndl < x < bndu (element-wise)

This minimization is performed using a projected gradient method with step lengths computed using the Barzilai-Borwein method. See 10.1007/s00211-004-0569-y for a description.

Arguments
A (n x n) NumPy array appearing in cost function B (n) NumPy array appearing in cost function bndl (n) NumPy array giving lower boundaries for the variables bndu (n) NumPy array giving upper boundaries for the variables x0 (n) NumPy array providing an initial guess
Optional Arguments
threshold Criterion to consider the iterations converged status Return also the number of iterations performed
quickff.tools.set_ffatypes(system, level)[source]

A method to guess atom types. This will overwrite ffatypes that are already defined in the system.

Arguments:

system
A yaff system instance
level

If level is a string containing comma’s, it is assumed to be an ordered list containing the atom type of each atom in the system. Otherwise, level is assumed to be a string defining how to guess the atom types from the local topology. Possible levels are:

  • low: based on atomic number
  • medium: based on atomic number and number of neighbors
  • high: based on atomic number, number of neighbors and atomic number of neighbors
  • highest: based on index in the molecule
quickff.tools.term_sort_atypes(ffatypes, indexes, kind)[source]

Routine to sort the atoms defined in indexes to give consistent term names. This routine returns the sorted atom indexes as well as the corresponding atom types.

quickff.tools.get_multiplicity(n1, n2)[source]

Routine to estimate m from local topology

quickff.tools.get_restvalue(values, m, thresshold=0.3490658503988659, mode=1)[source]

Get a rest value of 0.0, 360/(2*m) or None depending on the given equilbrium values.

For mode=0, the rest value is:

0, if all ‘values modulo per’ are in the interval [0,thresshold] U [per-thresshold,per]

per/2, if all ‘values modulo per’ are in the interval [per/2-thresshold,per/2+thresshold]

None, in all other cases

For mode=1, the rest value is determined as follows:

first the values are folded in the interval [0,per/2] by first taking the module with per and then mirroring values in [per/2,per] on [0,per/2]. Next, the mean and std of the folded values are computed. If the std is larger then the thresshold, the values are considered to be too scattered and no rest value can be computed (None is returned). If the std is small enough, the rest value will be determined based on the mean. If the mean is close enough to 0, the rest value will be 0. If the mean is close enough to per/2, the rest value will be per/2. In all other cases no rest value will be computed and None is returned.
quickff.tools.get_ei_radii(numbers)[source]

Routine to return atomic radii for use in the Gaussian charge distribution. These radii are computed according to the procedure of Chen and Slater:

First the Slater exponent is computed from the hardness using the formula of Rappe and Goddard (hardness of Pearson and Parr is used)

Next the gaussian exponent alpha is fitted by minimizing the L2-difference between the between the homonuclear Coulomb integral over Slater orbitals and over Gaussian orbitals.

quickff.tools.digits(x, n)[source]

returns a string representation of x with exactly n digits if possible.

quickff.tools.average(data, ffatypes, fmt='full', verbose=False)[source]

Average the atomic parameters stored in data over atoms of the same atom type.

Arguments

data
a list or numpy array containing the data
ffatypes
a listor numpy array containing the atom types. Should have equal length as data

Keywork arguments

fmt

Should be either full, dict or sort. In case of full, the result will be returned as an numpy array of equal length as data and ffatypes. In case of dict, the result will be returned as a dictionairy of the following format:

{ffatype0: value0, ffatype1: value1, …}

in which value0, … is the mean value of the given ffatype. In case of sort, a dictionairy of the following format will be returned:

{ffatype0: values0, ffatype1: values1, …}

in which values0, … is a list of the values for the given ffatype.

quickff.valence – ValenceFF class storing all valence ff terms

class quickff.valence.ValenceFF(system, settings)[source]

Class to collect all valence terms in the force field for which parameters need to be estimated.

Arguments

system
an instance of the Yaff System class containing all system properties
settings
an instance of Settings containing all QuickFF settings
add_term(pot, ics, basename, tasks, units, diag_term_indexes=[])[source]

Adds new term both to the Yaff vlist object and a new QuickFF list (self.terms) which holds all information about the term for easy access later in QuickFF.

Arguments

pot
an instance of ValenceTerm from yaff.pes.vlist.py representing the potential chosen for the term
ics
list of InternalCoordinate instances from yaff.pes.iclist.py
basename
base name for the term, identical for master and each slave
tasks
List of strings defining all tasks that have to be performed for this term. Possible tasks are: PT_ALL, HC_FC_DIAG, HC_FC_CROSS, EQ_RV.
units
Units for all parameters in this term (ordered in the same way as parameters are stored in yaff.vlist.vtab)

Optional arguments

diag_term_indexes
Indexes of the diagonal terms that correspond to the ics composing an off-diagonal term. Empty if the term is diagonal.
apply_constraints(constraints)[source]

Routine to apply equality constraints in the force field fitting

Arguments

contraints
A dictionairy containing (master: slaves) definitions in which master is a string defining the master basename and slaves is a list of strings defining the slave basenames.
check_params(term, labels)[source]

Check whether the given term has all given pars defined in labels.

Arguments

term
An instance of the Term class defining the term that has to be checked
labels
A list of strings defining which parameters should be checked. only arguments of the only option of Valence.get_params are allowed.
get_hessian_contrib(index, fc=None)[source]

Get the contribution to the covalent hessian of term with given index (and its slaves). If fc is given, set the fc of the master and its slave to the given fc.

init_bend_terms(thresshold=0.17453292519943295)[source]

Initialize all bend terms in the system based on the bends attribute of the system instance. All bend terms are given harmonic potentials either in the angle or the cosine of the angle.

init_bond_terms()[source]

Initialize all bond terms in the system based on the bonds attribute of the system instance. All bond terms are given harmonic potentials.

init_cross_angle_terms()[source]

Initialize cross terms between bonds and bends.

init_cross_dihed_terms()[source]

Initialize cross terms between diheds and bonds,bends.

init_dihedral_terms(thresshold=0.3490658503988659)[source]

Initialize the dihedral potentials from the local topology. The dihedral potential will be one of the two following possibilities:

The multiplicity m is determined from the local topology, i.e. the number of neighbors of the central two atoms in the dihedral

If the equilibrium value of all instances of the torsion are within thresshold of 0 deg or per/2 with per = 180deg/m, the following potential will be chosen:

0.5*K*(1-cos(m*psi-m*psi0)) with psi0 = 0 or 360/(2*m)
init_oop_terms(thresshold_zero=0.09448630669606262)[source]

Initialize all out-of-plane terms in the system based on the oops attribute of the system instance. All oops are given harmonic potentials.

is_negligible(term_index)[source]

Return True if the given term can be neglected (e.g. in parameter files).

iter_masters(label=None, use_re=False)[source]

Iterate over all master terms in the valence force field. If label is given, only iterate of the terms with the label in its name.

iter_terms(label=None, use_re=False)[source]

Iterate over all terms in the valence force field. If label is given, only iterate over terms that matches label.

modify_term(term_index, pot, ics, basename, tasks, units)[source]

Modify the term with given index to a new valence term.

Scripts

An overview is provided of the scripts in QuickFF as well as their basic usage and options. The information below can also be obtained by running the script with the --help option.

usage: qff.py [-h] [--version] [-s] [-v] [-V] [-l LOGFILE] [--scoop]
              [-c CONFIG_FILE]
              [-m {MakeTrajectories,PlotTrajectories,DeriveFF}]
              [--fn-traj FN_TRAJ] [--only-traj ONLY_TRAJ] [-p PLOT_TRAJ] [-x]
              [--suffix SUFFIX] [--ei EI] [--ei-rcut EI_RCUT] [--vdw VDW]
              [--vdw-rcut VDW_RCUT] [--covres COVRES]
              [--ffatypes {None,list_of_atypes,low,medium,high,highest}]
              fn [fn ...]

This script will apply QuickFF to derive a covalent force field for the given
system from the ab initio input given in the input files.

positional arguments:
  fn                    Input file name that specify the system and ab initio
                        reference data. Multiple file names are allowed, but
                        at least one should be given. Files later in the list
                        overwrite information from earlier files. Allowed file
                        formats are MolMod checkpoint files (file.chk),
                        Gaussian formatted checkpoint files (file.fchk) and
                        VASP xml files (file.xml).

optional arguments:
  -h, --help            show this help message and exit
  --version             show program's version number and exit
  -s, --silent          Swith of all logging completely, overwrites all other
                        verbosity options.
  -v, --verbose         Increases verbosity, is overwriten if SILENT or
                        VERY_VERBOSE is switched on.
  -V, --very-verbose    Increases verbosity to highest level, is overwriten if
                        SILENT is switched on.
  -l LOGFILE, --logfile LOGFILE
                        Redirect logger output to a file with name LOGFILE.
  --scoop               Flag to enable parallelisation using SCOOP. With
                        SCOOP, the command to run QuickFF is slightly
                        different, the absolute path to quickff.py should be
                        used. For example, to run on 4 cores: python -m scoop
                        -n4 /path/to/qff.py --scoop [options] fns

General QuickFF specifications:
  -c CONFIG_FILE, --config-file CONFIG_FILE
                        Specify a configuration file to read all QuickFF
                        settings from.
  -m {MakeTrajectories,PlotTrajectories,DeriveFF}, --program-mode {MakeTrajectories,PlotTrajectories,DeriveFF}
                        Specify the program mode which defines the set of
                        instructions that will be executed.
  --fn-traj FN_TRAJ     Read/write the perturbation trajectories from/to
                        FN_TRAJ. If the given file exists, the trajectories
                        are read from the file. Otherwise, the trajectories
                        are written to the given file.
  --only-traj ONLY_TRAJ
                        Construct the perturbation trajectory only for the
                        terms with the given basenames. This options is only
                        applied in the MakeTrajectories program.
  -p PLOT_TRAJ, --plot-traj PLOT_TRAJ
                        If set to final, plots the various energy
                        contributions along the perturbation trajectories to
                        using the final force field. If set to all, plots the
                        contributions along the trajectories using all
                        intermediate force fields (given suffixes _Apt1,
                        _Bhc1, _Cpt2 and _Dhc2) as well as the final force
                        field (given the suffix _Ehc3).
  -x, --xyz-traj        Write the perturbation trajectories in XYZ format.
  --suffix SUFFIX       Suffix that will be added to all output files.

Options related to the definition and derivation of the force field:
  --ei EI               A Yaff parameters file defining the electrostatic
                        contribution of the force field.
  --ei-rcut EI_RCUT     The real space cut off for the electrostatic
                        interactions. If the system is periodic, the ewald
                        parameters will be adjusted to this cut off.
  --vdw VDW             A Yaff parameters file defining the van der Waals
                        contribution of the force field.
  --vdw-rcut VDW_RCUT   The real space cut off for the van der Waals
                        interactions.
  --covres COVRES       A Yaff parameters file defining a residual
                        contribution to the covalent part of the force field.

Options related to the definition of the system:
  --ffatypes {None,list_of_atypes,low,medium,high,highest}
                        Assign atom types in the system by parsing an ordered
                        list ofatom types as argument or through the automatic
                        built-in detection (see documentation). By default (or
                        if None is given), the atom types are assumed to be
                        defined in the input files. [default=None]
usage: qff-input-ei.py [-h] [-v]
                       [--ffatypes {None,list_of_atypes,low,medium,high,highest}]
                       [--gaussian] [--bci]
                       [--bci-constraints BCI_CONSTRAINTS]
                       fn_sys charges [fn_out]

This script reads atomic charges from an input file and makes a Yaff
parameters file suitable for the QuickFF option --ei.

positional arguments:
  fn_sys                Any file from which the system can be extracted
                        (MolMod CHK, Gaussian FCHK, XYZ, ...).
  charges               The atomic charges to be used. This argument has the
                        form fn_charges:path, where fn_charges is a filename
                        that contains the charges and path refers to a
                        location within that file where the charges can be
                        found. If fn_charges is an HDF5 file, path is the
                        location of a dataset containing the charges, e.g.
                        '/charges'. If fn_charges is a MolMod CHK file, path
                        is the label of the dataset that contains the atomic
                        charges.
  fn_out                Name of the Yaff file to write the parameters to.
                        [default=pars_ei.txt]

optional arguments:
  -h, --help            show this help message and exit
  -v, --verbose         Increase verbosity of the script [default=False].
  --ffatypes {None,list_of_atypes,low,medium,high,highest}
                        Assign atom types in the system by parsing an ordered
                        list ofatom types as argument or through the automatic
                        built-in detection (see documentation). By default (or
                        if None is given), the atom types are assumed to be
                        defined in the input files. [default=None]
  --gaussian            Use gaussian smeared charges. The radii are taken from
                        the input file fn_in (from dataset /path/radii for
                        HDF5 file or from label `radii` for CHK file) if the
                        data is present, otherwise the radii are estimated
                        according to the procedure of Chen et al. See
                        ``quickff.tools.get_ei_radii`` for more info.
  --bci                 Convert averaged atomic charges to bond charge
                        increments, i.e. charge transfers along the chemical
                        bonds in the system. In this way, one is certain of
                        ending up with a globally neutral system even if bcis
                        from different systems are combined. This option
                        requires the definition of bonds, if the bonds cannot
                        be read from the system file, they will be estimated
                        from the interatomic distances using the detect_bonds
                        routine in the Yaff System class. [default=False]
  --bci-constraints BCI_CONSTRAINTS
                        A file containing constraints for the charge to bci
                        fit in a master: slave0,slave1,...: sign format. A new
                        line should be used for each master and the format is
                        insensitive towards spaces.Sign should be 1.0 or -1.0
                        indicating wheter or not a sign switch should be
                        introduced when mapping the slaves to the master.