cipsi

CIPSI algorithm.

The run_stochastic_cipsi() and run_cipsi() subroutines start with a single determinant, or with the wave function in the EZFIO database if determinants read_wf is true.

The run_cipsi() subroutine iteratively:

  • Selects the most important determinants from the external space and adds them to the internal space
  • If determinants s2_eig is true, it adds all the necessary determinants to allow the eigenstates of \(\hat H\) to be eigenstates of \(\widehat{S^2}\)
  • Diagonalizes \(\hat H\) in the enlarged internal space
  • Computes the PT2 contribution to the energy stochastically [13] or deterministically, depending on perturbation do_pt2
  • Extrapolates the variational energy by fitting \(E=E_\text{FCI} - \alpha\, E_\text{PT2}\)

The difference between run_stochastic_cipsi() and run_cipsi() is that run_stochastic_cipsi() selects the determinants on the fly with the computation of the stochastic PT2 [13]. Hence, it is a semi-stochastic selection. It

  • Selects the most important determinants from the external space and adds them to the internal space, on the fly with the computation of the PT2 with the stochastic algorithm presented in [13].
  • If determinants s2_eig is true, it adds all the necessary determinants to allow the eigenstates of \(\hat H\) to be eigenstates of \(\widehat{S^2}\)
  • Extrapolates the variational energy by fitting \(E=E_\text{FCI} - \alpha\, E_\text{PT2}\)
  • Diagonalizes \(\hat H\) in the enlarged internal space

The number of selected determinants at each iteration will be such that the size of the wave function will double at every iteration. If determinants s2_eig is true, then the number of selected determinants will be 1.5x the current number, and then all the additional determinants will be added.

By default, the program will stop when more than one million determinants have been selected, or when the PT2 energy is below \(10^{-4}\).

The variational and PT2 energies of the iterations are stored in the EZFIO database, in the iterations module.

Computation of the PT2 energy

At each iteration, the PT2 energy is computed considering the Epstein-Nesbet zeroth-order Hamiltonian:

\[E_{\text{PT2}} = \sum_{ \alpha } \frac{|\langle \Psi_S | \hat{H} | \alpha \rangle|^2} {E - \langle \alpha | \hat{H} | \alpha \rangle}\]

where the \(|\alpha \rangle\) determinants are generated by applying all the single and double excitation operators to all the determinants of the wave function \(\Psi_G\).

When the hybrid-deterministic/stochastic algorithm is chosen (default), \(Psi_G = \Psi_S = \Psi\), the full wavefunction expanded in the internal space. When the deterministic algorithm is chosen (perturbation do_pt2 is set to false), \(Psi_G\) is a truncation of \(|\Psi \rangle\) using determinants threshold_generators, and \(Psi_S\) is a truncation of \(|\Psi \rangle\) using determinants threshold_selectors, and re-weighted by \(1/\langle \Psi_s | \Psi_s \rangle\).

At every iteration, while computing the PT2, the variance of the wave function is also computed:

\[\begin{split}\sigma^2 & = \langle \Psi | \hat{H}^2 | \Psi \rangle - \langle \Psi | \hat{H} | \Psi \rangle^2 \\ & = \sum_{i \in \text{FCI}} \langle \Psi | \hat{H} | i \rangle \langle i | \hat{H} | \Psi \rangle - \langle \Psi | \hat{H} | \Psi \rangle^2 \\ & = \sum_{ \alpha } \langle |\Psi | \hat{H} | \alpha \rangle|^2.\end{split}\]

The expression of the variance is the same as the expression of the PT2, with a denominator of 1. It measures how far the wave function is from the FCI solution. Note that the absence of denominator in the Heat-Bath selected CI method is selection method by minimization of the variance, whereas CIPSI is a selection method by minimization of the energy.

If perturbation do_pt2 is set to false, then the stochastic PT2 is not computed, and an approximate value is obtained from the CIPSI selection. The calculation is faster, but the extrapolated FCI value is less accurate. This way of running the code should be used when the only goal is to generate a wave function, as for using CIPSI wave functions as trial wave functions of QMC calculations for example.

The PT2 program reads the wave function of the EZFIO database and computes the energy and the PT2 contribution.

State-averaging

Extrapolated FCI energy

An estimate of the FCI energy is computed by extrapolating

\[E=E_\text{FCI} - \alpha\, E_\text{PT2}\]

This extrapolation is done for all the requested states, and excitation energies are printed as energy differences between the extrapolated energies of the excited states and the extrapolated energy of the ground state.

