kp_fciqmc_read_inp Subroutine

public subroutine kp_fciqmc_read_inp(file_reader, kp)

Arguments

Type IntentOptional Attributes Name
class(FileReader_t), intent(inout) :: file_reader
type(kp_fciqmc_data), intent(inout) :: kp

Contents

Source Code


Source Code

    subroutine kp_fciqmc_read_inp(file_reader, kp)

        use CalcData, only: tWritePopsNorm, iPopsFileNoRead, tPairedReplicas
        use FciMCData, only: alloc_popsfile_dets
        use LoggingData, only: tIncrementPops
        use perturbations, only: init_perturbation_creation, init_perturbation_annihilation
        use SystemData, only: tAllSymSectors, tSpn

        type(kp_fciqmc_data), intent(inout) :: kp
        class(FileReader_t), intent(inout) :: file_reader
        logical :: eof
        type(TokenIterator_t) :: tokens
        character(len=100) :: w
        integer :: i, j, niters_temp, nvecs_temp, nreports_temp, npert, nexcit
        integer :: ras_size_1, ras_size_2, ras_size_3, ras_min_1, ras_max_3
        character(len=*), parameter :: t_r = "kp_fciqmc_read_inp", this_routine = t_r

        ! Default values.
        tExcitedStateKP = .false.
        kp%nconfigs = 1
        kp%nreports = 0
        kp%nrepeats = 1
        kp%nvecs = 0
        niters_temp = 1
        memory_factor_kp = 1.0_dp

        tFiniteTemp = .false.
        tExcitedInitState = .false.

        tCalcSpin = .false.

        nwalkers_per_site_init = 1.0_dp
        av_mc_excits_kp = 0.0_dp
        kp_hamil_exact_frac = 1.0_dp

        tMultiplePopStart = .false.
        tExactHamil = .false.
        tExactHamilSpawning = .false.
        tFullyStochasticHamil = .false.
        tInitCorrectNWalkers = .false.
        tOccDetermInit = .false.
        vary_niters = .false.
        tUseInitConfigSeeds = .false.
        tAllSymSectors = .false.
        tOutputAverageKPMatrices = .false.
        tOverlapPert = .false.
        tScalePopulation = .false.
        scaling_factor = 1.0_dp

        tPairedReplicas = .true.
#if defined(PROG_NUMRUNS_)
        nreplicas = 2
#endif
        tOrthogKPReplicas = .false.
        orthog_kp_iter = 0

        read_inp: do while (file_reader%nextline(tokens, skip_empty=.true.))
            w = to_upper(tokens%next())
            select case (w)
            case ("END-KP-FCIQMC")
                exit read_inp
            case ("EXCITED-STATE-KP")
#ifdef PROG_NUMRUNS_
                tExcitedStateKP = .true.
                kp%nvecs = to_int(tokens%next())
#else
                call stop_all(t_r, "NECI must be compiled with multiple replicas (mneci) to use the excited-state-kp option.")
