Previous topic

Gamess(US) tutorial

This Page

API Reference

Calculators

Base Calculator

class Calculator(exevar=None, executable=None, runopts=None, scratch=None)[source]

Abstract class that should be subclassed when adding a new code interace.

accomplished(fname)[source]

Return True if the job completed without errors

property executable

Set the executable

abstract parse(fname, objective, regexp=None)[source]

Parse the value of the objective from the file fname

abstract run(*args)[source]

run a single job

abstract run_multiple(*args)[source]

Run multiple jobs

property scratch

Return scratch

abstract write_input()[source]

write the input file in the format of the code used

Dalton

class Dalton(name='Dalton', **kwargs)[source]

Wrapper for running Dalton program

parse(fname, objective, regularexp=None)[source]

Parse a value from the output file fname based on the objective.

If the value of the objective is regexp then the regularexp will be used to parse the file.

run(fname)[source]

Run a single job

Args:
fnamedict

A dictionary with keys mol and dal and their respective file name strings as values

Returns:
outstr

Name of the dalton output file

run_multiple(fnames)[source]

Spawn two single jobs as paralell processes

write_input(fname, template, basis, mol, core)[source]

Write dalton input files: fname.dal and system.mol

Args:
fnamestr

Name of the input file .dal

templatedict

Dictionary with templates for the dal and mol with those strings as keys and actual templates as values

basisdict

An instance of BasisSet class or a dictionary of BasisSet objects with element symbols as keys

molchemtools.molecule.Molecule

Molecule object with the system geometry

corestr

Core definition

DMFT

class Dmft(logfile)[source]

dmft class

parse_gamess()[source]

Parse gamess-us log file to get the neccessary data to write dmft input and set defaults.

run()[source]

Run a single dmft job.

write_input(functional=None, a1=None, b1=None, a2=None, b2=None, printlevel=None, analyze=None)[source]

Write dmft input based on information in the gamess log file.

GAMESS-US

class GamessUS(name='GamessUS', version='00', runopts=None, **kwargs)[source]

Container object for Gamess-us jobs.

accomplished(outfile)[source]

Return True if Gamess(US) job finished without errors.

get_command(inpfile)[source]

Return the command to execute

parse(output, objective, regexp=None)[source]

Parser GAMESS(US) output file to get the objective.

remove_dat(inpfile)[source]

Remove the gamess dat file if it exists in the scratch directory.

run(inpfile, logfile=None, remove_dat=True)[source]

Run a single gamess job interactively - without submitting to the queue.

run_multiple(inputs)[source]

Run multiple jobs

property runopts

Return the runopts

property version

Return the version.

write_input(fname, template=None, mol=None, basis=None, core=None)[source]

Write the molpro input to “fname” file based on the information from the keyword arguments.

Args:
molchemtools.molecule.Molecule

Molecule object instance

basisdict or BasisSet

An instance of BasisSet class or a dictionary of BasisSet objects with element symbols as keys

corelist of ints

Molpro core specification

templatestr

Template of the input file

fnamestr

Name of the input file to be used

Molpro

class Molpro(name='Molpro', **kwargs)[source]

Wrapper for the Molpro program.

accomplished(fname)[source]

Return True if the job completed without errors

parse(fname, objective, regularexp=None)[source]

Parse a value from the output file fname based on the objective.

If the value of the objective is regexp then the regularexp will be used to parse the file.

run(inpfile)[source]

Run a single molpro job interactively - without submitting to the queue.

run_multiple(inputs)[source]

Run a single molpro job interactively - without submitting to the queue.

write_input(fname=None, template=None, mol=None, basis=None, core=None)[source]

Write the molpro input to “fname” file based on the information from the keyword arguments.

Args:
molchemtools.molecule.Molecule

Molecule object instance

basisdict or BasisSet

An instance of BasisSet class or a dictionary of BasisSet objects with element symbols as keys

corelist of ints

Molpro core specification

templatestr

Template of the input file

fnamestr

Name of the input file to be used

PSI4

class Psi4(name='Psi4', **kwargs)[source]

Wrapper for the Psi4 program.

parse(fname, objective, regularexp=None)[source]

Parse a value from the output file fname based on the objective.

If the value of the objective is regexp then the regularexp will be used to parse the file.

run(inpfile)[source]

Run a single Psi4 job interactively - without submitting to the queue.

run_multiple(inputs)[source]

Run a single Psi4 job interactively - without submitting to the queue.

write_input(fname, template, mol=None, basis=None, core=None)[source]

Write the Psi4 input to “fname” file based on the information from the keyword arguments.

Args:
molchemtools.molecule.Molecule

Molecule object instance

basisdict or BasisSet

An instance of BasisSet class or a dictionary of BasisSet objects with element symbols as keys

corelist of ints

