state_to_state_cross_sections_mod Module

This module provides functions and subroutines for calculating and analyzing state-to-state cross-sections. It is divided into three parts:

(1) Calculating cross-sections: Functions for computing state-to-state cross-sections based on quantum states, S-matrix, and scattering parameters ("calculate_state_to_state_cross_section", "compute_individual_cross_section", "get_block_indices", "sum_cross_section_contributions", "compute_real_component", "compute_imag_component") (2) Printing cross-sections: Subroutines to output the largest partial cross-sections, providing both basic and detailed information. ("print_largest_partial_cross_sections", "print_basic_cross_section_info", "print_detailed_cross_section_info") (3) Threshold checking: Subroutine to check if the computed cross-sections meet specified convergence conditions ("check_cross_section_thresholds")


Uses

  • module~~state_to_state_cross_sections_mod~~UsesGraph module~state_to_state_cross_sections_mod state_to_state_cross_sections_mod module~physics_utilities_mod physics_utilities_mod module~state_to_state_cross_sections_mod->module~physics_utilities_mod iso_fortran_env iso_fortran_env module~state_to_state_cross_sections_mod->iso_fortran_env module~utility_functions_mod utility_functions_mod module~state_to_state_cross_sections_mod->module~utility_functions_mod module~global_variables_mod global_variables_mod module~state_to_state_cross_sections_mod->module~global_variables_mod module~array_operations_mod array_operations_mod module~state_to_state_cross_sections_mod->module~array_operations_mod module~physics_utilities_mod->iso_fortran_env module~physics_utilities_mod->module~utility_functions_mod module~physics_utilities_mod->module~global_variables_mod module~physics_utilities_mod->module~array_operations_mod module~utility_functions_mod->iso_fortran_env module~global_variables_mod->iso_fortran_env module~array_operations_mod->iso_fortran_env

Used by

  • module~~state_to_state_cross_sections_mod~~UsedByGraph module~state_to_state_cross_sections_mod state_to_state_cross_sections_mod program~scattering SCATTERING program~scattering->module~state_to_state_cross_sections_mod

Contents


Functions

private function compute_imag_component(initial_index_, final_index_, s_matrix_imag_) result(imag_contribution_)

Computes the imaginary part of the cross-section contribution for given indices.

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: initial_index_

index pointing to the basis element involving initial state

integer(kind=int32), intent(in) :: final_index_

index pointing to the basis element involving final state

real(kind=dp), intent(in) :: s_matrix_imag_(:,:)

imaginary part of the S-matrix

Return Value real(kind=dp)

(output) contribution to the cross-section from imaginary part of the S-matrix

private function compute_individual_cross_section(initial_state_, final_state_, open_basis_levels_, basis_wavevectors_, s_matrix_real_, s_matrix_imag_, channel_indices_, channel_l_values_, total_angular_momentum_) result(cross_section_)

Calculates cross-section for a given initial and final state.

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: initial_state_

index pointing to the intial state in basis arrays

integer(kind=int32), intent(in) :: final_state_

index pointing to the final state in basis arrays

integer(kind=int32), intent(in) :: open_basis_levels_(:)

holds indices to the basis arrays which correspond to open channels

real(kind=dp), intent(in) :: basis_wavevectors_(:)

holds wavevectors k_{i}

real(kind=dp), intent(in) :: s_matrix_real_(:,:)

real and imaginary parts of the S-matrix

real(kind=dp), intent(in) :: s_matrix_imag_(:,:)

real and imaginary parts of the S-matrix

integer(kind=int32), intent(in) :: channel_indices_(:)

holds the indices pointing to the basis arrays

integer(kind=int32), intent(in) :: channel_l_values_(:)

holds all values of l

integer(kind=int32), intent(in) :: total_angular_momentum_

total angular momentum

Return Value real(kind=dp)

(output) cross-section

private function compute_real_component(initial_index_, final_index_, l_initial_, l_final_, s_matrix_real_) result(real_contribution_)

Computes the real part of the cross-section contribution for given indices.

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: initial_index_

index pointing to the basis element involving initial state

integer(kind=int32), intent(in) :: final_index_

index pointing to the basis element involving final state

integer(kind=int32), intent(in) :: l_initial_

pre-collisional \(l\) value

integer(kind=int32), intent(in) :: l_final_

post-collisional \(l\) value

real(kind=dp), intent(in) :: s_matrix_real_(:,:)

real part of the S-matrix

Return Value real(kind=dp)

(output) contribution to the cross-section from real part of the S-matrix

