richmol.rot package
|
Describes rigid molecule |
|
Solves rotational eigenvalue problem |
|
Represents rotational solutions for different values of J and symmetry |
|
Declares new molecular-frame Cartesian tensor which is dynamically rotated to a chosen molecular frame |
|
Represents rotational matrix elements of laboratory-frame Cartesian tensor operator |
Rigid molecule description
- class richmol.rot.Molecule
Describes rigid molecule
- XYZ()
Cartesian coordinates and masses of atoms
- Parameters:
arg –
Atomic labels, Cartesian coordinates of atoms, and units in a tuple, for example,
water.XYZ = ("bohr", "O", 0.00000000, 0.00000000, 0.12395915, "H", 0.00000000, -1.43102686, -0.98366080, "H", 0.00000000, 1.43102686, -0.98366080)
To define a non-standard isotoplogue, add the corresponding number next to the atom label, e.g., “O18”, “H2”. Cartesian coordinates can also be read from an XYZ file, in this case arg should contain the name of XYZ file.
- Returns:
- structured numpy array
XYZ[‘xyz’] - Cartesian coordinates of atoms in Angstrom, XYZ[‘mass’] - atomic masses, XYZ[‘label’] - atomic labels
- frame()
Molecular frame embedding
- Parameters:
arg –
String defining molecular frame
water.frame = "diag(inertia)" # rotates to a frame where the inertia tensor # becomes diagonal water.frame = "zxy" # permutes the x, y, and z axes water.pol = [[9.1369, 0, 0], [0, 9.8701, 0], [0, 0, 9.4486]] water.frame = "pol" # rotates frame with water.pol matrix water.frame = "diag(pol)" # rotates to a frame where water.pol tensor # becomes diagonal water.frame = "None" # or None, resets frame to the one defined # by the input molecular geometry
- Returns:
- array (3,3)
Frame rotation matrix
- B()
Molecular rotational constants
- Parameters:
arg – Tuple with Bx, By, and Bz user-defined rotational constants (in \(cm^{-1}\)).
- Returns:
Tuple with Bx, By, and Bz user-defined rotational constants. If constants were not defined, returns constants computed from the input molecular geometry.
- B_geom()
Returns Bx, By, and Bz rotational constants (in units of \(cm^{-1}\)), calculated from molecular geometry input
- sym()
Molecular point symmetry group
- Parameters:
arg – str Molecular symmetry label, e.g., “C1”, “D2”, “C2v”
- Returns:
richmol.rot.symmetry.SymtopSymmetry
class.
- store_xyz(file_name, comment='')
Stores Cartesian coordinates of atoms into XYZ file
- Parameters:
file_name – str Name of XYZ file
comment – str Optional comment line
- imom()
Computes inertia tensor
- gmat()
Computes rotational kinetic energy matrix (in units of \(cm^{-1}\))
- linear()
Returns True/False if molecule is linear/non-linear
- kappa()
Returns rotational asymmetry parameter kappa = (2*B-A-C)/(A-C)
- store(filename, name=None, comment=None, replace=False)
Stores object in HDF5 file
- Parameters:
filename – str Name of HDF5 file
name – str Name of the data group, by default name of the variable is used
comment – str User comment
replace – bool If True, the existing data set will be replaced
- read(filename, name=None)
Reads object from HDF5 file
- Parameters:
filename – str Name of HDF5 file
name – str Name of the data group, if None, first group with matching “__class_name__” attribute will be loaded
- richmol.rot.mol_tensor(val)
Declares new molecular-frame Cartesian tensor which is dynamically rotated to a chosen molecular frame
Molecule.frame()
whenever the latter is changed- Parameters:
val – array Cartesian tensor
- Returns:
object, needs to be assigned to an attribute of class
Molecule
Rotational solutions
- richmol.rot.solve(mol, Jmin=0, Jmax=10, verbose=False, filter=<function <lambda>>)
Solves rotational eigenvalue problem
- Parameters:
mol – Molecule Molecular parameters.
Jmin – int Min value of J quantum number.
Jmax – int Max value of J quantum number.
verbose – bool If True, some log will be printed.
filter –
function(**kw) State filter, takes as arguments state quantum numbers J and m and symmetry sym, and returns True or False depending if the corresponding state needs to be included (True) or excluded (False) form the basis. By default, all states spanned by J=Jmin..Jmax will be included. The following keyword arguments are passed into the filter function:
- Jint
J quantum number
- symstr
Symmetry
- mint
m quantum number
- Returns:
Solution
Nested dictionary containing rotational wave functions
richmol.rot.SymtopBasis
for different values of J=Jmin..Jmax and different symmetries
- class richmol.rot.Solution(val=None)
Represents rotational solutions for different values of J and symmetry
This is a subclass of
collections.UserDict
class.- Parameters:
val – dict Nested dictionary containing rotational wave functions
richmol.rot.SymtopBasis
for different values of J=Jmin..Jmax and different symmetries
- store(filename, name=None, comment=None, replace=False)
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 is used
comment – str User comment
replace – bool If True, the existing data set will be replaced
- read(filename, name=None)
Reads object from HDF5 file
- Parameters:
filename – str Name of HDF5 file
name – str Name of the data group, if None, first group with matching “__class_name__” attribute will be loaded
Matrix elements
- class richmol.rot.LabTensor(arg, basis, thresh=None, **kwargs)
Represents rotational matrix elements of laboratory-frame Cartesian tensor operator
This is a subclass of
richmol.field.CarTens
class.- Attrs:
- Usnumpy.complex128 2D array
Cartesian-to-spherical tensor transformation matrix.
- Uxnumpy.complex128 2D array
Spherical-to-Cartesian tensor transformation matrix, Ux = np.conj(Us.T)
- tens_flat1D array
Flattened molecular-frame Cartesian tensor.
- molecule
richmol.rot.Molecule
Molecular parameters.
- Parameters:
arg – numpy.ndarray, list or tuple,
richmol.rot.Molecule
, function, or str Depending on the type of argument, initializes lab-frame molecular tensor operator (arg is numpy.ndarray, list, or tuple), field-free Hamiltonian operator (arg isrichmol.rot.Molecule
), arbitrary function of spherical coordinates theta and phi (arg is function(theta, phi)), or matrix elements of some named operators (arg is str). Named operators: “cos2theta” for cos²(theta), “costheta” for cos(theta).basis –
richmol.rot.Solution
Rotational wave functions.thresh – float Threshold for neglecting matrix elements.
- Kwargs:
- 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 if the corresponding state needs to be included (True) or excluded (False) form the basis. By default, all states spanned by basis 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 state’s rotational or ro-vibrational quanta and energy combined in a string
- wig_jmaxint
Maximal value of J quantum number used in the expansion of input user-defined function of spherical coordinates in terms of Wigner D-functions (used when arg is function).
- leb_degint (deprecated, use leb_ind)
Degree of angular Lebedev quadrature used for computing the Wigner expansion coefficients (used when arg is function)
- leb_indint
Index (0..32) of the angular Lebedev quadrature used for computing the Wigner expansion coefficients (used when arg is a function)