InitLanczosCalc Subroutine

public subroutine InitLanczosCalc(this, det_list, print_info, hamil_type, n_states, max_lanczos_vecs, t_store_subspace_basis, t_orthogonalise, max_restarts, energy_precision, ritz_overlap_precision)

Arguments

Type IntentOptional Attributes Name
type(LanczosCalcType), intent(out) :: this
integer, intent(in) :: det_list(:,:)
logical, intent(in) :: print_info
integer, intent(in) :: hamil_type
integer, intent(in) :: n_states
integer, intent(in) :: max_lanczos_vecs
logical, intent(in) :: t_store_subspace_basis
logical, intent(in) :: t_orthogonalise
integer, intent(in) :: max_restarts
integer, intent(in) :: energy_precision
integer, intent(in) :: ritz_overlap_precision

Contents

Source Code


Source Code

    subroutine InitLanczosCalc(this, det_list, print_info, hamil_type, n_states, max_lanczos_vecs, &
                               t_store_subspace_basis, t_orthogonalise, max_restarts, energy_precision, ritz_overlap_precision)
        type(LanczosCalcType), intent(out) :: this
        integer, intent(in) :: det_list(:, :), hamil_type, n_states, max_lanczos_vecs, max_restarts, &
                               energy_precision, ritz_overlap_precision
        logical, intent(in) :: t_store_subspace_basis, t_orthogonalise
        character(len=*), parameter :: t_r = "init_lanczos"
        logical, intent(in) :: print_info

        integer :: i, HFindex, mem_reqd, residual_mem_reqd, ierr
        integer(MPIArg) :: mpi_temp
        real(dp), allocatable :: hamil_diag_temp(:)
        real(dp), allocatable :: lowest_energies(:)
        integer, allocatable :: lowest_energy_det_indices(:)

        call InitHamiltonianCalc(this%super, print_info, hamil_type, max_lanczos_vecs, t_store_subspace_basis, t_orthogonalise)
        associate ( &
            space_size => this%super%space_size, &
            max_subspace_size => this%super%max_subspace_size &
            )

            if (n_states > space_size) then
                call stop_all(t_r, "Not enough determinants in the space to produce the required number of approximate eigenpairs")
            end if

#if(0)
            ! TODO: given a target multiplicity for the many body ground state, an initial
            ! lanczos vector guess is chosen to have a good overlap with each of the
            ! ground_state_multiplicity states in the ground state WF

            ! verify the given target multiplicity for the ground state
            if (mod(nel * ground_state_multiplicity) == 0) then
                ! odd target multiplicities correspond to singlet, triplet, ... states
                ! these are only possible in closed shell systems
                ! even target multiplicities correspond to doublet, quadruplet, ... states
                ! these are only possible in open shell systems
                call stop_all(t_r, "Invalid target multiplicites for Lanczos initialisation")
            end if
#endif

            this%n_states = n_states
            this%max_restarts = max_restarts
            this%convergence_error = 10**(-real(energy_precision, dp))
            this%orthog_tolerance = 10**(-real(ritz_overlap_precision, dp))
            ! these will be larger than neccessary until the final iteration:
            safe_malloc(this%ritz_values, (max_subspace_size))
            safe_calloc(this%eigenvalues, (n_states), 0.0_dp)
            safe_calloc(this%t_states_converged, (n_states), .false.)
            safe_malloc(this%ritz_values_old, (n_states))
            safe_malloc(this%T_eigenvectors, (max_subspace_size, max_subspace_size))

            safe_calloc(this%lanczos_vector, (space_size), 0.0_dp)

            if (.not. t_store_subspace_basis) then
                safe_malloc(this%first_v, (max_subspace_size))
                safe_malloc(this%current_v, (max_subspace_size))
                safe_malloc(this%old_v, (max_subspace_size))
            end if

            ! check for previous allocation of the eigenvector estimates
            if (allocated(this%ritz_vectors)) then
                deallocate(this%ritz_vectors, stat=ierr)
                call logmemdealloc(t_r, lanczosTag, ierr)
            end if

            ! nstates columns, holding approx eigenvectors of length space_size
            safe_calloc_e(this%ritz_vectors, (space_size, n_states), 0.0_dp, ierr)
            call logmemalloc("ritz_vectors", space_size * n_states, 8, t_r, lanczosTag, ierr)

            safe_calloc_e(this%eigenvectors, (space_size, n_states), 0.0_dp, ierr)
            call logmemalloc("eigenvectors", space_size * n_states, 8, t_r, lanczosTag, ierr)

            if (iProcIndex == root) then
                safe_malloc(hamil_diag_temp, (size(hamil_diag)))
                hamil_diag_temp(:) = -hamil_diag(:)
                safe_calloc(lowest_energies, (n_states), 0.0_dp)
                safe_calloc(lowest_energy_det_indices, (n_states), 0)
                lowest_energy_det_indices(1) = maxloc(hamil_diag_temp, dim=1)
                lowest_energies(1) = hamil_diag(lowest_energy_det_indices(1))
                do i = 2, n_states
                    lowest_energy_det_indices(i) = maxloc(hamil_diag_temp, dim=1, mask=hamil_diag_temp < -lowest_energies(i - 1))
                    lowest_energies(i) = hamil_diag(lowest_energy_det_indices(i))
                end do
                safe_free(hamil_diag_temp)
            end if

            ! If there is only one determinant per state in the space being diagonalised:
            if (space_size == n_states) then
                if (iProcIndex == root) then
                    this%eigenvalues(:) = lowest_energies(:)
                    do i = 1, n_states
                        this%eigenvectors(lowest_energy_det_indices(i), i) = 1.0_dp
                    end do
                end if
                this%super%skip_calc = .true.
                return
            end if

            if (iProcIndex == root) then
                ! with the memory management done, now we can start setting up the pre-iteration
                ! numerics. Set the first basis vector to be an equal, normalised superposition
                ! of the nstates lowest lying determinants

                do i = 1, n_states
                    write(stdout, *) det_list(:, lowest_energy_det_indices(i))
                    this%lanczos_vector(lowest_energy_det_indices(i)) = 1.0_dp / sqrt(real(n_states, dp))
                end do

                this%lanczos_vector(:) = this%lanczos_vector(:) / euclidean_norm(this%lanczos_vector(:))

                safe_free(lowest_energies)
                safe_free(lowest_energy_det_indices)
            end if
        end associate

    end subroutine InitlanczosCalc