TABLE OF CONTENTS


ABINIT/get_all_gkr [ Functions ]

[ Top ] [ Functions ]

NAME

 get_all_gkr

FUNCTION

 This routine determines what to do with the rspace
   matrix elements of the el phon coupling (to disk or in memory),
   then reads those given in the gkq file and Fourier Transforms them

COPYRIGHT

 Copyright (C) 2004-2018 ABINIT group (MVer)
 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 .

INPUTS

   elph_ds = elphon datastructure with data and dimensions
   gprim = reciprocal space lattice vectors
   natom = number of atoms
   nrpt = number of real-space points used for FT
   onegkksize = size of one record of the new gkk output file, in bytes
   rpt = positions of real-space points for FT
   qpt_full = qpoint coordinates
   wghatm = weights for real-space rpt in FT

OUTPUT

   elph_ds%gkr = real space elphon matrix elements.

PARENTS

      elphon

CHILDREN

      ftgkk

SOURCE

 40 #if defined HAVE_CONFIG_H
 41 #include "config.h"
 42 #endif
 43 
 44 #include "abi_common.h"
 45 
 46 
 47 subroutine get_all_gkr (elph_ds,gprim,natom,nrpt,onegkksize,rpt,qpt_full,wghatm)
 48 
 49  use defs_basis
 50  use defs_elphon
 51  use m_profiling_abi
 52  use m_errors
 53 
 54 !This section has been created automatically by the script Abilint (TD).
 55 !Do not modify the following lines by hand.
 56 #undef ABI_FUNC
 57 #define ABI_FUNC 'get_all_gkr'
 58  use interfaces_77_ddb, except_this_one => get_all_gkr
 59 !End of the abilint section
 60 
 61  implicit none
 62 
 63 !Arguments ------------------------------------
 64 !scalars
 65  integer,intent(in) :: natom,nrpt,onegkksize
 66  type(elph_type),intent(inout) :: elph_ds
 67 !arrays
 68  real(dp),intent(in) :: gprim(3,3),rpt(3,nrpt),qpt_full(3,elph_ds%nqpt_full)
 69  real(dp),intent(in) :: wghatm(natom,natom,nrpt)
 70 
 71 !Local variables-------------------------------
 72 !scalars
 73  integer :: ikpt_phon0,iost,qtor,sz2,sz3,sz4,sz5
 74 
 75 ! *************************************************************************
 76 
 77 !
 78 !WARNING : disk file used for large arrays gkk_rpt and
 79 !(eventually) gkk2
 80 !
 81 !allocate (gkk_rpt(2,elph_ds%nbranch,elph_ds%nFSband,elph_ds%nFSband,&
 82 !&  elph_ds%k_phon%nkpt,nrpt))
 83  elph_ds%unit_gkk_rpt = 36
 84 !see if the gkk_rpt should be written to a file (only available option now)
 85  if (elph_ds%gkk_rptwrite == 1) then
 86 !  file is not present : we need to do the FT
 87    open (unit=elph_ds%unit_gkk_rpt,file='gkk_rpt_file',access='direct',&
 88 &   recl=onegkksize,form='unformatted',&
 89 &   status='new',iostat=iost)
 90    if (iost /= 0) then
 91      MSG_ERROR('get_all_gkr : error opening gkk_rpt_file as new')
 92    end if
 93    write(std_out,*) ' get_all_gkr : will write real space gkk to a disk file.'
 94    write(std_out,*) ' size = ', 4.0*dble(onegkksize)*dble(nrpt)/&
 95 &   1024.0_dp/1024.0_dp, ' Mb'
 96 
 97 !  else if (elph_ds%gkk_rptwrite  == 0) then
 98  else
 99    write(std_out,*) ' get_all_gkr : will keep real space gkk in memory.'
100    write(std_out,*) ' size = ', 4.0*dble(onegkksize)*dble(nrpt)/&
101 &   1024.0_dp/1024.0_dp, ' Mb'
102    sz2=elph_ds%ngkkband*elph_ds%ngkkband
103    sz3=elph_ds%nbranch*elph_ds%nbranch
104    sz4=elph_ds%k_phon%nkpt
105    sz5=elph_ds%nsppol
106    ABI_ALLOCATE(elph_ds%gkk_rpt,(2,sz2,sz3,sz4,sz5,nrpt))
107 !  write(std_out,*) ' get_all_gkr: invalid value for gkk_rptwrite'
108 !  stop
109  end if
110  write(std_out,*) '    about to FT the recip space gkk to real space '
111  qtor = 1
112 
113 !
114 !NOTE: should be very easy to parallelize!
115 !
116  ikpt_phon0 = 1
117  call ftgkk (wghatm,elph_ds%gkk_qpt,elph_ds%gkk_rpt,&
118 & elph_ds%gkqwrite,elph_ds%gkk_rptwrite,gprim,1,natom,&
119 & elph_ds%k_phon%nkpt,elph_ds%ngkkband,elph_ds%k_phon%nkpt,elph_ds%nqpt_full,&
120 & nrpt,elph_ds%nsppol,qtor,rpt,qpt_full,elph_ds%unit_gkk_rpt,elph_ds%unitgkq)
121 
122 !call ftgkk (elph_ds,gprim,ikpt_phon0,natom,nrpt,qtor,rpt,qpt_full,wghatm)
123  write(std_out,*) ' get_all_gkr : done with FT of gkk to real space'
124 
125 !No longer need the gkk_qpt?
126 !if (elph_ds%gkqwrite == 0) deallocate (elph_ds%gkk_qpt)
127 
128 !!DEBUG
129 !Test the FT of the gkk elements.
130 !call test_ftgkk(elph_ds,gprim,natom,nrpt,rpt,qpt_full,wghatm)
131 !!ENDDEBUG
132 
133 !DEBUG
134 !do irpt=1,nrpt
135 !do ipert1=1,elph_ds%nbranch
136 !write(std_out,'(6(F16.5,1x))') elph_ds%gkk_rpt(:,ipert1,1,1,1,irpt)
137 !end do
138 !end do
139 !ENDDEBUG
140 
141 
142 end subroutine get_all_gkr