forum@abinit.org
Subject: The ABINIT Users Mailing List ( CLOSED )
List archive
- From: Pascal Boulet <pascal.boulet@univ-provence.fr>
 - To: forum@abinit.org
 - Subject: Re: [abinit-forum] Can I generate "prtgeo 1" _GEO output files afterwards?
 - Date: Mon, 25 May 2009 14:36:59 +0200
 
| 
Hello Rob, I am not sure this post will answer your question. Here is a piece of basic code I wrote to parse the ABINIT output file. Basically, it gets the coordinates out of the ouput file and dump them into a xyz "trajectory" suitable for VMD. It should also work during the course of the calculation or if the job stopped for some reasons. Only the coordinates are read (no energies, no gradients). I would be glad to know if it works for you. Hope this help, Best Regards Pascal Rob wrote: 255164.18295.qm@web33308.mail.mud.yahoo.com" type="cite"> -- Dr. Pascal Boulet, Computational Chemist University Aix-Marseille I Laboratoire Chimie Provence, UMR6264 Centre Saint-Jerome, Bat. MADIREL F-13397 MARSEILLE Cedex 20, France Tel. +33 (0) 491 63 71 17 Fax. +33 (0) 491 63 71 11 courriel: pascal.boulet@univ-provence.fr http://www.lc-provence.fr http://allos.up.univ-mrs.fr/boulet  | 
program abinit2xmol
 implicit none
 real(kind=8),parameter :: r0=0.52917720859d0
 character(len=2),dimension(54),parameter :: sym=(/'H ','He',          &
 'Li','Be','B ','C ','N ','O ','F ','Ne',                              &
 'Na','Mg','Al','Si','P ','S ','Cl','Ar',                              &
 'K ','Ca','Sc','Ti','V ','Cr','Mn','Fe','Co','Ni','Cu','Zn',          &
 'Ga','Ge','As','Se','Br','Kr',                                        &
 'Rb','Sr','Y ','Zr','Nb','Mo','Tc','Ru','Rh','Pd','Ag','Cd',          &
 'In','Sn','Sb','Te','I ','Xe'/)
 real(kind=8),dimension(:,:),allocatable :: x
 character(len=2),dimension(:),allocatable :: elem
 real(kind=8),dimension(3) :: a 
 integer,dimension(:),allocatable :: typatom
 integer,dimension(:),allocatable :: znuc
 integer :: i,j,k,l,n,ierr
 integer(kind=8) :: natom
 integer(kind=8) :: ntypatom
 character(len=80) :: line
 character(len=25) :: abinitfn
 
 
  print*,'abinit output filename?'
  read*,abinitfn
  open(10,file=abinitfn,status='old')
  open(11,file='abinit.xmol',status='unknown')
 
  line=' '
  do while (index(line,'acell').eq.0)
    read(10,'(A80)') line
  enddo
 
  read(line(12:),*) (a(i),i=1,3)
  a=a*r0
 
  do while (index(line,'natom').eq.0)
    read(10,'(A80)') line
  enddo
 
  read(line(12:),*) natom
 
  allocate(x(natom,3))
  allocate(typatom(natom))
  allocate(elem(natom))
 
  do while (index(line,'ntypat').eq.0)
    read(10,'(A80)') line
  enddo
 
  read(line(12:),*) ntypatom
 
  allocate(znuc(ntypatom))
  
  line=' '
  do while (index(line,'typat').eq.0)
    read(10,'(A80)') line
  enddo
  
  i=int(natom/20)
  j=mod(natom,20)
  if (i.eq.0) then
    read(line(12:),*) (typatom(k),k=1,natom)
  else if (i.eq.1) then
    read(line(12:),*) (typatom(k),k=1,20)
    n=k
    read(10,'(A80)') line
    read(line(12:),*) (typatom(k),k=n,n+j-1)
  else if (i.gt.1) then
    read(line(12:),*) (typatom(k),k=1,20)
    n=k
    do l=1,i-1
      read(10,'(A80)') line
      read(line(12:),*) (typatom(k),k=n,n+19)
      n=k
    enddo
    read(10,'(A80)') line
    read(line(12:),*) (typatom(k),k=n,n+j-1)
  endif
  
  line=' '
  do while (index(line,'znucl').eq.0)
    read(10,'(A80)') line
  enddo
  
  read(line(12:),*) (znuc(i),i=1,ntypatom)
  
  l1:do
    do while (index(line,'Cartesian coordinates (bohr)').eq.0)
      read(10,'(A80)',iostat=ierr) line
      if (ierr.ne.0) exit l1
    enddo
    
    do i=1,natom
      read(10,'(A80)',iostat=ierr) line
      if (ierr.ne.0) exit l1
      read(line,*,iostat=ierr) (x(i,j),j=1,3)
      if (ierr.ne.0) exit l1
      elem(i)=sym(znuc(typatom(i)))
    enddo
    
    write(11,'(I5)') natom
    write(11,'(A)') ''
    x=x*r0
    do i=1,natom
      write(11,'(A2,X,3(F15.8))') elem(i),(x(i,j),j=1,3)
    enddo
      
  enddo l1
  
  close(10)
  close(11)
end program abinit2xmol
- [abinit-forum] Can I generate "prtgeo 1" _GEO output files afterwards?, Rob, 05/25/2009
- Re: [abinit-forum] Can I generate "prtgeo 1" _GEO output files afterwards?, matthieu verstraete, 05/25/2009
 
- <Possible follow-up(s)>
 - Re: [abinit-forum] Can I generate "prtgeo 1" _GEO output files afterwards?, Rob, 05/25/2009
- Re: [abinit-forum] Can I generate "prtgeo 1" _GEO output files afterwards?, Pascal Boulet, 05/25/2009
 
 
 
Archive powered by MHonArc 2.6.15.