TABLE OF CONTENTS
ABINIT/m_atprj [ Modules ]
NAME
m_atprj
FUNCTION
Module to output atomic projections of phonon modes
COPYRIGHT
Copyright (C) 2011-2024 ABINIT group (MJV) This file is distributed under the terms of the GNU General Public Licence, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
SOURCE
18 #if defined HAVE_CONFIG_H 19 #include "config.h" 20 #endif 21 22 #include "abi_common.h" 23 24 module m_atprj 25 26 use defs_basis 27 use m_abicore 28 use m_errors 29 30 use m_io_tools, only : get_unit, open_file 31 use m_fstrings, only : int2char4 32 33 implicit none 34 35 private
m_atprj/atprj_destroy [ Functions ]
[ Top ] [ m_atprj ] [ Functions ]
NAME
atprj_destroy
FUNCTION
deallocate atomic projection datastructure and close files INPUT
OUTPUT
t_atprj = container object for atomic projections
SOURCE
202 subroutine atprj_destroy(t_atprj) 203 204 type(atprj_type), intent(inout) :: t_atprj 205 206 if (allocated(t_atprj%iatprj_bs)) then 207 ABI_FREE(t_atprj%iatprj_bs) 208 end if 209 210 if (allocated(t_atprj%filename)) then 211 ABI_FREE(t_atprj%filename) 212 end if 213 214 end subroutine atprj_destroy 215 216 end module m_atprj
m_atprj/atprj_init [ Functions ]
[ Top ] [ m_atprj ] [ Functions ]
NAME
atprj_init
FUNCTION
initialize atprj datastructure INPUT natom = number of atoms natprj_bs = number of atoms to project on iatprj_bs = indices of atoms to project on outfile_radix = base file name for output files
OUTPUT
t_atprj = container object for atomic projections
SOURCE
83 subroutine atprj_init(t_atprj, natom, natprj_bs, iatprj_bs, outfile_radix) 84 85 type(atprj_type), intent(out) :: t_atprj 86 integer, intent(in) :: natom 87 integer, intent(in) :: natprj_bs 88 integer, intent(in) :: iatprj_bs(natprj_bs) 89 character(len=*), intent(in) :: outfile_radix 90 91 !Local variables------------------------------- 92 !scalars 93 integer :: iatom, imode, iunit 94 character(len=10) :: imodestring, iatomstring 95 character(len=500) :: msg 96 ! ************************************************************************* 97 98 t_atprj%natprj_bs = natprj_bs 99 t_atprj%natom = natom 100 101 ABI_MALLOC(t_atprj%iatprj_bs,(natprj_bs)) 102 t_atprj%iatprj_bs = iatprj_bs 103 104 ! for each phonon mode and atom for projection, open a file 105 ABI_MALLOC(t_atprj%filename ,(3*natom,natprj_bs)) 106 iunit = get_unit() 107 do imode = 1, 3*natom 108 call int2char4(imode, imodestring) 109 ABI_CHECK((imodestring(1:1)/='#'),'Bug: string length too short!') 110 do iatom = 1, natprj_bs 111 call int2char4(iatom, iatomstring) 112 ABI_CHECK((iatomstring(1:1)/='#'),'Bug: string length too short!') 113 t_atprj%filename(imode,iatom) = trim(outfile_radix)//"_mod"//trim(imodestring)//"_iat"//trim(iatomstring) 114 if (open_file(t_atprj%filename(imode,iatom), msg, newunit=iunit, form="formatted", action="write") /= 0) then 115 ABI_ERROR(msg) 116 end if 117 ! print header 118 write (unit=iunit, fmt='(a)') '##' 119 write (unit=iunit, fmt='(a,I6,a)') '## This file contains abinit phonon frequencies for mode number ', & 120 & imode, ' along a path in reciprocal space,' 121 write (unit=iunit, fmt='(a,I6)') '## the third column is the projection along atom number ',& 122 & t_atprj%iatprj_bs(iatom) 123 write (unit=iunit, fmt='(a)') '##' 124 125 close (iunit) 126 end do 127 end do 128 129 end subroutine atprj_init
m_atprj/atprj_print [ Functions ]
[ Top ] [ m_atprj ] [ Functions ]
NAME
atprj_print
FUNCTION
print out 1 line per atomic projection file INPUT t_atprj = container object for atomic projections phfrq = phonon frequencies for present q point eigvec = eigenvectors for present q point
OUTPUT
writes to files
SOURCE
150 subroutine atprj_print(t_atprj, iq, phfrq, eigvec) 151 152 !arguments 153 integer, intent(in) :: iq 154 type(atprj_type), intent(in) :: t_atprj 155 real(dp), intent(in) :: phfrq(3*t_atprj%natom) 156 real(dp), intent(in) :: eigvec(2,3,t_atprj%natom,3,t_atprj%natom) 157 158 !local variables 159 integer :: jatom, idir, iatom, imode, iunit, jdir 160 real(dp) :: proj 161 162 ! ************************************************************************* 163 164 iunit = get_unit() 165 do iatom = 1, t_atprj%natom 166 do idir = 1, 3 167 imode = idir + 3*(iatom-1) 168 do jatom = 1, t_atprj%natprj_bs 169 open (unit=iunit, file=t_atprj%filename(imode,jatom), position='append') 170 write (unit=iunit, fmt='(a,I4,a)') '# atom ', jatom, ' sum of directions' 171 proj = sum(eigvec(:,:,jatom,idir,iatom)**2) 172 write (unit=iunit, fmt='(I6,2E20.10)') iq, phfrq(imode), proj 173 174 do jdir = 1, 3 175 write (unit=iunit, fmt='(2a,I4,a,I4)') ch10, '# atom ', jatom, ' directions ', jdir 176 proj = sum(eigvec(:,jdir,jatom,idir,iatom)**2) 177 write (unit=iunit, fmt='(I6,2E20.10)') iq, phfrq(imode), proj 178 end do 179 close (iunit) 180 end do 181 end do 182 end do 183 184 end subroutine atprj_print
m_atprj/atprj_type [ Types ]
NAME
atprj_type
FUNCTION
Container for atomic projection file data
SOURCE
47 type, public :: atprj_type 48 49 integer :: natprj_bs 50 integer :: natom 51 52 integer, allocatable :: iatprj_bs(:) 53 character(len=fnlen), allocatable :: filename(:,:) 54 55 end type atprj_type 56 57 public :: atprj_init 58 public :: atprj_print 59 public :: atprj_destroy 60 61 contains