TABLE OF CONTENTS


ABINIT/append_xyz [ Functions ]

[ Top ] [ Functions ]

NAME

 append_xyz

FUNCTION

 Translate the data from a xyz file (xyz_fname),
 and add it at the end of the usual ABINIT input data string (string),
 taking into account the dtset (dtset_char)

COPYRIGHT

 Copyright (C) 2002-2017 ABINIT group (MJV).
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  dtset_char*2=possible dtset label
  xyz_fname = name of the xyz file
  strln=maximal number of characters of string, as declared in the calling routine

OUTPUT

SIDE EFFECTS

  lenstr=actual number of characters in string
  string*(strln)=string of characters  (upper case) to which the xyz data are appended

PARENTS

      importxyz

CHILDREN

      atomdata_from_symbol,wrtout

SOURCE

 37 #if defined HAVE_CONFIG_H
 38 #include "config.h"
 39 #endif
 40 
 41 #include "abi_common.h"
 42 
 43 subroutine append_xyz(dtset_char,lenstr,string,xyz_fname,strln)
 44 
 45  use defs_basis
 46  use m_errors
 47  use m_profiling_abi
 48  use m_atomdata
 49 
 50  use m_io_tools,       only : open_file
 51 
 52 !This section has been created automatically by the script Abilint (TD).
 53 !Do not modify the following lines by hand.
 54 #undef ABI_FUNC
 55 #define ABI_FUNC 'append_xyz'
 56  use interfaces_14_hidewrite
 57 !End of the abilint section
 58 
 59  implicit none
 60 
 61 !Arguments ------------------------------------
 62 !scalars
 63  integer,intent(in) :: strln
 64  integer,intent(inout) :: lenstr
 65  character(len=2),intent(in) :: dtset_char
 66  character(len=fnlen),intent(in) :: xyz_fname
 67  character(len=strln),intent(inout) :: string
 68 
 69 !Local variables-------------------------------
 70  character :: blank=' '
 71 !scalars
 72  integer :: unitxyz, iatom, natom, mu
 73  integer :: lenstr_new
 74  integer :: lenstr_old
 75  integer :: ntypat
 76  real(dp) :: znucl
 77  character(len=5) :: string5
 78  character(len=20) :: string20
 79  character(len=500) :: message
 80  type(atomdata_t) :: atom
 81 !arrays
 82  real(dp),allocatable :: xangst(:,:)
 83  integer, save :: atomspecies(200) = 0
 84  character(len=500), save :: znuclstring = ""
 85  character(len=2),allocatable :: elementtype(:)
 86 
 87 !************************************************************************
 88 
 89  lenstr_new=lenstr
 90 
 91  if (dtset_char == "-1") then
 92 !  write znucl
 93    lenstr_old=lenstr_new
 94    lenstr_new=lenstr_new+7+len_trim(znuclstring)+1
 95    string(lenstr_old+1:lenstr_new)=" ZNUCL"//blank//trim(znuclstring)//blank
 96 
 97 !  write ntypat
 98    ntypat = sum(atomspecies)
 99    write(string20,'(i10)') ntypat
100    lenstr_old=lenstr_new
101    lenstr_new=lenstr_new+8+len_trim(string20)+1
102    string(lenstr_old+1:lenstr_new)=" NTYPAT"//blank//trim(string20)//blank
103 
104    return
105  end if
106 
107 !open file with xyz data
108  if (open_file(xyz_fname, message, newunit=unitxyz, status="unknown") /= 0) then
109    MSG_ERROR(message)
110  end if
111  write(message, '(3a)')' importxyz : Opened file ',trim(xyz_fname),'; content stored in string_xyz'
112  call wrtout(std_out,message,'COLL')
113 
114 !check number of atoms is correct
115  read(unitxyz,*) natom
116  
117  write(string5,'(i5)')natom
118  lenstr_old=lenstr_new
119  lenstr_new=lenstr_new+7+len_trim(dtset_char)+1+5
120  string(lenstr_old+1:lenstr_new)=" _NATOM"//trim(dtset_char)//blank//string5
121 
122  ABI_ALLOCATE(xangst,(3,natom))
123  ABI_ALLOCATE(elementtype,(natom))
124 
125 !read dummy line
126  read(unitxyz,*)
127 
128 !read atomic types and positions
129  do iatom = 1, natom
130    read(unitxyz,*) elementtype(iatom), xangst(:,iatom)
131 !  extract znucl for each atom type
132    call atomdata_from_symbol(atom,elementtype(iatom))
133    znucl = atom%znucl
134    if (znucl > 200) then
135      write (message,'(5a)')&
136 &     'Error: found element beyond Z=200 ', ch10,&
137 &     'Solution: increase size of atomspecies in append_xyz', ch10
138      MSG_ERROR(message)
139    end if
140 !  found a new atom type 
141    if (atomspecies(int(znucl)) == 0) then
142      write(string20,'(f10.2)') znucl
143      znuclstring = trim(znuclstring) // " " // trim(string20) // " "
144    end if
145    atomspecies(int(znucl)) = 1
146  end do
147  close (unitxyz)
148 
149 
150 !Write the element types
151  lenstr_old=lenstr_new
152  lenstr_new=lenstr_new+7+len_trim(dtset_char)+1
153  string(lenstr_old+1:lenstr_new)=" _TYPAX"//trim(dtset_char)//blank
154  do iatom=1,natom
155    lenstr_old=lenstr_new
156    lenstr_new=lenstr_new+3
157    string(lenstr_old+1:lenstr_new)=elementtype(iatom)//blank
158  end do
159  lenstr_old=lenstr_new
160  lenstr_new=lenstr_new+3
161  string(lenstr_old+1:lenstr_new)="XX " ! end card for TYPAX
162 
163 !Write the coordinates
164  lenstr_old=lenstr_new
165  lenstr_new=lenstr_new+8+len_trim(dtset_char)+1
166  string(lenstr_old+1:lenstr_new)=" _XANGST"//trim(dtset_char)//blank
167 
168  do iatom=1,natom
169    do mu=1,3
170      write(string20,'(f20.12)')xangst(mu,iatom)
171      lenstr_old=lenstr_new
172      lenstr_new=lenstr_new+20
173      string(lenstr_old+1:lenstr_new)=string20
174    end do
175  end do
176 
177  ABI_DEALLOCATE(elementtype)
178  ABI_DEALLOCATE(xangst)
179 
180 !Check the length of the string
181  if(lenstr_new>strln)then
182    write(message,'(3a)')&
183 &   'The maximal size of the input variable string has been exceeded.',ch10,&
184 &   'The use of a xyz file is more character-consuming than the usual input file. Sorry.'
185    MSG_BUG(message)
186  end if
187 
188 !Update the length of the string
189  lenstr=lenstr_new
190 
191 end subroutine append_xyz