Code Documentation

diffiqult.Task

class diffiqult.Task.Tasks(mol, name, verbose=False)[source]

This class manage the implemented tasks over a system included in DiffiQult.

Attributes:

sys : System_mol object
Contains the basis functions, geometry and information of the a molecular system.
name : string
An id of the task, it is used as prefix for outputs.
verbose: bool

It defines the output and screen options.

True it prints in screen all details and a file “name.out”.

status: bool

Keep track of the success of the different task.

True the SCF and/or the optimization converged.

__init__(mol, name, verbose=False)[source]

Initialize task object that contains the inital parameters.

Parameters:

mol : System_mol object
Contains the basis functions, geometry and information of the a molecular system.
name : string
An id of the task, it is used as prefix for outputs.
Options:
verbose : bool

It defines the output and screen options

True it prints in screen all details and a file “name.out”

status: bool

Keep track of the success of the different task.

True the SCF and/or the optimization converged.

energy(max_scf=30, max_d=300, printguess=None, output=False, **kwargs)[source]

This function handdles the single point calculations

Options:

max_scf : integer
Maximum number of scf steps, default 30.

printguess : str

File path if it is requiered to prepare an inital guess for the molecular orbital coefficients.

output : bool

True if it will print a molden file in case of success.
optimization(max_scf=100, log=True, scf=True, name='Output', readguess=None, output=False, argnum=[0], **kwargs)[source]

This function handdles the optimization procedure

Options:

argnum : list of integers

Parameter to optimize

0 widths

1 contraction coefficients

2 Gaussian centers

e.g. [0,1] to optimized widhts and contraction coefficients

max_scf : integer
Maximum number of scf steps, default 30
log : bool
If we are not optimizing the log of exponent, we highly recoment leave it as True, the default.
name : str
Output file name default Output
readguess : str
File path to a npy file in case on predefined initial guess of the density matrix
output : bool
True if it will print a molden file in case of success
runtask(task, **kwargs)[source]

This method run a given task and if it has success, it uptates system with the most recent energy value and basis function

Parameters:

task : str
If defines the task:

‘Energy’ is a single point calculation.

‘Opt’ an optimization of a given parameter.

Options:
Check documentation for each task

Returns:

success : bool
True if task ended successfully.

diffiqult.Molecule

class diffiqult.Molecule.Getbasis(molecule, basis, shifted=False)[source]

This class loads the basis function for a given molecule Attributes:

alpha : array
Gaussian withds
coef : array
Contraction coefficients
xyz : array
Gaussian centers
list_contr : list
List of integers with contractions, eg. [3,3] are two AO with three primitives each
l : list of tuples
List of tuple of integers with angular momentums eg. [(0,0,0),(0,1,0)] one s orbital and one p_x
__init__(molecule, basis, shifted=False)[source]
Parameters:
molecule : list

It contains spects of geometry [( atomic_number_atom_1,(x,y,z),

atomic_number_atom_1,(x,y,z), atomic_number_atom_1,(x,y,z)]

eg. water

(8,(0.0, 0.0, 0.091685801102911746)), (1,(1.4229678834888837, 0.0, -0.98120954931681137)), (1,(-1.4229678834888837, 0.0, -0.98120954931681137))]
basis : dict
{atomic_number: [(‘type’,[(exp,coef),
(exp,coef)]),
(‘type’,[(exp,coef),
(exp,coef)])],
atomic_number: [(‘type’,[(exp,coef),
(exp,coef)]),
(‘type’,[(exp,coef),
(exp,coef)])]}
class diffiqult.Molecule.Getgeom(molecule)[source]
This class loads the geometry in xyz coordinates
Parameters:
molecule : list

