10
1
mirror of https://gitlab.com/scemama/QCaml.git synced 2024-06-02 11:25:19 +02:00
QCaml/CI/CI.ml

383 lines
11 KiB
OCaml
Raw Normal View History

2019-02-22 00:18:32 +01:00
open Lacaml.D
2019-03-04 22:56:03 +01:00
module Ds = DeterminantSpace
2019-02-22 00:18:32 +01:00
type t =
{
2019-02-22 19:19:11 +01:00
det_space : Ds.t ;
2019-02-28 12:30:20 +01:00
m_H : Matrix.t lazy_t ;
m_S2 : Matrix.t lazy_t ;
2019-02-22 00:18:32 +01:00
eigensystem : (Mat.t * Vec.t) lazy_t;
2019-02-28 16:55:50 +01:00
n_states : int;
2019-02-22 00:18:32 +01:00
}
let det_space t = t.det_space
2019-02-28 16:55:50 +01:00
let n_states t = t.n_states
2019-02-26 11:58:53 +01:00
let m_H t = Lazy.force t.m_H
2019-02-22 00:18:32 +01:00
2019-02-26 11:58:53 +01:00
let m_S2 t = Lazy.force t.m_S2
2019-02-22 19:19:11 +01:00
2019-02-22 00:18:32 +01:00
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
2019-02-22 19:19:11 +01:00
let h_integrals mo_basis =
2019-02-22 00:18:32 +01:00
let one_e_ints = MOBasis.one_e_ints mo_basis
and two_e_ints = MOBasis.two_e_ints mo_basis
in
2019-02-22 19:19:11 +01:00
( (fun i j _ -> one_e_ints.{i,j}),
(fun i j k l s s' ->
if s' = Spin.other s then
ERI.get_phys two_e_ints i j k l
else
(ERI.get_phys two_e_ints i j k l) -.
(ERI.get_phys two_e_ints i j l k)
) )
let h_ij mo_basis ki kj =
let integrals =
List.map (fun f -> f mo_basis)
[ h_integrals ]
in
CIMatrixElement.make integrals ki kj
|> List.hd
2019-02-22 00:18:32 +01:00
2019-02-28 19:37:54 +01:00
let create_matrix_arbitrary f det_space =
lazy (
2019-03-04 22:56:03 +01:00
let ndet = Ds.size det_space in
let data =
2019-02-28 19:37:54 +01:00
match Ds.determinants det_space with
2019-03-04 22:56:03 +01:00
| Ds.Arbitrary a -> a
2019-02-28 19:37:54 +01:00
| _ -> assert false
in
2019-03-04 22:56:03 +01:00
let det_alfa = data.Ds.det_alfa
and det_beta = data.Ds.det_beta
and det = data.Ds.det
and index_start = data.Ds.index_start
in
2019-02-22 00:18:32 +01:00
2019-02-28 16:55:50 +01:00
2019-03-04 22:56:03 +01:00
(** Array of (list of singles, list of doubles) in the beta spin *)
let degree_bb =
Array.map (fun det_i ->
let deg = Spindeterminant.degree det_i in
let doubles =
Array.mapi (fun i det_j ->
let d = deg det_j in
if d < 3 then
Some (i,d,det_j)
else
None
) det_beta
|> Array.to_list
|> Util.list_some
in
let singles =
List.filter (fun (i,d,det_j) -> d < 2) doubles
|> List.map (fun (i,_,det_j) -> (i,det_j))
in
let doubles =
List.map (fun (i,_,det_j) -> (i,det_j)) doubles
in
(singles, doubles)
) det_beta
in
let task (i,i_dets) =
let shift = index_start.(i) in
2019-02-28 19:37:54 +01:00
2019-03-04 22:56:03 +01:00
let result =
Array.init (index_start.(i+1) - shift)
(fun _ -> [])
in
(** Update function when ki and kj are connected *)
let update i j ki kj =
let x = f ki kj in
if abs_float x > Constants.epsilon then
result.(i - shift) <- (j, x) :: result.(i - shift) ;
in
let i_alfa = det_alfa.(i) in
let deg_a = Spindeterminant.degree i_alfa in
Array.iteri (fun j j_dets ->
let j_alfa = det_alfa.(j) in
let degree_a = deg_a j_alfa in
begin
match degree_a with
| 2 ->
Array.iteri (fun i' i_b ->
try
Array.iteri (fun j' j_b ->
if j_b >= i_b then
( if j_b = i_b then
( let i_beta = det_beta.(i_b) in
let ki = Determinant.of_spindeterminants i_alfa i_beta in
let kj = Determinant.of_spindeterminants j_alfa i_beta in
update (index_start.(i) + i')
(index_start.(j) + j' + 1) ki kj);
raise Exit)
) j_dets
with Exit -> ()
) i_dets
| 1 ->
Array.iteri (fun i' i_b ->
let i_beta = det_beta.(i_b) in
let ki = Determinant.of_spindeterminants i_alfa i_beta in
let singles, _ = degree_bb.(i_b) in
let rec aux singles j' =
match singles with
| [] -> ()
| (js, j_beta)::r_singles ->
begin
match compare js j_dets.(j') with
| -1 -> aux r_singles j'
| 0 ->
let kj =
Determinant.of_spindeterminants j_alfa j_beta
in (update
(index_start.(i) + i') (index_start.(j) + j' + 1)
ki kj;
aux r_singles (j'+1);)
| 1 -> if (j' < Array.length j_dets) then aux singles (j'+1)
| _ -> assert false
end
in aux singles 0
) i_dets
| 0 ->
Array.iteri (fun i' i_b ->
let i_beta = det_beta.(i_b) in
let ki = Determinant.of_spindeterminants i_alfa i_beta in
let _, doubles = degree_bb.(i_b) in
let rec aux doubles j' =
match doubles with
| [] -> ()
| (js, j_beta)::r_doubles ->
begin
match compare js j_dets.(j') with
| -1 -> aux r_doubles j'
| 0 ->
let kj =
Determinant.of_spindeterminants j_alfa j_beta
in (update
(index_start.(i) + i') (index_start.(j) + j' + 1)
ki kj;
aux r_doubles (j'+1);)
| 1 -> if (j' < Array.length j_dets) then aux doubles (j'+1)
| _ -> assert false
end
in aux doubles 0
) i_dets
| _ -> ();
end
) det;
let r =
Array.map (fun l ->
List.rev l
|> Vector.sparse_of_assoc_list ndet
) result
in (i,r)
in
let result =
if Parallel.master then
Array.init ndet (fun _ -> Vector.sparse_of_assoc_list ndet [])
else
Array.init ndet (fun _ -> Vector.sparse_of_assoc_list ndet [])
in
let n_det_alfa =
Array.length det_alfa
in
Array.mapi (fun i i_dets -> i, i_dets) det
|> Array.to_list
|> Stream.of_list
|> Farm.run ~ordered:false ~f:task
|> Stream.iter (fun (k, r) ->
let shift = index_start.(k) in
let det_k = det.(k) in
Array.iteri (fun j r_j -> result.(shift+det_k.(j)) <- r_j) r;
Printf.eprintf "%8d / %8d\r%!" (k+1) n_det_alfa;
) ;
Matrix.sparse_of_vector_array result
2019-02-28 19:37:54 +01:00
)
(* Create a matrix using the fact that the determinant space is made of
the outer product of spindeterminants. *)
let create_matrix_spin f det_space =
lazy (
let ndet = Ds.size det_space in
let a, b =
match Ds.determinants det_space with
| Ds.Spin (a,b) -> (a,b)
| _ -> assert false
in
let n_beta = Array.length b in
2019-02-28 16:55:50 +01:00
2019-02-28 19:37:54 +01:00
(** Array of (list of singles, list of doubles) in the beta spin *)
let degree_bb =
Array.map (fun det_i ->
2019-02-28 18:18:26 +01:00
let deg = Spindeterminant.degree det_i in
let doubles =
Array.mapi (fun i det_j ->
2019-02-28 19:37:54 +01:00
let d = deg det_j in
if d < 3 then
Some (i,d,det_j)
else
None
2019-02-28 18:18:26 +01:00
) b
|> Array.to_list
|> Util.list_some
in
let singles =
List.filter (fun (i,d,det_j) -> d < 2) doubles
|> List.map (fun (i,_,det_j) -> (i,det_j))
in
let doubles =
List.map (fun (i,_,det_j) -> (i,det_j)) doubles
in
(singles, doubles)
) b
2019-02-28 19:37:54 +01:00
in
let a = Array.to_list a
and b = Array.to_list b
in
2019-03-02 18:50:12 +01:00
let task (i,i_alfa) =
2019-03-04 22:56:03 +01:00
let result =
Array.init n_beta (fun _ -> [])
in
2019-03-02 17:05:15 +01:00
(** Update function when ki and kj are connected *)
let update i j ki kj =
let x = f ki kj in
if abs_float x > Constants.epsilon then
result.(i) <- (j, x) :: result.(i) ;
in
2019-03-04 22:56:03 +01:00
2019-02-28 19:37:54 +01:00
let j = ref 1 in
let deg_a = Spindeterminant.degree i_alfa in
List.iter (fun j_alfa ->
let degree_a = deg_a j_alfa in
begin
match degree_a with
| 2 ->
2019-03-02 17:05:15 +01:00
let i' = ref 0 in
2019-02-28 19:37:54 +01:00
List.iteri (fun ib i_beta ->
let ki = Determinant.of_spindeterminants i_alfa i_beta in
let kj = Determinant.of_spindeterminants j_alfa i_beta in
update !i' (ib + !j) ki kj;
incr i';
) b;
| 1 ->
2019-03-02 17:05:15 +01:00
let i' = ref 0 in
2019-02-28 19:37:54 +01:00
List.iteri (fun ib i_beta ->
let ki = Determinant.of_spindeterminants i_alfa i_beta in
let singles, _ = degree_bb.(ib) in
List.iter (fun (j', j_beta) ->
let kj = Determinant.of_spindeterminants j_alfa j_beta in
update !i' (j' + !j) ki kj
) singles;
incr i';
) b;
| 0 ->
2019-03-02 17:05:15 +01:00
let i' = ref 0 in
2019-02-28 19:37:54 +01:00
List.iteri (fun ib i_beta ->
let ki = Determinant.of_spindeterminants i_alfa i_beta in
let _singles, doubles = degree_bb.(ib) in
List.iter (fun (j', j_beta) ->
let kj = Determinant.of_spindeterminants j_alfa j_beta in
update !i' (j' + !j) ki kj
) doubles;
incr i';
) b;
| _ -> ();
end;
j := !j + n_beta
) a;
2019-03-02 18:50:12 +01:00
let r =
Array.map (fun l ->
List.rev l
|> Vector.sparse_of_assoc_list ndet
) result
in (i,r)
2019-03-02 17:05:15 +01:00
in
2019-03-02 18:50:12 +01:00
2019-03-02 17:05:15 +01:00
let result =
2019-03-02 18:50:12 +01:00
if Parallel.master then
Array.init ndet (fun _ -> Vector.sparse_of_assoc_list ndet [])
else
Array.init ndet (fun _ -> Vector.sparse_of_assoc_list ndet [])
2019-03-02 17:05:15 +01:00
in
2019-03-02 18:50:12 +01:00
List.mapi (fun i i_alfa -> i*n_beta, i_alfa) a
|> Stream.of_list
|> Farm.run ~ordered:false ~f:task
|> Stream.iter (fun (k, r) ->
Array.iteri (fun j r_j -> result.(k+j) <- r_j) r;
Printf.eprintf "%8d / %8d\r%!" (k+1) ndet;
) ;
2019-03-02 17:05:15 +01:00
Matrix.sparse_of_vector_array result
2019-02-28 19:37:54 +01:00
)
let make ?(n_states=1) det_space =
2019-02-28 16:55:50 +01:00
2019-02-28 19:37:54 +01:00
let m_H =
2019-03-02 18:50:12 +01:00
2019-02-28 19:37:54 +01:00
let mo_basis = Ds.mo_basis det_space in
2019-03-02 18:50:12 +01:00
(* While in a sequential region, initiate the parallel
2019-03-04 22:56:03 +01:00
4-idx transformation to avoid nested parallel jobs
*)
2019-03-02 18:50:12 +01:00
ignore @@ MOBasis.two_e_ints mo_basis;
2019-02-28 19:37:54 +01:00
let f =
match Ds.determinants det_space with
| Ds.Arbitrary _ -> create_matrix_arbitrary
| Ds.Spin _ -> create_matrix_spin
in
f (fun ki kj -> h_ij mo_basis ki kj) det_space
2019-02-28 12:50:42 +01:00
in
2019-02-28 16:55:50 +01:00
2019-02-28 19:37:54 +01:00
let m_S2 =
let f =
match Ds.determinants det_space with
| Ds.Arbitrary _ -> create_matrix_arbitrary
| Ds.Spin _ -> create_matrix_spin
in
f (fun ki kj -> CIMatrixElement.make_s2 ki kj) det_space
2019-02-22 00:18:32 +01:00
in
2019-02-28 16:55:50 +01:00
2019-02-22 00:18:32 +01:00
let eigensystem = lazy (
2019-02-28 16:55:50 +01:00
let m_H =
Lazy.force m_H
in
2019-02-28 12:30:20 +01:00
let diagonal =
Vec.init (Matrix.dim1 m_H) (fun i -> Matrix.get m_H i i)
in
let matrix_prod psi =
Matrix.mm ~transa:`T m_H psi
in
2019-02-28 16:55:50 +01:00
Davidson.make ~n_states diagonal matrix_prod
2019-02-22 00:18:32 +01:00
)
in
2019-02-28 16:55:50 +01:00
{ det_space ; m_H ; m_S2 ; eigensystem ; n_states }
2019-02-22 00:18:32 +01:00
2019-02-27 14:56:59 +01:00