first_real_time_spawn Subroutine

public subroutine first_real_time_spawn(err)

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: err

Contents

Source Code


Source Code

    subroutine first_real_time_spawn(err)
        ! routine which first loops over the CurrentDets array and creates the
        ! first spawning list k1 and combines it to y(n) + k1/2
!         integer, intent(inout) :: n_determ_states
        implicit none
        integer, intent(out) :: err
        character(*), parameter :: this_routine = "first_real_time_spawn"
        integer :: idet, parent_flags, nI_parent(nel), unused_flags, ex_level_to_ref, &
                   ireplica, nspawn, ispawn, nI_child(nel), ic, ex(2, maxExcit), &
                   ex_level_to_hf
        integer(n_int) :: ilut_child(0:niftot)
        real(dp) :: parent_sign(lenof_sign), parent_hdiag, prob, child_sign(lenof_sign), &
                    unused_sign(lenof_sign), prefactor, unused_rdm_real
        logical :: tParentIsDeterm, tParentUnoccupied, tParity, break
        HElement_t(dp) :: HelGen, parent_hoffdiag
        integer :: TotWalkersNew, run, determ_index, MaxIndex
        real(dp) :: unused_fac, unused_avEx

        ! declare this is the first runge kutta step
        runge_kutta_step = 1
        ! use part of nicks code, and remove the parts, that dont matter
        prefactor = 1.0_dp
        determ_index = 1

        do idet = 1, int(TotWalkers)

            ! The 'parent' determinant from which spawning is to be attempted.
            parent_flags = 0_n_int

            ! Indicate that the scratch storage used for excitation generation from the
            ! same walker has not been filled (it is filled when we excite from the first
            ! particle on a determinant).
            fcimc_excit_gen_store%tFilled = .false.

            call extract_bit_rep(CurrentDets(:, idet), nI_parent, parent_sign, unused_flags, &
                                 idet, fcimc_excit_gen_store)

            ex_level_to_ref = FindBitExcitLevel(iLutRef, CurrentDets(:, idet), max_calc_ex_level)
            if (tRef_Not_HF) then
                ex_level_to_hf = FindBitExcitLevel(iLutHF_true, CurrentDets(:, idet), &
                                                   max_calc_ex_level)
            else
                ex_level_to_hf = ex_level_to_ref
            end if

            tParentIsDeterm = check_determ_flag(CurrentDets(:, idet), core_run)
            tParentUnoccupied = IsUnoccDet(parent_sign)

            if (tParentIsDeterm) then
                associate(rep => cs_replicas(core_run))
                    rep%indices_of_determ_states(determ_index) = idet
                    rep%partial_determ_vecs(:, determ_index) = parent_sign
                end associate
                determ_index = determ_index + 1
                if (tParentUnoccupied) cycle
            end if

            ! If this slot is unoccupied (and also not a core determinant) then add it to
            ! the list of free slots and cycle.
            if (tParentUnoccupied) then
                iEndFreeSlot = iEndFreeSlot + 1
                FreeSlot(iEndFreeSlot) = idet
                temp_iendfreeslot = temp_iendfreeslot + 1
                temp_freeslot(temp_iendfreeslot) = idet
                cycle
            end if

            ! get new population of observed orbitals
            if (numSnapShotOrbs > 0) call makePopSnapshot(idet)

            ! The current diagonal matrix element is stored persistently.
            parent_hdiag = det_diagH(idet)
            parent_hoffdiag = det_offdiagH(idet)

            if (tTruncInitiator) call CalcParentFlag(idet, parent_flags)

            ! do i need to calc. the energy contributions in the rt-fciqmc?
            ! leave it for now.. and figure out later..
            call SumEContrib(nI_parent, ex_level_to_ref, parent_sign, CurrentDets(:, idet), &
                             parent_hdiag, parent_hoffdiag, 1.0_dp, tPairedReplicas, idet)

            ! If we're on the Hartree-Fock, and all singles and
            ! doubles are in the core space, then there will be
            ! no stochastic spawning from this determinant, so
            ! we can the rest of this loop.
            if (ss_space_in%tDoubles .and. ex_level_to_hf == 0 .and. tDetermHFSpawning) cycle

            ! If this condition is not met (if all electrons have spin up or all have spin down)
            ! then there will be no determinants to spawn to, so don't attempt spawning.
            ! thats a really specific condition.. shouldnt be checked each
            ! cycle.. since this is only input dependent..
            do ireplica = 1, lenof_sign

                call decide_num_to_spawn(parent_sign(ireplica), AvMCExcits, nspawn)
                !call merge_spawn(nspawn,prefactor)
                do ispawn = 1, nspawn

                    ! Zero the bit representation, to ensure no extraneous data gets through.
                    ilut_child = 0_n_int

                    call generate_excitation(nI_parent, CurrentDets(:, idet), nI_child, &
                                             ilut_child, exFlag, ic, ex, tParity, prob, &
                                             HElGen, fcimc_excit_gen_store)

                    ! If a valid excitation.
                    if (.not. IsNullDet(nI_child)) then

                        call encode_child(CurrentDets(:, idet), ilut_child, ic, ex)
                        ilut_child(IlutBits%ind_flag) = 0_n_int

                        if (tSemiStochastic) then
                            break = check_semistoch_flags(ilut_child, nI_child, part_type_to_run(ireplica), tParentIsDeterm)
                            if (break) cycle
                        end if

                        child_sign = attempt_create(nI_parent, CurrentDets(:, idet), parent_sign, &
                                                    nI_child, ilut_child, prob, HElGen, ic, ex, tParity, &
                                                    ex_level_to_ref, ireplica, unused_sign, AvMCExcits, &
                                                    unused_rdm_real, 1.0_dp)
                        child_sign = child_sign * prefactor
                    else
                        child_sign = 0.0_dp
                    end if

                    ! If any (valid) children have been spawned.
                    if ((any(.not. near_zero(child_sign))) .and. (ic /= 0) .and. (ic <= 2)) then

                        ! not quite sure about how to collect the child
                        ! stats in the new rt-fciqmc ..
                        call new_child_stats_normal(iter_data_fciqmc, CurrentDets(:, idet), &
                                             nI_child, ilut_child, ic, ex_level_to_ref, &
                                             child_sign, parent_flags, ireplica)

                        call create_particle_with_hash_table(nI_child, ilut_child, child_sign, &
                                                             ireplica, CurrentDets(:, idet), &
                                                             iter_data_fciqmc, err)

                    end if ! If a child was spawned.

                end do ! Over mulitple particles on same determinant.

            end do ! Over the replicas on the same determinant.

            ! If this is a core-space determinant then the death step is done in
            ! determ_projection.

            ! in the 2nd RK loop of the real-time fciqmc the death step
            ! should not act on the currently looped over walker list y(n)+k1/2
            ! but should be added to the spawned list k2, to then apply it to
            ! the original y(n) + k2
            ! the iEndFreeSlot variable also gets influenced in the
            ! death-step (duh) -> so keep count of the temp iEndFreeSlot above
            if (.not. tParentIsDeterm) then
                call walker_death_realtime(iter_data_fciqmc, nI_parent, &
                                           CurrentDets(:, idet), parent_hdiag, parent_sign, idet, ex_level_to_ref)
            end if

        end do ! Over all determinants.

        if (tSemiStochastic) call real_time_determ_projection()

        ! also update the temp. variables to reuse in y(n) + k2 comb.
        ! this should be done before the annihilaiton step, as there these
        ! values get changed!
        ! have to do that above in the loop as it gets influenced in the
        ! death-step too,

        ! this is the original number of dets.
        TotWalkersNew = int(TotWalkers)

        ! have to call end_iter_stats to get correct acceptance rate
!         call end_iter_stats(TotWalkersNew)
        ! but end iter stats for me is only uses to get SumWalkersCyc ..

        ! the number TotWalkersNew changes below in annihilation routine
        ! Annihilation is done after loop over walkers
        call communicate_and_merge_spawns(MaxIndex, iter_data_fciqmc, .false.)
        call DirectAnnihilation(TotWalkersNew, MaxIndex, iter_data_fciqmc, err)

        TotWalkers = int(TotWalkersNew)

    end subroutine first_real_time_spawn