TABLE OF CONTENTS


ABINIT/get_all_gkk2 [ Functions ]

[ Top ] [ Functions ]

NAME

 get_all_gkk2

FUNCTION

 This routine determines where to store gkk2 matrix elements (disk or RAM)
   and calls interpolate_gkk to calculate them.
   This is the most time consuming step.

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

   acell = lengths of unit cell vectors
   amu = masses of atoms
   atmfrc = atomic force constants
   dielt = dielectric tensor
   dipdip = dipole-dipole contribution flag
   dyewq0 =
   elph_ds = datastructure for elphon data and dimensions
   kptirr_phon = irreducible set of fermi-surface kpoints
   kpt_phon = full set of fermi-surface kpoints
   ftwghtgkk = weights for FT of matrix elements
   gmet = metric in reciprocal space
   indsym = indirect mapping of atoms under symops
   mpert = maximum number of perturbations
   msym = maximum number of symmetries (usually nsym)
   nsym = number of symmetries
   ntypat = number of types of atoms
   onegkksize = size of one gkk record, in bytes
   rmet = real-space metric
   rprim = unit cell lattice vectors (dimensionless)
   rprimd = real-space unit-cell lattice vectors
   rpt = points in real space for FT, in canonical coordinates
   symrel = symmetry operations in reduced real space
   trans = Atomic translations : xred = rcan + trans
   typat = array of types of atoms
   ucvol = unit cell volume
   xred = reduced coordinates of atoms
   zeff = Born effective charges

OUTPUT

   elph_ds = calculated |gkk|^2 are in elph_ds%gkk2

PARENTS

      elphon

CHILDREN

      interpolate_gkk

SOURCE

 59 #if defined HAVE_CONFIG_H
 60 #include "config.h"
 61 #endif
 62 
 63 #include "abi_common.h"
 64 
 65 
 66 subroutine get_all_gkk2(crystal,ifc,elph_ds,kptirr_phon,kpt_phon)
 67 
 68  use defs_basis
 69  use defs_elphon
 70  use m_errors
 71  use m_profiling_abi
 72 
 73  use m_crystal,   only : crystal_t
 74  use m_ifc,       only : ifc_type
 75 
 76 !This section has been created automatically by the script Abilint (TD).
 77 !Do not modify the following lines by hand.
 78 #undef ABI_FUNC
 79 #define ABI_FUNC 'get_all_gkk2'
 80  use interfaces_77_ddb, except_this_one => get_all_gkk2
 81 !End of the abilint section
 82 
 83  implicit none
 84 
 85 !Arguments ------------------------------------
 86 !scalars
 87  type(crystal_t),intent(in) :: crystal
 88  type(ifc_type),intent(in) :: ifc
 89  type(elph_type),intent(inout) :: elph_ds
 90 !arrays
 91  real(dp),intent(in) :: kpt_phon(3,elph_ds%k_phon%nkpt)
 92  real(dp),intent(in) :: kptirr_phon(3,elph_ds%k_phon%nkptirr)
 93 
 94 !Local variables-------------------------------
 95 !scalars
 96  integer :: iost,onediaggkksize,sz1,sz2,sz3,sz4
 97  real(dp) :: realdp_ex
 98  !character(len=500) :: msg
 99 
100 ! *************************************************************************
101 
102  if (elph_ds%nsppol /= 1) then
103    MSG_ERROR('get_all_gkk2: nsppol>1 not coded yet!')
104  end if
105 
106  onediaggkksize = elph_ds%nbranch*elph_ds%k_phon%nkpt*kind(realdp_ex)
107 
108  elph_ds%unit_gkk2 = 37
109  if (elph_ds%gkk2write == 0) then
110    write(std_out,*) 'get_all_gkk2 : keep gkk2 in memory. Size = ',&
111 &   4.0*dble(elph_ds%k_phon%nkpt)*dble(onediaggkksize)/&
112 &   1024.0_dp/1024.0_dp, " Mb"
113    sz1=elph_ds%nbranch
114    sz2=elph_ds%ngkkband
115    sz3=elph_ds%ngkkband
116    sz4=elph_ds%k_phon%nkpt
117    ABI_ALLOCATE(elph_ds%gkk2,(sz1,sz2,sz3,sz4,elph_ds%k_phon%nkpt,1))
118    elph_ds%gkk2(:,:,:,:,:,:) = zero
119 
120  else if (elph_ds%gkk2write == 1) then
121    write(std_out,*) 'get_all_gkk2 : About to open gkk2 file : '
122    write(std_out,*) elph_ds%unit_gkk2,onediaggkksize
123    open (unit=elph_ds%unit_gkk2,file='gkk2file',access='direct',&
124 &   recl=onediaggkksize,form='unformatted',status='new',iostat=iost)
125    if (iost /= 0) then
126      MSG_ERROR('error opening gkk2file as new')
127    end if
128 !  rewind (elph_ds%unit_gkk2)
129    write(std_out,*) 'get_all_gkk2 : disk file with gkk^2 created'
130    write(std_out,*) '  calculate from real space gkk and phonon modes'
131    write(std_out,*) '  gkk2write = 1 is forced: can take a lot of time! '
132    write(std_out,*) ' size = ', 4.0*dble(onediaggkksize)*dble(elph_ds%k_phon%nkpt)/&
133 &   1024.0_dp/1024.0_dp, ' Mb'
134  else
135    MSG_ERROR('bad value of gkk2write')
136  end if
137 
138 !here do the actual calculation of |g_kk|^2
139  MSG_ERROR("MGNOTE: interpolate_gkk is broken")
140  ABI_UNUSED(kptirr_phon(1,1))
141  call interpolate_gkk (crystal,ifc,elph_ds,kpt_phon)
142 
143  !MG: This was the old coding in version 7.6.2:
144 
145 ! call interpolate_gkk (elph_ds,kptirr_phon,kpt_phon,natom,nrpt,phon_ds,rcan,wghatm)
146 !
147 ! and interpolate_gkk had the prototype:
148 !
149 !subroutine interpolate_gkk(elph_ds,kpt_phon,gprim,natom,nrpt,phon_ds,rpt,wghatm)
150 
151 ! hence we were associating kpt_phon to gprim!
152 
153 end subroutine get_all_gkk2