~ubuntu-branches/debian/sid/abinit/sid

1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
!{\src2tex{textfont=tt}}
!!****f* ABINIT/prtxcfermsurf
!! NAME
!! prtxcfermsurf
!!
!! FUNCTION
!!  Print the Fermi surface in format XCrysDen
!!
!! COPYRIGHT
!! Copyright (C) 2003-2007 ABINIT group (MVerstraete)
!! This file is distributed under the terms of the
!! GNU General Public License, see ~abinit/COPYING
!! or http://www.gnu.org/copyleft/gpl.txt .
!!
!! INPUTS
!!  eigen(mband,nkpt,nsppol)=eigenvalues in hartree
!!  fermie=Fermi energy (Hartree)
!!  gprimd(3,3)=dimensional primitive translations for reciprocal space(bohr^-1)
!!  klatt(3,3)=reciprocal of lattice vectors for full kpoint grid
!!  indkpt(nkpt_fullbz)=indexes of irred kpoints equivalent to kpt_fullbz
!!  kpt_fullbz(3,nkpt_fullbz)=kpoints in full brillouin zone
!!  mband=maximum number of bands
!!  nkpt=number of k points.
!!  nkpt_fullbz=number of kpoints in full brillouin zone
!!  nsppol=1 for unpolarized, 2 for spin-polarized
!!  shiftk(3)=shift vector for k point grid
!!
!! OUTPUT
!!  (only output)
!!
!! PARENTS
!!
!! CHILDREN
!!
!! SOURCE

#if defined HAVE_CONFIG_H
#include "config.h"
#endif

subroutine prtxcfermsurf(eigen,fermie,gprimd,klatt,indkpt,kpt_fullbz,&
&        mband,nkpt,nkpt_fullbz,nsppol,shiftk)

 use defs_basis
  use defs_datatypes

 implicit none

!Arguments ------------------------------------
!scalars
 integer,intent(in) :: mband,nkpt,nkpt_fullbz,nsppol
 real(dp),intent(in) :: fermie
!arrays
 integer,intent(in) :: indkpt(nkpt_fullbz)
 real(dp),intent(in) :: eigen(mband,nkpt,nsppol),gprimd(3,3),klatt(3,3)
 real(dp),intent(in) :: kpt_fullbz(3,nkpt_fullbz),shiftk(3)

!Local variables-------------------------------
!scalars
 integer :: formeig,iband,ikpt,isppol,maxband,minband,nk1,nk2,nk3,option
 real(dp) :: mkval
 character(len=fnlen) :: filename

! *************************************************************************

! Error if klatt is no simple orthogonal lattice (in red space)
!  for generalization to MP grids, need new version of XCrysDen
  if (abs(klatt(1,2))+abs(klatt(1,3))+abs(klatt(2,3)) > tol6) then
    write(*,*) 'prtxcfermsurf Warning! : klatt is not orthogonal'
    write(*,*) ' klatt is not orthogonal'
    write(*,*) ' action: correct kptrlatt or ngpkpt...'
  end if

!GS form of eigenvalues only
  formeig = 0
  option=1
  minband = mband
  maxband = 0
!for all bands and 2 sppols (set as different bands)
  do isppol=1,nsppol
! if Fermi surface crosses band, include in plot of FSurf
   do iband=1,mband
    if(minval(eigen(iband,:,isppol))-fermie < -tol6) then
     minband = iband
    end if
   end do
   do iband=mband,1,-1
    if(maxval(eigen(iband,:,isppol))-fermie > tol6) then
     maxband = iband
    end if
   end do
  end do
! end sppol

  filename = 'toto.bxsf'

  open(unit=tmp_unit,file=filename,status='replace',form='formatted')
  rewind(tmp_unit)
!print out header

  mkval = maxval(abs(klatt(:,1)))
  nk1=int(1.0/mkval)
  mkval = maxval(abs(klatt(:,2)))
  nk2=int(1.0/mkval)
  mkval = maxval(abs(klatt(:,3)))
  nk3=int(1.0/mkval)

  write(tmp_unit,*) '    BEGIN_INFO'
  write(tmp_unit,*) '      #'
  write(tmp_unit,*) '      # this is a Band-XCRYSDEN-Structure-File for'
  write(tmp_unit,*) '      # Visualization of Fermi Surface'
  write(tmp_unit,*) '      # generated by cut3d in ABINIT package'
  write(tmp_unit,*) '      #'
  write(tmp_unit,*) '      # Launch as: xcrysden --bxsf toto.bxsf.gz'
  write(tmp_unit,*) '      #'
  write(tmp_unit,'(a,F12.4)') '      Fermi Energy: ', fermie
  write(tmp_unit,*) '    END_INFO'
  write(tmp_unit,*)
  write(tmp_unit,*) '    BEGIN_BLOCK_BANDGRID_3D'
  write(tmp_unit,*) '    band_energies'
  write(tmp_unit,*) '    BANDGRID_3D_BANDS'
  write(tmp_unit,*) '    ', maxband-minband+1
  write(tmp_unit,*) '    ', nk1,nk2,nk3
  write(tmp_unit,*) '    ', shiftk(:)
  write(tmp_unit,*) '    ', gprimd(:,1)
  write(tmp_unit,*) '    ', gprimd(:,2)
  write(tmp_unit,*) '    ', gprimd(:,3)


!print out data for all relevant bands and full kpt grid (redundant, yes)
!  for each kpt in full zone, find equivalent irred kpt and print
!  eigenval
  do isppol=1,nsppol
  do iband=minband,maxband
  write(tmp_unit,*) '    BAND: ', iband
  write(tmp_unit,'(6(F12.6,1x))') (eigen(iband,indkpt(ikpt),isppol), &
  &                            ikpt=1,nkpt_fullbz)
  end do
  end do
  write(tmp_unit,*) '    END_BANDGRID_3D'
  write(tmp_unit,*) '    END_BLOCK_BANDGRID_3D'

end subroutine prtxcfermsurf
!!***