TABLE OF CONTENTS
ABINIT/dfpt_prtph [ Functions ]
NAME
dfpt_prtph
FUNCTION
Print the phonon frequencies, on unit 6 as well as the printing unit (except if the associated number -iout- is negative), and for the latter, in Hartree, meV, Thz, Kelvin or cm-1. If eivec==1,2, also print the eigenmodes : displacements in cartesian coordinates. If eivec==4, generate output files for band2eps (drawing tool for the phonon band structure
COPYRIGHT
Copyright (C) 1999-2018 ABINIT group (XG) 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
displ(2,3*natom,3*natom)= contains the displacements of atoms in cartesian coordinates. The first index means either the real or the imaginary part, The second index runs on the direction and the atoms displaced The third index runs on the modes. eivec=(if eivec==0, the eigendisplacements are not printed, if eivec==1,2, the eigendisplacements are printed, if eivec==4, files for band2eps enunit=units for output of the phonon frequencies : 0=> Hartree and cm-1, 1=> eV and Thz, other=> Ha,Thz,eV,cm-1 and K iout= unit for long print (if negative, the routine only print on unit 6, and in Hartree only). natom= number of atom phfreq(3*natom)= phonon frequencies in Hartree qphnrm=phonon wavevector normalisation factor qphon(3)=phonon wavevector
OUTPUT
Only printing
NOTES
called by one processor only
PARENTS
anaddb,m_effective_potential_file,m_ifc,m_phonons,respfn
CHILDREN
wrtout
SOURCE
51 #if defined HAVE_CONFIG_H 52 #include "config.h" 53 #endif 54 55 #include "abi_common.h" 56 57 58 subroutine dfpt_prtph(displ,eivec,enunit,iout,natom,phfrq,qphnrm,qphon) 59 60 use defs_basis 61 use m_profiling_abi 62 use m_errors 63 64 !This section has been created automatically by the script Abilint (TD). 65 !Do not modify the following lines by hand. 66 #undef ABI_FUNC 67 #define ABI_FUNC 'dfpt_prtph' 68 use interfaces_14_hidewrite 69 !End of the abilint section 70 71 implicit none 72 73 !Arguments ------------------------------- 74 !scalars 75 integer,intent(in) :: eivec,enunit,iout,natom 76 real(dp),intent(in) :: qphnrm 77 !arrays 78 real(dp),intent(in) :: displ(2,3*natom,3*natom),phfrq(3*natom),qphon(3) 79 80 !Local variables ------------------------- 81 !scalars 82 integer :: i,idir,ii,imode,jj 83 real(dp) :: tolerance 84 logical :: t_degenerate 85 character(len=500) :: message 86 !arrays 87 real(dp) :: vecti(3),vectr(3) 88 character(len=1) :: metacharacter(3*natom) 89 90 ! ********************************************************************* 91 92 !Check the value of eivec 93 if (all(eivec /= [0,1,2,4])) then 94 write(message, '(a,i0,a,a)' )& 95 & 'In the calling subroutine, eivec is',eivec,ch10,& 96 & 'but allowed values are between 0 and 4.' 97 MSG_BUG(message) 98 end if 99 100 !write the phonon frequencies on unit std_out 101 write(message,'(4a)' )' ',ch10,& 102 & ' phonon wavelength (reduced coordinates) , ','norm, and energies in hartree' 103 call wrtout(std_out,message,'COLL') 104 105 !The next format should be rewritten 106 write(message,'(a,4f5.2)' )' ',(qphon(i),i=1,3),qphnrm 107 call wrtout(std_out,message,'COLL') 108 do jj=1,3*natom,5 109 if (3*natom-jj<5) then 110 write(message,'(5es17.9)') (phfrq(ii),ii=jj,3*natom) 111 else 112 write(message,'(5es17.9)') (phfrq(ii),ii=jj,jj+4) 113 end if 114 call wrtout(std_out,message,'COLL') 115 end do 116 write(message,'(a,es17.9)')' Zero Point Motion energy (sum of freqs/2)=',sum(phfrq(1:3*natom))/2 117 call wrtout(std_out,message,'COLL') 118 119 !Put the wavevector in nice format 120 if(iout>=0)then 121 call wrtout(iout,' ','COLL') 122 if(qphnrm/=0.0_dp)then 123 write(message, '(a,3f9.5)' )& 124 & ' Phonon wavevector (reduced coordinates) :',(qphon(i)/qphnrm+tol10,i=1,3) 125 else 126 write(message, '(3a,3f9.5)' )& 127 & ' Phonon at Gamma, with non-analyticity in the',ch10,& 128 & ' direction (cartesian coordinates)',qphon(1:3)+tol10 129 end if 130 call wrtout(iout,message,'COLL') 131 132 ! Write it, in different units. 133 if(enunit/=1)then 134 write(iout, '(a)' )' Phonon energies in Hartree :' 135 do jj=1,3*natom,5 136 if (3*natom-jj<5) then 137 write(message, '(1x,5es14.6)') (phfrq(ii),ii=jj,3*natom) 138 else 139 write(message, '(1x,5es14.6)') (phfrq(ii),ii=jj,jj+4) 140 end if 141 call wrtout(iout,message,'COLL') 142 end do 143 end if 144 if(enunit/=0)then 145 write(iout, '(a)' )' Phonon energies in meV :' 146 do jj=1,3*natom,5 147 if (3*natom-jj<5) then 148 write(message, '("-",5es14.6)') (phfrq(ii)*Ha_eV*1.0d3,ii=jj,3*natom) 149 else 150 write(message, '("-",5es14.6)') (phfrq(ii)*Ha_eV*1.0d3,ii=jj,jj+4) 151 end if 152 call wrtout(iout,message,'COLL') 153 end do 154 end if 155 if(enunit/=1)then 156 write(iout, '(a)' )' Phonon frequencies in cm-1 :' 157 do jj=1,3*natom,5 158 if (3*natom-jj<5) then 159 write(message, '("-",5es14.6)') (phfrq(ii)*Ha_cmm1,ii=jj,3*natom) 160 else 161 write(message, '("-",5es14.6)') (phfrq(ii)*Ha_cmm1,ii=jj,jj+4) 162 end if 163 call wrtout(iout,message,'COLL') 164 end do 165 end if 166 if(enunit/=0)then 167 write(iout, '(a)' )' Phonon frequencies in Thz :' 168 do jj=1,3*natom,5 169 if (3*natom-jj<5) then 170 write(message, '("-",5es14.6)') (phfrq(ii)*Ha_THz,ii=jj,3*natom) 171 else 172 write(message, '("-",5es14.6)') (phfrq(ii)*Ha_THz,ii=jj,jj+4) 173 end if 174 call wrtout(iout,message,'COLL') 175 end do 176 end if 177 if(enunit/=0.and.enunit/=1)then 178 write(iout, '(a)' )' Phonon energies in Kelvin :' 179 do jj=1,3*natom,5 180 if (3*natom-jj<5) then 181 write(message, '("-",5es14.6)') (phfrq(ii)/kb_HaK,ii=jj,3*natom) 182 else 183 write(message, '("-",5es14.6)') (phfrq(ii)/kb_HaK,ii=jj,jj+4) 184 end if 185 call wrtout(iout,message,'COLL') 186 end do 187 end if 188 end if 189 190 !Take care of the eigendisplacements 191 if(eivec==1 .or. eivec==2)then 192 write(message, '(a,a,a,a,a,a,a,a)' ) ch10,& 193 & ' Eigendisplacements ',ch10,& 194 & ' (will be given, for each mode : in cartesian coordinates',ch10,& 195 & ' for each atom the real part of the displacement vector,',ch10,& 196 & ' then the imaginary part of the displacement vector - absolute values smaller than 1.0d-7 are set to zero)' 197 call wrtout(std_out,message,'COLL') 198 if(iout>=0) then 199 call wrtout(iout,message,'COLL') 200 end if 201 202 ! Examine the degeneracy of each mode. The portability of the echo of the eigendisplacements 203 ! is very hard to obtain, and has not been attempted. 204 do imode=1,3*natom 205 ! The degenerate modes are not portable 206 t_degenerate=.false. 207 if(imode>1)then 208 if(phfrq(imode)-phfrq(imode-1)<tol6)t_degenerate=.true. 209 end if 210 if(imode<3*natom)then 211 if(phfrq(imode+1)-phfrq(imode)<tol6)t_degenerate=.true. 212 end if 213 metacharacter(imode)=';'; if(t_degenerate)metacharacter(imode)='-' 214 end do 215 216 do imode=1,3*natom 217 write(message,'(a,i4,a,es16.6)' )' Mode number ',imode,' Energy',phfrq(imode) 218 call wrtout(std_out,message,'COLL') 219 if(iout>=0)then 220 write(message, '(a,i4,a,es16.6)' )' Mode number ',imode,' Energy',phfrq(imode) 221 call wrtout(iout,message,'COLL') 222 end if 223 tolerance=1.0d-7 224 if(abs(phfrq(imode))<1.0d-5)tolerance=2.0d-7 225 if(phfrq(imode)<1.0d-5)then 226 write(message,'(3a)' )' Attention : low frequency mode.',ch10,& 227 & ' (Could be unstable or acoustic mode)' 228 call wrtout(std_out,message,'COLL') 229 if(iout>=0)then 230 write(iout, '(3a)' )' Attention : low frequency mode.',ch10,& 231 & ' (Could be unstable or acoustic mode)' 232 end if 233 end if 234 do ii=1,natom 235 do idir=1,3 236 vectr(idir)=displ(1,idir+(ii-1)*3,imode) 237 if(abs(vectr(idir))<tolerance)vectr(idir)=0.0_dp 238 vecti(idir)=displ(2,idir+(ii-1)*3,imode) 239 if(abs(vecti(idir))<tolerance)vecti(idir)=0.0_dp 240 end do 241 write(message,'(i4,3es16.8,a,4x,3es16.8)' ) ii,vectr(:),ch10,vecti(:) 242 call wrtout(std_out,message,'COLL') 243 if(iout>=0)then 244 write(message,'(a,i3,3es16.8,2a,3x,3es16.8)') metacharacter(imode),ii,vectr(:),ch10,& 245 & metacharacter(imode), vecti(:) 246 call wrtout(iout,message,'COLL') 247 end if 248 end do 249 end do 250 end if 251 252 end subroutine dfpt_prtph