save_open_basis_levels Subroutine

public subroutine save_open_basis_levels(number_of_open_basis_levels, open_basis_levels, basis_wavevectors)

saves indices to open levels in the basis and corresponding wavevectors (in A^2)

Arguments

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

number of energetically accessible levels in the basis

integer(kind=int32), intent(inout), allocatable :: open_basis_levels(:)

array holding indices to energetically accessible levels in the basis

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

array holding wavevectors calculated w.r.t energetically accessible levels in the basis


Calls

proc~~save_open_basis_levels~~CallsGraph proc~save_open_basis_levels save_open_basis_levels interface~allocate_1d allocate_1d proc~save_open_basis_levels->interface~allocate_1d proc~is_open is_open proc~save_open_basis_levels->proc~is_open proc~wavevector_squared_from_energy wavevector_squared_from_energy proc~save_open_basis_levels->proc~wavevector_squared_from_energy proc~total_energy total_energy proc~is_open->proc~total_energy proc~wavevector_squared_from_energy->proc~total_energy proc~write_error write_error proc~wavevector_squared_from_energy->proc~write_error proc~write_message write_message proc~write_error->proc~write_message

Called by

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

Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), public :: count_
integer(kind=int32), public :: level_index_1_

Source Code

      subroutine save_open_basis_levels(number_of_open_basis_levels,           &
         open_basis_levels, basis_wavevectors)
         !! saves indices to open levels in the basis and corresponding
         !! wavevectors (in A^2)
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: number_of_open_basis_levels
            !! number of energetically accessible levels in the basis
         integer(int32), intent(inout), allocatable :: open_basis_levels(:)
            !! array holding indices to energetically accessible levels
            !! in the basis
         real(dp), intent(inout), allocatable :: basis_wavevectors(:)
            !! array holding wavevectors calculated w.r.t energetically
            !! accessible levels in the basis
         !---------------------------------------------------------------------!
         integer(int32) :: count_, level_index_1_
         !---------------------------------------------------------------------!
         call allocate_1d(open_basis_levels, number_of_open_basis_levels)
         call allocate_1d(basis_wavevectors, number_of_open_basis_levels)
         !---------------------------------------------------------------------!
         count_ = 0
         do level_index_1_ = 1, number_of_basis_levels
            if (is_open(internal_energies(level_index_1_))) then
               count_ = count_ + 1
               open_basis_levels(count_) = level_index_1_
               basis_wavevectors(count_) = sqrt(                               &
                  wavevector_squared_from_energy(internal_energies(            &
                  level_index_1_)) ) / bohr_to_angstrom 
            endif
         enddo
         !---------------------------------------------------------------------!
      end subroutine save_open_basis_levels