TABLE OF CONTENTS


ABINIT/inpgkk [ Functions ]

[ Top ] [ Functions ]

NAME

 inpgkk

FUNCTION

 read in gkk file and return eigenvalue matrix
 Only works for a single gkk matrix (1 perturbation and qpoint) in the file
 like the files produced by outgkk

COPYRIGHT

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

  filegkk= filename

OUTPUT

  eigen1 = response function 1st order eigenvalue matrix

PARENTS

      read_el_veloc

CHILDREN

      hdr_fort_read,hdr_free,wrtout

SOURCE

 33 #if defined HAVE_CONFIG_H
 34 #include "config.h"
 35 #endif
 36 
 37 #include "abi_common.h"
 38 
 39 
 40 subroutine inpgkk(eigen1,filegkk,hdr1)
 41 
 42  use defs_basis
 43  use defs_abitypes
 44  use m_profiling_abi
 45  use m_errors
 46  use m_hdr
 47 
 48  use m_io_tools,        only : open_file
 49 
 50 !This section has been created automatically by the script Abilint (TD).
 51 !Do not modify the following lines by hand.
 52 #undef ABI_FUNC
 53 #define ABI_FUNC 'inpgkk'
 54  use interfaces_14_hidewrite
 55 !End of the abilint section
 56 
 57  implicit none
 58 
 59 !Arguments ------------------------------------
 60 !scalars
 61  character(len=fnlen),intent(in) :: filegkk
 62  type(hdr_type), intent(out) :: hdr1
 63 !arrays
 64  real(dp),allocatable,intent(out) :: eigen1(:)
 65 
 66 !Local variables-------------------------------
 67 !scalars
 68  integer :: bantot1
 69  integer :: isppol, ikpt, mband, ikb
 70  integer :: unitgkk, fform, ierr, n1wf, i1wf
 71  type(hdr_type) :: hdr0
 72  real(dp), allocatable :: eigen(:)
 73  character(len=500) :: message
 74 
 75 ! *************************************************************************
 76 
 77  if (open_file(filegkk,message,newunit=unitgkk,form='unformatted',status='old') /= 0) then
 78    MSG_ERROR(message)
 79  end if
 80 
 81 
 82 !read in header of GS file and eigenvalues
 83  call hdr_fort_read(hdr0, unitgkk, fform)
 84  ABI_CHECK(fform /= 0, "hdr_fort_read returned fform == 0")
 85  
 86  mband = maxval(hdr0%nband(:))
 87  ABI_ALLOCATE(eigen,(mband))
 88  call wrtout(std_out,'inpgkk : try to reread GS eigenvalues','COLL')
 89 
 90  do isppol=1,hdr0%nsppol
 91    do ikpt=1,hdr0%nkpt
 92      read (unitgkk,IOSTAT=ierr) eigen(1:hdr0%nband(ikpt))
 93      ABI_CHECK(ierr==0,'reading eigen from gkk file')
 94    end do
 95  end do
 96 
 97  read(unitgkk,IOSTAT=ierr) n1wf
 98  ABI_CHECK(ierr==0,"reading n1wf from gkk file")
 99 
100  ABI_DEALLOCATE(eigen)
101  call hdr_free(hdr0)
102 
103  if (n1wf > 1) then
104    write(message,'(3a)')&
105 &   'several 1wf records were found in the file,',ch10, &
106 &   'which is not allowed for reading with this routine'
107    MSG_ERROR(message)
108  end if
109 
110 !read in header of 1WF file
111  call hdr_fort_read(hdr1, unitgkk, fform)
112  if (fform == 0) then
113    write(message,'(a,i0,a)')' 1WF header number ',i1wf,' was mis-read. fform == 0'
114    MSG_ERROR(message)
115  end if
116 
117  bantot1 = 2*hdr1%nsppol*hdr1%nkpt*mband**2
118  ABI_ALLOCATE(eigen1, (bantot1))
119 
120 
121 !retrieve 1WF <psi_k+q | H | psi_k> from gkk file and echo to output
122  ikb = 0
123  do isppol=1,hdr1%nsppol
124    do ikpt=1,hdr1%nkpt
125      read (unitgkk,IOSTAT=ierr) eigen1(ikb+1:ikb+2*hdr1%nband(ikpt)**2)
126      ikb = ikb + 2*hdr1%nband(ikpt)**2
127      if (ierr /= 0) then
128        write(message,'(a,2i0)')'reading eigen1 from gkk file, spin, kpt_idx',isppol,ikpt
129        MSG_ERROR(message)
130      end if
131    end do
132  end do
133 
134  close(unitgkk)
135 
136 end subroutine inpgkk