generate_initiator_space Subroutine

public subroutine generate_initiator_space(space_in)

Arguments

Type IntentOptional Attributes Name
type(subspace_in) :: space_in

Contents


Source Code

    subroutine generate_initiator_space(space_in)

        ! A wrapper to call the correct generating routine.

        use bit_rep_data, only: flag_initiator
        use bit_reps, only: set_flag, encode_sign
        use FciMCData, only: SpawnedParts
        use searching, only: remove_repeated_states
        use SystemData, only: tAllSymSectors

        type(subspace_in) :: space_in

        integer :: ndets_this_proc, run
        !real(dp) :: zero_sign(lenof_sign)
        character(len=*), parameter :: t_r = "generate_initiator_space"

        ndets_this_proc = 0

        ! Call the requested generating routines.
        if (space_in%tHF) call add_state_to_space(ilutHF, SpawnedParts, ndets_this_proc)
        if (space_in%tPops) call generate_space_most_populated(space_in%npops, &
                           space_in%tApproxSpace, space_in%nApproxSpace, SpawnedParts, ndets_this_proc, GLOBAL_RUN, CurrentDets, TotWalkers)
        if (space_in%tRead) call generate_space_from_file(space_in%read_filename, SpawnedParts, ndets_this_proc)
        if (space_in%tDoubles) call generate_sing_doub_determinants(SpawnedParts, ndets_this_proc, space_in%tHFConn)
        if (space_in%tCAS) call generate_cas(space_in%occ_cas, space_in%virt_cas, SpawnedParts, ndets_this_proc)
        if (space_in%tRAS) call generate_ras(space_in%ras, SpawnedParts, ndets_this_proc)
        if (space_in%tOptimised) call generate_optimised_space(space_in%opt_data, space_in%tLimitSpace, &
                                                               SpawnedParts, ndets_this_proc, space_in%max_size)
        if (space_in%tMP1) call generate_using_mp1_criterion(space_in%mp1_ndets, SpawnedParts, ndets_this_proc)
        if (space_in%tFCI) then
            if (tAllSymSectors) then
                call gndts_all_sym_this_proc(SpawnedParts, .false., ndets_this_proc)
            else
                call generate_fci_core(SpawnedParts, ndets_this_proc)
            end if
            !else if (space_in%tHeisenbergFCI) then
            !call generate_heisenberg_fci(SpawnedParts, ndets_this_proc)
        end if

        ! If two different spaces have been called then there may be some
        ! repeated states. We don't want repeats, so remove them and update
        ! ndets_this_proc accordingly.
        call remove_repeated_states(SpawnedParts, ndets_this_proc)

        !zero_sign = 0.0_dp
        !do i = 1, ndets_this_proc
        !    call encode_sign(SpawnedParts(:,i), zero_sign)

        !    if (tTruncInitiator) then
        !        do run = 1, inum_runs
        !            call set_flag(SpawnedParts(:,i), get_initiator_flag_by_run(run))
        !        end do
        !    end if
        !end do

        ! Set the initiator space size for this process.
        initiator_sizes(iProcIndex) = int(ndets_this_proc, MPIArg)

        ! If requested, remove high energy orbitals so that the space size is below some max.
        if (space_in%tLimitSpace) then
            call remove_high_energy_orbs(SpawnedParts(0:NIfTot, 1:ndets_this_proc), ndets_this_proc, &
                                         space_in%max_size, .true.)
            initiator_sizes(iProcIndex) = int(ndets_this_proc, MPIArg)
        end if

    end subroutine generate_initiator_space