1
0
mirror of https://github.com/TREX-CoE/fparser.git synced 2024-10-02 06:21:08 +02:00

determinant file read successfully

This commit is contained in:
Ravindra Shinde 2021-03-03 15:53:47 +01:00
parent ff575f7a93
commit 61f7711ffd
4 changed files with 50 additions and 54 deletions

View File

@ -35,6 +35,7 @@ PROGRAM iochamp
real(dp) :: float_value
character(len=20) :: real_format = '(A, T20, F14.8)'
character(len=20) :: int_format = '(A, T20, I8)'
character(len=80) :: string_format = '(A, T40, A)'
! for determinants sections
integer :: nelectrons, nexcitation, iostat
@ -52,26 +53,29 @@ PROGRAM iochamp
! Get the directory where the pooled data is kept
path_pool = fdf_string('pool', './')
write(6,'(A)') 'pool directory location :: ', path_pool
write(6,fmt=string_format) 'pool directory location :: ', path_pool
! Get all the filenames from which the data is to be read
file_molecule = fdf_load_filename('molecule', 'default.xyz')
write(6,fmt=string_format) 'filename molecule :: ', trim(file_molecule)
file_pseudo = fdf_load_filename('pseudopot', 'default.psp')
write(6,'(A)') 'filename pseuodpotential :: ', file_pseudo
write(6,fmt=string_format) 'filename pseuodpotential :: ', trim(file_pseudo)
file_basis = fdf_load_filename('basis', 'default.bas')
write(6,'(A)') 'filename basis :: ', file_basis
write(6,fmt=string_format) 'filename basis :: ', trim(file_basis)
file_determinants = fdf_load_filename('determinants', 'default.det')
write(6,'(A)') 'filename determinants :: ', file_determinants
write(6,fmt=string_format) 'filename determinants :: ', trim(file_determinants)
file_orbitals = fdf_load_filename('orbitals', 'default.orb')
write(6,'(A)') 'filename orbitals :: ', file_orbitals
write(6,fmt=string_format) 'filename orbitals :: ', trim(file_orbitals)
file_jastrow = fdf_load_filename('jastrow', 'default.jas')
write(6,'(A)') 'filename jastrow :: ', file_jastrow
write(6,fmt=string_format) 'filename jastrow :: ',trim(file_jastrow)
file_jastrow_deriv = fdf_load_filename('jastrow_deriv', 'default.jasder')
write(6,'(A)') 'filename jastrow derivatives :: ', file_jastrow_deriv
write(6,fmt=string_format) 'filename jastrow derivatives :: ', trim(file_jastrow_deriv)
! &optwf ioptwf 1 ioptci 1 ioptjas 1 ioptorb 1
@ -403,59 +407,50 @@ PROGRAM iochamp
write(6,*) '------------------------------------------------------'
if (fdf_block('determinants', bfdf)) then
if (.not. fdf_block('determinants', bfdf)) then
! External file reading
write(6,*) 'Beginning of external file determinant block '
write(6,*) 'Reading the determinants block from an external file '
ia = 1
! call io_status()
! call fdf_printfdf()
print*, "printing label ", bfdf%label , trim(bfdf%mark%pline%line)
! print*, "printing label ", bfdf%label , trim(bfdf%mark%pline%line)
! print*, "pline obtained", (fdf_bline(bfdf, pline))
print*, "pline obtained", (fdf_bline(bfdf, pline))
open (unit=11,file='TZ_1M_500.det', iostat=iostat, action='read' )
open (unit=11,file=file_determinants, iostat=iostat, action='read' )
read(11,*) temp1, temp2, nelectrons, temp3, nalpha
write(*,'(a,1x,i3,1x,i3)') "write after read1", nelectrons, nalpha
read(11,*) temp1, ndeterminants, nexcitation
allocate(det_coeff(ndeterminants))
write(*,'(a,1x,i3, 1x, i3)') "write after read2", ndeterminants, nexcitation
read(11,*) temp1, ndeterminants, nexcitation
if (.not. allocated(det_coeff)) allocate(det_coeff(ndeterminants))
read(11,*) (det_coeff(i), i=1,ndeterminants)
write(fmt,*) '(', ndeterminants, '(f11.8,1x))'
write(*,fmt) (det_coeff(i), i=1,ndeterminants)
! write(*,'(<ndeterminants>(f11.8, 1x))') (det_coeff(i), i=1,ndeterminants) ! for Intel Fortran
close(11)
! write(*,'(<ndeterminants>(f11.8, 1x))') (det_coeff(i), i=1,ndeterminants) ! for Intel Fortran
nbeta = nelectrons - nalpha
! allocate the orbital mapping array
if (.not. allocated(iworbd)) allocate(iworbd(nelectrons, ndeterminants))
write(*,*) "total number of electrons ", nelectrons
write(*,*) " number of alpha electrons ", nalpha
write(*,*) " number of beta electrons ", nbeta
do i = 1, ndeterminants
read(11,*) (iworbd(j,i), j=1,nelectrons)
enddo
write(fmt,*) '(i4,1x)'
do i = 1, ndeterminants
write(*,'(<nelectrons>(i4, 1x))') (iworbd(j,i), j=1,nelectrons)
enddo
read(11,*) temp1
if (temp1 == "end" ) write(*,*) "Determinant File read successfully "
if(fdf_bsearch(pline, "&electrons")) then
nelectrons = integers(pline, 1)
nalpha = integers(pline, 2)
nbeta = nelectrons - nalpha
write(*,*) "total number of electrons ", nelectrons
write(*,*) " number of alpha electrons ", nalpha
write(*,*) " number of beta electrons ", nalpha
if (.not. allocated(det_alpha)) allocate(det_alpha(nalpha))
if (.not. allocated(det_beta)) allocate(det_beta(nbeta))
endif
close(11)
if(fdf_bsearch(pline, "determinants")) then
ndeterminants = fdf_bintegers(pline, 1)
nexcitation = fdf_bintegers(pline, 2)
write(*,*) "total number of determinants ", ndeterminants
write(*,*) " number of excitations ", nexcitation
if (.not. allocated(det_coeff)) allocate(det_coeff(ndeterminants))
endif
na = nintegers(pline)
write(*,'(2(a,i0),a)') 'number of integers: ', na, ' integers'
write(*,'(tr5,a,<nalpha>(tr1,i0))') 'list: ', det_alpha(1:nalpha)
! endif
! enddo
endif
write(6,'(A)')

