TABLE OF CONTENTS
ABINIT/write_eig [ Functions ]
NAME
write_eig
FUNCTION
Write the eigenvalues band by band and k point by k point in a NetCDF file format
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR) 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
filname = Filename of the file where the history will be stored
OUTPUT
(only writing)
PARENTS
clnup1
CHILDREN
ab_define_var
SOURCE
33 #if defined HAVE_CONFIG_H 34 #include "config.h" 35 #endif 36 37 #include "abi_common.h" 38 39 40 subroutine write_eig(eigen,filename,kptns,mband,nband,nkpt,nsppol) 41 42 use defs_basis 43 use m_profiling_abi 44 use m_errors 45 #ifdef HAVE_NETCDF 46 use netcdf 47 use m_nctk, only : ab_define_var 48 #endif 49 50 !This section has been created automatically by the script Abilint (TD). 51 !Do not modify the following lines by hand. 52 #undef ABI_FUNC 53 #define ABI_FUNC 'write_eig' 54 !End of the abilint section 55 56 implicit none 57 58 !Arguments ------------------------------------ 59 !scalars 60 character(len=fnlen),intent(in) :: filename 61 integer,intent(in) :: nkpt,nsppol,mband 62 !arrays 63 integer,intent(in) :: nband(nkpt*nsppol) 64 real(dp),intent(in) :: eigen(mband*nkpt*nsppol) 65 real(dp),intent(in) :: kptns(3,nkpt) 66 67 !Local variables------------------------------- 68 !scalars 69 integer :: ncerr,ncid,ii 70 integer :: xyz_id,nkpt_id,mband_id,nsppol_id 71 integer :: eig_id,kpt_id,nbk_id,nbk 72 integer :: ikpt,isppol,nband_k,band_index 73 real(dp):: convrt 74 !arrays 75 integer :: dimEIG(3),dimKPT(2),dimNBK(2) 76 integer :: count2(2),start2(2) 77 integer :: count3(3),start3(3) 78 real(dp):: band(mband) 79 80 ! ********************************************************************* 81 82 #if defined HAVE_NETCDF 83 84 convrt=1.0_dp 85 86 !1. Create netCDF file 87 ncerr = nf90_create(path=trim(filename),cmode=NF90_CLOBBER, ncid=ncid) 88 NCF_CHECK_MSG(ncerr," create netcdf EIG file") 89 90 !2. Define dimensions 91 ncerr = nf90_def_dim(ncid,"xyz",3,xyz_id) 92 NCF_CHECK_MSG(ncerr," define dimension xyz") 93 94 ncerr = nf90_def_dim(ncid,"mband",mband,mband_id) 95 NCF_CHECK_MSG(ncerr," define dimension mband") 96 97 ncerr = nf90_def_dim(ncid,"nkpt",nkpt,nkpt_id) 98 NCF_CHECK_MSG(ncerr," define dimension nkpt") 99 100 ncerr = nf90_def_dim(ncid,"nsppol",nsppol,nsppol_id) 101 NCF_CHECK_MSG(ncerr," define dimension nsppol") 102 103 !Dimensions for EIGENVALUES 104 dimEIG = (/ mband_id, nkpt_id, nsppol_id /) 105 !Dimensions for kpoint positions 106 dimKPT = (/ xyz_id, nkpt_id /) 107 !Dimensions for number kpoints per band and spin 108 !dimNBK = (/ nkpt_id, nsppol_id /) 109 dimNBK = (/ nkpt_id, nsppol_id /) 110 111 !3. Define variables 112 113 call ab_define_var(ncid, dimEIG, eig_id, NF90_DOUBLE,& 114 & "Eigenvalues",& 115 & "Values of eigenvalues",& 116 & "Hartree") 117 call ab_define_var(ncid, dimKPT, kpt_id, NF90_DOUBLE,"Kptns",& 118 & "Positions of K-points in reciprocal space",& 119 & "Dimensionless") 120 call ab_define_var(ncid, dimNBK, nbk_id, NF90_INT,"NBandK",& 121 & "Number of bands per kpoint and Spin",& 122 & "Dimensionless") 123 124 !4. End define mode 125 ncerr = nf90_enddef(ncid) 126 NCF_CHECK_MSG(ncerr," end define mode") 127 128 !5 Write kpoint positions 129 do ikpt=1,nkpt 130 start2 = (/ 1, ikpt /) 131 count2 = (/ 3, 1 /) 132 ncerr = nf90_put_var(ncid, kpt_id,& 133 & kptns(1:3,ikpt),& 134 & start = start2,& 135 & count = count2) 136 NCF_CHECK_MSG(ncerr," write variable kptns") 137 end do 138 139 140 !6 Write eigenvalues 141 band_index=0 142 do isppol=1,nsppol 143 do ikpt=1,nkpt 144 nband_k=nband(ikpt+(isppol-1)*nkpt) 145 start3 = (/ 1, ikpt, isppol /) 146 count3 = (/ mband, 1, 1 /) 147 band(:)=zero 148 do ii=1,nband_k 149 band(ii)=eigen(band_index+ii) 150 end do 151 ncerr = nf90_put_var(ncid, eig_id,& 152 & band,& 153 & start = start3,& 154 & count = count3) 155 NCF_CHECK_MSG(ncerr," write variable band") 156 157 band_index=band_index+nband_k 158 end do 159 end do 160 161 !6 Write Number of bands per kpoint and Spin 162 163 do isppol=1,nsppol 164 do ikpt=1,nkpt 165 start2 = (/ ikpt, 1 /) 166 count2 = (/ 1, 1 /) 167 nbk=nband(ikpt+(isppol-1)*nkpt) 168 ncerr = nf90_put_var(ncid, nbk_id,& 169 & nbk,& 170 & start = start2) 171 NCF_CHECK_MSG(ncerr," write variable nband") 172 end do 173 end do 174 175 !7 Close file 176 177 ncerr = nf90_close(ncid) 178 NCF_CHECK_MSG(ncerr," close netcdf EIG file") 179 #endif 180 181 end subroutine write_eig