TABLE OF CONTENTS


ABINIT/prt_cif [ Functions ]

[ Top ] [ Functions ]

NAME

 prt_cif

FUNCTION

   print out CIF format file, from abinit internal data

COPYRIGHT

 Copyright (C) 2010-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

OUTPUT

NOTES

PARENTS

      outscfcv

CHILDREN

SOURCE

 29 #if defined HAVE_CONFIG_H
 30 #include "config.h"
 31 #endif
 32 
 33 #include "abi_common.h"
 34 
 35 
 36 subroutine prt_cif(brvltt, ciffname, natom, nsym, ntypat, rprimd, &
 37 &   spgaxor, spgroup, spgorig, symrel, tnon, typat, xred, znucl)
 38 
 39  use defs_basis
 40  use m_profiling_abi
 41  use m_errors
 42  use m_atomdata
 43 
 44  use m_io_tools,       only : open_file
 45  use m_fstrings,       only : int2char10
 46 
 47 !This section has been created automatically by the script Abilint (TD).
 48 !Do not modify the following lines by hand.
 49 #undef ABI_FUNC
 50 #define ABI_FUNC 'prt_cif'
 51  use interfaces_41_geometry, except_this_one => prt_cif
 52 !End of the abilint section
 53 
 54  implicit none
 55 
 56 !Arguments ------------------------------------
 57 !scalars
 58  integer,intent(in) :: natom, ntypat, nsym
 59  integer, intent(in) :: brvltt, spgaxor, spgroup, spgorig
 60 !arrays
 61  integer, intent(in) :: typat(natom)
 62  integer, intent(in) :: symrel(3,3,nsym)
 63  character(len=*), intent(in) :: ciffname
 64  real(dp), intent(in) :: tnon(3,nsym)
 65  real(dp), intent(in) :: rprimd(3,3)
 66  real(dp), intent(in) :: xred(3,natom)
 67  real(dp), intent(in) :: znucl(ntypat)
 68 
 69 !Local variables -------------------------------
 70 !scalars
 71  integer :: unitcif
 72  integer :: iatom, isym
 73  integer :: sporder
 74  integer :: itypat, nat_this_type
 75  real(dp) :: ucvol
 76  type(atomdata_t) :: atom
 77 
 78 !arrays
 79  character(len=80) :: tmpstring
 80 
 81  character(len=1) :: brvsb
 82  character(len=15) :: intsb,ptintsb,ptschsb,schsb
 83  character(len=35) :: intsbl
 84 
 85  character(len=10) :: str_nat_type
 86  character(len=100) :: chemformula
 87  character(len=500) :: msg
 88 
 89  real(dp) :: angle(3)
 90  real(dp) :: gprimd(3,3)
 91  real(dp) :: rmet(3,3), gmet(3,3)
 92 
 93 !*************************************************************************
 94 
 95 !open file in append mode xlf and other compilers refuse append mode
 96  if (open_file(ciffname,msg,newunit=unitcif) /=0) then
 97    MSG_WARNING(msg)
 98    return
 99  end if
