Public API: richmol package
Subpackages
Unit conversion (convert_units
)
We use scipy.constants for the reference values of all physcial constants
|
Converts MHz to \(cm^{-1}\) |
|
Converts Joule to \(cm^{-1}\) |
|
Converts dipole moment form \(Debye\) to international system of units |
|
Converts dipole moment from \(Debye\) to atomic units |
|
Converts dipole moment from \(Debye\) to \(erg^{1/2}*cm^{3/2}\) |
|
Converts product of dipole moment (in \(Debye\)) with field (in \(Volts/meter\)) to \(cm^{-1}\) |
|
Converts product of dipole moment (in \(Debye\)) with field (in \(Volts/meter\)) to \(MHz\) |
|
Converts quadrupole moment form \(Buckingham\) to international system of units |
|
Converts quadrupole moment from \(Buckingham\) to atomic units |
|
Converts quadrupole moment from \(Buckinghom\) to \(erg^{1/2}*cm^{5/2}\) |
|
Converts product of dipole moment (in atomic units) with field (in \(Volts/meter\)) to \(cm^{-1}\) |
|
Converts product of polarizability (in atomic units) with field (in \(Volts^2/meter^2\)) to \(cm^{-1}\) |
Molecule-field interaction potential (field
)
|
General class for laboratory-frame Cartesian tensor operator |
|
Applies state selection filters to input Cartesian tensor |
- class richmol.field.CarTens(filename=None, name=None, **kwargs)
General class for laboratory-frame Cartesian tensor operator
- Parameters:
filename – str Name of the HDF5 file from which tensor data is loaded.
name – str Name of the data group in HDF5 file.
- Kwargs:
- threshfloat
Threshold for neglecting matrix elements when reading from file
- bra, ketfunction(**kw)
State filters for bra and ket basis sets, take as arguments state quantum numbers, symmetry, and energy, i.e., J, m, k, sym, and enr, and return True or False depending on if the corresponding state needs to be included or excluded form the basis. By default, all states stored in the file are included. The following keyword arguments are passed into the bra and ket functions:
- Jfloat (round to first decimal)
Value of J (or F) quantum number
- symstr
State symmetry
- enrfloat
State energy
- mstr
State assignment in the M subspace, usually just a value of the m quantum number
- kstr
State assignment in the K subspace, which are the rotational or ro-vibrational quanta joined in a string
- Attrs:
- rankint
Rank of tensor operator.
- cartlist of str
Contains string labels of tensor Cartesian components (e.g, ‘x’, ‘y’, ‘z’, ‘xx’, ‘xy’, ‘xz’, …)
- os[ (omega, sigma) for omega in range(nirrep)
for sigma in range(-omega, omega + 1) ]
List of spherical-tensor indices (omega, sigma), here nirrep is the number of tensor irreducible representations.
- Jlist1, Jlist2list
Lists of J quanta, spanned by the bra and ket basis sets, respectively.
- symlist1, symlist2dict
List of symmetries, for each J, spanned by the bra and ket basis sets, respectively. Example:
[sym for J in self.Jlist1 for sym in symlist1[J]]
- dim1, dim2nested dict
Matrix dimensions of tensor for different J and symmetry, spanned by the bra and ket basis sets, respectively. Example:
[ dim for J in self.Jlist1 for sym in symlist1[J] for dim in dim1[J][sym] ]
- dim_m1, dim_k1nested dict
Dimensions of M and K tensors for different J and symmetry, spanned by the bra basis. Example:
[ dim for J in self.Jlist1 for sym in symlist1[J] for dim in dim_m1[J][sym] ]
- dim_m2, dim_k2nested dict
Same as dim_m1 and dim_m2 but for the ket basis.
- quanta_m1, quanta_k1nested dict
M and ro-vibrational quantum numbers, respectively, for different J and symmetry, spanned by the bra basis set. The elements of quanta_m1[J][sym] list represent the m quantum number, while the elements of quanta_k1[J][sym] list are tuples (q, enr), where q is the string of ro-vibrational quantum numbers and enr is the ro-vibrational energy (or None). Example:
[ int(m) for J in self.Jlist1 for sym in symlist1[J] for m in quanta_m1[J][sym] ] [ (k, energy) for J in self.Jlist1 for sym in symlist1[J] for m in quanta_k1[J][sym] ]
- quanta_m2, quanta_k2nested dict
Same as quanta_m1 and quanta_k1 but for the ket basis.
- kmatnested dict
K-tensor matrix elements (in CSR format) for different pairs of bra and ket J quanta, symmetries, and irreducible spherical-tensor components. Example:
for (J1, J2), kmat_J in kmat.items(): for (sym1, sym2), kmat_sym in kmat_J.items(): for irrep, kmat_irrep in kmat_sym.items(): # K-subspace matrix elements print(type(kmat_irrep)) # scipy.sparse.spmatrix
- mmatnested dict
M-tensor matrix elements (in CSR format) for different pairs of bra and ket J quanta, symmetries, irreducible spherical-tensor and Cartesian components. Example:
for (J1, J2), mmat_J in mmat.items() for (sym1, sym2), mmat_sym in mmat_J.items(): for irrep, mmat_irrep in mmat_sym.items(): for cart, mmat_cart in mmat_irrep.items(): # M-subspace matrix elements print(type(mmat_cart)) # scipy.sparse.spmatrix
- mfmatnested dict
M-tensor matrix elements contracted with field. Produced after multiplication of tensor with a vector of X, Y, and Z field values (see
field()
). Has the same structure askmat
.- symtop_basis: nested dict
Symmetric-top basis representation of wave functions. Present if Cartesian tensor represents field-free solutions.
- eigvecnested dict
Eigenvectors (in dense matrix format) for different J quanta and different symmetries.
NOTE: state selection filters (
filter()
) are not implemented, and thus have no effect, for eigvec attribute.Example:
for J in eigvec.items(): for sym in eigvec[J].items(): print(type(eigvec[J][sym]))
- filter(thresh=None, bra=<function CarTens.<lambda>>, ket=<function CarTens.<lambda>>)
Applies state selection filters to tensor matrix elements
- tomat(form='block', repres='csr_matrix', thresh=None, cart=None)
Returns full matrix representation of tensor
- Parameters:
form – str For form = ‘block’, the matrix representation is build as dictionary containing matrix blocks for different pairs of J and symmetries. For form = ‘full’, full 2D matrix is constructed.
repres – str Defines representation for matrix blocks or full 2D matrix. Can be set to the name of one of
scipy.sparse
matrix classes, e.g., repres = ‘coo_matrix’. Alternatively it can be set to ‘dense’.thresh – float Threshold for neglecting matrix elements when converting into the sparse form.
cart – str Desired Cartesian component of tensor, e.g., cart = ‘z’ or cart = ‘xx’. If set to None (default), the function will attempt to return a matrix representation of the corresponding potential (i.e., product of tensor and field)
- Returns:
- nested dict or 2D array
For form = ‘block’, returns dictionary containing matrix blocks for different pairs of J and symmetries. For form = ‘full’, returns 2D matrix.
- assign(form='block')
Returns assignments of bra and ket basis states
- Parameters:
form – str Form of the assignment output (see form argument of
CarTens.tomat()
function)- Returns:
- dict
Assignments of bra and ket states, respectively. For form = ‘block’, assign[J][sym][‘m’] and assign[J][sym][‘k’] contain list of m and ro-vibrational quantum numbers, respectively, for states with given values of J quantum number and symmetry sym. For form = ‘full’, assign[‘m’], assign[‘k’], assign[‘sym’], and assign[‘J’] contain list of m quanta, ro-vibrational quanta, symmetries, and J values for all states in the basis. The ordering of elements in assign1 and assign2 lists corresponds to the ordering of rows and columns in a matrix returned by
CarTens.tomat()
function.
- Return type:
assign1, assign2
- full_form(mat, repres='csr_matrix', thresh=None)
Converts block representation of tensor matrix into 2D matrix form
- Parameters:
mat – nested dict Block representation of matrix for different values of bra and ket J quanta and different symmetries
repres – str Set to the name of one of
scipy.sparse
matrix classes to output the matrix in sparse format, e.g., sparse = ‘coo_matrix’ or ‘csr_matrix’. Alternatively set to ‘dense’.thresh – float Threshold for neglecting matrix elements when converting into sparse form.
- Returns:
- array
2D matrix representation.
- block_form(mat)
Converts 2D tensor matrix into a block form
- Parameters:
mat – array (numpy.ndarray or scipy.sparse.spmatrix) 2D tensor matrix.
- Returns:
- nested dict
Block form of tensor matrix, split for different values of bra and ket J quanta and different symmetries.
- mul(arg)
In-place multiplication of tensor with a scalar arg
- field(field, thresh=None)
In-place multiplication of tensor with field
The result is stored in the attribute
mfmat
.- Parameters:
field –
- array (3)
Contains field’s X, Y, Z components
- threshfloat
Threshold for neglecting field’s components and their products, as well as the elements of the resulting product tensor.
- vec(vec, matvec_lib='scipy')
Computes product of tensor with vector
- Parameters:
vec –
nested dict Dictionary containing vector elements for different J and symmetries.
for J in vec.keys(): for sym in vec[J].keys(): print( vec[J][sym].shape == self.dim2[J][sym] ) # must be True
matvec_lib – str The python library to use for SpMV. By default, this is set to ‘scipy’. Alternatively this can be set to ‘numba’ or ‘cupy’.
- Returns:
- nested dict
Resulting vector, has same structure as input vec.
- store(filename, name=None, comment=None, replace=False, replace_k=False, replace_m=False, store_k=True, store_m=True, thresh=None)
Stores object into HDF5 file
- Parameters:
filename – str Name of HDF5 file.
name – str Name of the data group, by default, name of the variable will be used.
comment – str User comment.
replace – bool If True, the existing in file complete tensor data group will be replaced.
replace_k – bool If True, the existing in file K-tensor data sets will be replaced.
replace_m – bool If True, the existing in file M-tensor data sets will be replaced.
store_k – bool If False, writing K-tensor into file will be skipped.
store_m – bool If False, writing M-tensor into file will be skipped.
thresh – float Threshold for neglecting matrix elements (M and K tensors) when writing into file.
- read(filename, name=None, thresh=None, **kwargs)
Reads object from HDF5 file
- Parameters:
filename – str Name of HDF5 file.
name – str Name of the data group, by default the name of the variable will be used.
thresh – float Threshold for neglecting matrix elements when reading from file.
- richmol.field.filter(obj, bra, ket, thresh=None)
- Applies state selection filters to input Cartesian tensor
Rotational density (rotdens
)
|
Computes rotational densities |Ψ*Ψ'| on a grid of Euler angles (α, β, γ). |
|
Computes wave functions on a grid of Euler angles (α, β, γ). |
- richmol.rotdens.dens_grid(h: CarTens, alpha: ndarray, beta: ndarray, gamma: ndarray, diag_only: bool = False, form: str = 'block') Union[Dict[Any, Any], ndarray]
Computes rotational densities |Ψ*Ψ'| on a grid of Euler angles (α, β, γ).
- Parameters:
h (CarTens) – Cartesian tensor representing rotational or rovibrational field-free solutions.
alpha (np.ndarray) – Values of the Euler α angle in radians.
beta (np.ndarray) – Values of the Euler β angle in radians.
gamma (np.ndarray) – Values of the Euler γ angle in radians.
form (str, optional) –
Determines the format of the output. - ‘block’: Returns a dictionary with densities grouped by J and symmetry. - ‘full’: Returns densities concatenated across different J and symmetry.
Default is ‘block’.
diag_only (bool, optional) – If set to True, only diagonal densities |Ψ*Ψ| are calculated. If False, all diagonal and non-diagonal elements |Ψ*Ψ'| are calculated. Default is False.
- Returns:
Densitties on the grid of Euler angles. If form is ‘block’, it returns a dictionary. For ‘full’, it returns a 5D array if diag_only is False, and a 4D array if diag_only is True. See examples below.
Example for ‘block’:
dens = dens_grid(h, alpha, beta, gamma, form='block') assign1, assign2 = h.assign(form='block') # bra and ket state assignments for (J1, J2), dens_J in dens.items(): for (sym1, sym2), dens_sym in dens_J.items(): assgn1 = assign1[J1][sym1] assgn2 = assign2[J2][sym2] for istate in range(dens_sym.shape[-2]): for jstate in range(dens_sym.shape[-1]): print(J1, sym1, istate, assgn1['m'][istate], assgn1['k'][istate], " | ", J2, sym2, jstate, assgn2['m'][jstate], assgn2['k'][jstate], " ", dens_sym[ia, ib, ig, istate, jstate]) # where `ia`, `ib`, and `ig` are indices of Euler angles α, β, γ # for diagonal-only contributions dens = dens_grid(h, alpha, beta, gamma, form='block', diag_only=True) assign1, assign2 = h.assign(form='block') for (J1, J2), dens_J in dens.items(): for (sym1, sym2), dens_sym in dens_J.items(): assgn1 = assign1[J1][sym1] assgn2 = assign2[J2][sym2] for istate in range(dens_sym.shape[-1]): print(J1, sym1, istate, assgn1['m'][istate], assgn1['k'][istate], " | ", J2, sym2, istate, assgn2['m'][istate], assgn2['k'][istate], " ", dens_sym[ia, ib, ig, istate]) # where `ia`, `ib`, and `ig` are indices of Euler angles α, β, γ
Example for ‘full’:
assign, _ = h.assign(form='full') J, sym, m, k = [assign[key] for key in ('J', 'sym', 'm', 'k')] dens = dens_grid(h, alpha, beta, gamma, form='full') for istate in range(dens.shape[-2]): for jstate in range(dens.shape[-1]): print(istate, J[istate], sym[istate], m[istate], k[istate], " | ", jstate, J[jstate], sym[jstate], m[jstate], k[jstate], " ", dens[ia, ib, ig, istate, jstate]) # where `ia`, `ib`, and `ig` are indices of Euler angles α, β, γ # for diagonal-only contributions dens = dens_grid(h, alpha, beta, gamma, form='full', diag_only=True) for istate in range(dens.shape[-1]): print(istate, J[istate], sym[istate], m[istate], k[istate], dens[ia, ib, ig, istate]) # where `ia`, `ib`, and `ig` are indices of Euler angles α, β, γ
- Return type:
Union[Dict[Any, Any], np.ndarray]
- richmol.rotdens.psi_grid(h: CarTens, alpha: ndarray, beta: ndarray, gamma: ndarray, form: str = 'block', sum_over_vib: bool = True) Union[Dict[Any, Any], ndarray]
Computes wave functions on a grid of Euler angles (α, β, γ).
- Parameters:
h (CarTens) – Cartesian tensor representing rotational or rovibrational field-free solutions.
alpha (np.ndarray) – Values of the Euler α angle in radians.
beta (np.ndarray) – Values of the Euler β angle in radians.
gamma (np.ndarray) – Values of the Euler γ angle in radians.
form (str, optional) –
Determines the format of the output. - ‘block’: Returns a dictionary with wave functions grouped by J and symmetry. - ‘full’: Returns wave functions concatenated across different J and symmetry.
Default is ‘block’.
sum_over_vib (bool, optional) – If set to True, the returned function sums over basis set contributions with different vibrational quanta. If False, these contributions are kept separately in the leading dimension of the output array. Default is True.
- Returns:
Wave functions on the grid of Euler angles. If form is ‘block’, it returns a dictionary. For ‘full’, it returns a 5D array, if sum_over_vib is False, and a 4D array, if sum_over_vib is True. See examples below.
Example for ‘block’:
psi = psi_grid(h, alpha, beta, gamma, form='block', sum_over_vib=False) assign, _ = h.assign(form='block') # state assignment for J, psi_J in psi.items(): for sym, (psi_sym, vib_ind) in psi_J.items(): assgn = assign[J][sym] for istate in range(psi_sym.shape[-1]): print(J, sym, istate, assgn['m'][istate], assgn['k'][istate], psi_sym[iv, ia, ib, ig, istate]) # where `ia`, `ib`, and `ig` are indices of Euler angles # α, β, γ, and `iv` is vibrational basis index
Example for ‘full’:
psi = psi_grid(h, alpha, beta, gamma, form='full') assign, _ = h.assign(form='full') # state assignment J, sym, m, k = [assign[key] for key in ('J', 'sym', 'm', 'k')] for istate in range(psi.shape[-1]): print(istate, J[istate], sym[istate], m[istate], k[istate], psi[ia, ib, ig, istate]) # where `ia`, `ib`, and `ig` are indices of Euler angles α, β, γ
- Return type:
Union[Dict[Any, Any], np.ndarray]
TROVE bindings (trove
)
|
To initialise Cartesian tensor operator from TROVE output files. |
- class richmol.trove.CarTensTrove(filename=None, matelem=None, coefs=None, **kwargs)
To initialise Cartesian tensor operator from TROVE output files.
TROVE format includes ‘states file’, containing information about rovibrational energies and assignments, and ‘matrix elements files’, containing rovibrational matrix elements of an operator. Matrix elements for different pairs of angular momentum J (or F) quanta are stored in different files.
This is a subclass of
richmol.field.CarTens
class.- Parameters:
filename – str Name of states file.
matelem – str In TROVE output format, matrix elements of Cartesian tensors for different values of bra and ket J (or F) quanta are stored in separate files. The argument matelem provides a template for generating the names of these files. For example, for matelem = ‘me<j1>_<j2>’, the following files will be searched: ‘me0_0’, ‘me0_1”, “me1_0”, ‘me2_0’, ‘me2_1’, and so on, i.e. <j1> and <j2> will be replaced by integer values of bra and ket J quanta, respectively. For half-integer values of the F quantum number, replace <j1> and <j2> in the template by <f1> and <f2>, these will then be substituted by the floating point values of bra and ket F quanta rounded to the first decimal.
coefs – str Name of file containing rovibrational coefficients, used to for plotting rotational densities and wavepackets.
- read_states(filename, coef_file=None, **kwargs)
Reads states file (rovibrational state information) and optionally file with rovibrational coefficients, which is used for estimations of rotational density function (mainly for plotting)
- Parameters:
filename – str Name of states file.
coef_file – str Name of coefficients file. This file is used to estimate rotational probability density functions.
- read_trans(filename, thresh=None, **kwargs)
Reads matrix elements of Cartesian tensor
NOTE: call
read_states()
beforeread_trans()
to load rovibrational states information stored in states file.- Parameters:
filename – str Matrix elements of Cartesian tensors for different values of bra and ket J (or F) quanta are stored in separate files. The parameter filename provides a template for generating the names of these files. For example, for filename = ‘me<j1>_<j2>’, the following files will be searched: ‘me0_0’, ‘me0_1’, ‘me1_0’, ‘me2_0’, ‘me2_1’, and so on, i.e. <j1> and <j2> will be replaced by integer values of bra and ket J quanta, respectively. For half-integer values of the F quantum number, replace <j1> and <j2> in the template by <f1> and <f2>, these will then be substituted by the floating point values of bra and ket F quanta rounded to the first decimal.
thresh – float Threshold for neglecting matrix elements when reading from file.
Time-dependent solutions (tdse
)
|
Class for time-propagator |
- class richmol.tdse.TDSE(**kwargs)
Class for time-propagator
- Attrs:
- _t_startint, float
Starting time of time-propagation
- _t_endint, float
Terminal time of time-propagation
- _dtint, float
Time-step of propagation
- _t_to_sfloat
Factor to convert time into units of second
- _enr_to_Jfloat
Factor to convert energy into units of Joule
- _time_gridtuple()
Lower borders, upper borders and centers of time-intervals to propagate over
- _exp_fac_H0numpy.ndarray
Matrix exponential of the field-free part used during split operator approach
- __init__(**kwargs)
Initializes TDSE object
- time_grid(grid_type='equidistant')
Generates time-grid to propagate on
- init_state(H, **kwargs)
Generates initial state vectors
- update(H, vecs, **kwargs)
Propagates vectors by one time-step
- time_grid(grid_type='equidistant')
Generates time-grid to propagate on
- Parameters:
grid_type – str type of time grid
- Returns:
- numpy.ndarray
times at which to evaluate Hamiltonian
- Return type:
t_c
- init_state(H, **kwargs)
Generates initial state vectors
- Initial state vectors are eigenfunctions of Hamiltonian H. If
temp is None, all eigenfunctions are returned unweighted. If temp is 0, only lowest energy eigenfunction is returned. If temp is > 0, all eigenfunctions are returned Boltzmann- -weighted.
- Parameters:
H –
field.CarTens
Hamiltonian operator
- Kwargs:
- tempfloat
Temperature in Kelvin
- threshfloat
Collective threshold below which to neglect higher energy states
- zpefloat
Zero-point energy, by default lowest eigenvalue of H is taken
- Returns:
- numpy.ndarray
Initial state vectors (each row represents an initial state)
- Return type:
vecs
Spectra (spectrum
)
|
Class for field-free spectrum |
- class richmol.spectrum.FieldFreeSpec(filename, matelem=None, names=None, **kwargs)
Class for field-free spectrum
- Attrs:
- filenamestr
Name of the HDF5 file from which tensor data is loaded. Alternatively, one can load tensor from the old-format ASCII files, by providing in filename the name of the richmol states file and in matelem a template for generating the names of the richmol matrix elements files.
- matelemstr
In the old-format, matrix elements of Cartesian tensors for different values of bra and ket J (or F) quanta are stored in separate files.
- names: list
Names (str) of two data groups in HDF5 file. First name represents group name of field-free data, while second name represents group name of matrix element data.
- out_fstr
Name of the HDF5 file to which results are written.
- typestr
The type of coupling moment (‘elec’, ‘magn’). Set to ‘elec’ by default.
- orderstr
The order of coupling moment. If matelem is given, order will be deduced from matelem. Else must be specified here (‘dip’, ‘quad’).
- j_minint, float
The minimum value of J (or F). Set to 0.0 (for integer) or 0.5 (for half-integer) by default.
- j_maxint, float
The maximum value of J (or F). Set to 1.0 (for ‘dip’) or 2.0 (for ‘quad’) for integer by default. Set to 0.5 (for ‘dip’) or 1.5 (for ‘quad’) for half-integer by default.
- e_maxint, float
The maximum value of energy. Set to 1.0e9 by default.
- unitsstr
The units of coupling moment. Set to ‘Debye’ by default.
- linestr_threshint, float
The linestrength threshold below which to neglect transitions. Set to 0 by default.
- abunint, float
The natural terrestrial isotopic abundance. Set to 1.0 by default.
- tempint, float
The temperature. Set to 296.0 by default.
- part_sumint, float
The total internal partition sum. Set to 1.0 by default.
- filters: list
User defined filtering functions for transitions. The folowing keyword arguments are passed into the filter functions:
- symtuple
Symmetries of bra and ket state.
- qstrtuple
Qstrs of bra and ket state.
Set to [def all_(**kwargs): return True]
- abs_intens_threshint, float
The absorption intensity threshold below which to neglect transitions. Set to 0 by default.
- __init__(filename, matelem=None, **kwargs)
Initializes field-free spectrum object, computes assignments.
- linestr(**kwargs)
Computes linestrengths.
- abs_intens(**kwargs)
Computes, filters and plots absorption intensities.
- totxt(**kwargs)
Writes results into ASCII file.
- linestr(**kwargs)
Wrapper for computing assignments and linestrengths.
- Kwargs:
- threshfloat
The linestrength threshold below which to neglect transitions. Set to 0 by default.
- abs_intens(temp, part_sum, **kwargs)
- Wrapper for computing, filtering and plotting absorption
intensities.
- Parameters:
temp – int, float The temperature.
part_sum – int, float The total internal partition sum.
- Kwargs:
- abunfloat
Natural terrestrail isotopic abundance Set to 1.0 by default
- filterslist
User defined filtering functions for transitions (see filters in kwargs of
field_free_spectrum
).- threshint, float
The absorption intensity threshold below which to neglect transitions. Set to 0 by default.
- totxt(form='default')
Writes assignments and spectrum into ASCII files.
- Parameters:
form – str The file format of the ASCII file. Set to ‘default’.