Psi4 core specification

templatestr

Template of the input file

fnamestr

Name of the input file to be used

Basis Set Tools

BasisSet module

class BasisSet(name, element, family=None, kind=None, functions=None, info=None)[source]

Basis set module supporting basic operation on basis sets.

Args:
namestr

Name of the basis set

elementstr

Symbol of the element

Kwargs:
kindstr

Classification of the basis set functions, diffuse, tight

familystr

basis set family

functionsdict

Dict of functions with s, p, d, f, … as keys

paramslist of dicts

Parameters for generating the functions according to the model

append(other)[source]

Append functions from another BasisSet object

Args:
otherBasisSet

BasisSet object whose functions will be added to the existing ones

completeness_profile(zetas)[source]

Calculate the completeness profile of each shell of the basis set

Args:
zetasnumpy.array

Scaning exponents

Returns:
outnumpy.array

numpy array with values of the profile (shells, zetas)

contraction_matrix(shell)[source]

Return the contraction coefficients for a given shell in a matrix form with size ne * nc, where ne is the number of exponents and nc is the number of contracted functions

Args:
shellstr

shell label, s, p, d, …

Returns:
out2D numpy.array

2D array with contraction coefficients

contraction_scheme()[source]

Return a string describing the contraction scheme.

contraction_type()[source]

Try to determine the contraction type: segmented, general, uncontracted, unknown.

contractions_per_shell()[source]

Calculate how many contracted functions are in each shell.

Returns:

out : list of ints

classmethod from_file(fname=None, fmt=None, name=None)[source]

Read and parse a basis set from file and return a BasisSet object

Args:
fnamestr

File name

fmtstr

Format of the basis set in the file (molpro, gamessus)

namestr

Name of the basis set

Returns:
outBasisSet or dict

Basisset object parsed from file or dictionary of BasisSet objects

classmethod from_json(jsonstring, **kwargs)[source]

Instantiate the BasisSet object from a JSON string

Args:
jsonstring: str

A JSON serialized string with the basis set

classmethod from_optpars(x0, funs=None, name=None, element=None, explogs=False)[source]

Return a basis set object generated from a sequence based on the specified arguments.

Args:
x0list or numpy.array

Parameters to generate the basis set as a continuous list or array

funslist of tuples

A list of tuple specifying the shell type, number of functions and parameters, e.g. [(‘s’, ‘et’, 4, (0.5, 2.0)), (‘p’, ‘et’, 3, (1.0, 3.0))]

namestr

Name of the basis set

elementstr

Chemical symbol of the element

Returns:

out : BasisSet

static from_pickle(fname, **kwargs)[source]

Read a pickled BasisSet object from a pickle file

Args:
fnamestr

File name containing the BasisSet

kwargsdict

Extra arguments for the pickle.load method

Raises:
UnicodeDecodeError

When you try to read a python2 pickle using python3, to fix that use encoding=’latin1’ option

classmethod from_sequence(funs=None, name=None, element=None)[source]

Return a basis set object generated from a sequence based on the specified arguments.

Args:
funslist of tuples

A list of tuple specifying the shell type, number of functions and parameters, e.g. [(‘s’, ‘et’, 4, (0.5, 2.0)), (‘p’, ‘et’, 3, (1.0, 3.0))]

namestr

Name of the basis set

elementstr

Chemical symbol of the element

Returns:

out : BasisSet

classmethod from_str(string, fmt=None, name=None)[source]

Parse a basis set from string

Args:
stringstr

A string with the basis set

fmtstr

Format of the basis set in the file: molpro, gamessus

namestr

Name of the basis set

Returns:
outBasisSet or dict

Basisset object parsed from string or dictionary of BasisSet objects

get_exponents(asdict=True)[source]

Return the exponents of a given shell or if the shell isn’t specified return all of the available exponents

Args:
asdict (bool): if True a dict with exps per shell is

returned

nf(spherical=True)[source]

Calculate the number of basis functions

Args:
sphericalbool

flag indicating if spherical or cartesian functions should be used, default: True

Returns:
resint

number of basis functions

normalization()[source]

For each function (contracted) calculate the norm and return a list of tuples containing the shell, function index and the norm respectively.

normalize()[source]

Normalize contraction coefficients for each contracted functions based on the primitive overlaps so that the norm is equal to 1.

nprimitive(spherical=True)[source]

Return the number of primitive functions assuming sphrical or cartesian gaussians.

Args:
sphericalbool

A flag to select either spherical or cartesian gaussians

Returns:
outint

Number of primitive function in the basis set

partial_wave_expand()[source]

From a given basis set with shells spdf… return a list of basis sets that are subsets of the entered basis set with increasing angular momentum functions included [s, sp, spd, spdf, …]

primitives_per_contraction()[source]

