Derivatives of the energy (qcip_tools.derivatives)

Tools to help the representation of derivatives of the energy.

Theory

Supposing a function \(\mathbf{f}:\mathbb{R}^{n}\rightarrow \mathbb{R}^{m}\), a function which takes a vector \(\mathbf{x}\in\mathbb{R}^n\) as input and produce a vector \(\mathbf{f}(\mathbf{x})\in\mathbb{R}^m\), one can write its Maclaurin series expansion as

(1)\[\mathbf{M}(\mathbf{q}) = \sum_{r=0}^{r_{max}} \frac{1}{r!}\frac{\partial^r\mathbf{f}(0)}{\partial\mathbf{q}^r}\,\mathbf{q}^r,\]

where \(r_{max}=\infty\) if the series is not truncated.

In this expansion, one recognise, for example, \(\frac{\partial\mathbf{f}(0)}{\partial\mathbf{q}}\), which is the jacobian of the function (and a linear application \(\mathbb{R}^n\rightarrow\mathbb{R}^m\)) evaluated at \(\mathbf{q}=\mathbf{0}\), which correspond to a \(n\times m\) matrix \(\alpha\), for which the elements \(i\) and \(j\) are given by

\[\alpha_{ij} = \left.\frac{\partial \mathbf{f}_i}{\partial \mathbf{q}_j}\right|_{\mathbf{q}=\mathbf{0}}\]

where \(\mathbf{f}_i\) is the application that associate to \(\mathbf{x}\in\mathbb{R}^n\) the \(i^\text{th}\) coordinates of \(\mathbf{f}(\mathbf{x})\in\mathbb{R}^m\). By facility, it will be referred as the term “tensor of order 2”. For any \(r\) in (1), there is a corresponding tensor of order \(r+1\) (corresponding to a multilinear application to \(\mathbb{R}^m\) and wich can be flatten into a \(n^r\, m\) dimension vector).

Due to the Shwarz’s theorem (referred as “Kleinman symmetry” in the field of nonlinear optics), some component of those tensors are equals to each others.

Three kind of derivatives are considered for the moment:

  • Electrical, with respect to a static (F) or dynamic (D or d) electric field (though the first one is a special case of the second). The input space is \(\mathbf{f}\in\mathbb{R}^3\), and

    \[\begin{split}\begin{align} E(\mathbf{f}) = E_0 &+ \sum_i \mu_i\,\mathbf{f}_i+ \frac{1}{2!}\,\sum_{ij}\alpha_{ij}\mathbf{f}_i\mathbf{f}_j + \frac{1}{3!}\sum_{ijk}\beta_{ijk}\mathbf{f}_i\mathbf{f}_j\mathbf{f}_k \\ &+ \frac{1}{4!}\,\sum_{ijkl}\gamma_{ijkl}\mathbf{f}_i\mathbf{f}_j\mathbf{f}_k\mathbf{f}_l + \ldots \end{align}\end{split}\]

    where \(\mu\) is the electric dipole moment, \(\alpha\) is the polarizability, \(\beta\) and \(\gamma\) are the first and second hyperpolarizabilites. Note that in the case of dynamic electric fields, \(\mathbf{f}_i(t)=\mathbf{f}^0_i\,(e^{-i\omega t}+c.c.)\), where \(\omega\) is the pulsation (\(\omega=2\pi\nu\)).

  • Geometrical, with respect to cartesian (G) or normal mode (N) coordinates. The input space is \(\mathbf{q}\in\mathbf{R}^{3N}\) where \(3N\) is the number of degrees of freedom of \(N\) atoms expressed in cartesian coordinates, and

    \[V(\mathbf{q}) = V_0 + \sum_i F_i\mathbf{q}_i + \frac{1}{2!}\,\sum_{ij}H_{ij}\mathbf{q}_i\mathbf{q}_j + \frac{1}{3!}\,\sum_{ij}F_{ijk}\mathbf{q}_i\mathbf{q}_j\mathbf{q}_k+ \ldots\]

    where \(F_i=\left.\frac{\partial V}{\partial\mathbf{q}_i}\right|_{\mathbf{q}=0}\) is the gradient (so the forces, wich are close to 0 if the system is in equilibrium), \(H_{ij}=\left.\frac{\partial^2 V}{\partial\mathbf{q}_i\mathbf{q}_j}\right|_{\mathbf{q}=0}\) is the hessian (force constant matrix) and \(F_{ijk}\) is the cubic force constants matrix (and so on).

  • Excitation energies: not exactly derivatives themselves, but obtained by single/double/… residue on the response functions. Excited state are represented by #, while excitation are represented by !. Therefore, #F means “excited state dipole moment” while !F means “(ground to excited) transition dipole moment”, and #!F means “excited to excited transition dipole moment”.

