mirror of
https://github.com/LCPQ/quantum_package
synced 2024-12-31 16:45:54 +01:00
Fixed Pseudo and dummy atoms
This commit is contained in:
parent
5e1b077576
commit
2dcb4eba0d
244
ocaml/Pseudo.ml
244
ocaml/Pseudo.ml
@ -124,23 +124,27 @@ let to_string t =
|
|||||||
let find in_channel element =
|
let find in_channel element =
|
||||||
In_channel.seek in_channel 0L;
|
In_channel.seek in_channel 0L;
|
||||||
|
|
||||||
let element_read, old_pos =
|
let loop, element_read, old_pos =
|
||||||
ref Element.X,
|
ref true,
|
||||||
|
ref None,
|
||||||
ref (In_channel.pos in_channel)
|
ref (In_channel.pos in_channel)
|
||||||
in
|
in
|
||||||
while !element_read <> element
|
|
||||||
|
while !loop
|
||||||
do
|
do
|
||||||
let buffer =
|
|
||||||
old_pos := In_channel.pos in_channel;
|
|
||||||
match In_channel.input_line in_channel with
|
|
||||||
| Some line -> String.split ~on:' ' line
|
|
||||||
|> List.hd_exn
|
|
||||||
| None -> ""
|
|
||||||
in
|
|
||||||
try
|
try
|
||||||
element_read := Element.of_string buffer
|
let buffer =
|
||||||
|
old_pos := In_channel.pos in_channel;
|
||||||
|
match In_channel.input_line in_channel with
|
||||||
|
| Some line -> String.split ~on:' ' line
|
||||||
|
|> List.hd_exn
|
||||||
|
| None -> raise End_of_file
|
||||||
|
in
|
||||||
|
element_read := Some (Element.of_string buffer);
|
||||||
|
loop := !element_read <> (Some element)
|
||||||
with
|
with
|
||||||
| Element.ElementError _ -> ()
|
| Element.ElementError _ -> ()
|
||||||
|
| End_of_file -> loop := false
|
||||||
done ;
|
done ;
|
||||||
In_channel.seek in_channel !old_pos;
|
In_channel.seek in_channel !old_pos;
|
||||||
!element_read
|
!element_read
|
||||||
@ -148,124 +152,126 @@ let find in_channel element =
|
|||||||
|
|
||||||
(** Read the Pseudopotential in GAMESS format *)
|
(** Read the Pseudopotential in GAMESS format *)
|
||||||
let read_element in_channel element =
|
let read_element in_channel element =
|
||||||
ignore (find in_channel element);
|
match find in_channel element with
|
||||||
|
| Some e when e = element ->
|
||||||
let rec read result =
|
|
||||||
match In_channel.input_line in_channel with
|
|
||||||
| None -> result
|
|
||||||
| Some line ->
|
|
||||||
if (String.strip line = "") then
|
|
||||||
result
|
|
||||||
else
|
|
||||||
read (line::result)
|
|
||||||
in
|
|
||||||
|
|
||||||
let data =
|
|
||||||
read []
|
|
||||||
|> List.rev
|
|
||||||
in
|
|
||||||
|
|
||||||
let debug_data =
|
|
||||||
String.concat ~sep:"\n" data
|
|
||||||
in
|
|
||||||
|
|
||||||
let decode_first_line = function
|
|
||||||
| first_line :: rest ->
|
|
||||||
begin
|
begin
|
||||||
let first_line_split =
|
let rec read result =
|
||||||
String.split first_line ~on:' '
|
match In_channel.input_line in_channel with
|
||||||
|> List.filter ~f:(fun x -> (String.strip x) <> "")
|
| None -> result
|
||||||
|
| Some line ->
|
||||||
|
if (String.strip line = "") then
|
||||||
|
result
|
||||||
|
else
|
||||||
|
read (line::result)
|
||||||
in
|
in
|
||||||
match first_line_split with
|
|
||||||
| e :: "GEN" :: n :: p ->
|
|
||||||
{ element = Element.of_string e ;
|
|
||||||
n_elec = Int.of_string n |> Positive_int.of_int ;
|
|
||||||
local = [] ;
|
|
||||||
non_local = []
|
|
||||||
}, rest
|
|
||||||
| _ -> failwith (
|
|
||||||
Printf.sprintf "Unable to read Pseudopotential : \n%s\n"
|
|
||||||
debug_data )
|
|
||||||
end
|
|
||||||
| _ -> failwith ("Error reading pseudopotential\n"^debug_data)
|
|
||||||
in
|
|
||||||
|
|
||||||
let rec loop create_primitive accu = function
|
let data =
|
||||||
| (0,rest) -> List.rev accu, rest
|
read []
|
||||||
| (n,line::rest) ->
|
|> List.rev
|
||||||
begin
|
in
|
||||||
match
|
|
||||||
String.split line ~on:' '
|
let debug_data =
|
||||||
|> List.filter ~f:(fun x -> String.strip x <> "")
|
String.concat ~sep:"\n" data
|
||||||
with
|
in
|
||||||
| c :: i :: e :: [] ->
|
|
||||||
let i =
|
let decode_first_line = function
|
||||||
Int.of_string i
|
| first_line :: rest ->
|
||||||
in
|
begin
|
||||||
let elem =
|
let first_line_split =
|
||||||
( create_primitive
|
String.split first_line ~on:' '
|
||||||
(Float.of_string e |> AO_expo.of_float)
|
|> List.filter ~f:(fun x -> (String.strip x) <> "")
|
||||||
(i-2 |> R_power.of_int),
|
in
|
||||||
Float.of_string c |> AO_coef.of_float
|
match first_line_split with
|
||||||
)
|
| e :: "GEN" :: n :: p ->
|
||||||
in
|
{ element = Element.of_string e ;
|
||||||
loop create_primitive (elem::accu) (n-1, rest)
|
n_elec = Int.of_string n |> Positive_int.of_int ;
|
||||||
|
local = [] ;
|
||||||
|
non_local = []
|
||||||
|
}, rest
|
||||||
|
| _ -> failwith (
|
||||||
|
Printf.sprintf "Unable to read Pseudopotential : \n%s\n"
|
||||||
|
debug_data )
|
||||||
|
end
|
||||||
| _ -> failwith ("Error reading pseudopotential\n"^debug_data)
|
| _ -> failwith ("Error reading pseudopotential\n"^debug_data)
|
||||||
end
|
in
|
||||||
| _ -> failwith ("Error reading pseudopotential\n"^debug_data)
|
|
||||||
in
|
|
||||||
|
|
||||||
let decode_local (pseudo,data) =
|
let rec loop create_primitive accu = function
|
||||||
let decode_local_n n rest =
|
| (0,rest) -> List.rev accu, rest
|
||||||
let result, rest =
|
| (n,line::rest) ->
|
||||||
loop Primitive_local.of_expo_r_power [] (Positive_int.to_int n,rest)
|
begin
|
||||||
|
match
|
||||||
|
String.split line ~on:' '
|
||||||
|
|> List.filter ~f:(fun x -> String.strip x <> "")
|
||||||
|
with
|
||||||
|
| c :: i :: e :: [] ->
|
||||||
|
let i =
|
||||||
|
Int.of_string i
|
||||||
|
in
|
||||||
|
let elem =
|
||||||
|
( create_primitive
|
||||||
|
(Float.of_string e |> AO_expo.of_float)
|
||||||
|
(i-2 |> R_power.of_int),
|
||||||
|
Float.of_string c |> AO_coef.of_float
|
||||||
|
)
|
||||||
|
in
|
||||||
|
loop create_primitive (elem::accu) (n-1, rest)
|
||||||
|
| _ -> failwith ("Error reading pseudopotential\n"^debug_data)
|
||||||
|
end
|
||||||
|
| _ -> failwith ("Error reading pseudopotential\n"^debug_data)
|
||||||
in
|
in
|
||||||
{ pseudo with local = result }, rest
|
|
||||||
in
|
let decode_local (pseudo,data) =
|
||||||
match data with
|
let decode_local_n n rest =
|
||||||
| n :: rest ->
|
let result, rest =
|
||||||
let n =
|
loop Primitive_local.of_expo_r_power [] (Positive_int.to_int n,rest)
|
||||||
String.strip n
|
in
|
||||||
|> Int.of_string
|
{ pseudo with local = result }, rest
|
||||||
|> Positive_int.of_int
|
|
||||||
in
|
in
|
||||||
decode_local_n n rest
|
match data with
|
||||||
| _ -> failwith ("Unable to read (non-)local pseudopotential\n"^debug_data)
|
| n :: rest ->
|
||||||
in
|
let n =
|
||||||
|
String.strip n
|
||||||
let decode_non_local (pseudo,data) =
|
|> Int.of_string
|
||||||
let decode_non_local_n proj n (pseudo,data) =
|
|> Positive_int.of_int
|
||||||
let result, rest =
|
in
|
||||||
loop (Primitive_non_local.of_proj_expo_r_power proj)
|
decode_local_n n rest
|
||||||
[] (Positive_int.to_int n, data)
|
| _ -> failwith ("Unable to read (non-)local pseudopotential\n"^debug_data)
|
||||||
in
|
in
|
||||||
{ pseudo with non_local = pseudo.non_local @ result }, rest
|
|
||||||
in
|
|
||||||
let rec new_proj (pseudo,data) proj =
|
|
||||||
match data with
|
|
||||||
| n :: rest ->
|
|
||||||
let n =
|
|
||||||
String.strip n
|
|
||||||
|> Int.of_string
|
|
||||||
|> Positive_int.of_int
|
|
||||||
in
|
|
||||||
let result =
|
|
||||||
decode_non_local_n proj n (pseudo,rest)
|
|
||||||
and proj_next =
|
|
||||||
(Positive_int.to_int proj)+1
|
|
||||||
|> Positive_int.of_int
|
|
||||||
in
|
|
||||||
new_proj result proj_next
|
|
||||||
| _ -> pseudo
|
|
||||||
in
|
|
||||||
new_proj (pseudo,data) (Positive_int.of_int 0)
|
|
||||||
in
|
|
||||||
|
|
||||||
decode_first_line data
|
let decode_non_local (pseudo,data) =
|
||||||
|> decode_local
|
let decode_non_local_n proj n (pseudo,data) =
|
||||||
|> decode_non_local
|
let result, rest =
|
||||||
|
loop (Primitive_non_local.of_proj_expo_r_power proj)
|
||||||
|
[] (Positive_int.to_int n, data)
|
||||||
|
in
|
||||||
|
{ pseudo with non_local = pseudo.non_local @ result }, rest
|
||||||
|
in
|
||||||
|
let rec new_proj (pseudo,data) proj =
|
||||||
|
match data with
|
||||||
|
| n :: rest ->
|
||||||
|
let n =
|
||||||
|
String.strip n
|
||||||
|
|> Int.of_string
|
||||||
|
|> Positive_int.of_int
|
||||||
|
in
|
||||||
|
let result =
|
||||||
|
decode_non_local_n proj n (pseudo,rest)
|
||||||
|
and proj_next =
|
||||||
|
(Positive_int.to_int proj)+1
|
||||||
|
|> Positive_int.of_int
|
||||||
|
in
|
||||||
|
new_proj result proj_next
|
||||||
|
| _ -> pseudo
|
||||||
|
in
|
||||||
|
new_proj (pseudo,data) (Positive_int.of_int 0)
|
||||||
|
in
|
||||||
|
|
||||||
|
decode_first_line data
|
||||||
|
|> decode_local
|
||||||
|
|> decode_non_local
|
||||||
|
end
|
||||||
|
| _ -> empty element
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
include To_md5
|
include To_md5
|
||||||
|
@ -1,7 +1,7 @@
|
|||||||
program e_curve
|
program e_curve
|
||||||
use bitmasks
|
use bitmasks
|
||||||
implicit none
|
implicit none
|
||||||
integer :: i,j,k, nab, m, l
|
integer :: i,j,k, kk, nab, m, l
|
||||||
double precision :: norm, E, hij, num, ci, cj
|
double precision :: norm, E, hij, num, ci, cj
|
||||||
integer, allocatable :: iorder(:)
|
integer, allocatable :: iorder(:)
|
||||||
double precision , allocatable :: norm_sort(:)
|
double precision , allocatable :: norm_sort(:)
|
||||||
@ -60,7 +60,7 @@ program e_curve
|
|||||||
num = 0.d0
|
num = 0.d0
|
||||||
norm = 0.d0
|
norm = 0.d0
|
||||||
m = 0
|
m = 0
|
||||||
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(k,l,det_i,det_j,ci,cj,hij) REDUCTION(+:norm,m,num)
|
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(k,kk,l,det_i,det_j,ci,cj,hij) REDUCTION(+:norm,m,num)
|
||||||
allocate( det_i(N_int,2), det_j(N_int,2))
|
allocate( det_i(N_int,2), det_j(N_int,2))
|
||||||
!$OMP DO SCHEDULE(guided)
|
!$OMP DO SCHEDULE(guided)
|
||||||
do k=1,n_det
|
do k=1,n_det
|
||||||
@ -68,15 +68,19 @@ program e_curve
|
|||||||
cycle
|
cycle
|
||||||
endif
|
endif
|
||||||
ci = psi_bilinear_matrix_values(k,1)
|
ci = psi_bilinear_matrix_values(k,1)
|
||||||
det_i(:,1) = psi_det_alpha_unique(:,psi_bilinear_matrix_rows(k))
|
do kk=1,N_int
|
||||||
det_i(:,2) = psi_det_beta_unique(:,psi_bilinear_matrix_columns(k))
|
det_i(kk,1) = psi_det_alpha_unique(kk,psi_bilinear_matrix_rows(k))
|
||||||
|
det_i(kk,2) = psi_det_beta_unique(kk,psi_bilinear_matrix_columns(k))
|
||||||
|
enddo
|
||||||
do l=1,n_det
|
do l=1,n_det
|
||||||
if (psi_bilinear_matrix_values(l,1) == 0.d0) then
|
if (psi_bilinear_matrix_values(l,1) == 0.d0) then
|
||||||
cycle
|
cycle
|
||||||
endif
|
endif
|
||||||
cj = psi_bilinear_matrix_values(l,1)
|
cj = psi_bilinear_matrix_values(l,1)
|
||||||
det_j(:,1) = psi_det_alpha_unique(:,psi_bilinear_matrix_rows(l))
|
do kk=1,N_int
|
||||||
det_j(:,2) = psi_det_beta_unique(:,psi_bilinear_matrix_columns(l))
|
det_j(kk,1) = psi_det_alpha_unique(kk,psi_bilinear_matrix_rows(l))
|
||||||
|
det_j(kk,2) = psi_det_beta_unique(kk,psi_bilinear_matrix_columns(l))
|
||||||
|
enddo
|
||||||
call i_h_j(det_i, det_j, N_int, hij)
|
call i_h_j(det_i, det_j, N_int, hij)
|
||||||
num = num + ci*cj*hij
|
num = num + ci*cj*hij
|
||||||
enddo
|
enddo
|
||||||
|
Loading…
Reference in New Issue
Block a user