1
!{\src2tex{textfont=tt}}
2
!!****f* ABINIT/prtxcfermsurf
7
!! Print the Fermi surface in format XCrysDen
10
!! Copyright (C) 2003-2007 ABINIT group (MVerstraete)
11
!! This file is distributed under the terms of the
12
!! GNU General Public License, see ~abinit/COPYING
13
!! or http://www.gnu.org/copyleft/gpl.txt .
16
!! eigen(mband,nkpt,nsppol)=eigenvalues in hartree
17
!! fermie=Fermi energy (Hartree)
18
!! gprimd(3,3)=dimensional primitive translations for reciprocal space(bohr^-1)
19
!! klatt(3,3)=reciprocal of lattice vectors for full kpoint grid
20
!! indkpt(nkpt_fullbz)=indexes of irred kpoints equivalent to kpt_fullbz
21
!! kpt_fullbz(3,nkpt_fullbz)=kpoints in full brillouin zone
22
!! mband=maximum number of bands
23
!! nkpt=number of k points.
24
!! nkpt_fullbz=number of kpoints in full brillouin zone
25
!! nsppol=1 for unpolarized, 2 for spin-polarized
26
!! shiftk(3)=shift vector for k point grid
37
#if defined HAVE_CONFIG_H
41
subroutine prtxcfermsurf(eigen,fermie,gprimd,klatt,indkpt,kpt_fullbz,&
42
& mband,nkpt,nkpt_fullbz,nsppol,shiftk)
49
!Arguments ------------------------------------
51
integer,intent(in) :: mband,nkpt,nkpt_fullbz,nsppol
52
real(dp),intent(in) :: fermie
54
integer,intent(in) :: indkpt(nkpt_fullbz)
55
real(dp),intent(in) :: eigen(mband,nkpt,nsppol),gprimd(3,3),klatt(3,3)
56
real(dp),intent(in) :: kpt_fullbz(3,nkpt_fullbz),shiftk(3)
58
!Local variables-------------------------------
60
integer :: formeig,iband,ikpt,isppol,maxband,minband,nk1,nk2,nk3,option
62
character(len=fnlen) :: filename
64
! *************************************************************************
66
! Error if klatt is no simple orthogonal lattice (in red space)
67
! for generalization to MP grids, need new version of XCrysDen
68
if (abs(klatt(1,2))+abs(klatt(1,3))+abs(klatt(2,3)) > tol6) then
69
write(*,*) 'prtxcfermsurf Warning! : klatt is not orthogonal'
70
write(*,*) ' klatt is not orthogonal'
71
write(*,*) ' action: correct kptrlatt or ngpkpt...'
74
!GS form of eigenvalues only
79
!for all bands and 2 sppols (set as different bands)
81
! if Fermi surface crosses band, include in plot of FSurf
83
if(minval(eigen(iband,:,isppol))-fermie < -tol6) then
88
if(maxval(eigen(iband,:,isppol))-fermie > tol6) then
95
filename = 'toto.bxsf'
97
open(unit=tmp_unit,file=filename,status='replace',form='formatted')
101
mkval = maxval(abs(klatt(:,1)))
103
mkval = maxval(abs(klatt(:,2)))
105
mkval = maxval(abs(klatt(:,3)))
108
write(tmp_unit,*) ' BEGIN_INFO'
109
write(tmp_unit,*) ' #'
110
write(tmp_unit,*) ' # this is a Band-XCRYSDEN-Structure-File for'
111
write(tmp_unit,*) ' # Visualization of Fermi Surface'
112
write(tmp_unit,*) ' # generated by cut3d in ABINIT package'
113
write(tmp_unit,*) ' #'
114
write(tmp_unit,*) ' # Launch as: xcrysden --bxsf toto.bxsf.gz'
115
write(tmp_unit,*) ' #'
116
write(tmp_unit,'(a,F12.4)') ' Fermi Energy: ', fermie
117
write(tmp_unit,*) ' END_INFO'
119
write(tmp_unit,*) ' BEGIN_BLOCK_BANDGRID_3D'
120
write(tmp_unit,*) ' band_energies'
121
write(tmp_unit,*) ' BANDGRID_3D_BANDS'
122
write(tmp_unit,*) ' ', maxband-minband+1
123
write(tmp_unit,*) ' ', nk1,nk2,nk3
124
write(tmp_unit,*) ' ', shiftk(:)
125
write(tmp_unit,*) ' ', gprimd(:,1)
126
write(tmp_unit,*) ' ', gprimd(:,2)
127
write(tmp_unit,*) ' ', gprimd(:,3)
130
!print out data for all relevant bands and full kpt grid (redundant, yes)
131
! for each kpt in full zone, find equivalent irred kpt and print
134
do iband=minband,maxband
135
write(tmp_unit,*) ' BAND: ', iband
136
write(tmp_unit,'(6(F12.6,1x))') (eigen(iband,indkpt(ikpt),isppol), &
137
& ikpt=1,nkpt_fullbz)
140
write(tmp_unit,*) ' END_BANDGRID_3D'
141
write(tmp_unit,*) ' END_BLOCK_BANDGRID_3D'
143
end subroutine prtxcfermsurf