generate_entire_ras_space Subroutine

public subroutine generate_entire_ras_space(ras, classes, space_size, ilut_list, parities)

Arguments

Type IntentOptional Attributes Name
type(ras_parameters), intent(in) :: ras
type(ras_class_data), intent(in) :: classes(:)
integer, intent(in) :: space_size
integer(kind=n_int), intent(out) :: ilut_list(0:NIfTot,space_size)
integer, intent(out), optional :: parities(space_size)

Contents


Source Code

    subroutine generate_entire_ras_space(ras, classes, space_size, ilut_list, parities)

        type(ras_parameters), intent(in) :: ras
        type(ras_class_data), intent(in) :: classes(:)
        integer, intent(in) :: space_size
        integer(n_int), intent(out) :: ilut_list(0:NIfTot, space_size)
        integer, intent(out), optional :: parities(space_size)
        integer :: nI(nel)
        integer :: string(tot_nelec)
        integer, allocatable, dimension(:, :) :: string_list
        integer, allocatable, dimension(:, :) :: min_indices
        integer :: temp_class
        integer :: string_address, block_address
        integer :: i, j, k, l, m, n, o, counter
        logical :: none_left
        integer :: par

        allocate(string_list(tot_nelec, ras%num_strings))

        ! min_indices(i,j) stores the index of the first state in class i with symmetry label j.
        allocate(min_indices(ras%num_classes, 0:7))

        ! Loop over all classes, and for each, generate all states in that class. Store all these
        ! states in blocks, one after another, in a 1d array (sting_list). Within each block
        ! find_address and classes(i)%address_map are used to to find an address *relative to
        ! the first state* in that class.
        block_address = 0
        do i = 1, ras%num_classes
            ! Address of the last string in the last class.
            associate(ind => i)
            if (i > 1) block_address = block_address + classes(ind - 1)%class_size
            end associate

            do j = 0, 7
                min_indices(i, j) = block_address + 1
                do k = 0, j - 1
                    min_indices(i, j) = min_indices(i, j) + classes(i)%num_sym(k)
                end do
            end do

            call generate_first_full_string(string, ras, classes(i))
            string_address = classes(i)%address_map(get_address(classes(i), string)) + &
                             block_address
            string_list(:, string_address) = string

            do j = 2, classes(i)%class_size
                call generate_next_string(string, ras, classes(i), none_left)
                string_address = classes(i)%address_map(get_address(classes(i), string)) + &
                                 block_address
                string_list(:, string_address) = string
            end do
        end do

        ! Now combine the alpha and beta strings together to create the full list of states.

        counter = 0
        ilut_list = 0
        ! Loop over all classes, call it class_i.
        do i = 1, ras%num_classes
            ! For class_i, loop over all classes which can be combined with this one.
            do j = 1, classes(i)%num_comb
                temp_class = classes(i)%allowed_combns(j)
                ! Loop over all symmetries.
                do k = 0, 7
                    l = ieor(HFSym_ras, k)
                    ! Add all combinations of states in these symmetry blocks to ilut_list.
                    do m = min_indices(i, k), min_indices(i, k) + &
                        classes(i)%num_sym(k) - 1
                        do n = min_indices(temp_class, l), min_indices(temp_class, l) + &
                            classes(temp_class)%num_sym(l) - 1

                            ! Beta string.
                            nI(1:tot_nelec) = string_list(:, m) * 2 - 1
                            ! Alpha string.
                            nI(tot_nelec + 1:nel) = string_list(:, n) * 2

                            ! Replace all orbital numbers, orb, with the true orbital
                            ! numbers, BRR(orb). Also, sort this list.
                            do o = 1, nel
                                nI(o) = BRR(nI(o))
                            end do

                            call sort_orbitals(nI, par)

                            ! Find bitstring representation.
                            counter = counter + 1
                            call EncodeBitDet(int(nI), ilut_list(:, counter))
                            if (present(parities)) then
                                parities(counter) = par
                            end if
                        end do
                    end do
                end do
            end do
        end do

        if (counter /= space_size) call stop_all("generate_entire_ras_space", "Wrong number of &
                                                 &states found in generating loop.")

        deallocate(string_list)
        deallocate(min_indices)

    end subroutine generate_entire_ras_space