lib_bound: bound state calculations
Potential energy matrix elements
Potential matrix elements: atom + open-shell diatom, body fixed
- qmsV = qms_V_bf_DLMK2(sjmkf, sjmki, sLMK, qnames)
Compute matrix elements of the complex conjugate of a Wigner rotation matrix in a symmetric top basis.
Input arguments
- sjmkf:
A
q_set
data structure (seelib_qmat
) containing the set of bra quantum numbers: \(\{j_f m_f k_f; f=1,2,\ldots\}\)Or:
A Scilab
matrix
with two columns, withsjmkf(f,:) = [j_f m_f k_f]
- sjmki:
A
q_set
data structure (seelib_qmat
) containing the set of ket quantum numbers \(\{j_i m_i k_i; i=1,2,\ldots\}\)Or:
A Scilab
matrix
with two columns, withsjmki(i,:) = [j_i m_i k_i]
- sLMK:
A
q_set
data structure (seelib_qmat
) containing \(LMK\) of the Wigner D-matrix element \(D^{(L)\ast}_{MK}\).Or:
A Scilab \(n \times 3\)
matrix
Q
with elements \(Q(i,:)=[L_i M_i K_i]\)- qnames:
A
vector
ofstrings
with the names of the quantum numbers. By defaultqnames = ['j' 'm' 'k' 'L' 'M' 'K']
Output argument
- qmsV:
A
q_mats
data structure (seelib_qmat
) containing a list of matrices with Gaunt coefficients for all values of \(L\)qmsV.s1 = sjmkf qmsV.s2 = sjmki qmsV.sa = sLMK qmsV.As(i) = V_n // list of matrices V_n(f, i) = <sjmkf(f,:)|D*L_MK|sjmki(i,:)>Example:
Here is an example of how to use this function:
// Example usage of the function qms_V_bf_PL2 lib_load(['lib_bound' 'lib_special' 'lib_tenop']) jmk = [0 0 0 1 0 0 1 0 1 2 0 0 2 0 1 2 0 2] LMKs = [0 0 0 1 0 0 1 0 1 2 0 0 2 0 1 2 0 2] qmsV1 = qms_V_bf_DLMK2(jmk, jmk, LMKs, ['j' 'm' 'k' 'L' 'M' 'K']) // |J(3), K(3)>A(6;6)<(same)| [L(5); 5] // Or, create q_sets first s1 = q_set(jmk, ['j' 'm' 'k']) sL = q_set(LMKs, ['L' 'M' 'K']) // LMKs has 3 columns qmsV2 = qms_V_bf_DLMK2(s1, s1, sL, [s1.qnames sL.qnames]) if qmsV1==qmsV2 mprintf('The same!\n') else mprintf('qmsV1 <> qmsV2\n') end
Gaunt coefficients
- qmsV = qms_V_bf_PL2(sjkf, sjki, sL, qnames)
This function computes matrices with Gaunt coefficients:
where \(Y_{jk}\) are spherical harmonics with the Condon and Shortley phase convention and \(P_L\) are Legendre polynomials.
Input arguments
- sjkf:
A
q_set
data structure (seelib_qmat
) containing the set of bra quantum numbers: \(\{j_f k_f; f=1,2,\ldots\}\)Or:
A Scilab
matrix
with two columns, withsjkf(f,:) = [j_f k_f]
- sjki:
A
q_set
data structure (seelib_qmat
) containing the set of ket quantum numbers \(\{j_i k_i; i=1,2,\ldots\}\)Or:
A Scilab
matrix
with two columns, withsjki(i,:) = [j_i k_i]
- sL:
A
q_set
data structure (seelib_qmat
) containing the orders \(L\) or the Legendre polynomials \(P_L\)Or:
A Scilab
vector
with the values of \(L\)- qnames:
A
vector
ofstrings
with the names of the quantum numbers. By defaultqnames = ['j' 'k' 'L']
Output argument
- qmsV:
A
q_mats
data structure (seelib_qmat
) containing a list of matrices with Gaunt coefficients for all values of \(L\)qmsV.s1 = sjkf qmsV.s2 = sjki qmsV.sa = sL qmsV.As(i) = VL // list of Gaunt coefficient matrices for: L = sL(i, 'L') // VL(f, i) = <sjkf(f, :)|P_L|sjki(i,:)>Example:
Here is an example of how to use this function:
// Example usage of the function qms_V_bf_PL2 lib_load(['lib_bound' 'lib_special' 'lib_tenop']) jk = [0 0 1 0 1 1 2 0 2 1 2 2] Ls = 0:4 qmsV1 = qms_V_bf_PL2(jk, jk, Ls, ['j' 'k' 'L']) // |J(3), K(3)>A(6;6)<(same)| [L(5); 5] // Or, create q_sets first s1 = q_set(jk, ['j' 'k']) sL = q_set(Ls(:), 'L') // q_set requires a column vector here qmsV2 = qms_V_bf_PL2(s1, s1, sL, [s1.qnames 'L']) if qmsV1==qmsV2 mprintf('The same!\n') else mprintf('qmsV1 <> qmsV2\n') end
G. Dhont, W. B. Zeimen, G. C. Groenenboom, and A. van der Avoird. Theoretical study of the He–HF⁺ complex; ii. rovibronic states from coupled diabatic potential energy surfaces. J. Chem. Phys., 120:103–116, 2004. doi:10.1063/1.1629672.
Taha Selim, Ad van der Avoird, and Gerrit C. Groenenboom. State-to-state rovibrational transition rates for CO₂ in the bend mode in collisions with he atoms. J. Chem. Phys., 159:164310, 2023. arXiv:2309.03781, doi:10.1063/5.0174787.