private function get_block_indices(v_, j_, channel_indices_) result(block_indices_)

Returns indices from channel_indices_ that match the specified quantum state.

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: v_

vibrational quantum number

integer(kind=int32), intent(in) :: j_

rotational quantum number

integer(kind=int32), intent(in) :: channel_indices_(:)

holds the indices pointing to the basis arrays

Return Value integer(kind=int32), allocatable, (:)

(output) indices that match the specified quantum state

private function sum_cross_section_contributions(init_indices_, final_indices_, s_matrix_real_, s_matrix_imag_, channel_l_values_) result(sum_contributions_)

Sum the contributions to the cross-section from S-matrix components

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: init_indices_(:)

indices pointing to basis element involving initial state

integer(kind=int32), intent(in) :: final_indices_(:)

indices pointing to basis element involving final state

real(kind=dp), intent(in) :: s_matrix_real_(:,:)

real and imaginary parts of the S-matrix

real(kind=dp), intent(in) :: s_matrix_imag_(:,:)

real and imaginary parts of the S-matrix

integer(kind=int32), intent(in) :: channel_l_values_(:)

holds all values of l

Return Value real(kind=dp)

(output) contribution to the cross-section from S-matrix components


Subroutines

public subroutine add_cross_sections(size_, partial_, accumulated_)

Add partial cross-sections to accumulated cross-sections

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: size_

number of open basis levels (defines the size of both arrays)

real(kind=dp), intent(in) :: partial_(size_*size_)

array holding partial cross-sections

real(kind=dp), intent(inout) :: accumulated_(size_*size_)

array holding accumulated cross-sections

public subroutine calculate_state_to_state_cross_section(total_angular_momentum_, open_basis_levels_, basis_wavevectors_, s_matrix_real_, s_matrix_imag_, channel_indices_, channel_l_values_, cross_section_array_)

Calculates all state-to-state cross-sections.

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: total_angular_momentum_

total angular momentum

integer(kind=int32), intent(in) :: open_basis_levels_(:)

holds indices to the basis arrays which correspond to open channels

real(kind=dp), intent(in) :: basis_wavevectors_(:)

holds wavevectors k_{i}

real(kind=dp), intent(in) :: s_matrix_real_(:,:)

real and imaginary parts of the S-matrix

real(kind=dp), intent(in) :: s_matrix_imag_(:,:)

real and imaginary parts of the S-matrix

integer(kind=int32), intent(in) :: channel_indices_(:)

holds the indices pointing to the basis arrays

integer(kind=int32), intent(in) :: channel_l_values_(:)

holds all values of l

real(kind=dp), intent(inout) :: cross_section_array_(:)

array holding all XSs

public subroutine check_cross_section_thresholds(largest_elastic_xs_, largest_inelastic_xs_, consecutive_elastic_, consecutive_inelastic_, terminate_)

Checks if the elastic_xs_threshold (threshold for elastic XS) and inelastic_xs_threshold (threshold for inelastic XS) conditions are already fulfilled.

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: largest_elastic_xs_

largest elastic XS in the block

real(kind=dp), intent(in) :: largest_inelastic_xs_

largest inelastic XS in the block

integer(kind=int32), intent(inout) :: consecutive_elastic_

number of consecutive blocks meeting condition on elastic XS

integer(kind=int32), intent(inout) :: consecutive_inelastic_

number of consecutive blocks meeting condition on inelastic XS

logical, intent(inout) :: terminate_

flag to indicate termination of loop based on thresholds

public subroutine determine_largest_cross_sections(open_basis_levels_, cross_sections_, max_elastic_cross_section_, max_inelastic_cross_section_, max_elastic_index_, max_inelastic_index_1_, max_inelastic_index_2_)

Determine the largest partial elastic and inelastic cross-sections in a given set of cross-sections.

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: open_basis_levels_(:)

array holding indices to open basis levels

real(kind=dp), intent(in) :: cross_sections_(:)

array holding partial cross-sections

real(kind=dp), intent(out) :: max_elastic_cross_section_

largest elastic cross-section

real(kind=dp), intent(out) :: max_inelastic_cross_section_

largest inelastic cross-section

integer(kind=int32), intent(out) :: max_elastic_index_

index pointing to the largest elastic cross-section

integer(kind=int32), intent(out) :: max_inelastic_index_1_

first index pointing to the largest inelastic cross-section

integer(kind=int32), intent(out) :: max_inelastic_index_2_

second index pointing to the largest inelastic cross-section

