calculates the non-zero algebraic coefficients \( g_{{\lambda},\gamma,\gamma'}^{Jp} \) for a single matrix element - see Eq. (1) in the "Coupling matrix" section; algebraic coefficients are saved to nonzero_algebraic_coefficients array; corresponding indices to legendre_indices are saved to nonzero_legendre_indices array
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | j_ |
pre-collisional rotational angular momentum |
||
integer(kind=int32), | intent(in) | :: | j_prime_ |
post-collisional rotational angular momentum |
||
integer(kind=int32), | intent(in) | :: | omega_ |
\(\bar{\Omega}\) |
||
integer(kind=int32), | intent(inout) | :: | count_nonzero_algebraic_coefficients |
running index counting non-zero coupling coefficients, \( g_{{\lambda},\gamma,\gamma'}^{Jp} \) in the whole matrix; incremented in this subroutine |
||
integer(kind=int32), | intent(inout) | :: | count_nonzero_legendre_terms |
number of non-zero terms in Legendre expansion for a single matrix element of the interaction potential |
||
integer(kind=int32), | intent(inout) | :: | nonzero_legendre_indices(:) |
holds indices pointing to legendre_indices, which correspond to the non-vanishing elements of the sum over \(\lambda\) for each non-zero element of the PES matrix; |
||
real(kind=dp), | intent(inout) | :: | nonzero_algebraic_coefficients(:) |
holds values of the non-zero algebraic coefficients |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer(kind=int32), | private | :: | lambda_ | ||||
integer(kind=int32), | private | :: | legendre_term_index_ | ||||
real(kind=dp), | private | :: | pscoeff |
subroutine process_single_matrix_element(j_, j_prime_, omega_, &
count_nonzero_algebraic_coefficients, count_nonzero_legendre_terms, &
nonzero_legendre_indices, nonzero_algebraic_coefficients)
!! calculates the non-zero algebraic coefficients
!! \\( g\_{{\lambda},\gamma,\gamma'}^{Jp} \\) for a single matrix
!! element - see Eq. (1) in the "Coupling matrix" section;
!! algebraic coefficients are saved to nonzero_algebraic_coefficients
!! array; corresponding indices to legendre_indices are saved to
!! nonzero_legendre_indices array
!---------------------------------------------------------------------!
integer(int32), intent(in) :: j_
!! pre-collisional rotational angular momentum
integer(int32), intent(in) :: j_prime_
!! post-collisional rotational angular momentum
integer(int32), intent(in) :: omega_
!! \\(\bar{\Omega}\\)
integer(int32), intent(inout) :: count_nonzero_algebraic_coefficients
!! running index counting non-zero coupling coefficients,
!! \\( g\_{{\lambda},\gamma,\gamma'}^{Jp} \\) in the whole matrix;
!! incremented in this subroutine
integer(int32), intent(inout) :: count_nonzero_legendre_terms
!! number of non-zero terms in Legendre expansion for a single
!! matrix element of the interaction potential
integer(int32), intent(inout) :: nonzero_legendre_indices(:)
!! holds indices pointing to legendre_indices, which correspond to
!! the non-vanishing elements of the sum over \\(\lambda\\)
!! for each non-zero element of the PES matrix;
real(dp), intent(inout) :: nonzero_algebraic_coefficients(:)
!! holds values of the non-zero algebraic coefficients
!---------------------------------------------------------------------!
integer(int32) :: legendre_term_index_, lambda_
real(dp) :: pscoeff
!---------------------------------------------------------------------!
count_nonzero_legendre_terms = 0
do legendre_term_index_ = 1, number_of_legendre_indices
lambda_ = legendre_indices(legendre_term_index_)
!------------------------------------------------------------------!
! check the condition on 3-j symbol with zero projections
!------------------------------------------------------------------!
if (.not. zero_projections_3j_condition(j_, j_prime_, lambda_)) cycle
!------------------------------------------------------------------!
count_nonzero_algebraic_coefficients = &
count_nonzero_algebraic_coefficients + 1
!------------------------------------------------------------------!
! calculate the Percival-Seaton coefficient
!------------------------------------------------------------------!
pscoeff = percival_seaton_coefficient(j_, j_prime_, lambda_, omega_)
!------------------------------------------------------------------!
! save the Percival-Seaton coefficient
!------------------------------------------------------------------!
nonzero_algebraic_coefficients( &
count_nonzero_algebraic_coefficients) = pscoeff
!------------------------------------------------------------------!
! save indices to legendre_indices corresponding to \\(\lambda\\)
!------------------------------------------------------------------!
nonzero_legendre_indices(count_nonzero_algebraic_coefficients)&
= legendre_term_index_
!------------------------------------------------------------------!
count_nonzero_legendre_terms = count_nonzero_legendre_terms + 1
enddo
!---------------------------------------------------------------------!
end subroutine process_single_matrix_element