TABLE OF CONTENTS
ABINIT/mknormpath [ Functions ]
NAME
mknormpath
FUNCTION
Please do not use this routine, use make_normpath instead. mknormpath should be removed This simple routine generates a normalized path that can be used to plot a band structures in an easy way. For normalized path we mean a path where the number of division on each segment is proportional to the length of the segment itself. To generate the above mentioned path, the subroutine must be called twice. The first call reports the total number of divisions in the normalized path, dimension that is required to correctly allocate the array. The second call calculates the reduced coordinates of the circuit.
COPYRIGHT
Copyright (C) 2007-2018 ABINIT group (MG) 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
nbounds=number of points defining the path ndiv_small=number of points to be used to sample the smallest segment defined by bounds(:,1:nbounds) bounds(3,nbounds)=points defining the path gmet(3,3)=metric
OUTPUT
ndiv(nbounds-1)= number of divisions for each segment npt_tot=total number of points sampled along the circuit path(3,npt_tot)= normalized path in reciprocal space
TODO
Do not use this routine, it is obsolete and should be replaced by make_path in m_bz_mesh.
PARENTS
inkpts
CHILDREN
wrtout
SOURCE
47 #if defined HAVE_CONFIG_H 48 #include "config.h" 49 #endif 50 51 #include "abi_common.h" 52 53 54 subroutine mknormpath(nbounds,bounds,gmet,ndiv_small,ndiv,npt_tot,path) 55 56 use defs_basis 57 use m_profiling_abi 58 use m_errors 59 60 !This section has been created automatically by the script Abilint (TD). 61 !Do not modify the following lines by hand. 62 #undef ABI_FUNC 63 #define ABI_FUNC 'mknormpath' 64 use interfaces_14_hidewrite 65 !End of the abilint section 66 67 implicit none 68 69 !Arguments ------------------------------------ 70 !F95 construct, interface required but we can call mknormpath once 71 !real(dp),pointer :: path(:,:) 72 !scalars 73 integer,intent(in) :: nbounds,ndiv_small 74 integer,intent(inout) :: npt_tot 75 !arrays 76 integer,intent(inout) :: ndiv(nbounds-1) 77 real(dp),intent(in) :: bounds(3,nbounds),gmet(3,3) 78 real(dp),intent(out),optional :: path(3,npt_tot) 79 80 !Local variables------------------------------- 81 !scalars 82 integer :: idx,ii,jp 83 real(dp) :: fct 84 character(len=500) :: message 85 !arrays 86 real(dp) :: dd(3),lng(nbounds-1) 87 88 ! ************************************************************************* 89 90 if (ndiv_small<=0) then 91 write(message,'(3a,i0)')& 92 & 'The argument ndiv_small should be a positive number,',ch10,& 93 & 'however, ndiv_small=',ndiv_small 94 MSG_ERROR(message) 95 end if 96 97 do ii=1,nbounds-1 98 dd(:)=bounds(:,ii+1)-bounds(:,ii) 99 lng(ii)= sqrt( dd(1)*gmet(1,1)*dd(1)+ & 100 & dd(2)*gmet(2,2)*dd(2)+ & 101 & dd(3)*gmet(3,3)*dd(3)+ & 102 & 2.0d0*(dd(1)*gmet(1,2)*dd(2)+ & 103 & dd(1)*gmet(1,3)*dd(3)+ & 104 & dd(2)*gmet(2,3)*dd(3)) & 105 & ) 106 end do 107 write(std_out,*)lng 108 fct=minval(lng) 109 110 !Avoid division by zero if k(:,i+1)=k(:,i) 111 if (abs(fct)<tol6) then 112 write(message,'(3a)')& 113 & 'found two consecutive points in the path which are equal',ch10,& 114 & 'This is not allowed, please modify the path in your input file' 115 MSG_ERROR(message) 116 end if 117 118 fct=fct/ndiv_small 119 ndiv(:)=nint(lng(:)/fct) 120 !The 1 stand for the first point 121 npt_tot=sum(ndiv)+1 122 123 !allocate(path(3,npt_tot) 124 if (.not.present(path)) then 125 write(message,'(2a,i8)')ch10,& 126 & ' mknormpath : total number of points on the path : ',npt_tot 127 call wrtout(std_out,message,'COLL') 128 write(message,'(2a)')ch10,' Number of divisions for each segment of the normalized path : ' 129 call wrtout(std_out,message,'COLL') 130 do ii=1,nbounds-1 131 write(message,'(2(3f8.5,a),i5,a)')& 132 bounds(:,ii),' ==> ',bounds(:,ii+1),' ( ndiv : ',ndiv(ii),' )' 133 call wrtout(std_out,message,'COLL') 134 end do 135 write(message,'(a)')ch10 136 call wrtout(std_out,message,'COLL') 137 else 138 write(message,'(2a)')ch10,' Normalized Path : ' 139 call wrtout(std_out,message,'COLL') 140 idx=1 141 do ii=1,nbounds-1 142 do jp=1,ndiv(ii) 143 path(:,idx)=bounds(:,ii)+(jp-1)*(path(:,ii+1)-path(:,ii))/ndiv(ii) 144 write(message,'(i4,4x,3(f8.5,1x))')idx,path(:,idx) 145 call wrtout(std_out,message,'COLL') 146 idx=idx+1 147 end do 148 end do 149 end if 150 151 end subroutine mknormpath