public subroutine initialize_cross_section_arrays(size_, partial_block_, partial_jtot_, accumulated_)

allocate arrays keeping accumulated and partial cross-sections in each jtot and parity block

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: size_

number of open basis levels (defines the size of both arrays)

real(kind=dp), intent(inout), allocatable :: partial_block_(:)

array holding partial cross-sections in a parity block

real(kind=dp), intent(inout), allocatable :: partial_jtot_(:)

array holding partial cross-sections in a jtot block

real(kind=dp), intent(inout), allocatable :: accumulated_(:)

array holding accumulated cross-sections

public subroutine print_cross_sections_for_jtot(total_angular_momentum_, open_basis_levels_, cross_sections_)

Prints information about cross-sections at the end of each total angular momentum (jtot) block

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: total_angular_momentum_

total angular momentum

integer(kind=int32), intent(in) :: open_basis_levels_(:)

array holding indices to open basis levels

real(kind=dp), intent(in) :: cross_sections_(:)

holds values of the cross-sections

public subroutine print_final_cross_sections(open_basis_levels_, cross_sections_)

Prints information about cross-sections at the end of the program

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: open_basis_levels_(:)

array holding indices to open basis levels

real(kind=dp), intent(in) :: cross_sections_(:)

holds values of the cross-sections

public subroutine print_largest_partial_cross_sections(total_angular_momentum_, largest_elastic_xs_, largest_inelastic_xs_, elastic_index_, inelastic_index_1_, inelastic_index_2_, open_basis_levels_)

Print the largest partial elastic and inelastic state-to-state cross-sections in a given block.

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: total_angular_momentum_

total angular momentum

real(kind=dp), intent(in) :: largest_elastic_xs_

the largest partial elastic state-to-state XS in the block

real(kind=dp), intent(in) :: largest_inelastic_xs_

the largest partial inelastic state-to-state XS in the block

integer(kind=int32), intent(in) :: elastic_index_

index pointing indirectly to quantum numbers associated with the largest partial elastic state-to-state XS in the block

integer(kind=int32), intent(in) :: inelastic_index_1_

indices pointing indirectly to quantum numbers associated with the largest partial inelastic state-to-state XS in the block

integer(kind=int32), intent(in) :: inelastic_index_2_

indices pointing indirectly to quantum numbers associated with the largest partial inelastic state-to-state XS in the block

integer(kind=int32), intent(in) :: open_basis_levels_(:)

array holding indices to open basis levels

public subroutine save_partial_xs_file_header()

save "header" of the partial cross-sections file: -- label, "itype", number of levels in the basis, reduced mass of the system -- vibrational and rotational quantum numbers -- rovibrational energies -- index pointing to the initial level and the kinetic/total energy

Arguments

None

public subroutine save_partial_xs_single_block(jtot_, block_number_, open_basis_levels_, xs_block)

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: jtot_
integer(kind=int32), intent(in) :: block_number_
integer(kind=int32), intent(in) :: open_basis_levels_(:)
real(kind=dp), intent(in) :: xs_block(:)

private subroutine print_all_cross_sections(open_basis_levels_, cross_sections_)

Prints information about cross-sections from provided "cross_sections_" array

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: open_basis_levels_(:)

array holding indices to open basis levels

real(kind=dp), intent(in) :: cross_sections_(:)

holds values of the cross-sections

private subroutine print_basic_cross_section_info(total_angular_momentum_, cross_section_value_, type_label_)

Prints basic information about the largest elastic and inelastic state-to-state xs (print_level <= 2)

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: total_angular_momentum_

total angular momentum

real(kind=dp), intent(in) :: cross_section_value_

value of the cross-section

character(len=*), intent(in) :: type_label_

"elastic" or "inelastic"

private subroutine print_detailed_cross_section_info(total_angular_momentum_, cross_section_value_, index_1_, index_2_, open_basis_levels_, type_label_)

Prints detailed information about the largest elastic and inelastic state-to-state xs (print_level >= 3)

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: total_angular_momentum_

total angular momentum

real(kind=dp), intent(in) :: cross_section_value_

value of the cross-section

integer(kind=int32), intent(in) :: index_1_

indices pointing indirectly to quantum numbers associated with the largest partial XS in the block

integer(kind=int32), intent(in) :: index_2_

indices pointing indirectly to quantum numbers associated with the largest partial XS in the block

integer(kind=int32), intent(in) :: open_basis_levels_(:)

array holding indices to open basis levels

character(len=*), intent(in) :: type_label_

"elastic" or "inelastic"