read_input_file Subroutine

public subroutine read_input_file()

reads the input file prepared by the user using NAMELIST feature the code uses 3 namelists: "input", "basis" and "potential"

Arguments

None

Calls

proc~~read_input_file~~CallsGraph proc~read_input_file read_input_file proc~file_io_status file_io_status proc~read_input_file->proc~file_io_status interface~allocate_3d allocate_3d proc~read_input_file->interface~allocate_3d proc~to_lowercase to_lowercase proc~read_input_file->proc~to_lowercase interface~allocate_1d allocate_1d proc~read_input_file->interface~allocate_1d proc~check_namelist_input check_namelist_input proc~read_input_file->proc~check_namelist_input proc~write_message write_message proc~read_input_file->proc~write_message proc~check_namelist_basis check_namelist_basis proc~read_input_file->proc~check_namelist_basis proc~check_namelist_potential check_namelist_potential proc~read_input_file->proc~check_namelist_potential proc~input_summary input_summary proc~read_input_file->proc~input_summary proc~integer_to_character integer_to_character proc~file_io_status->proc~integer_to_character proc~write_error write_error proc~file_io_status->proc~write_error proc~char_to_lowercase char_to_lowercase proc~to_lowercase->proc~char_to_lowercase proc~check_namelist_input->proc~to_lowercase proc~check_namelist_input->proc~write_message proc~check_namelist_input->proc~integer_to_character proc~check_namelist_input->proc~write_error interface~incorrect_value incorrect_value proc~check_namelist_input->interface~incorrect_value proc~check_namelist_basis->proc~integer_to_character proc~check_namelist_basis->interface~incorrect_value proc~check_namelist_potential->proc~integer_to_character proc~check_namelist_potential->interface~incorrect_value proc~input_summary->proc~write_message proc~float_to_character float_to_character proc~input_summary->proc~float_to_character proc~input_summary->proc~integer_to_character proc~total_energy total_energy proc~input_summary->proc~total_energy proc~write_error->proc~write_message proc~incorrect_value_ch incorrect_value_ch interface~incorrect_value->proc~incorrect_value_ch proc~incorrect_value_dp incorrect_value_dp interface~incorrect_value->proc~incorrect_value_dp proc~incorrect_value_int32 incorrect_value_int32 interface~incorrect_value->proc~incorrect_value_int32 proc~incorrect_value_sp incorrect_value_sp interface~incorrect_value->proc~incorrect_value_sp proc~incorrect_value_ch->proc~write_error proc~incorrect_value_dp->proc~write_error proc~incorrect_value_int32->proc~write_error proc~incorrect_value_sp->proc~write_error

Called by

proc~~read_input_file~~CalledByGraph proc~read_input_file read_input_file program~scattering SCATTERING program~scattering->proc~read_input_file

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: coupling_index_
character(len=200), private :: err_message
integer(kind=int32), private :: io_status
integer(kind=int32), private :: level_index_1_
integer(kind=int32), private :: level_index_2_