View File

@ -86,6 +86,7 @@ MODULE keywords
public :: nctype, ntypes_atom ! number of atom/center types
public :: ncent, natoms, ncenters, ncentres ! number of atoms/centers
public :: iwctype ! specify atom-type for each atom
public :: znuc, znuclear, atomic_number ! nuclear charge
public :: cent !atom_coords ! atom positions
public :: ndet, ndeterminants ! number of determinants in wavefunction
@ -152,6 +153,7 @@ MODULE keywords
character(len=132) :: file_input, file_output
character(len=132) :: file_basis
character(len=132) :: file_pseudo
character(len=132) :: file_molecule
character(len=132) :: file_orbitals
character(len=132) :: file_determinants
character(len=132) :: file_jastrow
@ -250,6 +252,7 @@ MODULE keywords
integer, pointer :: natoms => ncent, ncenters => ncent, ncentres => ncent
integer :: iwctype
integer, allocatable :: iworbd(:,:) ! to store orbital mapping in determinants
integer, target :: znuc
integer, pointer :: znuclear => znuc, atomic_number => znuc
@ -268,7 +271,6 @@ MODULE keywords
real(dp), allocatable :: cdet(:)
! integer, pointer :: det_coeffs => cdet ! issues in performance
integer :: iworbd
integer, target :: ianalyt_lap
integer, pointer :: analytic_laplacian => ianalyt_lap

View File

@ -7,10 +7,10 @@ basis BFD-T-normf0
%include global.fdf
%block molecule < caffeine.xyz
%block determinants < sample.det
load determinants TZ_1M_500.det
load molecule caffeine.xyz
load basis BFD-T-normf0.bas
load determinants "TZ_1M_500.det"
optimize_wavefunction 1

View File

@ -1,7 +1,8 @@
title A Sample QMC input file parsed by libfdf interfaced to CHAMP
pool ./pool
molecule caffeine.xyz
pseudopot default.psp # default value
basis default.bas # default value
basis BFD-T-normf0.bas
determinants TZ_1M_500.det
orbitals default.orb # default value
jastrow default.jas # default value
@ -57,5 +58,3 @@ etrial -408.1744362 eV
H 3.2249 1.5791 0.8255
H 2.6431 2.5130 -0.5793
%endblock molecule
%block determinants
&electrons nelec 22 nup 11