transfer_to_block_form Subroutine

public subroutine transfer_to_block_form(ras, classes, full_vec, ras_vec)

Arguments

Type IntentOptional Attributes Name
type(ras_parameters), intent(in) :: ras
type(ras_class_data), intent(in) :: classes(ras%num_classes)
real(kind=dp), intent(in) :: full_vec(:)
type(ras_vector), intent(inout) :: ras_vec(ras%num_classes,ras%num_classes,0:7)

Contents


Source Code

    subroutine transfer_to_block_form(ras, classes, full_vec, ras_vec)

        ! Take a standard vector (as a 1d array) and transform it to a RAS vector, stored
        ! in matrix form. This can then be input into the perform_multiplication routine.

        ! Important: Because orbitals are ordered differently when using alpha and beta
        ! strings compared to the rest of NECI, some components of the vector will differ
        ! by a factor or -1. However, this routine does *not* apply these minus signs,
        ! so one should not take a vector from somewhere in NECI, use this routine and
        ! then the perform_multiplication routine. The extra minuses should be applied
        ! separately (although there is no routine to do this currently...)

        type(ras_parameters), intent(in) :: ras
        type(ras_class_data), intent(in) :: classes(ras%num_classes)
        real(dp), intent(in) :: full_vec(:)
        type(ras_vector), intent(inout) :: ras_vec(ras%num_classes, ras%num_classes, 0:7)
        integer :: class_i, class_j, j, sym_i, sym_j, ind_i, ind_j
        integer :: counter

        counter = 0

        ! Over all RAS classes.
        do class_i = 1, ras%num_classes
            ! Over all classes allowed with class_i.
            do j = 1, classes(class_i)%num_comb
                class_j = classes(class_i)%allowed_combns(j)
                ! Over all symmetry labels.
                do sym_i = 0, 7
                    ! The symmetry label of string_j is fixed by that of string_i.
                    sym_j = ieor(HFSym_ras, sym_i)
                    ! If there are no strings with the correct symmetries in these classes.
                    if (classes(class_i)%num_sym(sym_i) == 0) cycle
                    if (classes(class_j)%num_sym(sym_j) == 0) cycle
                    ! Over all strings with symmetry sym_i and class class_i.
                    do ind_i = 1, classes(class_i)%num_sym(sym_i)
                        ! Over all strings with symmetry sym_j and class class_j.
                        do ind_j = 1, classes(class_j)%num_sym(sym_j)
                            counter = counter + 1
                            ! Copy the amplitude across.
                            ras_vec(class_i, class_j, sym_i)%elements(ind_i, ind_j) = full_vec(counter)
                        end do
                    end do
                end do
            end do
        end do

    end subroutine transfer_to_block_form