The extrapolations are given considering the 2 last points, the 3 last points, …, the 7 last points. The extrapolated value should be chosen such that the extrpolated value is stable with the number of points.

Providers

initialize_pt2_e0_denominator

File : cipsi/energy.irp.f

logical :: initialize_pt2_e0_denominator

If true, initialize pt2_E0_denominator

Needed by:

nthreads_pt2

File : cipsi/environment.irp.f

integer :: nthreads_pt2

Number of threads for Davidson

Needs:

pt2_collector:()

File : cipsi/pt2_stoch_routines.irp.f

subroutine pt2_collector(zmq_socket_pull, E, relative_error, pt2, error,  &

variance, norm, b, N_)

Needs:

Called by:

  • zmq_pt2()

Calls:

  • add_to_selection_buffer()
  • check_mem()
  • create_selection_buffer()
  • delete_selection_buffer()
  • end_zmq_to_qp_run_socket()
  • pull_pt2_results()
  • sleep()
  • sort_selection_buffer()
  • wall_time()
pt2_cw

File : cipsi/pt2_stoch_routines.irp.f

double precision, allocatable   :: pt2_w        (N_det_generators)
double precision, allocatable   :: pt2_cw       (0:N_det_generators)
double precision        :: pt2_w_t
double precision        :: pt2_u_0
integer, allocatable    :: pt2_n_0      (pt2_N_teeth+1)

Needs:

Needed by:

pt2_e0_denominator

File : cipsi/energy.irp.f

double precision, allocatable   :: pt2_e0_denominator   (N_states)

E0 in the denominator of the PT2

Needs:

pt2_f

File : cipsi/pt2_stoch_routines.irp.f

integer, allocatable    :: pt2_f        (N_det_generators)
integer :: pt2_n_tasks_max

Needs:

pt2_j

File : cipsi/pt2_stoch_routines.irp.f

integer, allocatable    :: pt2_j        (N_det_generators)
integer, allocatable    :: pt2_r        (N_det_generators)

Needs:

pt2_mindetinfirstteeth

File : cipsi/pt2_stoch_routines.irp.f

integer :: pt2_n_teeth
integer :: pt2_mindetinfirstteeth

Needs:

Needed by:

pt2_n_0

File : cipsi/pt2_stoch_routines.irp.f

double precision, allocatable   :: pt2_w        (N_det_generators)
double precision, allocatable   :: pt2_cw       (0:N_det_generators)
double precision        :: pt2_w_t
double precision        :: pt2_u_0
integer, allocatable    :: pt2_n_0      (pt2_N_teeth+1)

Needs:

Needed by:

pt2_n_tasks

File : cipsi/pt2_stoch_routines.irp.f

integer :: pt2_n_tasks

Number of parallel tasks for the Monte Carlo

Needs:

Needed by:

pt2_n_tasks_max

File : cipsi/pt2_stoch_routines.irp.f

integer, allocatable    :: pt2_f        (N_det_generators)
integer :: pt2_n_tasks_max

Needs:

pt2_n_teeth

File : cipsi/pt2_stoch_routines.irp.f

integer :: pt2_n_teeth
integer :: pt2_mindetinfirstteeth

Needs:

Needed by:

pt2_r

File : cipsi/pt2_stoch_routines.irp.f

integer, allocatable    :: pt2_j        (N_det_generators)
integer, allocatable    :: pt2_r        (N_det_generators)

Needs:

pt2_stoch_istate

File : cipsi/pt2_stoch_routines.irp.f

integer :: pt2_stoch_istate

State for stochatsic PT2

Needed by:

pt2_u

File : cipsi/pt2_stoch_routines.irp.f

double precision, allocatable   :: pt2_u        (N_det_generators)

Needs:

Needed by:

pt2_u_0

File : cipsi/pt2_stoch_routines.irp.f

double precision, allocatable   :: pt2_w        (N_det_generators)
double precision, allocatable   :: pt2_cw       (0:N_det_generators)
double precision        :: pt2_w_t
double precision        :: pt2_u_0
integer, allocatable    :: pt2_n_0      (pt2_N_teeth+1)

Needs:

Needed by:

pt2_w

File : cipsi/pt2_stoch_routines.irp.f

double precision, allocatable   :: pt2_w        (N_det_generators)
double precision, allocatable   :: pt2_cw       (0:N_det_generators)
double precision        :: pt2_w_t
double precision        :: pt2_u_0
integer, allocatable    :: pt2_n_0      (pt2_N_teeth+1)

Needs:

Needed by:

pt2_w_t

