TABLE OF CONTENTS


ABINIT/get_all_gkq [ Functions ]

[ Top ] [ Functions ]

NAME

 get_all_gkq

FUNCTION

 This routine determines what to do with the initial qspace
   matrix elements of the electron phonon coupling (to disk or in memory),
   then reads those given in the gkk file and completes them
   (for kpts, then perturbations)
   01/2010: removed completion on qpoints here (MJV)

COPYRIGHT

 Copyright (C) 2004-2018 ABINIT group (MVer, MG)
 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
   Cryst<crystal_t>=Info on the unit cell and on its symmetries.
   Ifc<ifc_type>=Object containing the interatomic force constants.
   Bst<ebands_t>=GS energies, occupancies and Fermi level.
   FSfullpqtofull = mapping of k+q to another k
   kphon_full2full = mapping of FS kpoints under symops
   kpt_phon = fermi surface kpoints
   %k_phon%wtk = integration weights for bands and kpoints near the FS
   gkk_flag = flag to
   nband = number of bands
   n1wf = number of file headers from perturbation calculations
      which are present in the initial gkk input file.
   onegkksize = size of one record of the new gkk output file, in bytes
   qpttoqpt = mapping of qpoints onto each other under symmetries
   unitgkk = fortran unit for initial gkk input file
   xred = reduced coordinates of atoms

OUTPUT

   elph_ds%gkq = recip space elphon matrix elements.

PARENTS

      elphon

CHILDREN

      complete_gkk,int2char4,read_gkk,wrtout

SOURCE

 50 #if defined HAVE_CONFIG_H
 51 #include "config.h"
 52 #endif
 53 
 54 #include "abi_common.h"
 55 
 56 subroutine get_all_gkq (elph_ds,Cryst,ifc,Bst,FSfullpqtofull,nband,n1wf,onegkksize,&
 57 &    qpttoqpt,ep_prt_yambo,unitgkk,ifltransport)
 58 
 59  use defs_basis
 60  use defs_datatypes
 61  use defs_elphon
 62  use m_profiling_abi
 63  use m_errors
 64  use m_io_tools
 65  use m_xmpi
 66 
 67  use m_fstrings, only : int2char4
 68  use m_crystal,  only : crystal_t
 69  use m_ifc,      only : ifc_type
 70 
 71 !This section has been created automatically by the script Abilint (TD).
 72 !Do not modify the following lines by hand.
 73 #undef ABI_FUNC
 74 #define ABI_FUNC 'get_all_gkq'
 75  use interfaces_14_hidewrite
 76  use interfaces_77_ddb, except_this_one => get_all_gkq
 77 !End of the abilint section
 78 
 79  implicit none
 80 
 81 !Arguments ------------------------------------
 82 !scalars
 83  integer,intent(in) :: n1wf,nband,onegkksize,unitgkk,ep_prt_yambo,ifltransport
 84  type(crystal_t),intent(in) :: Cryst
 85  type(ifc_type),intent(in) :: ifc
 86  type(ebands_t),intent(in) :: Bst
 87  type(elph_type),intent(inout) :: elph_ds
 88 !arrays
 89  integer,intent(in) :: FSfullpqtofull(elph_ds%k_phon%nkpt,elph_ds%nqpt_full)
 90  integer,intent(in) :: qpttoqpt(2,Cryst%nsym,elph_ds%nqpt_full)
 91 
 92 !Local variables-------------------------------
 93 !scalars
 94  integer :: iost,ierr,me,sz2,sz3,sz4,sz5,sz6
 95  character(len=10) :: procnum
 96  character(len=500) :: message
 97  character(len=fnlen) :: fname
 98 !arrays
 99  integer,allocatable :: gkk_flag(:,:,:,:,:) 