All derivatives are written with respect to the energy (assuming a corresponding Maclaurrin expansion). For example, GG corresponds to the hessian expressed in cartesian coordinates (second order derivative of the energy), while F is the electric dipole moment. Therefore, FG (or GF) is the first-order derivative of the dipole moment with respect to the cartesian coordinates.

Note

The name of a derivative is expressed as its electrical-derivative counterpart with respect to the geometrical parameters, though the opposite is valid (GF is also the derivative of the gradient with respect to electric field). The first is just more used in the literature (for IR and Raman spectra, for example).

Implementation

The object Derivative handle the manipulation of derivative on a symbolic level : you can define a derivative and differenciate it. It computes the shape, dimension and order of the corresponding derivative. Note that you must provide spacial_dof, the number of spacial degrees of freeedom of the corresponding system if the derivative contains a geometrical one (G or N), and nstate if the derivative contains an excitation (! or #).

The object Tensor allows to store this derivative in the components array (but without the specific operation that you can perform on certain tensors, e.g. the isotropic value for the polarizability and so all). It is provided with a Derivative (self.representation) for efficiency, so the components array is shaped accordingly. You must also provide frequency if the derivatives contains a dynamic electric field (D or d).

Warning

Always use the representation() function to get the actual representation of the derivatives. The representation of the derivative is ordered so that geometrical derivatives come before electric one. Therefore, the content of the components array in Tensor must be given with the indices in that order.

API documentation

qcip_tools.derivatives.ALLOWED_DERIVATIVES = ('!', '#', 'G', 'N', 'F', 'D', 'd', 'X')

In term of derivative of the energy G = geometrical derivative, N = normal mode derivative, F = static electric field derivative, D = dynamic electric field derivative (which can be static), d = inverse of the dynamic electric field (-w). X = anything but -w and w ! = excited state (for example !F stands for “excited state dipole moment”) # = transition (for example, #F stands for “transition dipole moment”) Note: do not change the order, except if you know what your are doing

qcip_tools.derivatives.COORDINATES = {0: 'x', 1: 'y', 2: 'z'}

spacial 3D coordinates

class qcip_tools.derivatives.Derivative(from_representation=None, basis=None, spacial_dof=None, nstates=None)

Represent a quantity, which is derivable

Parameters:
  • from_representation (str) – representation of the derivative

  • basis (Derivative) – basis for the representation, if any

  • spacial_dof (int) – spacial degrees of freedom (3N)

  • nstates (int) – number of excited states

classmethod apply_permutation(iterable, permutations)

Apply the permutations (list of tuple (from, to)) over iterable

Parameters:
  • iterable (list) – an iterable that you can copy (a list)

  • permutations (list|iterator) – list of permutation

Returns:

new iterable

Return type:

list

components_to_flatten_component(components)

Given a component, give the index if the array would be flatten

Parameters:

components (list|tuple) – components

Return type:

int

classmethod correct_components(ideal_representation, representation, list_of_components)

Correct the order of the components

Parameters:
  • ideal_representation (str) – representation when the different type of derivatives nicely follows each others

  • representation (str) – real representation

  • list_of_components (list) – list of components to correct

differentiate(derivatives_representation, spacial_dof=None, nstate=None)

Create a new derivative from the differentiation of of the current one.

Parameters:
  • derivatives_representation (str) – the representation of the derivatives

  • spacial_dof (int) – spacial degrees of freedom

Returns:

a new derivative

Return type:

Derivative

dimension(raw=False)

Return the dimension of the (full) flatten tensor

Parameters:

raw (bool) – return the shape of the raw representation

Returns:

the size

Return type:

int

classmethod expend_list(iterable, sub_iterable)

For each element of iterable, create a new list by adding each element of sub_iterable.

Parameters:
  • iterable (list) – main list

  • sub_iterable (list) – sub list

Return type:

list

flatten_component_to_components(i)

From the index if the array would be flatten, gives the components

Parameters:

i (int) – flatten index

Return type:

list

inverse_smart_iterator(element, as_flatten=False)

Back-iterate over all the other components : from a coordinates, give all the other ones that are equivalents

Parameters:
  • element (int|tuple) – the coordinates, either a flatten index (if as_flatten=True) or a tuple

  • as_flatten (bool) – yield full components, not flatten ones

order()

Get the order of derivation with respect to energy

Return type:

int

raw_representation(exclude=None)

Get raw representation (simply the parent representation + obj representation).

Warning

Dangerous to use, prefer representation().

Parameters:

exclude (tuple|list|str) – exclude some element from the representation

Return type:

str

representation()

Get the full representation (mix basis and current)

Note

First come the excitations, then the geometrical derivatives (G or N), then the electrical ones (F, D or d).

Return type:

str

shape(raw=False)

Return the shape of the (full) tensor

Parameters:

raw (bool) – return the shape of the raw representation

Returns:

the shape

Return type:

list

smart_iterator(as_flatten=False)

Apply the Shwarz’s theorem and only return as subset of independant coordinates. Order normally guaranteed.

Note

Ideally, the different type derivatives follows each other. This is not always the case for electrical ones, so there is a stage of re-ordering.

Parameters:

as_flatten (bool) – yield full components, not flatten ones

qcip_tools.derivatives.ORDERS = {1: 'first', 2: 'second', 3: 'third', 4: 'fourth', 5: 'fifth'}

number to x*th*.

exception qcip_tools.derivatives.RepresentationError
class qcip_tools.derivatives.Tensor(representation, components=None, spacial_dof=None, frequency=None, nstates=None, name='')

Create a tensor to store a given derivative.

Parameters:
  • representation (str|Derivative) – representation of the derivative

  • spacial_dof (int) – number of spacial degrees of freedom

  • frequency (str|float) – frequency if dynamic electric field

  • nstates (int) – number of excited states

project_over_normal_modes(mwh)

Project the tensor over normal modes (get its normal mode equivalent)

Parameters:

mwh (qcip_tools.derivatives_g.MassWeightedHessian) – mass weighted hessian

Return type:

qcip_tools.derivatives.Tensor

to_string(threshold=1e-05, columns_per_line=6, molecule=None, skip_trans_plus_rot_dof=0, **kwargs)

Print the tensor in a more or less textual version.

Parameters:
  • threshold (float) – show 0 instead of the value if lower than threshold

  • columns_per_line (int) – number of columns per “line” of the tensor

  • molecule (qcip_tools.molecule.Molecule) – use molecule to gives the atom instead of Gxx

  • skip_trans_plus_rot_dof (int) – skip the n first N (normal) modes

Returns:

representation

Return type:

str

qcip_tools.derivatives.compute_numerical_derivative_of_tensor(basis, derivative_repr, tensor_func, k_max, min_field, ratio, frequency=None, dry_run=False, accuracy_level=1, force_choice=None, **kwargs)

Compute numerical derivative of tensor (taking advantage of the so-called smart iterator).

Note

The parameter accuracy_level is provided for retro-compatibility, but should remain to “1”.

Parameters:
  • basis (qcip_tools.derivatives.Derivative) – basis of differentiation (representation for the base tensors)

  • derivative_repr (qcip_tools.derivatives.Derivative) – representation of the further derivatives of the tensor

  • tensor_func – function to access to the tensor

  • dry_run (bool) – do not fill the tensor or perform romberg analysis

  • k_max (int) – Romberg triangle k maximum

  • min_field (float) – minimal value

  • ratio (float) – ratio (a)

  • frequency (float|str) – frequency if electrical derivative

  • accuracy_level – level of accuracy. Default value (1) does not change the behavior, but “0” and lower use forward derivatives for order of differentiation == 3 (retro-compatibility behavior)

  • force_choice (tuple) – force the choice in the Romberg triangle

  • kwargs (dict) – args passed to tensor_func

Returns:

tensor and romberg triangles

Return type:

qcip_tools.derivatives.Tensor, collection.OrderedDict

qcip_tools.derivatives.is_electrical(derivative)

Return if the derivatives contains an electrical one

Return type:

bool

qcip_tools.derivatives.is_excitation(derivative)

Return if the derivatives contains an excitation

Return type:

bool

qcip_tools.derivatives.is_geometrical(derivative)

Return if the derivatives contains a geometrical one

Return type:

bool

qcip_tools.derivatives.representation_to_operator(representation, component=None, molecule=None)

Use the representation to find a translation for the operator

Parameters:
  • representation (str) – the representation

  • component (int) – (flatten) component

  • molecule (qcip_tools.molecule.Molecule) – the molecule (print the symbol of the atoms instead of a number if given)

Return type:

str