TABLE OF CONTENTS
ABINIT/prtposcar [ Functions ]
NAME
prtposcar
FUNCTION
output VASP style POSCAR and FORCES files for use with frozen phonon codes, like PHON from Dario Alfe' or frophon IMPORTANT: the order of atoms is fixed such that typat is re-grouped. First typat=1 then typat=2, etc...
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 .
INPUTS
fcart = forces on atoms in cartesian coordinates natom = number of atoms ntypat = number of types of atoms rprimd = lattice vectors for the primitive cell typat = type for each of the natom atoms ucvol = unit cell volume xred = reduced positions of the atoms znucl = nuclear charge of each atomic type OUTPUTS Only files written
PARENTS
afterscfloop
CHILDREN
atomdata_from_znucl
SOURCE
39 #if defined HAVE_CONFIG_H 40 #include "config.h" 41 #endif 42 43 #include "abi_common.h" 44 45 46 subroutine prtposcar(fcart, fnameradix, natom, ntypat, rprimd, typat, ucvol, xred, znucl) 47 48 use defs_basis 49 use m_profiling_abi 50 use m_atomdata 51 use m_errors 52 53 use m_io_tools, only : open_file 54 55 !This section has been created automatically by the script Abilint (TD). 56 !Do not modify the following lines by hand. 57 #undef ABI_FUNC 58 #define ABI_FUNC 'prtposcar' 59 !End of the abilint section 60 61 implicit none 62 63 !Arguments ------------------------------------ 64 !scalars 65 integer, intent(in) :: natom, ntypat 66 real(dp), intent(in) :: ucvol 67 !arrays 68 integer, intent(in) :: typat(natom) 69 real(dp), intent(in) :: fcart(3,natom) 70 real(dp), intent(in) :: rprimd(3,3) 71 real(dp), intent(in) :: xred(3,natom) 72 real(dp), intent(in) :: znucl(ntypat) 73 character(len=fnlen), intent(in) :: fnameradix 74 75 !Local variables------------------------------- 76 !scalars 77 integer :: iatom, itypat, iout 78 type(atomdata_t) :: atom 79 ! arrays 80 integer :: natoms_this_type(ntypat) 81 character(len=2) :: symbol 82 character(len=7) :: natoms_this_type_str 83 character(len=100) :: chem_formula, natoms_all_types 84 character(len=500) :: msg 85 character(len=fnlen) :: fname 86 87 !************************************************************************ 88 89 !output POSCAR file for positions, atom types etc 90 fname = trim(fnameradix)//"_POSCAR" 91 if (open_file(fname,msg,newunit=iout) /= 0) then 92 MSG_ERROR(msg) 93 end if 94 95 natoms_this_type = 0 96 do itypat=1, ntypat 97 do iatom=1,natom 98 if (typat(iatom) == itypat) natoms_this_type(itypat) = natoms_this_type(itypat)+1 99 end do 100 end do 101 102 chem_formula = "" 103 do itypat=1, ntypat 104 call atomdata_from_znucl(atom,znucl(itypat)) 105 symbol = atom%symbol 106 if (natoms_this_type(itypat) < 10) then 107 write (natoms_this_type_str, '(I1)') natoms_this_type(itypat) 108 else if (natoms_this_type(itypat) < 100) then 109 write (natoms_this_type_str, '(I2)') natoms_this_type(itypat) 110 else if (natoms_this_type(itypat) < 1000) then 111 write (natoms_this_type_str, '(I3)') natoms_this_type(itypat) 112 end if 113 chem_formula = trim(chem_formula) // symbol // trim(natoms_this_type_str) 114 end do 115 write (iout,'(2a)') "ABINIT generated POSCAR file. Title string - should be chemical formula... ",& 116 & trim(chem_formula) 117 118 write (iout,'(E24.14)') -ucvol*Bohr_Ang*Bohr_Ang*Bohr_Ang 119 write (iout,'(3E24.14,1x)') Bohr_Ang*rprimd(:,1) ! (angstr? bohr?) 120 write (iout,'(3E24.14,1x)') Bohr_Ang*rprimd(:,2) 121 write (iout,'(3E24.14,1x)') Bohr_Ang*rprimd(:,3) 122 natoms_all_types = " " 123 do itypat=1, ntypat 124 write (natoms_this_type_str, '(I7)') natoms_this_type(itypat) 125 natoms_all_types = trim(natoms_all_types) // " " // trim(natoms_this_type_str) 126 end do 127 write (iout,'(a)') trim(natoms_all_types) 128 write (iout,'(a)') "Direct" 129 do itypat=1, ntypat 130 do iatom=1,natom 131 if (typat(iatom) /= itypat) cycle 132 write (iout,'(3(E24.14,1x))') xred(:,iatom) 133 end do 134 end do 135 close (iout) 136 137 138 !output FORCES file for forces in same order as positions above 139 fname = trim(fnameradix)//"_FORCES" 140 if (open_file(fname,msg,newunit=iout) /= 0 ) then 141 MSG_ERROR(msg) 142 end if 143 !ndisplacements 144 !iatom_displaced displacement_red_coord(3) 145 !forces_cart_ev_Angstr(3) 146 !... 147 !<repeat for other displaced atoms> 148 write (iout,'(I7)') 1 149 write (iout,'(a)') '1 0 0 0 ! TO BE FILLED IN ' 150 do itypat=1, ntypat 151 do iatom=1,natom 152 if (typat(iatom) /= itypat) cycle 153 write (iout,'(3(E24.14,1x))') Ha_eV/Bohr_Ang*fcart(:,iatom) 154 end do 155 end do 156 close (iout) 157 158 end subroutine prtposcar