#endif
            case ("FINITE-TEMPERATURE")
                tFiniteTemp = .true.

            case ("CALC-SPIN")
                tCalcSpin = .true.

            case ("MULTIPLE-POPS")
                tMultiplePopStart = .true.
                tIncrementPops = .true.
                iPopsFileNoRead = -1
            case ("NUM-INIT-CONFIGS")
                kp%nconfigs = to_int(tokens%next())
            case ("NUM-REPEATS-PER-INIT-CONFIG")
                kp%nrepeats = to_int(tokens%next())
            case ("NUM-KRYLOV-VECS")
                nvecs_temp = to_int(tokens%next())
                if (kp%nvecs /= 0) then
                    if (kp%nvecs /= nvecs_temp) call stop_all(t_r, 'The number of values specified for the number of iterations &
                        &between Krylov vectors is not consistent with the number of Krylov vectors requested.')
                else
                    kp%nvecs = nvecs_temp
                    allocate(kp%niters(kp%nvecs))
                    kp%niters(kp%nvecs) = 0
                end if
            case ("NREPORTS")
                nreports_temp = to_int(tokens%next())
                if (kp%nreports /= 0) then
                    if (kp%nreports /= nreports_temp) call stop_all(t_r, 'The number of values specified for the number of &
                        &iterations between reports is not consistent with the number of reports requested.')
                else
                    kp%nreports = nreports_temp
                    allocate(kp%niters(kp%nreports))
                    kp%niters(kp%nreports) = 0
                end if
            case ("NUM-ITERS-BETWEEN-VECS")
                niters_temp = to_int(tokens%next())
            case ("NUM-ITERS-BETWEEN-REPORTS")
                niters_temp = to_int(tokens%next())
            case ("NUM-ITERS-BETWEEN-VECS-VARY")
                vary_niters = .true.
                if (kp%nvecs /= 0) then
                    if (kp%nvecs /= tokens%size()) call stop_all(t_r, 'The number of values specified for the number of iterations &
                        &between Krylov vectors is not consistent with the number of Krylov vectors requested.')
                end if
                kp%nvecs = tokens%size()
                if (.not. allocated(kp%niters)) allocate(kp%niters(kp%nvecs))
                do i = 1, kp%nvecs - 1
                    kp%niters(i) = to_int(tokens%next())
                end do
                kp%niters(kp%nvecs) = 0
            case ("NUM-ITERS-BETWEEN-REPORTS-VARY")
                vary_niters = .true.
                if (kp%nreports /= 0) then
                    if (kp%nreports /= tokens%size()) call stop_all(t_r, 'The number of values specified for the number of &
                        &iterations between reports is not consistent with the number of reports requested.')
                end if
                kp%nreports = tokens%size()
                if (.not. allocated(kp%niters)) allocate(kp%niters(kp%nreports))
                do i = 1, kp%nreports - 1
                    kp%niters(i) = to_int(tokens%next())
                end do
                kp%niters(kp%nreports) = 0
            case ("MEMORY-FACTOR")
                memory_factor_kp = to_realdp(tokens%next())
            case ("NUM-WALKERS-PER-SITE-INIT")
                nwalkers_per_site_init = to_realdp(tokens%next())
            case ("AVERAGEMCEXCITS-HAMIL")
                av_mc_excits_kp = to_realdp(tokens%next())
            case ("EXACT-HAMIL-SPAWNING")
                tExactHamilSpawning = .true.
            case ("EXACT-HAMIL-SPAWNING-FRAC")
                kp_hamil_exact_frac = to_realdp(tokens%next())
            case ("EXACT-HAMIL")
                tExactHamil = .true.
            case ("FULLY-STOCHASTIC-HAMIL")
                tFullyStochasticHamil = .true.
            case ("INIT-CORRECT-WALKER-POP")
                tInitCorrectNWalkers = .true.
            case ("OCCUPY-DETERM-INIT")
                tOccDetermInit = .true.
            case ("INIT-CONFIG-SEEDS")
                if (kp%nconfigs == 0) call stop_all(t_r, 'Please use the num-init-configs option to enter the number of &
                                                          &initial configurations before using the init-vec-seeds option.')
                tUseInitConfigSeeds = .true.
                allocate(init_config_seeds(kp%nconfigs))
                do i = 1, kp%nconfigs
                    init_config_seeds(i) = to_int(tokens%next())
                end do
            case ("ALL-SYM-SECTORS")
                tAllSymSectors = .true.
            case ("WRITE-AVERAGE-MATRICES")
                tOutputAverageKPMatrices = .true.
            case ("SCALE-POPULATION")
                tScalePopulation = .true.
            case ("EXCITED-INIT-STATE")
                tExcitedInitState = .true.
                ! Read in the number of excited states to be used.
                nexcit = to_int(tokens%next())
                allocate(kpfciqmc_ex_labels(nexcit))
                allocate(kpfciqmc_ex_weights(nexcit))

                ! Read in which excited states to use.
                if (file_reader%nextline(tokens, skip_empty=.false.)) then
                    do j = 1, nexcit
                        kpfciqmc_ex_labels(j) = to_int(tokens%next())
                    end do
                else
                    call stop_all(this_routine, 'Unexpected EOF reached.')
                end if

                ! Read in the relative weights for the trial excited states in
                ! the initial state.
                if (file_reader%nextline(tokens, skip_empty=.false.)) then
                    do j = 1, nexcit
                        kpfciqmc_ex_weights(j) = to_realdp(tokens%next())
                    end do
                else
                    call stop_all(this_routine, 'Unexpected EOF reached.')
                end if

                ! Normalise weights so that they sum to 1.
                kpfciqmc_ex_weights = kpfciqmc_ex_weights / sum(kpfciqmc_ex_weights)

            case ("OVERLAP-PERTURB-ANNIHILATE")
                alloc_popsfile_dets = .true.
                tOverlapPert = .true.
                tWritePopsNorm = .true.

                ! Read in the number of perturbation operators which are about
                ! to be read in.
                npert = to_int(tokens%next())
                if (.not. allocated(overlap_pert)) then
                    allocate(overlap_pert(npert))
                else
                    if (npert /= size(overlap_pert)) then
                        call stop_all(t_r, "A different number of creation and annihilation perturbation have been requested.")
                    end if
                end if

                do i = 1, npert
                    if (file_reader%nextline(tokens, skip_empty=.false.)) then
                        overlap_pert(i)%nannihilate = tokens%remaining_items()
                        allocate(overlap_pert(i)%ann_orbs(tokens%remaining_items()))
                        do j = 1, size(overlap_pert(i)%ann_orbs)
                            overlap_pert(i)%ann_orbs(j) = to_int(tokens%next())
                        end do
                        ! Create the rest of the annihilation-related
                        ! components of the pops_pert object.
                        call init_perturbation_annihilation(overlap_pert(i))
                    else
                        call stop_all(this_routine, 'Unexpected EOF reached.')
                    end if
                end do

            case ("OVERLAP-PERTURB-CREATION")
                alloc_popsfile_dets = .true.
                tOverlapPert = .true.
                tWritePopsNorm = .true.

                ! Read in the number of perturbation operators which are about
                ! to be read in.
                npert = to_int(tokens%next())
                if (.not. allocated(overlap_pert)) then
                    allocate(overlap_pert(npert))
                else
                    if (npert /= size(overlap_pert)) then
                        call stop_all(t_r, "A different number of creation and annihilation perturbation have been requested.")
                    end if
                end if

                do i = 1, npert
                    if (file_reader%nextline(tokens, .false.)) then
                        overlap_pert(i)%ncreate = tokens%remaining_items()
                        allocate(overlap_pert(i)%crtn_orbs(tokens%remaining_items()))
                        do j = 1, size(overlap_pert(i)%crtn_orbs)
                            overlap_pert(i)%crtn_orbs(j) = to_int(tokens%next())
                        end do
                        ! Create the rest of the creation-related
                        ! components of the pops_pert object.
                        call init_perturbation_creation(overlap_pert(i))
                    else
                        call stop_all(this_routine, 'Unexpected EOF reached.')
                    end if
                end do

            case ("DOUBLES-TRIAL")
                kp_trial_space_in%tDoubles = .true.
            case ("CAS-TRIAL")
                kp_trial_space_in%tCAS = .true.
                tSpn = .true.
                kp_trial_space_in%occ_cas = to_int(tokens%next())  ! Number of electrons in the CAS
                kp_trial_space_in%virt_cas = to_int(tokens%next()) ! Number of virtual spin-orbitals in the CAS
            case ("RAS-TRIAL")
                kp_trial_space_in%tRAS = .true.
                ras_size_1 = to_int(tokens%next())  ! Number of spatial orbitals in RAS1.
                ras_size_2 = to_int(tokens%next())  ! Number of spatial orbitals in RAS2.
                ras_size_3 = to_int(tokens%next())  ! Number of spatial orbitals in RAS3.
                ras_min_1 = to_int(tokens%next()) ! Min number of electrons (alpha and beta) in RAS1 orbs.
                ras_max_3 = to_int(tokens%next()) ! Max number of electrons (alpha and beta) in RAS3 orbs.
                kp_trial_space_in%ras%size_1 = int(ras_size_1, sp)
                kp_trial_space_in%ras%size_2 = int(ras_size_2, sp)
                kp_trial_space_in%ras%size_3 = int(ras_size_3, sp)
                kp_trial_space_in%ras%min_1 = int(ras_min_1, sp)
                kp_trial_space_in%ras%max_3 = int(ras_max_3, sp)
            case ("MP1-TRIAL")
                kp_trial_space_in%tMP1 = .true.
                kp_trial_space_in%mp1_ndets = to_int(tokens%next())
            case ("HF-TRIAL")
                kp_trial_space_in%tHF = .true.
            case ("POPS-TRIAL")
                kp_trial_space_in%tPops = .true.
                kp_trial_space_in%npops = to_int(tokens%next())
            case ("READ-TRIAL")
                kp_trial_space_in%tRead = .true.
            case ("FCI-TRIAL")
                kp_trial_space_in%tFCI = .false.
            case ("UNPAIRED-REPLICAS")
                tPairedReplicas = .false.
#if defined(PROG_NUMRUNS_)
                nreplicas = 1
#endif
            case ("ORTHOGONALISE-REPLICAS")
                tOrthogKPReplicas = .true.
                if (tokens%remaining_items() > 0) then
                    orthog_kp_iter = to_int(tokens%next())
                end if
            case default
                call stop_all(this_routine, "Keyword "//trim(w)//" not recognized in kp-fciqmc block")
            end select
        end do read_inp

        if ((.not. vary_niters)) then
            if (tExcitedStateKP) then
                if (.not. allocated(kp%niters)) allocate(kp%niters(kp%nreports))
                kp%niters = niters_temp
                kp%niters(kp%nreports) = 0
            else
                if (.not. allocated(kp%niters)) allocate(kp%niters(kp%nvecs))
                kp%niters = niters_temp
                kp%niters(kp%nvecs) = 0
            end if
        end if

    end subroutine kp_fciqmc_read_inp