TABLE OF CONTENTS


ABINIT/outgkk [ Functions ]

[ Top ] [ Functions ]

NAME

 outgkk

FUNCTION

 output gkk file for one perturbation (used for elphon calculations in anaddb)

COPYRIGHT

 Copyright (C) 2008-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

  bantot0 = total number of bands for all kpoints
  bantot1 = total number of matrix elements for 1st order eigenvalues
  eigen0 = GS eigenvalues
  eigen1 = response function 1st order eigenvalue matrix
  hdr0 = GS header
  hdr1 = RF header
  mpi_enreg=information about MPI parallelization

PARENTS

      dfpt_looppert

CHILDREN

      hdr_fort_write,wrtout

SOURCE

 33 #if defined HAVE_CONFIG_H
 34 #include "config.h"
 35 #endif
 36 
 37 #include "abi_common.h"
 38 
 39 subroutine outgkk(bantot0,bantot1,outfile,eigen0,eigen1,hdr0,hdr1,mpi_enreg,phasecg)
 40 
 41  use defs_basis
 42  use defs_abitypes
 43  use m_profiling_abi
 44  use m_errors
 45  use m_hdr
 46 
 47  use m_io_tools,    only : open_file
 48 
 49 !This section has been created automatically by the script Abilint (TD).
 50 !Do not modify the following lines by hand.
 51 #undef ABI_FUNC
 52 #define ABI_FUNC 'outgkk'
 53  use interfaces_14_hidewrite
 54 !End of the abilint section
 55 
 56  implicit none
 57 
 58 !Arguments ------------------------------------
 59 !scalars
 60  integer,intent(in) :: bantot0,bantot1
 61  character(len=fnlen),intent(in) :: outfile
 62  type(MPI_type),intent(in) :: mpi_enreg
 63  type(hdr_type),intent(inout) :: hdr0,hdr1
 64 !arrays
 65  real(dp),intent(in) :: eigen0(bantot0),eigen1(2*bantot1)
 66  real(dp),intent(in) :: phasecg(2,bantot1)
 67 
 68 !Local variables-------------------------------
 69 !scalars
 70  integer :: fform,iband,ikpt,isppol,me,ntot,unitout
 71  integer :: iband_off, mband, ierr
 72  character(len=500) :: msg
 73  real(dp), allocatable :: tmpeig(:)
 74 
 75 ! *************************************************************************
 76 
 77 !only master should be writing to disk
 78 !Init me
 79  me=mpi_enreg%me_kpt
 80  if (me /= 0) return
 81 
 82  call wrtout(std_out,' writing gkk file: '//outfile,"COLL")
 83 
 84 !initializations
 85  fform = 42
 86  ntot = 1
 87 
 88 !open gkk file
 89  if (open_file(outfile, msg, newunit=unitout, form='unformatted', status='unknown', action="write") /= 0) then
 90    MSG_ERROR(msg)
 91  end if
 92 
 93 !output GS header
 94  call hdr_fort_write(hdr0, unitout, fform, ierr)
 95  ABI_CHECK(ierr == 0 , "hdr_fort_write returned ierr != 0")
 96 
 97 !output GS eigenvalues
 98  iband=0
 99  do isppol=1,hdr0%nsppol
100    do ikpt=1,hdr0%nkpt
101      write (unitout) eigen0(iband+1:iband+hdr0%nband(ikpt))
102      iband=iband+hdr0%nband(ikpt)
103    end do
104  end do
105 
106 !output number of gkk in this file (1)
107  write (unitout) ntot
108 
109 !output RF header
110  call hdr_fort_write(hdr1, unitout, fform, ierr)
111  ABI_CHECK(ierr == 0 , "hdr_fort_write returned ierr != 0")
112 
113 !output RF eigenvalues
114  mband = maxval(hdr1%nband(:))
115  ABI_ALLOCATE(tmpeig,(2*mband**2))
116  iband_off = 0
117  tmpeig(1) = phasecg(1, 1)
118  do isppol = 1, hdr1%nsppol
119    do ikpt = 1, hdr1%nkpt
120      tmpeig = zero
121      do iband = 1, hdr1%nband(ikpt)**2
122        tmpeig (2*(iband-1)+1) = eigen1(2*(iband_off+iband-1)+1)
123        tmpeig (2*(iband-1)+2) = eigen1(2*(iband_off+iband-1)+2)
124      end do
125      write (unitout) tmpeig(1:2*hdr1%nband(ikpt)**2)
126      iband_off = iband_off + hdr1%nband(ikpt)**2
127    end do
128  end do
129  ABI_DEALLOCATE(tmpeig)
130 
131 !close gkk file
132  close (unitout)
133 
134 end subroutine outgkk