perform_spawn Subroutine

public subroutine perform_spawn(ilut_parent, nI, parent_sign, hdiag, tCoreDet)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilut_parent(0:niftot)
integer, intent(in) :: nI(nel)
real(kind=dp), intent(in) :: parent_sign(lenof_sign)
real(kind=dp), intent(in) :: hdiag
logical, intent(in) :: tCoreDet

Contents

Source Code


Source Code

    subroutine perform_spawn(iLut_parent, nI, parent_sign, hdiag, tCoreDet)
        real(dp), intent(in) :: parent_sign(lenof_sign), hdiag
        integer, intent(in) :: nI(nel)
        integer(n_int), intent(in) :: ilut_parent(0:niftot)
        logical, intent(in) :: tCoreDet
        integer :: part, nspawn, ispawn, nI_child(nel), ic, ex(2, maxExcit), unused_ex_level
        integer(n_int) :: ilut_child(0:niftot)!, ilut_parent_init(0:niftot)
        real(dp) :: prob, child_sign(lenof_sign), unused_rdm_real, unused_sign(nel)
        real(dp) :: unused_precond_fac, unused_avEx
        logical :: tParity, break
        integer :: err
        HElement_t(dp) :: HElGen

        unused_ex_level = 0
        do part = 1, lenof_sign
            call decide_num_to_spawn(parent_sign(part), AvMCExcits, nspawn)
            do ispawn = 1, nspawn
                ilut_child = 0_n_int

                call generate_excitation(nI, iLut_parent, nI_child, ilut_child, exFlag, ic, &
                                         ex, tParity, prob, HElGen, fcimc_excit_gen_store)

                if (.not. IsNullDet(nI_child)) then
                    call encode_child(ilut_parent, ilut_child, ic, ex)
                    ilut_child(IlutBits%ind_flag) = 0_n_int

                    ! treating semi-stochastic space
                    ! note that diagonal event either are not in the core space at all
                    ! or from the core space to the core space, i.e. they are covered fully
                    ! by this branching
                    if (tSemiStochastic) then
                        break = check_semistoch_flags(ilut_child, nI_child, &
                                                      part_type_to_run(part), tCoredet)
                        if (break) cycle
                    end if

                    child_sign = attempt_create(nI, iLut_parent, parent_sign, nI_child, iLut_child, &
                                                prob, HElGen, ic, ex, tParity, unused_ex_level, part, &
                                                unused_sign, unused_avEx, unused_rdm_real, unused_precond_fac)
                else
                    child_sign = 0.0_dp
                end if

                if ((any(abs(child_sign) > EPS)) .and. (ic /= 0) .and. (ic <= 2)) then
                    call create_particle_with_hash_table(nI_child, ilut_child, child_sign, &
                                                         part, ilut_parent, iter_data_fciqmc, err)
                    if (err /= 0) return
                end if ! If a child was spawned.

            end do ! Over mulitple particles on same determinant.
        end do
        if (.not. tCoreDet) then
            child_sign = -attempt_die_realtime(hdiag, parent_sign, unused_ex_level)
            if (any(abs(child_sign) > EPS)) then
                ! diagonal events are treated the same way as offdiagonal ones,
                ! except for an extra flag indicating the diagonal spawn
                !ilut_parent_init = ilut_parent
                !call set_flag(ilut_parent_init,get_initiator_flag(part_type))
                ! we have to think about whether diagonal spawns get special treatment
                ! I think they should be subject to the same rules as offdiagonal spawns
                call create_diagonal_with_hashtable(nI, iLut_parent, child_sign, err)
                if (err /= 0) return
            end if
        end if

    end subroutine perform_spawn