TABLE OF CONTENTS
- ABINIT/m_spectra
- m_spectra/dump_Qlist
- m_spectra/spectra_free
- m_spectra/spectra_init
- m_spectra/spectra_repr
- m_spectra/spectra_t
- m_spectra/spectra_write
ABINIT/m_spectra [ Modules ]
NAME
m_spectra
FUNCTION
This module contains the definition of the specta_type data type used to store results related to optical spectra with or without nonlocal field effects as well as the electron energy loss function for a given q. These quantities are obtained from the dielectric matrix as calculated in the GW part of ABINIT (screening.F90)
COPYRIGHT
Copyright (C) 2008-2024 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 .
SOURCE
20 #if defined HAVE_CONFIG_H 21 #include "config.h" 22 #endif 23 24 #include "abi_common.h" 25 26 MODULE m_spectra 27 28 use defs_basis 29 use m_errors 30 use m_abicore 31 32 use m_io_tools, only : open_file 33 use m_fstrings, only : strcat 34 35 implicit none 36 37 private
m_spectra/dump_Qlist [ Functions ]
[ Top ] [ m_spectra ] [ Functions ]
NAME
dump_Qlist
FUNCTION
Helper function used to write the list of q-points used for the long-wavelength limit.
INPUTS
Spectra=The Object containing the spectra
OUTPUT
Only writing.
SOURCE
272 subroutine dump_Qlist() 273 274 integer :: iqpt 275 write(unt,'(a,i3)')'# Q-point list, No. ',Spectra%nqpts 276 do iqpt=1,Spectra%nqpts 277 write(unt,'(a,i3,a,3f9.6,a)')'# ',iqpt,') [',Spectra%qpts(:,iqpt),'] r.l.u. ' 278 end do 279 end subroutine dump_Qlist 280 281 end subroutine spectra_write
m_spectra/spectra_free [ Functions ]
[ Top ] [ m_spectra ] [ Functions ]
NAME
spectra_free
FUNCTION
Deallocate all associated pointers defined in the structure.
INPUTS
OUTPUT
SOURCE
157 subroutine spectra_free(Spectra) 158 159 !Arguments ------------------------------------ 160 class(spectra_t),intent(inout) :: spectra 161 162 ! ********************************************************************* 163 164 ABI_SFREE(Spectra%omega) 165 ABI_SFREE(Spectra%qpts) 166 ABI_SFREE(Spectra%emacro_lf) 167 ABI_SFREE(Spectra%emacro_nlf) 168 ABI_SFREE(Spectra%eelf) 169 170 end subroutine spectra_free
m_spectra/spectra_init [ Functions ]
[ Top ] [ m_spectra ] [ Functions ]
NAME
spectra_init
FUNCTION
Initialize the object.
INPUTS
OUTPUT
SOURCE
115 subroutine spectra_init(Spectra,nomega,omega,nqpts,qpts) 116 117 !Arguments ------------------------------------ 118 !scalars 119 integer,intent(in) :: nomega,nqpts 120 !arrays 121 real(dp),intent(in) :: omega(nomega),qpts(3,nqpts) 122 type(spectra_t),intent(out) :: Spectra 123 124 ! ********************************************************************* 125 126 Spectra%nomega = nomega 127 Spectra%nqpts = nqpts 128 129 ABI_MALLOC(Spectra%qpts,(3,nqpts)) 130 Spectra%qpts = qpts 131 132 ABI_MALLOC(Spectra%omega,(nomega)) 133 Spectra%omega = omega 134 135 ABI_CALLOC(Spectra%emacro_lf,(nomega,nqpts)) 136 ABI_CALLOC(Spectra%emacro_nlf,(nomega,nqpts)) 137 ABI_CALLOC(Spectra%eelf,(nomega,nqpts)) 138 139 end subroutine spectra_init
m_spectra/spectra_repr [ Functions ]
[ Top ] [ m_spectra ] [ Functions ]
NAME
spectra_repr
FUNCTION
Returns a string reporting info on the calculated dielectric constant.
INPUTS
Spectra=The Object containing the spectra
OUTPUT
SOURCE
300 subroutine spectra_repr(Spectra,str) 301 302 !Arguments ------------------------------------ 303 !scalars 304 class(spectra_t),intent(in) :: Spectra 305 character(len=*),intent(out) :: str 306 307 !Local variables------------------------------- 308 !scalars 309 integer :: iqpt 310 real(dp) :: epsilon0,epsilon0_nlf 311 character(len=15), parameter :: format_f84 = '(1x,a, f8.4,a)' 312 character(len=15), parameter :: format_es135 = '(1x,a,es13.5,a)' 313 character(len=15) :: format_diel 314 character(len=500) :: msg 315 316 ! ********************************************************************* 317 318 !istatic = -1 319 !do io = Spectra%nomega 320 ! if (ABS(REAL(Spectra%omega(io)))<1.e-3.and.ABS(AIMAG(Spectra%omega(io)))<1.e-3) then 321 ! istatic = io 322 ! EXIT 323 ! end if 324 !end do 325 326 str = "" 327 do iqpt=1,Spectra%nqpts 328 epsilon0 = REAL(Spectra%emacro_lf (1,iqpt)) 329 epsilon0_nlf= REAL(Spectra%emacro_nlf(1,iqpt)) 330 write(msg,'(a,3f9.6,a)')' For q-point: ',Spectra%qpts(:,iqpt),ch10 331 format_diel=format_f84 332 if(abs(epsilon0)>1000.0d0 .or. abs(epsilon0_nlf)>1000.0d0) format_diel=format_es135 333 str = strcat(str,msg) 334 write(msg,format_diel)' dielectric constant = ',epsilon0,ch10 335 str = strcat(str,msg) 336 write(msg,format_diel)' dielectric constant without local fields = ',epsilon0_nlf,ch10 337 str = strcat(str,msg) 338 end do 339 340 end subroutine spectra_repr 341 342 !---------------------------------------------------------------------- 343 344 END MODULE m_spectra
m_spectra/spectra_t [ Types ]
[ Top ] [ m_spectra ] [ Types ]
NAME
spectra_t
FUNCTION
Object used to store optical spectra with or without non-local field effects.
SOURCE
49 type,public :: spectra_t 50 51 !scalars 52 integer :: nomega 53 ! number of frequencies 54 55 integer :: nqpts 56 ! number of q-points 57 58 !arrays 59 real(dp),allocatable :: omega(:) 60 ! omega(nomega) 61 ! Real frequency mesh for optical spectra. 62 63 real(dp),allocatable :: qpts(:,:) 64 ! qpts(3,nqpoints) 65 ! The list of q-points used for the spectra 66 67 real(dp),allocatable :: eelf(:,:) 68 ! eelf(nomega,nqpoints) 69 ! contains the Electron Energy Loss Function i.e. -\Im{ e^{-1}_{G1=0,G2=0}(q-->0,nomega)} 70 71 complex(dpc),allocatable :: emacro_lf(:,:) 72 ! emacro_lf(nomega,nqpoints) 73 ! contains 1/e^{-1}_{G1=0,G2=0}(q-->0,nomega) (with Local field effects) 74 75 complex(dpc),allocatable :: emacro_nlf(:,:) 76 ! emacro_nlf(nomega,nqpoints) 77 ! contains e_{G1=0,G2=0}(q-->0,nomega) (without Local field effects) 78 79 contains 80 81 procedure :: free => spectra_free 82 ! Free memory. 83 84 procedure :: write => spectra_write 85 ! Write results on file. 86 87 procedure :: repr => spectra_repr 88 ! Return info on Macroscopic diel. constant in form of a string. 89 90 end type spectra_t
m_spectra/spectra_write [ Functions ]
[ Top ] [ m_spectra ] [ Functions ]
NAME
spectra_write
FUNCTION
Write the optical spectra stored in the object on an external formatted file.
INPUTS
Spectra=The Object containing the spectra write_bits=Positive integer defining the quantities to be written (bit representation is used) fname=Name of the file to be written.
OUTPUT
SOURCE
191 subroutine spectra_write(Spectra,write_bits,fname) 192 193 !Arguments ------------------------------------ 194 !scalars 195 class(spectra_t),intent(in) :: Spectra 196 integer,intent(in) :: write_bits 197 character(len=*),intent(in) :: fname 198 199 !Local variables------------------------------- 200 !scalars 201 integer :: unt,io,iqpt 202 real(dp) :: mino,maxo 203 character(len=100) :: fmt 204 character(len=500) :: msg 205 206 ! ********************************************************************* 207 208 if (write_bits<0) RETURN 209 210 mino = MINVAL(Spectra%omega) 211 maxo = MAXVAL(Spectra%omega) 212 213 if (open_file(fname,msg,newunit=unt) /= 0) then 214 ABI_ERROR(msg) 215 end if 216 217 !write(unt,'(a,i5,2(a,f9.1),a)')'# nomega : ',Spectra%nomega,' from ',mino*Ha_eV,' up to ',maxo*Ha_eV,' [eV] ' 218 219 if ( IAND(write_bits,W_EM_NLF ) == W_EM_NLF ) then 220 write(unt,'(a)')'#' 221 write(unt,'(a)')'# Macroscopic Dielectric Function without local fields' 222 call dump_qlist() 223 write(unt,'(a)')'# Omega [eV] Re epsilon_M IM eps_M ' 224 write(fmt,'(a,i3,a)') '(1x,f7.3,7x,',Spectra%nqpts,'(2(e12.4),2x))' 225 do io=1,Spectra%nomega 226 write(unt,fmt) Spectra%omega(io)*Ha_eV, ( Spectra%emacro_nlf(io,iqpt), iqpt=1,Spectra%nqpts ) 227 end do 228 end if 229 ! 230 if ( IAND(write_bits, W_EM_LF ) == W_EM_LF ) then 231 write(unt,'(a)')'#' 232 write(unt,'(a)')'# Macroscopic Dielectric Function with local fields included' 233 call dump_qlist() 234 write(unt,'(a)')'# Omega [eV] Re epsilon_M Im eps_M ' 235 write(fmt,'(a,i3,a)') '(1x,f7.3,7x,',Spectra%nqpts,'(2(e12.4),2x))' 236 do io=1,Spectra%nomega 237 write(unt,fmt) Spectra%omega(io)*Ha_eV, ( Spectra%emacro_lf(io,iqpt), iqpt=1,Spectra%nqpts ) 238 end do 239 end if 240 ! 241 if ( IAND(write_bits,W_EELF) == W_EELF) then 242 write(unt,'(a)')'#' 243 write(unt,'(a)')'# Electron Energy Loss Function -Im(1/epsilon_M)' 244 call dump_qlist() 245 write(unt,'(a)')'# Omega [eV] -Im(1/epsilon_M)' 246 write(fmt,'(a,i3,a)') '(1x,f7.3,7x,',Spectra%nqpts,'(e12.4,2x))' 247 do io=1,Spectra%nomega 248 write(unt,fmt) Spectra%omega(io)*Ha_eV, ( Spectra%eelf(io,iqpt), iqpt=1,Spectra%nqpts ) ! -AIMAG(chi0(1,1,io)) 249 end do 250 end if 251 252 close(unt) 253 254 CONTAINS