TABLE OF CONTENTS


ABINIT/outbsd [ Functions ]

[ Top ] [ Functions ]

NAME

 outbsd

FUNCTION

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

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (XG, DRH, MB)
 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.

INPUTS

  bdeigrf=number of bands for which the derivatives of the eigenvalues have been computed
  dtset = dataset variable for run flags
  eig2nkq= second ordre eigenvalue (or electron lifetime) that must be printed out 
  mpert= maximum number of perturbations
  nkpt_rbz= number of k-points for perturbation
  unitout= writting unit of file

 OUTPUTS
  to file

PARENTS

      dfpt_looppert,eig2tot

CHILDREN

SOURCE

 34 #if defined HAVE_CONFIG_H
 35 #include "config.h"
 36 #endif
 37 
 38 #include "abi_common.h"
 39 
 40 
 41 subroutine outbsd(bdeigrf,dtset,eig2nkq,mpert,nkpt_rbz,unitout)
 42 
 43  use defs_basis
 44  use defs_abitypes
 45  use m_profiling_abi
 46 
 47 !This section has been created automatically by the script Abilint (TD).
 48 !Do not modify the following lines by hand.
 49 #undef ABI_FUNC
 50 #define ABI_FUNC 'outbsd'
 51 !End of the abilint section
 52 
 53  implicit none
 54 
 55 !Arguments ------------------------------------
 56 !scalars
 57  integer,intent(in) :: bdeigrf,mpert,nkpt_rbz,unitout
 58  type(dataset_type),intent(in) :: dtset
 59 !arrays
 60  real(dp),intent(in) :: eig2nkq(2,dtset%mband*dtset%nsppol,nkpt_rbz,3,mpert,3,mpert)
 61 
 62 !Local variables -------------------------
 63 !scalars
 64  integer :: bandtot_index,iband,idir1,idir2,ikpt,ipert1,ipert2,isppol
 65 
 66 ! *********************************************************************
 67 
 68 !DEBUG
 69 !write(std_out,*)' outbsd : enter'
 70 !write(std_out,*)' eig2nkq(1,1,1,1,1,1)=',eig2nkq(1,1,1,1,1,1,1)
 71 !ENDDEBUG
 72 
 73 !output information in this file
 74  write(unitout,*)
 75  write(unitout,'(a,i8)') ' 2nd eigenvalue derivatives   - # elements :', 9*dtset%natom**2
 76  write(unitout,'(a,3es16.8,a)') ' qpt', dtset%qptn(:), ' 1.0'
 77 
 78 !output RF eigenvalues
 79 
 80  do ikpt=1,nkpt_rbz
 81 !  bandtot_index differs from zero only in the spin-polarized case
 82    bandtot_index=0
 83    write (unitout,'(a,3es16.8)') ' K-point:', dtset%kptns(:,ikpt)
 84    do isppol=1,dtset%nsppol
 85      do iband=1,bdeigrf
 86        write (unitout,'(a,i5)') ' Band:', iband+bandtot_index
 87 !      write (unitout,*) 'ipert1     ','idir1     ','ipert2     ','idir2    ','Real    ','Im    '
 88        do ipert2=1,mpert
 89          do idir2=1,3
 90            do ipert1=1,mpert
 91              do idir1=1,3
 92                write (unitout,'(4i4,2d22.14)') idir1,ipert1,idir2,ipert2,&
 93 &               eig2nkq(1,iband+bandtot_index,ikpt,idir1,ipert1,idir2,ipert2),&
 94 &               eig2nkq(2,iband+bandtot_index,ikpt,idir1,ipert1,idir2,ipert2)
 95              end do !idir2
 96            end do !ipert2
 97          end do !idir1
 98        end do !ipert1
 99      end do !iband
100      bandtot_index = bandtot_index + dtset%mband
101    end do !isppol
102  end do !ikpt
103 
104 !close bsd file
105  close (unitout)
106 
107 end subroutine outbsd