calc_trial_states_qmc Subroutine

public subroutine calc_trial_states_qmc(space_in, nexcit, qmc_iluts, qmc_ht, paired_replicas, ndets_this_proc, trial_iluts, evecs_this_proc, space_sizes, space_displs)

Arguments

Type IntentOptional Attributes Name
type(subspace_in) :: space_in
integer, intent(in) :: nexcit
integer(kind=n_int), intent(in) :: qmc_iluts(0:,:)
type(ll_node), intent(inout), pointer :: qmc_ht(:)
logical, intent(in) :: paired_replicas
integer, intent(out) :: ndets_this_proc
integer(kind=n_int), intent(out) :: trial_iluts(0:,:)
real(kind=dp), intent(out), allocatable :: evecs_this_proc(:,:)
integer(kind=MPIArg), intent(out) :: space_sizes(0:nProcessors-1)
integer(kind=MPIArg), intent(out) :: space_displs(0:nProcessors-1)

Contents

Source Code


Source Code

    subroutine calc_trial_states_qmc(space_in, nexcit, qmc_iluts, qmc_ht, paired_replicas, ndets_this_proc, &
                                     trial_iluts, evecs_this_proc, space_sizes, space_displs)

        use CalcData, only: subspace_in
        use DetBitOps, only: ilut_lt, ilut_gt
        use FciMCData, only: ilutHF, CurrentDets, TotWalkers
        use Parallel_neci, only: nProcessors
        use semi_stoch_gen
        use sort_mod, only: sort
        use SystemData, only: tAllSymSectors

        type(subspace_in) :: space_in
        integer, intent(in) :: nexcit
        integer(n_int), intent(in) :: qmc_iluts(0:, :)
        type(ll_node), pointer, intent(inout) :: qmc_ht(:)
        logical, intent(in) :: paired_replicas
        integer, intent(out) :: ndets_this_proc
        integer(n_int), intent(out) :: trial_iluts(0:, :)
        real(dp), allocatable, intent(out) :: evecs_this_proc(:, :)
        integer(MPIArg), intent(out) :: space_sizes(0:nProcessors - 1), space_displs(0:nProcessors - 1)

        integer(n_int), allocatable :: ilut_list(:, :)
        integer, allocatable :: det_list(:, :)
        integer :: i, j, ierr
        integer(MPIArg) :: ndets_all_procs, ndets_this_proc_mpi
        character(len=*), parameter :: this_routine = "calc_trial_states_qmc"

        if (paired_replicas) then
            ASSERT(nexcit == (lenof_sign.div.2))
        else
            ASSERT(nexcit == lenof_sign)
        end if

        ndets_this_proc = 0
        trial_iluts = 0_n_int

        ! Choose the correct generating routine.
        if (space_in%tHF) call add_state_to_space(ilutHF, trial_iluts, ndets_this_proc)
        if (space_in%tPops) call generate_space_most_populated(space_in%npops, &
                                                     space_in%tApproxSpace, space_in%nApproxSpace, trial_iluts, ndets_this_proc, GLOBAL_RUN)
        if (space_in%tRead) call generate_space_from_file(space_in%read_filename, trial_iluts, ndets_this_proc)
        if (space_in%tDoubles) then
            if (tGUGA) then
                call generate_sing_doub_guga(trial_iluts, ndets_this_proc, .false.)
            else
                call generate_sing_doub_determinants(trial_iluts, ndets_this_proc, .false.)
            end if
        end if
        if (space_in%tCAS) call generate_cas(space_in%occ_cas, space_in%virt_cas, trial_iluts, ndets_this_proc)
        if (space_in%tRAS) call generate_ras(space_in%ras, trial_iluts, ndets_this_proc)
        if (space_in%tOptimised) call generate_optimised_space(space_in%opt_data, space_in%tLimitSpace, &
                                                               trial_iluts, ndets_this_proc, space_in%max_size)
        if (space_in%tMP1) call generate_using_mp1_criterion(space_in%mp1_ndets, trial_iluts, ndets_this_proc)
        if (space_in%tFCI) then
            if (tAllSymSectors) then
                call gndts_all_sym_this_proc(trial_iluts, .true., ndets_this_proc)
            else
                call generate_fci_core(trial_iluts, ndets_this_proc)
            end if
        end if

        if (.not. (space_in%tPops .or. space_in%tRead .or. space_in%tDoubles .or. space_in%tCAS .or. &
                   space_in%tRAS .or. space_in%tOptimised .or. space_in%tMP1 .or. space_in%tFCI)) then
            call stop_all(this_routine, "A space for the trial functions was not chosen.")
        end if

        ndets_this_proc_mpi = int(ndets_this_proc, MPIArg)
        call MPIAllGather(ndets_this_proc_mpi, space_sizes, ierr)
        ndets_all_procs = sum(space_sizes)

        space_displs(0) = 0_MPIArg
        do i = 1, nProcessors - 1
            space_displs(i) = space_displs(i - 1) + space_sizes(i - 1)
        end do

        call sort(trial_iluts(:, 1:ndets_this_proc), ilut_lt, ilut_gt)

        allocate(evecs_this_proc(nexcit, ndets_this_proc), stat=ierr)

        call get_qmc_trial_weights(trial_iluts, ndets_this_proc, qmc_iluts, qmc_ht, nexcit, paired_replicas, evecs_this_proc)

    end subroutine calc_trial_states_qmc