Calculate how many primities are used in each contracted function.

Returns:

out : list of ints

primitives_per_shell()[source]

Calculate how many primitive functions are in each shell.

Returns:

out : list of ints

print_functions(efmt='20.10f', cfmt='15.8f')[source]

Return a string with the basis set.

Args:
efmtstr

string describing output format for the exponents, default: “20.10f”

cfmtstr

string describing output format for the contraction coefficients, default: “15.8f”

Returns:
resstr

basis set string

shell_overlap(shell)[source]

Calculate the overlap integrals for a given shell

Args:
shellstr

Shell

Returns:
outnumpy.array

Overlap integral matrix

sort(reverse=False)[source]

Sort shells in the order of increasing angular momentum and for each shell sort the exponents.

Args:
reversebool

If False sort the exponents in each shell in the descending order (default), else sort exponents in ascending order

to_cfour(comment='', efmt='15.8f', cfmt='15.8f')[source]

Return a string with the basis set in (new) CFOUR format.

Args:
commentstr

comment string

efmtstr

string describing output format for the exponents, default: “20.10f”

cfmtstr

string describing output format for the contraction coefficients, default: “15.8f”

Returns:
resstr

basis set string in Cfour format

to_dalton(fmt='prec')[source]

Return a string with the basis set in DALTON format.

Args:
fmtstr
string describing output format for the exponents and coefficents
  • prec “20.10f”

  • default “10.4f”

  • or python format string e.g. “15.8f”

Returns:
resstr

basis set string in Dalton format

to_gamessus(efmt='20.10f', cfmt='15.8f')[source]

Return a string with the basis set in GAMESS(US) format.

Args:
efmtstr

string describing output format for the exponents, default: “20.10f”

cfmtstr

string describing output format for the contraction coefficients, default: “15.8f”

Returns:
resstr

basis set string in Gamess(US) format

to_gaussian(efmt='20.10f', cfmt='15.8f')[source]

Return a string with the basis set in Gaussian format.

Args:
efmtstr

string describing output format for the exponents, default: “20.10f”

cfmtstr

string describing output format for the contraction coefficients, default: “15.8f”

Returns:
resstr

basis set string in Gaussian format

to_json(fname=None, **kwargs)[source]

Serizalize the basis set object to JSON format

to_latex(efmt='20.10f', cfmt='15.8f')[source]

Return a string with the basis set as LaTeX table/

Args:
efmtstr

Output format for the exponents, default: “20.10f”

cfmtstr

Output format for the contraction coefficients, default: “15.8f”

Returns:
resstr

basis set string in LaTeX format

to_molpro(withpars=False, efmt='20.10f', cfmt='15.8f')[source]

Return a string with the basis set in MOLPRO format.

Args:
withparsbool

A flag to indicate whether to wrap the basis with basis={ } string

efmtstr

Output format for the exponents, default: “20.10f”

cfmtstr

Output format for the contraction coefficients, default: “15.8f”

Returns:
resstr

basis set string

to_nwchem(efmt='20.10f', cfmt='15.8f')[source]

Return a string with the basis set in NWChem format.

Args:
efmtstr

string describing output format for the exponents, default: “20.10f”

cfmtstr

string describing output format for the contraction coefficients, default: “15.8f”

Returns:
resstr

basis set string in NwChem format

to_pickle(fname=None)[source]

Save the basis set in pickle format under the filename fname

Args:
fnamestr

File name

uncontract(copy=False)[source]

Uncontract the basis set. This replaces the contraction coefficients in the current object.

Args:
copybool

If True return an uncontracted copy of the basis set rather than uncontracting in place, default is False.

has_consecutive_indices(shell)[source]

Check if all the contracted functions have consecutive indices

Args:
shelldict

Basis functions for a given shell asa dict with structure {'e' : np.array(), 'cf': [np.array(), np.array(), ...]}

reorder_shell_to_consecutive(shell)[source]

Reorder the exponents so that the indices of the contracted functions have consecutive inidices.

Args:
shelldict

Basis functions for a given shell asa dict with structure {'e' : np.array(), 'cf': [np.array(), np.array(), ...]}

Returns:
shelldict

Same shell as on input but with reordered exponents and relabelled contracted functions

merge(first, other)[source]

Merge functions from two BasisSet objects

Args:
firstBasisSet

First BasisSet object to merge

otherBasisSet

Second BasisSet object to merge

Returns:
ourBasisSet

BasisSet instance with functions from first and other merged

primitive_overlap(l, a, b)[source]

Calculate the overlap integrals for a given shell l and two sets of exponents

\[S(\zeta_i, \zeta_j) = \frac{2^{l + \frac{3}{2}} \left(\zeta_{1} \zeta_{2}\right)^{\frac{l}{2} + \frac{3}{4}}}{\left(\zeta_{1} + \zeta_{2}\right)^{l + \frac{3}{2}}}\]
Args:
lint