File : cipsi/pt2_stoch_routines.irp.f

double precision, allocatable   :: pt2_w        (N_det_generators)
double precision, allocatable   :: pt2_cw       (0:N_det_generators)
double precision        :: pt2_w_t
double precision        :: pt2_u_0
integer, allocatable    :: pt2_n_0      (pt2_N_teeth+1)

Needs:

Needed by:

selection_weight

File : cipsi/selection.irp.f

double precision, allocatable   :: selection_weight     (N_states)

Weights used in the selection criterion

Needs:

  • n_states

Subroutines / functions

add_to_selection_buffer:()

File : cipsi/selection_buffer.irp.f

subroutine add_to_selection_buffer(b, det, val)

Needs:

Called by:

  • fill_buffer_double()
  • pt2_collector()
  • selection_collector()

Calls:

  • sort_selection_buffer()
bitstring_to_list_in_selection:()

File : cipsi/selection.irp.f

subroutine bitstring_to_list_in_selection( string, list, n_elements, Nint)

Gives the inidices(+1) of the bits set to 1 in the bit string

Called by:

  • splash_pq()
  • spot_isinwf()
create_selection_buffer:()

File : cipsi/selection_buffer.irp.f

subroutine create_selection_buffer(N, siz_, res)

Needs:

Called by:

  • pt2_collector()
  • run_pt2_slave()
  • run_selection_slave()
  • selection_collector()
  • zmq_pt2()
  • zmq_selection()

Calls:

  • check_mem()
delete_selection_buffer:()

File : cipsi/selection_buffer.irp.f

subroutine delete_selection_buffer(b)

Called by:

  • pt2_collector()
  • run_pt2_slave()
  • run_selection_slave()
  • selection_collector()
  • zmq_pt2()
  • zmq_selection()
fill_buffer_double:()

File : cipsi/selection.irp.f

subroutine fill_buffer_double(i_generator, sp, h1, h2, bannedOrb, banned, fock_diag_tmp, E0, pt2, variance, norm, mat, buf)

Needs:

Called by:

  • select_singles_and_doubles()

Calls:

  • add_to_selection_buffer()
  • apply_holes()
  • apply_particles()
get_d0:()

File : cipsi/selection.irp.f

subroutine get_d0(gen, phasemask, bannedOrb, banned, mat, mask, h, p, sp, coefs)

Needs:

Called by:

  • splash_pq()

Calls:

  • apply_particles()
  • get_mo_two_e_integrals()
  • i_h_j()
get_d1:()

File : cipsi/selection.irp.f

subroutine get_d1(gen, phasemask, bannedOrb, banned, mat, mask, h, p, sp, coefs)

Needs:

Called by:

  • splash_pq()

Calls:

  • apply_particles()
  • get_mo_two_e_integrals()
  • i_h_j()
get_d2:()

File : cipsi/selection.irp.f

subroutine get_d2(gen, phasemask, bannedOrb, banned, mat, mask, h, p, sp, coefs)

Needs:

  • n_states

Called by:

  • splash_pq()
get_mask_phase:()

File : cipsi/selection.irp.f

subroutine get_mask_phase(det1, pm, Nint)

Called by:

  • splash_pq()
get_phase_bi:()

File : cipsi/selection.irp.f

double precision function get_phase_bi(phasemask, s1, s2, h1, p1, h2, p2, Nint)
make_selection_buffer_s2:()

File : cipsi/selection_buffer.irp.f

subroutine make_selection_buffer_s2(b)

Needs:

  • elec_alpha_num

Called by:

  • zmq_pt2()
  • zmq_selection()

Calls:

  • check_mem()
  • dsort()
  • i8sort()
  • occ_pattern_to_dets()
  • occ_pattern_to_dets_size()
merge_selection_buffers:()

File : cipsi/selection_buffer.irp.f

subroutine merge_selection_buffers(b1, b2)

Merges the selection buffers b1 and b2 into b2

Needs:

Called by:

  • run_pt2_slave()
  • run_selection_slave()

Calls:

  • check_mem()
past_d1:()

File : cipsi/selection.irp.f

subroutine past_d1(bannedOrb, p)

Needs:

Called by:

  • splash_pq()
past_d2:()

File : cipsi/selection.irp.f

subroutine past_d2(banned, p, sp)

Needs:

Called by:

  • splash_pq()
provide_everything:()

File : cipsi/slave_cipsi.irp.f

subroutine provide_everything

Needs:

Called by:

  • run_slave_cipsi()
pt2_find_sample:()

File : cipsi/pt2_stoch_routines.irp.f