It contains spects of geometry [( atomic_number_atom_1,(x,y,z),

atomic_number_atom_1,(x,y,z), atomic_number_atom_1,(x,y,z)]
class diffiqult.Molecule.System_mol(mol, basis_set, ne, mol_name='molecule', angs=False, shifted=False)[source]

This class contains all the information of the system extracted from mol and basis

Parameters:

alpha : array
Gaussian withds
coef : array
Contraction coefficients
xyz : array
Gaussian centers
list_contr : list
List of integers with contractions, eg. [3,3] are two AO with three primitives each
l : list of tuples
List of tuple of integers with angular momentums eg. [(0,0,0),(0,1,0)] one s orbital and one p_x
__init__(mol, basis_set, ne, mol_name='molecule', angs=False, shifted=False)[source]

This function starts an object System_mol

Parameters:
mol : list
It contains spects of geometry
[( atomic_number_atom_1,(x,y,z),
atomic_number_atom_1,(x,y,z), atomic_number_atom_1,(x,y,z)]
basis: dict
{atomic_number: [(‘type’,[(exp,coef),
(exp,coef)]),
(‘type’,[(exp,coef),
(exp,coef)])],
atomic_number: [(‘type’,[(exp,coef),
(exp,coef)]),
(‘type’,[(exp,coef),
(exp,coef)])]}
ne : int
Number of electron
mol_name : str
An id
shifted : bool

False in case standard basis functions

True in case costumized basis functions

angs : bool

False is atomic units

True if lenght units are angstroms

printcurrentgeombasis(tape)[source]

This method prints the current values of sys on tape Parameters:

tape : file obj
Output file

diffiqult.Energy

diffiqult.Energy.rhfenergy(alpha_old, coef2, xyz, l, charges, xyz_atom, natoms, nbasis, contr_list, ne, max_scf, max_d, log, eigen, printguess, readguess, name, write, dtype)[source]

This function returns the rhf function

Parameters:
alpha_old : array
Gaussian exponents
coef2 : array
Contraction coeffients
xyz : array 3N
Gaussian centers
l : array 3N
Angular momentum each entry is a vector eg. s orbital (0,0,0) or pz (1,0,0)
charges : array
Atom charges
nbasis : int
Number of basis
contr_list: list of integers
Specify the number of orbitals in each atom
ne : int
Number of electrons
max_scf : int
maximum number of scf cycles
log : bool
The exponents are given in log
printguess: str or None
File to print coeff matrix initial guess
readguess : str or None
File that contains coeff matrix initial guess
name : str
Output file name
write : bool
True if printing
dtype : type of output
This is the directive to know if algopy will be used or not np.float64(1.0) if it is a single point calculation otherwise, it specify the size of the UTMP, autodifferentiation
Returns:
energy : float
RHF energy

diffiqult.Optimize

This routine was orginally taken from scipy

This module was taken from scipy and modified to optimize several parameters at a time with the BFGS Minimization

diffiqult.Optimize.algopy_wrapper(function)[source]

Returnfunction to obtain the gradient of the function Parameters ———- f : callable

The function of which to determine the gradient (partial derivatives). Should take xk as first argument, other arguments to f can be supplied in *args. Should return a scalar, the value of the function at xk.
grad_function : function
The partial derivatives of f to xk.
diffiqult.Optimize.fmin_bfgs(f, x0, fprime=None, args=(), argnum=None, gtol=1e-05, norm=inf, epsilon=1.4901161193847656e-08, maxiter=None, full_output=0, disp=1, retall=0, callback=None)[source]

Minimize a function using the BFGS algorithm.

f : callable f(x,*args)
Objective function to be minimized.
x0 : ndarray
Initial guess.
fprime : callable f’(x,*args), optional
Gradient of f.
args : tuple, optional
Extra arguments passed to f and fprime.
gtol : float, optional
Gradient norm must be less than gtol before successful termination.
norm : float, optional
Order of norm (Inf is max, -Inf is min)
epsilon : int or ndarray, optional
If fprime is approximated, use this value for the step size.
callback : callable, optional
An optional user-supplied function to call after each iteration. Called as callback(xk), where xk is the current parameter vector.
maxiter : int, optional
Maximum number of iterations to perform.
full_output : bool, optional
If True,return fopt, func_calls, grad_calls, and warnflag in addition to xopt.
disp : bool, optional
Print convergence message if True.
retall : bool, optional
Return a list of results at each iteration if True.
xopt : ndarray
Parameters which minimize f, i.e. f(xopt) == fopt.
fopt : float
Minimum value.
gopt : ndarray
Value of gradient at minimum, f’(xopt), which should be near 0.
Bopt : ndarray
Value of 1/f’‘(xopt), i.e. the inverse hessian matrix.
func_calls : int
Number of function_calls made.
grad_calls : int
Number of gradient calls made.
warnflag : integer
1 : Maximum number of iterations exceeded. 2 : Gradient and/or function calls not changing.
allvecs : list
OptimizeResult at each iteration. Only returned if retall is True.
minimize: Interface to minimization algorithms for multivariate
functions. See the ‘BFGS’ method in particular.

Optimize the function, f, whose gradient is given by fprime using the quasi-Newton method of Broyden, Fletcher, Goldfarb, and Shanno (BFGS)

Wright, and Nocedal ‘Numerical Optimization’, 1999, pg. 198. This function was originally taken and modified from scipy