Angular momentum quantum number of the shell

a (M,)numpy.array

First vector of exponents

b (N,)numpy.array

Second vector of exponents

Returns:
out (M, N)numpy.array

Overlap integrals

nspherical(l)[source]

Calculate the number of spherical components of a function with a given angular momentum value l.

ncartesian(l)[source]

Calculate the number of cartesian components of a function with a given angular momentum value l.

zetas2legendre(zetas, kmax)[source]

From a set of exponents (zetas), using least square fit calculate the expansion coefficients into the legendre polynomials of the order kmax.

Args:
kmaxint

length of the legendre expansion

zetaslist

list of exponents (floats) to be fitted

Returns:
coeffnumpy.array

numpy array of legendre expansion coeffcients of length kmax

zetas2eventemp(zetas)[source]

From a set of exponents (zetas), using least square fit calculate the approximate $alpha$ and $beta$ parameters from the even tempered expansion.

Args:
zetaslist

list of exponents (floats) to be fitted

Returns:
coeffnumpy.array

numpy array of legendre expansion coeffcients of length kmax

eventemp(numexp, params)[source]

Generate a sequence of nf even tempered exponents accodring to the even tempered formula

\[\zeta_i = \alpha \cdot \beta^{i-1}\]
Args:
numexpint

Number fo exponents to generate

paramstuple of floats

Alpha and beta parameters

Returns:
resnumpy array

Array of generated exponents (floats)

welltemp(numexp, params)[source]

Generate a sequence of nf well tempered exponents accodring to the well tempered fromula

\[\zeta_i = \alpha \cdot \beta^{i-1} \cdot \left[1 + \gamma \cdot \left(\frac{i}{N}\right)^{\delta}\right]\]
Args:
numexpint

Number fo exponents to generate

paramstuple of floats

Alpha, beta, gamma and delta parameters

Returns:
resnumpy.array

Array of generated exponents (floats)

legendre(numexp, coeffs)[source]

Generate a sequence of nf exponents from expansion in the orthonormal legendre polynomials as described in: Peterson, G. A. et.al J. Chem. Phys., Vol. 118, No. 3 (2003), eq. (7).

\[\ln \zeta_i = \sum^{k_{\max}}_{k=0} A_k P_k \left(\frac{2j-2}{N-1}-1\right)\]
Args:
numexpint

Number fo exponents to generate

paramstuple of floats

Polynomial coefficients (expansion parameters)

Returns:
resnumpy array

Array of generated exponents (floats)

xyzlist(l)[source]

Generate an array of \(l_x\), \(l_y\), \(l_z\) components of cartesian gaussian with a given angular momentum value in canonical order.