Source Code

      subroutine read_input_file
         !! reads the input file prepared by the user using NAMELIST feature
         !! the code uses 3 namelists: "input", "basis" and "potential"
         !---------------------------------------------------------------------!
         character(len = 200):: err_message
         integer(int32) :: level_index_1_, level_index_2_, coupling_index_,    &
            io_status
         !---------------------------------------------------------------------!
         namelist / INPUT / label, reduced_mass, relative_energy_flag, energy, &
            jtot_min, jtot_max, jtot_step, r_min, r_max, r_step, steps,        &
            potential_depth, consecutive_blocks_threshold,                     &
            elastic_xs_threshold, inelastic_xs_threshold,                      &
            number_of_basis_levels, initial_level, number_of_r_points,         &
            number_of_legendre_indices, total_number_of_coupling_terms,        &
            n_skip_lines, coupling_terms_r_unit, coupling_terms_file_name,     &
            s_matrix_file_name, print_partial_cross_sections,                  &
            partial_xs_file_name, print_level
         namelist / BASIS / vib_levels, rot_levels, internal_energies
         namelist / POTENTIAL / legendre_indices, vib_couplings, rot_couplings,&
            vib_prime_couplings, rot_prime_couplings
         !---------------------------------------------------------------------!
         open(unit=input_unit, action='read', form='formatted',                &
            access='sequential', status = 'old', iostat = io_status,           &
            iomsg = err_message)
         call file_io_status(io_status, err_message, input_unit, 'o')
         !---------------------------------------------------------------------!
         call write_message("-- Reading input file...")
         !---------------------------------------------------------------------!
         ! Read the "input" namelist:
         !---------------------------------------------------------------------!
         read(unit = input_unit, nml = INPUT, iostat = io_status,              &
            iomsg = err_message)
         call file_io_status(io_status, err_message, input_unit, 'r')
         !---------------------------------------------------------------------!
         ! Check if the variables from "input" namelist are supplied correctly:
         !---------------------------------------------------------------------!
         call check_namelist_input
         !---------------------------------------------------------------------!
         ! Set values/sizes of related variables
         !---------------------------------------------------------------------!
         if (jtot_max.eq.-1) jtot_max = 999999
         !---------------------------------------------------------------------!
         call allocate_1d(vib_levels,number_of_basis_levels)
         call allocate_1d(rot_levels,number_of_basis_levels)
         call allocate_1d(internal_energies,number_of_basis_levels)
         !---------------------------------------------------------------------!
         call allocate_1d(legendre_indices,number_of_legendre_indices)
         call allocate_1d(vib_couplings,total_number_of_coupling_terms)
         call allocate_1d(vib_prime_couplings,total_number_of_coupling_terms)
         call allocate_1d(rot_couplings,total_number_of_coupling_terms)
         call allocate_1d(rot_prime_couplings,total_number_of_coupling_terms)
         !---------------------------------------------------------------------!
         select case(to_lowercase(trim(adjustl(coupling_terms_r_unit))))
            case("bohr")
               radial_term_distance_converter = 1.0_dp
            case("angstrom")
               radial_term_distance_converter = bohr_to_angstrom
         end select
         !---------------------------------------------------------------------!
         radial_term_energy_converter = 1.0_dp / hartree_to_cm
         !---------------------------------------------------------------------!
         ! Read the "basis" namelist & check the values
         !---------------------------------------------------------------------!
         read(unit = input_unit, nml = BASIS, iostat = io_status,              &
            iomsg = err_message)
         call file_io_status(io_status, err_message, input_unit, 'r')
         call check_namelist_basis
         !---------------------------------------------------------------------!
         ! The code reads all the total_number_of_coupling_terms coupling terms,
         ! but some of them will not be used in the calculations. Here, the
         ! code prepares the arrays of minimal_number_of_coupling_terms size,
         ! that will hold only the necessary terms           !
         !---------------------------------------------------------------------!
         minimal_number_of_coupling_terms                                      &
            = number_of_basis_levels * (number_of_basis_levels + 1) / 2
         !---------------------------------------------------------------------!
         call allocate_1d(reduced_rot_couplings,                               &
            minimal_number_of_coupling_terms)
         call allocate_1d(reduced_rot_prime_couplings,                         &
            minimal_number_of_coupling_terms)
         call allocate_1d(reduced_vib_couplings,                               &
            minimal_number_of_coupling_terms)
         call allocate_1d(reduced_vib_prime_couplings,                         &
            minimal_number_of_coupling_terms)
         !---------------------------------------------------------------------!
         coupling_index_ = 0
         do level_index_1_ = 1, number_of_basis_levels
            do level_index_2_ = level_index_1_, number_of_basis_levels
               coupling_index_ = coupling_index_ + 1
               reduced_vib_couplings(coupling_index_)                          &
                  = vib_levels(level_index_1_)
               reduced_rot_couplings(coupling_index_)                          &
                  = rot_levels(level_index_1_)
               reduced_vib_prime_couplings(coupling_index_)                    &
                  = vib_levels(level_index_2_)
               reduced_rot_prime_couplings(coupling_index_)                    &
                  = rot_levels(level_index_2_)
            enddo
         enddo
         !---------------------------------------------------------------------!
         ! Read the "potential" namelist & check the values
         !---------------------------------------------------------------------!
         read(unit = input_unit, nml=POTENTIAL, iostat = io_status,            &
            iomsg = err_message)
         call file_io_status(io_status, err_message, input_unit, 'r')
         call check_namelist_potential
         !---------------------------------------------------------------------!
         close(input_unit)
         !---------------------------------------------------------------------!
         ! Set values/sizes of related variables
         !---------------------------------------------------------------------!
         call allocate_1d(r_grid,number_of_r_points)
         call allocate_3d(tabulated_coupling_terms,number_of_r_points,         &
            number_of_legendre_indices, total_number_of_coupling_terms)
         call allocate_3d(coupling_terms,number_of_r_points,                   &
            number_of_legendre_indices, minimal_number_of_coupling_terms)
         call allocate_3d(coupling_terms_b_coeffs,number_of_r_points,          &
            number_of_legendre_indices, minimal_number_of_coupling_terms)
         call allocate_3d(coupling_terms_c_coeffs,number_of_r_points,          &
            number_of_legendre_indices, minimal_number_of_coupling_terms)
         call allocate_3d(coupling_terms_d_coeffs,number_of_r_points,          &
            number_of_legendre_indices, minimal_number_of_coupling_terms)
         !---------------------------------------------------------------------!
         ! Summarize the input parameters
         !---------------------------------------------------------------------!
         call input_summary
         !---------------------------------------------------------------------!
      end subroutine read_input_file