integer function pt2_find_sample(v, w)

Needs:

pt2_find_sample_lr:()

File : cipsi/pt2_stoch_routines.irp.f

integer function pt2_find_sample_lr(v, w, l_in, r_in)

Needs:

pt2_slave_inproc:()

File : cipsi/pt2_stoch_routines.irp.f

subroutine pt2_slave_inproc(i)

Needs:

Called by:

  • zmq_pt2()

Calls:

  • run_pt2_slave()
pull_pt2_results:()

File : cipsi/run_pt2_slave.irp.f

subroutine pull_pt2_results(zmq_socket_pull, index, pt2, variance, norm, task_id, n_tasks, b)

Needs:

  • n_states

Called by:

  • pt2_collector()
pull_selection_results:()

File : cipsi/run_selection_slave.irp.f

subroutine pull_selection_results(zmq_socket_pull, pt2, variance, norm, val, det, N, task_id, ntask)

Needs:

  • n_states

Called by:

  • selection_collector()
push_pt2_results:()

File : cipsi/run_pt2_slave.irp.f

subroutine push_pt2_results(zmq_socket_push, index, pt2, variance, norm, b, task_id, n_tasks)

Needs:

  • n_states

Called by:

  • run_pt2_slave()
push_selection_results:()

File : cipsi/run_selection_slave.irp.f

subroutine push_selection_results(zmq_socket_push, pt2, variance, norm, b, task_id, ntask)

Needs:

  • n_states

Called by:

  • run_selection_slave()
remove_duplicates_in_selection_buffer:()

File : cipsi/selection_buffer.irp.f

subroutine remove_duplicates_in_selection_buffer(b)

Needs:

Called by:

  • zmq_pt2()

Calls:

  • check_mem()
  • i8sort()
run_cipsi:()

File : cipsi/cipsi.irp.f

subroutine run_cipsi

Selected Full Configuration Interaction with deterministic selection and stochastic PT2.

Needs:

Called by:

  • fci()

Calls:

  • check_mem()
  • diagonalize_ci()
  • ezfio_get_hartree_fock_energy()
  • ezfio_has_hartree_fock_energy()
  • make_s2_eigenfunction()
  • print_extrapolated_energy()
  • print_summary()
  • save_energy()
  • save_iterations()
  • save_wavefunction()
  • write_double()
  • zmq_pt2()
  • zmq_selection()

Touches:

run_pt2_slave:()

File : cipsi/run_pt2_slave.irp.f

subroutine run_pt2_slave(thread,iproc,energy)

Needs:

Called by:

  • pt2_slave_inproc()
  • run_slave_main()

Calls:

  • check_mem()
  • create_selection_buffer()
  • delete_selection_buffer()
  • end_zmq_push_socket()
  • end_zmq_to_qp_run_socket()
  • merge_selection_buffers()
  • push_pt2_results()
  • select_connected()
  • sleep()
  • sort_selection_buffer()
  • wall_time()
run_selection_slave:()

File : cipsi/run_selection_slave.irp.f

subroutine run_selection_slave(thread,iproc,energy)

Needs:

Called by:

  • run_slave_main()
  • selection_slave_inproc()

Calls:

  • create_selection_buffer()
  • delete_selection_buffer()
  • end_zmq_push_socket()
  • end_zmq_to_qp_run_socket()
  • merge_selection_buffers()
  • push_selection_results()
  • select_connected()
  • sleep()
  • sort_selection_buffer()
run_slave_cipsi:()

File : cipsi/slave_cipsi.irp.f

subroutine run_slave_cipsi

Helper program for distributed parallelism

Needs:

  • read_wf
  • distributed_davidson

Called by:

  • fci()
  • pt2()

Calls:

  • omp_set_nested()
  • provide_everything()
  • run_slave_main()
  • switch_qp_run_to_master()

Touches:

run_slave_main:()

File : cipsi/slave_cipsi.irp.f

subroutine run_slave_main

Needs:

Called by:

  • run_slave_cipsi()

Calls:

  • davidson_slave_tcp()
  • mpi_print()
  • omp_set_nested()
  • run_pt2_slave()
  • run_selection_slave()
  • sleep()
  • wait_for_states()
  • wall_time()
  • write_double()

Touches:

  • threshold_generators
run_stochastic_cipsi:()

File : cipsi/stochastic_cipsi.irp.f

subroutine run_stochastic_cipsi

Selected Full Configuration Interaction with Stochastic selection and PT2.

Needs:

Called by:

  • fci()