100 
101 !print title for dataset
102  write (unitcif,'(a)') 'data_set'
103 
104 !print cell parameters a,b,c, angles, volume
105  call metric(gmet,gprimd,-1,rmet,rprimd,ucvol)
106  angle(1)=acos(rmet(2,3)/sqrt(rmet(2,2)*rmet(3,3)))/two_pi*360.0_dp
107  angle(2)=acos(rmet(1,3)/sqrt(rmet(1,1)*rmet(3,3)))/two_pi*360.0_dp
108  angle(3)=acos(rmet(1,2)/sqrt(rmet(1,1)*rmet(2,2)))/two_pi*360.0_dp
109 
110  write (unitcif,'(a,E20.10)') '_cell_length_a                     ', sqrt(rmet(1,1))*Bohr_Ang
111  write (unitcif,'(a,E20.10)') '_cell_length_b                     ', sqrt(rmet(2,2))*Bohr_Ang
112  write (unitcif,'(a,E20.10)') '_cell_length_c                     ', sqrt(rmet(3,3))*Bohr_Ang
113  write (unitcif,'(a,E20.10)') '_cell_angle_alpha                  ', angle(1)
114  write (unitcif,'(a,E20.10)') '_cell_angle_beta                   ', angle(2)
115  write (unitcif,'(a,E20.10)') '_cell_angle_gamma                  ', angle(3)
116  write (unitcif,'(a,E20.10)') '_cell_volume                       ', ucvol*(Bohr_Ang)**3
117 
118 !print reduced positions
119  write (unitcif,'(a)') 'loop_'
120  write (unitcif,'(a,E20.10)') '  _atom_site_label                   '
121  write (unitcif,'(a,E20.10)') '  _atom_site_fract_x                 '
122  write (unitcif,'(a,E20.10)') '  _atom_site_fract_y                 '
123  write (unitcif,'(a,E20.10)') '  _atom_site_fract_z                 ' 
124  do iatom = 1, natom
125    call atomdata_from_znucl(atom,znucl(typat(iatom)))
126    write (unitcif,'(2a,3E20.10)') '  ', atom%symbol, xred(:,iatom)
127  end do
128 
129 !
130 !other specs in CIF dictionary which may be useful:
131 !GEOM_BOND GEOM_ANGLE GEOM_TORSION
132 !
133 
134 !print chemical composition in simplest form
135  chemformula = "'"
136  do itypat = 1, ntypat
137    nat_this_type = 0
138    do iatom = 1, natom
139      if (typat(iatom) == itypat) nat_this_type = nat_this_type+1
140    end do
141    call atomdata_from_znucl(atom,znucl(itypat))
142    call int2char10(nat_this_type, str_nat_type)
143    chemformula = trim(chemformula) // atom%symbol // trim(str_nat_type) // "  " 
144  end do
145  chemformula = trim(chemformula) // "'"
146  write (unitcif,'(2a)') '_chemical_formula_analytical              ', chemformula
147 
148 !FIXME: check that brvltt is correctly used here - is it equal to bravais(1) in the invars routines?
149  if     (brvltt==1)then 
150    write (unitcif,'(a)') '_symmetry_cell_setting             triclinic'
151  else if(brvltt==2)then 
152    write (unitcif,'(a)') '_symmetry_cell_setting             monoclinic'
153  else if(brvltt==3)then 
154    write (unitcif,'(a)') '_symmetry_cell_setting             orthorhombic'
155  else if(brvltt==4)then 
156    write (unitcif,'(a)') '_symmetry_cell_setting             tetragonal'
157  else if(brvltt==5)then 
158    write (unitcif,'(a)') '_symmetry_cell_setting             rhombohedral'
159  else if(brvltt==6)then 
160    write (unitcif,'(a)') '_symmetry_cell_setting             hexagonal'
161  else if(brvltt==7)then 
162    write (unitcif,'(a)') '_symmetry_cell_setting             cubic'
163  end if
164 
165  call spgdata(brvsb,intsb,intsbl,ptintsb,ptschsb,&
166 & schsb,spgaxor,spgroup,sporder,spgorig)
167 
168 !print symmetry operations
169  write (unitcif,'(a,I6)') "_symmetry_Int_Tables_number          ", spgroup
170  write (unitcif,'(5a)') "_symmetry_space_group_name_H-M        '", brvsb, " ", trim(intsb), "'"
171  write (unitcif,'(a)') ''
172  write (unitcif,'(a)') 'loop_'
173  write (unitcif,'(a)') '  _symmetry_equiv_pos_as_xyz           '
174  do isym = 1, nsym
175    call  symrel2string(symrel(:,:,isym), tnon(:,isym), tmpstring)
176    write (unitcif,'(2a)') '  ', trim(tmpstring)
177  end do
178 
179  close (unitcif)
180 
181 end subroutine prt_cif

ABINIT/symrel2string [ Functions ]

[ Top ] [ Functions ]

NAME

 symrel2string

FUNCTION

COPYRIGHT

 Copyright (C) 2010-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

OUTPUT

NOTES

PARENTS

      prt_cif

CHILDREN

SOURCE

211 subroutine symrel2string(symrel1, tnon, string)
212 
213  use defs_basis
214  use m_profiling_abi
215 
216 !This section has been created automatically by the script Abilint (TD).
217 !Do not modify the following lines by hand.
218 #undef ABI_FUNC
219 #define ABI_FUNC 'symrel2string'
220 !End of the abilint section
221 
222  implicit none
223  
224 !Arguments ------------------------------------
225 !scalars
226 !arrays
227  integer, intent(in) :: symrel1(3,3)
228  real(dp), intent(in) :: tnon(3)
229  character(len=80), intent(out) :: string
230 
231 !Local variables -------------------------
232 !scalars
233  integer :: i1,i2
234  character(len=1) :: xyz(3)
235 
236 ! *********************************************************************
237 
238  xyz(1) = 'x'
239  xyz(2) = 'y'
240  xyz(3) = 'z'
241 
242  string = ''
243  do i1=1,3
244    if (abs(tnon(i1)) > tol10) then
245 !    find fraction 1/n for tnon, otherwise do not know what to print
246      if (abs(one-two*tnon(i1)) < tol10) string = trim(string)//'1/2'
247      if (abs(one+two*tnon(i1)) < tol10) string = trim(string)//'-1/2'
248 
249      if (abs(one-three*tnon(i1)) < tol10) string = trim(string)//'1/3'
250      if (abs(one+three*tnon(i1)) < tol10) string = trim(string)//'-1/3'
251      if (abs(two-three*tnon(i1)) < tol10) string = trim(string)//'2/3'
252      if (abs(two+three*tnon(i1)) < tol10) string = trim(string)//'-2/3'
253 
254      if (abs(one-six*tnon(i1)) < tol10) string = trim(string)//'1/6'
255      if (abs(one+six*tnon(i1)) < tol10) string = trim(string)//'-1/6'
256      if (abs(five-six*tnon(i1)) < tol10) string = trim(string)//'5/6'
257      if (abs(five+six*tnon(i1)) < tol10) string = trim(string)//'-5/6'
258    end if
259    do i2=1,3
260 !    FIXME: check if this is correct ordering for symrel(i1,i2) looks ok
261      if (symrel1(i1,i2) == 1)  string = trim(string)//'+'//xyz(i2)
262      if (symrel1(i1,i2) == -1) string = trim(string)//'-'//xyz(i2)
263    end do
264    if (i1 /= 3) string = trim(string)//','
265  end do
266 
267 end subroutine symrel2string