100 
101 ! *************************************************************************
102 
103 !attribute file unit number
104  elph_ds%unitgkq = get_unit()
105 
106 !============================================
107 !save gkk for all qpts in memory or to disk
108 !============================================
109 
110 !DEBUG
111 !write(std_out,*) ' 4 bytes / ??'
112 !write(std_out,*) ' kind(real) = ', kind(one)
113 !write(std_out,*) ' elph_ds%ngkkband = ', elph_ds%ngkkband, '^2'
114 !write(std_out,*) ' elph_ds%nbranch = ', elph_ds%nbranch, '^2'
115 !write(std_out,*) ' elph_ds%k_phon%nkpt = ', elph_ds%k_phon%nkpt
116 !write(std_out,*) ' elph_ds%nsppol = ', elph_ds%nsppol
117 !write(std_out,*) ' elph_ds%nqptirred ', elph_ds%nqptirred
118 !ENDDEBUG
119 
120  write(message,'(a,f14.4,a)')&
121 & ' get_all_gkq : gkq file/array size = ',&
122  4.0*dble(onegkksize)*dble(elph_ds%k_phon%my_nkpt)*dble(elph_ds%nqptirred)/1024.0_dp/1024.0_dp/1024.0_dp,' Gb'
123  call wrtout(std_out,message,'COLL')
124 
125  if (elph_ds%gkqwrite == 0) then !calculate gkk(q) keeping all in memory
126 
127    call wrtout(std_out,' get_all_gkq : keep gkk(q) in memory ','COLL')
128 
129    sz2=elph_ds%ngkkband*elph_ds%ngkkband
130    sz3=elph_ds%nbranch*elph_ds%nbranch
131    sz4=elph_ds%k_phon%my_nkpt
132    sz5=elph_ds%nsppol
133    if (ifltransport == 3) then
134      sz6=elph_ds%nqpt_full
135    else
136      sz6=elph_ds%nqptirred
137    end if
138    ABI_STAT_ALLOCATE(elph_ds%gkk_qpt,(2,sz2,sz3,sz4,sz5,sz6), ierr)
139    ABI_CHECK(ierr==0, 'Trying to allocate array elph_ds%gkk_qpt')
140 
141    elph_ds%gkk_qpt = zero
142 
143  else if (elph_ds%gkqwrite == 1) then !calculate gkk(q) and write to file
144    me = xmpi_comm_rank(xmpi_world)
145    call int2char4(me,procnum)
146    ABI_CHECK((procnum(1:1)/='#'),'Bug: string length too short!')
147    fname=trim(elph_ds%elph_base_name) // "_P" // trim(procnum) // '_GKKQ'
148 
149    iost=open_file(file=fname,iomsg=message,newunit=elph_ds%unitgkq,access='direct',&
150 &   recl=onegkksize,form='unformatted')
151    if (iost /= 0) then
152      write (message,'(2a)')' get_all_gkq : ERROR- opening file ',trim(fname)
153      MSG_ERROR(message)
154    end if
155 
156    write (message,'(5a)')&
157 &   ' get_all_gkq : gkq matrix elements  will be written to file : ',trim(fname),ch10,&
158 &   ' Nothing is in files yet',ch10
159    call wrtout(std_out,message,'COLL')
160 
161  else
162    write(message,'(a,i0)')' gkqwrite must be 0 or 1 while it is : ',elph_ds%gkqwrite
163    MSG_BUG(message)
164  end if !if gkqwrite
165 
166 !=====================================================
167 !read in g_kk matrix elements for all bands, kpoints,
168 !and calculated qpoints
169 !=====================================================
170  call wrtout(std_out,' get_all_gkq : calling read_gkk to read in the g_kk matrix elements',"COLL")
171 
172  sz2=elph_ds%nbranch;sz3=elph_ds%k_phon%my_nkpt
173  sz4=elph_ds%nsppol;sz5=elph_ds%nqpt_full
174  ABI_STAT_ALLOCATE(gkk_flag,(sz2,sz2,sz3,sz4,sz5), ierr)
175  ABI_CHECK(ierr==0, "allocating gkk_flag")
176 
177  call read_gkk(elph_ds,Cryst,ifc,Bst,FSfullpqtofull,gkk_flag,n1wf,nband,ep_prt_yambo,unitgkk)
178 
179 !if (elph_ds%symgkq ==1) then
180 !MJV 01/2010 removed the completion on qpt here: it should be done after FS integration
181 !so that everything is lighter in memory etc... (only irred qpt)
182 ! if (0==1) then
183  if (ifltransport == 3) then !  bxu, complete gkk is necessary
184 
185 !  ==============================================================
186 !  complete gkk matrices for other qpoints on the full grid qpt_full
187 !  inspired and cannibalized from symdm9.f
188 !  FIXME: should add the possibility to copy over to other qpoints,
189 !  without full symmetrization, for testing purposes.
190 !  ==============================================================
191 
192    write(message,'(4a)')ch10,&
193 &   ' get_all_gkq : calling complete_gkk to complete ',ch10,&
194 &   ' gkk matrices for other qpoints on the full grid'
195    call wrtout(std_out,message,'COLL')
196 
197    call complete_gkk(elph_ds,gkk_flag,Cryst%gprimd,Cryst%indsym,&
198 &   Cryst%natom,Cryst%nsym,qpttoqpt,Cryst%rprimd,Cryst%symrec,Cryst%symrel)
199 
200    call wrtout(std_out,' get_all_gkq : out of complete_gkk','COLL')
201 
202  end if !symgkq
203 
204 !TODO Do we need gkk_flag in elphon?
205  ABI_DEALLOCATE(gkk_flag)
206 
207 end subroutine get_all_gkq