For exampe:
  • l = 0 generates the result array([[0, 0, 0]])

  • l = 1 generates array([[1, 0, 0], [0, 1, 0], [0, 0, 1])

  • etc.

The functions are coded by triples of powers of x, y, z, namely [1, 2, 3] corresponds to \(xy^{2}z^{3}\).

Args:
lint

Angular momentum value

Returns:
out ((l+1)*(l+2)/2, 3)numpy.array

Array of monomial powers, where columns correspond to x, y, and z respectively and rows correspond to functions.

zlmtoxyz(l)[source]

Generates the expansion coefficients of the real spherical harmonics in terms of products of cartesian components. Method based on [Schelegel1995]

[Schelegel1995]

Schlegel, H. B., & Frisch, M. J. (1995). “Transformation between Cartesian and pure spherical harmonic Gaussians”. International Journal of Quantum Chemistry, 54(2), 83–87. doi:10.1002/qua.560540202

Args:
lint

Angular momentum value

Returns:
out \(((l+1)*(l+2)/2, 2*l + 1)\)numpy.array

Expansion coefficients of real spherical harmonics in terms of cartesian gaussians

basisparse module

Parsing basis set from different formats.

class NumpyEncoder(*, skipkeys=False, ensure_ascii=True, check_circular=True, allow_nan=True, sort_keys=False, indent=None, separators=None, default=None)[source]
default(obj)[source]

If input object is an ndarray it will be converted into a dict holding the data and dtype.

get_l(shell)[source]

Return the orbital angular momentum quantum number for a given subshell

merge_exponents(a, b)[source]

Concatenate the arrays a and b using only the unique items from both arrays

Args:

a : numpy.array b : numpy.array

Returns:
res3-tuple of numpy arrays
  • res[0] sorted union of a and b

  • res[1] indices of a items in res[0]

  • res[2] indices of b items in res[0]

parse_basis(string, fmt=None)[source]

A wrapper for parsing the basis sets in different formats.

Args:
stringstr

A string with the basis set

fmtstr

Format in which the basis set is specified

Returns:
outdict

A dictionary of parsed basis sets with element symbols as keys and basis set functions as values

parse_gamessus_basis(string)[source]

Parse the basis set into a list of dictionaries from a string in gamess format.

parse_gamessus_function(lines)[source]

Parse a basis set function information from list of strings into three lists containg: exponents, indices, coefficients.

Remeber that python doesn’t recognise the 1.0d-3 format where d or D is used to the regex subsitution has to take care of that.

parse_gaussian_basis(string)[source]

Parse the basis set into a list of dictionaries from a string in gaussian format.

parse_gaussian_function(lines)[source]

Parse a basis set function information from list of strings into three lists containg: exponents, indices, coefficients.

Remeber that python doesn’t recognise the 1.0d-3 format where d or D is used to the regex subsitution has to take care of that.

parse_molpro_basis(string)[source]

Parse basis set from a string in Molpro format.

parse_molpro_shell(expsline, coeffs)[source]

Parse functions of one shell in molpro format.

basisopt module

Optimization of basis set exponents and contraction coefficients.

basisopt is a module containing flexible methods for optimizing primitive exponents of basis sets

class BSOptimizer(objective=None, core=None, template=None, regexp=None, verbose=False, code=None, optalg=None, mol=None, fsopt=None, staticbs=None, fname=None, uselogs=True, runcore=False, penalize=None, penaltykwargs=None, logfile=None)[source]

Basis Set Optimizer class is a convenient wrapper for optimizing primitive exponents of Gaussian basis sets using different code interfaces for performing the actual electronic structure calculations.

Args:
objectivestr or callable

Name of the objective using which the optimization will be evaluated, if it is a callable/function it should take the output name as argument

codeCalculator

Subclass of the Calculator specifying the electronic structure program to use

mol:py:class:` Molecule <chemtools.molecule.Molecule>`

Molecule specifying the system

staticbsdict or BasisSet

Basis set or dict of basis set with basis sets whose exponents are not going to be optimized

fsoptdict

A dictionary specifying the basis set to be optimized, the keys should be element/atom symbol and the values should contain lists of 4-tuples composed of: shell, type, number of functions and parameters/exponents, for example,

>>> fs = {'H' : [('s', 'et', 10, (0.5, 2.0))]}

describes 10 s-type exponents for H generated from even tempered formula with parameters 0.5 and 2.0

optalgdict

A dictionary specifying the optimization algorithm and its options

fnamestr

Name of the job/input file for the single point calculator

uselogsbool

Use natural logarithms of exponents in optimization rather than their values, this option is only relevant if functions asre given as exp or exponents

regexpstr

Regular expression to use in search for the objective if objective is regexp

runcorebool

Flag to mark wheather to run separate single point jobs with different numbers of frozen core orbitals to calculate core energy

penalizebool

Flag enabling the use of penalty function, when two exponent in any shell are too close the objective is multiplied by a penalty factor calculated in get_penalty

penaltykwargsdict

Keyword arguments for the penalty function, default {'alpha' : 25.0, 'smallestonly' : True}

get_basis(name=None, element=None)[source]

Construct the BasisSet object from the result of the optimized exponents and function definition.

Args:
namestr

Name to be assigned to the basis set

elementstr

Element symbol for the basis set

Returns:
basischemtools.basisset.BasisSet

BasisSet object with the optimized functions

get_x0()[source]

Collect all the parameters in an array of consecutive elements.

header()[source]

Return the basic information about the execution environment and optimization settings.

jobinfo()[source]

Return the information on the optimization objects and parameters

run()[source]

Start the basis set optimization

Returns:
resOptimizeResult

An instance of the scipy.optimize.OptimizeResult class

get_basis_dict(bso, x0)[source]

Return a dictionary with BasisSet objects as values and element symbols as keys. The dictionary is composed based on the current parameters x0 and attributes of the BSOptimizer including the staticbs

get_penalty(bsdict, alpha=25.0, smallestonly=True)[source]

For a given dict of basis sets calculate the penalty for pairs of exponents being too close together.

Args:
bsdictdict

Dictionary of BasisSet objects

alphafloat

Parameter controlling the magnitude and range of the penalty

smallestonlybool

A flag to mark whether to use only the smallest ratio to calculate the penalty or all smallest ratios from each shell and calculate the penalty as a product of individual penalties

For each basis and shell within the basis ratios between pairs of sorted exponents are calculated. The minimal ratio (closest to 1.0) is taken to calculate the penalty according to the formula

\[P = 1 + \exp(-\alpha (r - 1))\]

where \(r\) is the ratio between the two closest exponents (>1).

opt_multishell(shells=None, nfps=None, guesses=None, max_params=5, opt_tol=0.0001, save=False, bsopt=None, **kwargs)[source]

Optimize a basis set by saturating the function space shell by shell

Kwargs:
shells (list of strings):

list of shells to be optimized, in the order the optimization should be performed,

nfps (list of lists of integers):

list specifying a set of function numbers to be scanned per each shell,

guesses (list of lists of floats):

list specifying a set of starting parameters per each shell,

max_params (int)

maximal number of parameters to be used in the legendre expansion, (length of the expansion)

opt_tol (float):

threshold controlling the termination of the shell optimization, if energy difference between two basis sets with subsequent number of functionsis larger than this threshold, another function is added to this shell and parameters are reoptimized

kwargs:

options for the basis set optimization driver, see driver function from the basisopt module

opt_shell_by_nf(shell=None, nfs=None, max_params=5, opt_tol=0.0001, save=False, bsopt=None, **kwargs)[source]

For a given shell optimize the functions until the convergence criterion is reached the energy difference for two consecutive function numbers is less than the threshold

Kwargs:
shellstring

string label for the shell to be optimized

nfslist of ints

list of integers representing the number of basis functions to be inceremented in the optimization,

max_paramsint

maximal number of parameters to be used in the legendre expansion, (length of the expansion)

opt_tolfloat

threshold controlling the termination of the shell optimization, if energy difference between two basis sets with subsequent number of functionsis larger than this threshold, another function is added to this shell and parameters are reoptimized,

savebool

a flag to trigger saving all the optimized basis functions for each shell,

kwargsdict

options for the basis set optimization driver, see driver function from the basisopt module

Returns:

BasisSet object instance with optimized functions for the specified shell

Raises:
ValueError:

when shell is different from s, p, d, f, g, h, i when number of parameters equals 1 and there are more functions when there are more parameters than functions to optimize

run_core_energy(x0, *args)[source]

Funtion for running two single point calculations and parsing the resulting energy (or property) as specified by the objective function, primarily designed to extract core energy.

Args:
x0: list or numpy.array

contains a list of parameters to be optimized, may be explicit exponents or parametrized exponents in terms of some polynomial

args: tuple of dicts

bsopt, bsnoopt, code, job, mol, opt, needed for writing input and parsing output

Returns:

parsed result of the single point calculation as speficied by the objective function in the “job” dictionary

run_total_energy(x0, *args)[source]

Funtion for running a single point calculation and parsing the resulting energy (or property) as specified by the objective function.

Args:
x0list or numpy.array

contains a list of parameters to be optimized, may be explicit exponents or parametrized exponents in terms of some polynomial

argstuple of dicts

bsopt, bsnoopt, code, job, mol, opt, needed for writing input and parsing output

Returns:

parsed result of the single point calculation as speficied by the objective function in the “job” dictionary

CBS module

Complete Basis Set (CBS) extrapolation techniques.

Module for Complete Basis Set (CBS) Extrapolations.

expo()[source]

CBS extrapolation formula by exponential Dunning-Feller relation.

\[E^{HF}(X) = E_{CBS} + a\cdot\exp(-bX)\]
Returns:

function object

exposqrt(twopoint=True)[source]

Three-point formula for extrapolating the HF reference energy [2].

\[E^{HF}(X) = E_{CBS} + a\cdot \exp(-b\sqrt{X})\]
Args:
twpointbool

A flag marking the use of two point extrapolation with b=9.0

Returns:

funtion object

exposum()[source]

Three point extrapolation through sum of exponentials expression

\[E(X) = E_{CBS} + a \cdot\exp(-(X-1)) + b\cdot\exp(-(X-1)^2)\]
extrapolate(x, energy, method, **kwargs)[source]

An interface for performing CBS extrapolations using various methods.

Args:
xnumpy.array

A vector of basis set cardinal numbers

energynumpy.array

A vector of corresponding energies

methodstr

Method/formula to use to perform the extrapolation

Kwargs:

Keyword arguments to be passed to the requested extrapolation function using the method argument

poly(p=0.0, z=3.0, twopoint=True)[source]

CBS extrapolation by polynomial relation.

\[E(X) = E_{CBS} + \sum_{i}a_{i}\cdot (X + P)^{-b_{i}}\]
Kwargs:
twpointbool

A flag for choosing the two point extrapolation

zfloat or list of floats

Order of the polynomial, default=3.0

pfloat

A parameter modifying the cardinal number, default=0.0

uste(method='CI')[source]

CBS extrapolation using uniform singlet and triplet pair extrapolation (USTE) scheme [Varandas2007].

[Varandas2007]

Varandas, A. J. C. (2007). “Extrapolating to the one-electron basis-set limit in electronic structure calculations. The Journal of Chemical Physics, 126(24), 244105. doi:10.1063/1.2741259

Args:
xint

Cardinal number of the basis set

e_cbsfloat

Approximation to the energy value at the CBS limit

afloat

Empirical A3 parameter

methodstr

One of: ci, cc

Returns:

function object

Molecule module

Module for handling atoms and molecules.

class Atom(identifier, xyz=(0.0, 0.0, 0.0), dummy=False, id=None)[source]

Basic atom class representing an atom.

move(x=0.0, y=0.0, z=0.0)[source]

Move atom to a set of new coordinates given in xyz

class Molecule(name='', atoms=None, unique=None, sym='', charge=0, multiplicity=1)[source]

Molecule class handling all the operations on single molecules

get_distance(atom1, atom2)[source]

Calcualte the distance between two atoms.

nele()[source]

Get the total number of electrons in a molecule.

unique()[source]

Get a list of unique atom specified by unique keyword

parsetools module

Module with convenience functions for parsing files

contains(string, query)[source]

Check if string contains query

getchunk(filename, startlno, endlno)[source]

Get a list of lines from a file between specified line numbers startlno and endlno.

Args:
filenamestr

Name of the file to process

startlnoint

Number of the first line to obtain

endlnoint

Number of the last line to obtain

Returns:
lineslist

A list of lines from the file filename between line numbers startlno and endlno

getlines(filename, tolocate)[source]

Return the lines from the files based on tolocate

Args:
filenamestr

Name of the file

tolocatelist of tuples

List of tuples with strings to find (queries) as first elements and integer offset values as second

Return:

locatelinenos(filename, tolocate)[source]

Given a file and a list of strings return a dict with string as keys and line numbers in which they appear a values.

Args:
filenamestr

Name of the file

tolocatelist of tuples

List of tuples with strings to find (queries) as first elements and integer offset values as second

Returns:
outdict

Dictionary whose keys are indices corresponding to item in input list and values are lists of line numbers in which those string appear

TODO:
  • add option to ignore the case of the strings to search

parsepairs(los, sep='=')[source]

Parse a given list of strings “los” into a dictionary based on separation by “sep” character and return the dictionary.

sliceafter(seq, item, num)[source]

Return “num” elements of a sequence “seq” present after the item “item”.

slicebetween(string, start, end)[source]

Return a slice of the string between phrases start and end.

take(seq, num)[source]

Iterate over a sequence seq num times and return the list of the elements iterated over.

submitgamess module

main()[source]

Script for submitting gamessus jobs to the queue.

remove_dat(path, datfile)[source]

Remove the dat file from the ASCII scratch.

submit_ll(args)[source]

Write the run script for LoadLeveller and submit it to the queue.

submit_pbs(args)[source]

Write the run script for PBS and submit it to the queue.

submit_slurm(args)[source]

Write the run script for SLURM and submit it to the queue.

submitmolpro module

main()[source]

Script for submitting molpro job to the queue.

set_defaults(args)[source]

Set some useful default values and add them to args

submit_pbs(args)[source]

Write the run script for PBS and submit it to the queue.

gamessorbitals module

Orbitals class

class Orbitals(*args: Any, **kwargs: Any)[source]

A convenience class for handling GAMESS(US) orbitals.

assign_lz_values(decimals=6, tolv=0.01)[source]

Determine the eigenvalues of the Lz operator for each nondegenerate or a combination of degenerate orbitals and assign them to lzvals column in the DataFrame

The Lz integrals over atomic orbitals are read from the dictionary file record No. 379.

Args:
decimalsint

Number of decimals to keep when comparing float eigenvalues

tolvfloat

Threshold for keeping wights of the eiegenvalues (squared eigenvector components)

WARNING:

currently this only works if the eigenvalues on input are sorted

fragment_populations()[source]

Calculate the fragment populations and assign each orbital to a fragment

classmethod from_files(name=None, logfile=None, dictfile=None, datfile=None)[source]

Initialize the Orbitals instance based on orbital information parsed from the logfile and read from the dictfile.

Args:
namestr

One of hf or ci

logfilestr

Name of the GAMESS(US) log file

dictfilestr

Name of the GAMESS(US) dictionary file .F10

datfilestr :

Name of the GAMESS(US) dat file

lzmos(ETOLLZ=1e-06, TOLV=0.01)[source]

A python rewrite of the GAMESS(US) LZMOS subroutine (from symorb.src) for analyzing the Lz composition of the orbtials

Args:

ETOLLZ : float TOLV : float

check_duplicates(a, decimals=6)[source]

This funciton assumes that the array a is sorted

http://stackoverflow.com/questions/25264798/checking-for-and-indexing-non-unique-duplicate-values-in-a-numpy-array http://stackoverflow.com/questions/5426908/find-unique-elements-of-floating-point-array-in-numpy-with-comparison-using-a-d

gamessreader module

GamessReader : reading gamess binary files.

class BinaryFile(filename, mode='r', order='fortran')[source]

Class representing a binary file (C or Fortran ordered).

file

The file handler.

order

The order for file (‘c’ or ‘fortran’).

read(dtype, shape=(1,))[source]

Read an array of dtype and shape from current position.

shape must be any tuple made of integers or even () for scalars.

The current position will be updated to point to the end of read data.

seek(offset, whence=0)[source]

Move to new file position.

Argument offset is a byte count. Optional argument whence defaults to 0 (offset from start of file, offset should be >= 0); other values are 1 (move relative to current position, positive or negative), and 2 (move relative to end of file, usually negative, although many platforms allow seeking beyond the end of a file). If the file is opened in text mode, only offsets returned by tell() are legal. Use of other offsets causes undefined behavior.

tell()[source]

Returns current file position, an integer (may be a long integer).

write(arr)[source]

Write an arr to current position.

The current position will be updated to point to the end of written data.

class DictionaryFile(filename, irecln=4090, int_size=8)[source]

Wrapper for reading GAMESS(US) dictionary file (*.F10).

read_record(nrec, dtype=None)[source]

Read a logical record ‘rec’ from the dictionary file and return a numpy array of type defined in the ‘records’ list, and size defined through ‘self.ifilen’ array.

class GamessFortranReader(log)[source]
Class for holding method for reading gamess binary files:

$JOB.F08 : two electron integrals over AO’s, $JOB.F09 : two electron integrals over MO’s, $JOB.F10 : the dictionary file with 1-e integrals, orbitals etc., $JOB.F15 : GUGA and ORMAS two-electron reduced density matrix,

get_onee_size(aos=True)[source]

Get the size of the vector holding upper (or lower) triangle of a square matrix of size naos or nmos.

get_twoe_size(aos=False)[source]

Get the size of the 1d vector holding upper (or lower) triangle of a supermatrix of size nmos (2RDM and two-electrons integrals).

read_rdm2(filename=None, nmo=None)[source]

Read the 2rdm from the gamess-us file

read_twoeao(filename=None, nmo=None)[source]

Read the two electron integrals from the gamess-us file

read_twoemo(filename=None, nmo=None)[source]

Read the two electron integrals from the gamess-us file

class GamessReader(log)[source]
Class for holding method for reading gamess binary files:
  • $JOB.F08 : two electron integrals over AO’s,

  • $JOB.F09 : two electron integrals over MO’s,

  • $JOB.F15 : GUGA and ORMAS two-electron reduced density matrix,

TODO:

CI coefficients, and CI hamiltonian matrix elements.

get_onee_size(aos=True)[source]

Get the size of the vector holding upper (or lower) triangle of a square matrix of size naos or nmos.

get_twoe_size()[source]

Get the size of the 1d vector holding upper (or lower) triangle of a supermatrix of size nmos (2RDM and two-electrons integrals).

read_rdm2(filename=None, nmo=None)[source]

Read the 2rdm from the gamess-us file

read_twoeao(filename=None, nmo=None)[source]

Read the two electron integrals from the gamess-us file

class SequentialFile(filename, logfile=None)[source]
get_index_buffsize(buff_size, int_size)[source]

Return the index buffer size for reading 2-electron integrals

ijkl(i, j, k, l)[source]

Based on the four orbital indices i,j,k,l return the address in the 1d vector.

read_ci_coeffs()[source]

Read CI coefficients from file NFT12 and return them as a numpy array of floats.

readseq(buff_size=15000, int_size=8, mos=False, skip_first=False)[source]

Read FORTRAN sequential unformatted file with two-electron quantities:

  • two electron integrals over AO’s: .F08 file

  • two electron integrals over MO’s: .F09 file

  • elements of the two particle density matrix: .F15 file

Args:
buff_sizeint

size of the buffer holding values to be read, in gamessus it is stored under NINTMX variable and in Linux version is equal to 15000 which is the default value

large_labelsbool

a flag indicating if large labels should were used, if largest label i (of the MO) is i<255 then large_labels should be False (case LABSIZ=1 in gamessus), otherwise set to True (case LABSIZ=2 in gamess(us),

skip_firstbool

skips the first record of the file is set to True,

Returns:

numpy 1D array holding the values

factor(i, j, k, l)[source]

Based on the orbitals indices return the factor that takes into account the index permutational symmetry.

ijkl(i, j, k, l)[source]

Based on the four orbital indices i, j, k, l return the address in the 1d vector.

print_twoe(twoe, nbf)[source]

Print the two-electron integrals.

rec

alias of record

tri2full(vector, sym=True)[source]

Convert a triagonal matrix whose elements are stored in the vector into a rectangular matrix of the shape given by shape tuple.