Calls:

  • check_mem()
  • copy_h_apply_buffer_to_wf()
  • diagonalize_ci()
  • ezfio_get_hartree_fock_energy()
  • ezfio_has_hartree_fock_energy()
  • make_s2_eigenfunction()
  • print_extrapolated_energy()
  • print_summary()
  • save_energy()
  • save_iterations()
  • save_wavefunction()
  • write_double()
  • zmq_pt2()

Touches:

select_connected:()

File : cipsi/selection.irp.f

subroutine select_connected(i_generator,E0,pt2,variance,norm,b,subset,csubset)

Needs:

Called by:

  • run_pt2_slave()
  • run_selection_slave()

Calls:

  • build_fock_tmp()
  • select_singles_and_doubles()
select_singles_and_doubles:()

File : cipsi/selection.irp.f

subroutine select_singles_and_doubles(i_generator,hole_mask,particle_mask,fock_diag_tmp,E0,pt2,variance,norm,buf,subset,csubset)

WARNING /!: It is assumed that the generators and selectors are psi_det_sorted

Needs:

Called by:

  • select_connected()

Calls:

  • apply_hole()
  • bitstring_to_list_ab()
  • check_mem()
  • fill_buffer_double()
  • get_excitation_degree_spin()
  • isort()
  • splash_pq()
  • spot_isinwf()
selection_collector:()

File : cipsi/zmq_selection.irp.f

subroutine selection_collector(zmq_socket_pull, b, N, pt2, variance, norm)

Needs:

  • n_states

Called by:

  • zmq_selection()

Calls:

  • add_to_selection_buffer()
  • check_mem()
  • create_selection_buffer()
  • delete_selection_buffer()
  • end_zmq_to_qp_run_socket()
  • pull_selection_results()
  • sort_selection_buffer()
selection_slave_inproc:()

File : cipsi/zmq_selection.irp.f

subroutine selection_slave_inproc(i)

Needs:

Called by:

  • zmq_selection()

Calls:

  • run_selection_slave()
sort_selection_buffer:()

File : cipsi/selection_buffer.irp.f

subroutine sort_selection_buffer(b)

Needs:

Called by:

  • add_to_selection_buffer()
  • pt2_collector()
  • run_pt2_slave()
  • run_selection_slave()
  • selection_collector()

Calls:

  • check_mem()
  • dsort()
splash_pq:()

File : cipsi/selection.irp.f

subroutine splash_pq(mask, sp, det, i_gen, N_sel, bannedOrb, banned, mat, interesting)

Needs:

Called by:

  • select_singles_and_doubles()

Calls:

  • bitstring_to_list_in_selection()
  • get_d0()
  • get_d1()
  • get_d2()
  • get_mask_phase()
  • past_d1()
  • past_d2()
spot_isinwf:()

File : cipsi/selection.irp.f

subroutine spot_isinwf(mask, det, i_gen, N, banned, fullMatch, interesting)

Needs:

Called by:

  • select_singles_and_doubles()

Calls:

  • bitstring_to_list_in_selection()
testteethbuilding:()

File : cipsi/pt2_stoch_routines.irp.f

logical function testTeethBuilding(minF, N)

Needs:

Calls:

  • check_mem()
zmq_pt2:()

File : cipsi/pt2_stoch_routines.irp.f

subroutine ZMQ_pt2(E, pt2,relative_error, error, variance, norm, N_in)

Needs:

Called by:

  • run_cipsi()
  • run_stochastic_cipsi()

Calls:

  • check_mem()
  • create_selection_buffer()
  • delete_selection_buffer()
  • end_parallel_job()
  • fill_h_apply_buffer_no_selection()
  • make_selection_buffer_s2()
  • new_parallel_job()
  • omp_set_nested()
  • pt2_collector()
  • pt2_slave_inproc()
  • remove_duplicates_in_selection_buffer()
  • resident_memory()
  • write_double()
  • write_int()
  • zmq_selection()

Touches:

zmq_selection:()

File : cipsi/zmq_selection.irp.f

subroutine ZMQ_selection(N_in, pt2, variance, norm)

Needs:

Called by:

  • run_cipsi()
  • zmq_pt2()

Calls:

  • copy_h_apply_buffer_to_wf()
  • create_selection_buffer()
  • delete_selection_buffer()
  • end_parallel_job()
  • fill_h_apply_buffer_no_selection()
  • make_selection_buffer_s2()
  • new_parallel_job()
  • save_wavefunction()
  • selection_collector()
  • selection_slave_inproc()
  • write_double()

Touches: