This commit is contained in:
Anthony Scemama 2019-02-22 00:18:32 +01:00
parent 47691d0a0a
commit f2794fa54e
10 changed files with 206 additions and 1 deletions

159
CI/CI.ml Normal file
View File

@ -0,0 +1,159 @@
open Lacaml.D
module De = Determinant
module Ex = Excitation
module Sp = Spindeterminant
type t =
{
det_space : Determinant_space.t ;
h_matrix : Mat.t lazy_t ;
eigensystem : (Mat.t * Vec.t) lazy_t;
}
let det_space t = t.det_space
let h_matrix t = Lazy.force t.h_matrix
let eigensystem t = Lazy.force t.eigensystem
let eigenvectors t =
let (x,_) = eigensystem t in x
let eigenvalues t =
let (_,x) = eigensystem t in x
let h_ij mo_basis ki kj =
let one_e_ints = MOBasis.one_e_ints mo_basis
and two_e_ints = MOBasis.two_e_ints mo_basis
and degree_a, degree_b = De.degrees ki kj
in
if degree_a+degree_b > 2 then
0.
else
begin
let kia = De.alfa ki and kib = De.beta ki
and kja = De.alfa kj and kjb = De.beta kj
in
let integral =
ERI.get_phys two_e_ints
in
let anti_integral i j k l =
integral i j k l -. integral i j l k
in
let single h p same opposite =
let same_spin =
Sp.to_list same
|> List.fold_left (fun accu i -> accu +. anti_integral h i p i) 0.
and opposite_spin =
Sp.to_list opposite
|> List.fold_left (fun accu i -> accu +. integral h i p i) 0.
and one_e =
one_e_ints.{h,p}
in one_e +. same_spin +. opposite_spin
in
let diag_element () =
let mo_a = Sp.to_list kia
and mo_b = Sp.to_list kib
in
let one_e =
List.concat [mo_a ; mo_b]
|> List.fold_left (fun accu i -> accu +. one_e_ints.{i,i}) 0.
and two_e =
let two_index i j = integral i j i j
and anti_two_index i j = anti_integral i j i j
in
let rec aux_same accu = function
| [] -> accu
| i :: rest ->
let new_accu =
List.fold_left (fun accu j -> accu +. anti_two_index i j) accu rest
in
aux_same new_accu rest
in
let rec aux_opposite accu other = function
| [] -> accu
| i :: rest ->
let new_accu =
List.fold_left (fun accu j -> accu +. two_index i j) accu other
in
aux_opposite new_accu other rest
in
(aux_same 0. mo_a) +. (aux_same 0. mo_b) +. (aux_opposite 0. mo_a mo_b)
in
one_e +. two_e
in
match degree_a, degree_b with
| 1, 1 -> (* alpha-beta double *)
begin
let ha, pa, phase_a = Ex.single_of_spindet kia kja in
let hb, pb, phase_b = Ex.single_of_spindet kib kjb in
match phase_a, phase_b with
| Phase.Pos, Phase.Pos
| Phase.Neg, Phase.Neg -> +. integral ha hb pa pb
| Phase.Neg, Phase.Pos
| Phase.Pos, Phase.Neg -> -. integral ha hb pa pb
end
| 2, 0 -> (* alpha double *)
begin
let h1, p1, h2, p2, phase = Ex.double_of_spindet kia kja in
match phase with
| Phase.Pos -> +. anti_integral h1 h2 p1 p2
| Phase.Neg -> -. anti_integral h1 h2 p1 p2
end
| 0, 2 -> (* beta double *)
begin
let h1, p1, h2, p2, phase = Ex.double_of_spindet kib kjb in
match phase with
| Phase.Pos -> +. anti_integral h1 h2 p1 p2
| Phase.Neg -> -. anti_integral h1 h2 p1 p2
end
| 1, 0 -> (* alpha single *)
begin
let h, p, phase = Ex.single_of_spindet kia kja in
match phase with
| Phase.Pos -> +. single h p kia kib
| Phase.Neg -> -. single h p kia kib
end
| 0, 1 -> (* beta single *)
begin
let h, p, phase = Ex.single_of_spindet kib kjb in
match phase with
| Phase.Pos -> +. single h p kib kia
| Phase.Neg -> -. single h p kib kia
end
| 0, 0 -> (* diagonal element *)
diag_element ()
| _ -> assert false
end
let make det_space =
let ndet = Determinant_space.size det_space in
let det = Determinant_space.determinants det_space in
let mo_basis = Determinant_space.mo_basis det_space in
let h_matrix = lazy (
Array.init ndet (fun i ->
let ki = det.(i) in
Array.init ndet (fun j ->
let kj = det.(j) in
h_ij mo_basis ki kj
)
)
|> Mat.of_array
)
in
let eigensystem = lazy (
Lazy.force h_matrix
|> Util.diagonalize_symm
)
in
{ det_space ; h_matrix ; eigensystem }

View File

@ -41,6 +41,15 @@ let negate_phase t =
{ t with alfa = Spindeterminant.negate_phase t.alfa }
let degrees t t' =
Spindeterminant.degree t.alfa t'.alfa,
Spindeterminant.degree t.beta t'.beta
let degree t t' =
let a, b = degrees t t' in a+b
let of_lists a b =
let alfa = Spindeterminant.of_list a
and beta = Spindeterminant.of_list b

View File

@ -44,6 +44,13 @@ val single_excitation : Spin.t -> hole -> particle -> t -> t
val double_excitation : Spin.t -> hole -> particle -> Spin.t -> hole -> particle -> t -> t
(** Double excitation operator {% $T_{hh'}^{pp'} = a^\dagger_p a^\dagger_{p'} a_{h'} a_h$ %}. *)
val degrees : t -> t -> int*int
(** Returns degrees of excitation between two determinants in {% $\alpha$ %} and
{% $\beta$ %} spins as a pair. *)
val degree : t -> t -> int
(** Returns degree of excitation between two determinants. *)
(** {1 Creators} *)

View File

@ -14,6 +14,7 @@ let n_beta t = t.n_beta
let mo_class t = t.mo_class
let mo_basis t = t.mo_basis
let determinants t = t.determinants
let size t = Array.length t.determinants
let fci_of_mo_basis ?(frozen_core=true) mo_basis =
let s = MOBasis.simulation mo_basis in

View File

@ -21,6 +21,10 @@ val mo_basis : t -> MOBasis.t
val determinants : t -> Determinant.t array
(** All the determinants belonging to the space. *)
val size : t -> int
(** Size of the determinant space *)
val fci_of_mo_basis : ?frozen_core:bool -> MOBasis.t -> t
(** Creates a space of all possible ways to put [n_alfa] electrons in the {% $\alpha$ %}
[Active] MOs and [n_beta] electrons in the {% $\beta$ %} [Active] MOs.

View File

@ -76,6 +76,10 @@ let multiple_of_spindet t t' =
in
(phase, List.map2 (fun hole particle -> (hole, particle)) holes (List.rev particles) )
let double_of_spindet t t' =
match multiple_of_spindet t t' with
| (phase, (h1,p1)::(h2,p2)::[]) -> (h1, p1, h2, p2, phase)
| _ -> invalid_arg "t and t' are not doubly excited"
let multiple_of_det t t' =
let pa, a =

View File

@ -91,6 +91,14 @@ let particles_of t t' =
Z.logand (bitstring t') (Z.logxor (bitstring t) (bitstring t'))
|> bits_to_list []
let holes_particles_of t t' =
let x = Z.logxor (bitstring t) (bitstring t') in
let holes = Z.logand (bitstring t) x |> bits_to_list []
and particles = Z.logand (bitstring t') x |> bits_to_list []
in
List.map2 (fun h p -> (h,p)) holes particles
let negate_phase = function
| Some t -> Some { t with phase = Phase.neg t.phase }
| None -> None

View File

@ -45,7 +45,7 @@ val double_excitation : hole -> particle -> hole -> particle -> t -> t
(** Double excitation operator {% $T_{hh'}^{pp'} = a^\dagger_p a^\dagger_{p'} a_{h'} a_h$ %}. *)
val degree : t -> t -> int
(** Returns degree of excitation between two determinants. *)
(** Returns degree of excitation between two spin-determinants. *)
val holes_of : t -> t -> int list
(** Returns the list of holes in the excitation from one determinant to another. *)
@ -53,6 +53,9 @@ val holes_of : t -> t -> int list
val particles_of : t -> t -> int list
(** Returns the list of particles in the excitation from one determinant to another. *)
val holes_particles_of : t -> t -> (int*int) list
(** Returns the list of pairs of holes/particles in the excitation from one determinant to
another. *)
(** {1 Creation} *)

View File

@ -35,6 +35,10 @@ let mo_coef t = t.mo_coef
let eN_ints t = Lazy.force t.eN_ints
let ee_ints t = Lazy.force t.ee_ints
let kin_ints t = Lazy.force t.kin_ints
let two_e_ints t = Lazy.force t.ee_ints
let one_e_ints t =
Mat.add (NucInt.matrix @@ Lazy.force t.eN_ints)
(KinInt.matrix @@ Lazy.force t.kin_ints)
let mo_matrix_of_ao_matrix ~mo_coef ao_matrix =
xt_o_x ~x:mo_coef ~o:ao_matrix

View File

@ -40,6 +40,12 @@ val ee_ints : t -> ERI.t
val kin_ints : t -> KinInt.t
(** Kinetic energy integrals *)
val one_e_ints : t -> Mat.t
(** One-electron integrals {% $\hat{T} + V$ %} *)
val two_e_ints : t -> ERI.t
(** Electron-electron repulsion integrals *)
val size : t -> int
(** Number of molecular orbitals in the basis *)