TABLE OF CONTENTS
- ABINIT/m_bs_defs
- m_bs_defs/bs_parameters_free
- m_bs_defs/bsp_calctype2str
- m_bs_defs/excfiles
- m_bs_defs/excparam
- m_bs_defs/init_transitions
- m_bs_defs/print_bs_files
- m_bs_defs/print_bs_parameters
- m_bs_defs/repr_1trans
- m_bs_defs/repr_2trans
- m_bs_defs/transition
ABINIT/m_bs_defs [ Modules ]
NAME
m_bs_defs
FUNCTION
This module defines basic structures used for Bethe-Salpeter calculations.
COPYRIGHT
Copyright (C) 1992-2024 ABINIT and EXC group (L.Reining, V.Olevano, F.Sottile, S.Albrecht, G.Onida, 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 . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
SOURCE
17 #if defined HAVE_CONFIG_H 18 #include "config.h" 19 #endif 20 21 #include "abi_common.h" 22 23 MODULE m_bs_defs 24 25 use defs_basis 26 use m_abicore 27 use m_errors 28 29 implicit none 30 31 private
m_bs_defs/bs_parameters_free [ Functions ]
[ Top ] [ m_bs_defs ] [ Functions ]
NAME
bs_parameters_free
FUNCTION
Free all memory allocated in a structure of type excparam
SIDE EFFECTS
Bsp<excparam>=All associated pointers are deallocated.
SOURCE
286 subroutine bs_parameters_free(BSp) 287 288 !Arguments ------------------------------------ 289 class(excparam),intent(inout) :: BSp 290 291 !************************************************************************ 292 293 !@excparam 294 ABI_SFREE(BSp%q) 295 ABI_FREE(Bsp%nreh) 296 ABI_SFREE(Bsp%vcks2t) 297 ABI_SFREE(Bsp%omega) 298 ABI_SFREE(Bsp%lomo_spin) 299 ABI_SFREE(Bsp%homo_spin) 300 ABI_SFREE(Bsp%lumo_spin) 301 ABI_SFREE(Bsp%humo_spin) 302 ABI_SFREE(Bsp%nbndv_spin) 303 ABI_SFREE(Bsp%nbndc_spin) 304 ABI_SFREE(Bsp%Trans) 305 ABI_SFREE(Bsp%nreh_interp) 306 ABI_SFREE(Bsp%vcks2t_interp) 307 ABI_SFREE(Bsp%Trans_interp) 308 309 end subroutine bs_parameters_free
m_bs_defs/bsp_calctype2str [ Functions ]
[ Top ] [ m_bs_defs ] [ Functions ]
NAME
bsp_calctype2str
FUNCTION
Returns a string with the calculation type.
SOURCE
494 subroutine bsp_calctype2str(BSp, str) 495 496 !Arguments ------------------------------------ 497 !scalars 498 class(excparam),intent(in) :: BSp 499 character(len=500),intent(out) :: str 500 501 !************************************************************************ 502 503 SELECT CASE (Bsp%calc_type) 504 CASE (BSE_HTYPE_RPA_KS) 505 str = " RPA L0 with KS energies and KS wavefunctions" 506 CASE (BSE_HTYPE_RPA_QPENE) 507 str = " RPA L0 with QP energies and KS wavefunctions" 508 CASE (BSE_HTYPE_RPA_QP) 509 str = " RPA L0 with QP energies and QP wavefunctions" 510 CASE DEFAULT 511 str = " Unknown" 512 END SELECT 513 514 end subroutine bsp_calctype2str
m_bs_defs/excfiles [ Types ]
[ Top ] [ m_bs_defs ] [ Types ]
NAME
excfiles
FUNCTION
The excfiles derived data type contains file names and unit numbers used to store temporary or final results of the Bethe-Salpeter calculation.
SOURCE
237 type,public :: excfiles 238 239 character(len=fnlen) :: in_hreso = BSE_NOFILE 240 ! Name of the input file with the resonant part of the Hamiltonian (Hermitian). 241 242 character(len=fnlen) :: out_hreso = BSE_NOFILE 243 ! Name of the output file with the resonant part of the Hamiltonian (Hermitian). 244 245 character(len=fnlen) :: in_hcoup = BSE_NOFILE 246 ! Name of the input file with the coupling part of the Hamiltonian (Symmetric). 247 248 character(len=fnlen) :: out_hcoup = BSE_NOFILE 249 ! Name of the output file with the coupling part of the Hamiltonian (Symmetric). 250 251 character(len=fnlen) :: in_eig = BSE_NOFILE 252 ! Name of the input file with the eigenvalues and the eigenvectors of the Hamiltonian. 253 254 character(len=fnlen) :: out_eig = BSE_NOFILE 255 ! Name of the output file with the eigenvalues and the eigenvectors of the Hamiltonian. 256 257 character(len=fnlen) :: in_haydock_basename = BSE_NOFILE 258 ! Name of the input file used to restart Haydock algorithm. 259 260 character(len=fnlen) :: out_basename = BSE_NOFILE 261 ! Prefix to be used for other output files. 262 263 end type excfiles 264 265 public :: print_bs_files ! Printout of the excfiles data type.
m_bs_defs/excparam [ Types ]
[ Top ] [ m_bs_defs ] [ Types ]
NAME
excparam
FUNCTION
The excparam derived data type contains the parameters controlling the BS calculation.
SOURCE
106 type,public :: excparam 107 108 !scalars 109 integer :: algorithm ! Algorithm used for computing the BS dielectric function. 110 integer :: calc_type ! Calculation type (see Dtset%bs_calc_type). 111 integer :: hayd_term ! Option for the terminator used in the Haydock solver. 112 integer :: use_coupling ! Include off-diagonal block coupling resonant and anti-resonant transitions. 113 integer :: exchange_term ! Include the exchange term in the BS Hamiltonian. 114 integer :: inclvkb ! Option for the inclusion of the commutator [Vnl, r] for the optical limit. 115 integer :: mdlf_type ! Model dielectric function type. 116 integer :: nline ! Number of line minimization used for CG minimization. 117 integer :: nbdbuf ! Number of states in the buffer that will be excluded from the convergence check (CG only) 118 integer :: nstates ! Number of states that will be considered in the CG minimization. 119 integer :: npweps ! No. of G in the Screening. 120 integer :: npwwfn ! No. of G for wave functions. 121 !$integer :: npwx ! No. of G for the exchange part. 122 integer :: npwvec ! MAX between npwwfn and npweps 123 integer :: nbnds ! Total number of bands considered. 124 125 !integer :: nbndv ! No. of valence states treated (homo-lomo+1) 126 !integer :: nbndc ! No. of conduction states (humo-lumo+1) 127 !integer :: lomo,homo ! Lowest and highest occupied orbital considered. 128 !integer :: lumo,humo ! Lowest and highest unoccupied orbital considered. 129 130 ! new for spin 131 ! DO I need these? 132 integer :: lomo_min,homo_max ! Lowest and highest occupied orbital considered. 133 integer :: lumo_min,humo_max ! Lowest and highest unoccupied orbital considered. 134 integer :: maxnbndv, maxnbndc 135 136 integer,allocatable :: lomo_spin(:) ! Lowest occupied orbital considered for the different spins. 137 integer,allocatable :: homo_spin(:) ! Highest occupied orbital considered for the different spins. 138 integer,allocatable :: lumo_spin(:) ! Lowest unoccupied orbital considered for the different spins 139 integer,allocatable :: humo_spin(:) ! Highest unoccupied orbital considered for the different spins 140 integer,allocatable :: nbndv_spin(:) ! No. of valence states treated (homo-lomo+1) 141 integer,allocatable :: nbndc_spin(:) ! No. of conduction states (humo-lumo+1) 142 ! end new 143 144 integer :: niter ! No. of iterations for (Haydock|CG). 145 integer :: nkibz, nkbz ! No. of k-points in the IBZ and BZ (resp.) 146 integer :: nomega ! No. of frequencies for epsilon. 147 integer :: nq ! Number of "small" q for optical limit. 148 integer :: nsppol ! Number of independent spin polarizations. 149 integer :: wtype ! Option used for dealing with W (see BSE_WTYPE_) flags 150 151 !Interp@BSE 152 integer :: interp_mode ! Mode of interpolation 153 integer :: interp_method ! Method of interpolation 154 integer :: nkibz_interp,nkbz_interp ! Number of points in the interpolated kmesh 155 integer :: nstates_interp ! Number of states of interpolation 156 integer :: rl_nb ! Index of the nb in Rohlfing and Louie technique 157 158 real(dp) :: ecutwfn ! Cutoff energy for wavefunctions. 159 real(dp) :: ecuteps ! Cutoff energy for W. 160 real(dp) :: eps_inf ! Electronic dielectric constant used for the model dielectric function. 161 real(dp) :: mbpt_sciss ! Scissors energy (used if it absolute value is > tol6) 162 real(dp) :: omegai ! First omega for epsilon. 163 real(dp) :: omegae ! Last omega for epsilon (defaults to 10eV) 164 real(dp) :: domega ! Step of the frequency mesh. 165 real(dp) :: broad ! Lorentzian Broadening. 166 real(dp) :: ircut ! Infrared cutoff for transitions 167 real(dp) :: uvcut ! Ultraviolet cutoff for transitions. 168 real(dp) :: haydock_tol(2) ! Tolerance for stopping the Haydock algorithm. 169 real(dp) :: cg_tolwfr ! Tolerance for stopping the CG algorithm 170 171 !Interp@BSE 172 real(dp) :: interp_m3_width ! Width of the interpolated M3 method along the diagonal 173 174 logical :: use_diagonal_Wgg ! Use diagonal approximation for Wgg. 175 logical :: use_coulomb_term ! Include W term in the BS Hamiltonian. 176 logical :: have_complex_ene ! .TRUE. if energies have a non-zero imaginary part. 177 178 !Interp@BSE 179 logical :: use_interp ! .TRUE. if we use interpolation technique 180 logical :: prep_interp ! .TRUE. if we prepare interpolation with ABC 181 logical :: sum_overlaps ! .TRUE. if making the sum of the overlaps to 1 182 logical :: prt_ncham ! .TRUE. if we dump the hamiltonian in NetCDF 183 184 logical :: do_ep_renorm ! .TRUE. for electron-phonon renormalization of the spectrum 185 logical :: do_lifetime ! .TRUE. if using elphon lifetime (not yet implemented) 186 187 !arrays 188 integer :: mg0(3) ! For each reduced direction gives the max G0 component 189 ! to account for umklapp processes 190 191 integer,allocatable :: nreh(:) 192 ! nreh(nsppol) 193 ! Number of resonant electron-hole transitions for each spin. 194 195 integer,allocatable :: vcks2t(:,:,:,:) 196 ! vcks2t(v,c,ik_bz,spin) gives the transition index associated to (v,c,kbz,spin) 197 198 !Interp@BSE 199 integer :: interp_kmult(3) ! Factor to subdivide kmesh sampling 200 201 integer,allocatable :: nreh_interp(:) 202 ! nreh_interp(nsppol) 203 ! Number of transitions for the interpolated mesh 204 205 integer,allocatable :: vcks2t_interp(:,:,:,:) 206 ! vcks2t(v,c,ik_bz_dense,spin) : Transition index for the dense kmesh 207 208 real(dp),allocatable :: q(:,:) ! Q-points for optical limit (reduced coordinates). 209 210 complex(dpc),allocatable :: omega(:) 211 ! omega(nomega) 212 ! Frequency mesh for epsilon (including the complex imaginary shift) 213 214 type(transition),allocatable :: Trans(:,:) 215 ! Trans(max_nreh,nsppol) 216 217 type(transition),allocatable :: Trans_interp(:,:) 218 ! Transitions for interpolated mesh 219 220 end type excparam 221 222 public :: bs_parameters_free 223 public :: print_bs_parameters 224 public :: bsp_calctype2str
m_bs_defs/init_transitions [ Functions ]
[ Top ] [ m_bs_defs ] [ Functions ]
NAME
init_transitions
FUNCTION
Main creation method for the transition structured datatype.
INPUTS
lomo_spin(nsppol) humo_spin(nsppol) ir_cut,uv_cut nkbz=Number of k-points in the BZ. nbnds=Maximum number of bands. nkibz=Number of k-points in the IBZ. nsppol=Number of spins. nspinor=Number of spinor components. gw_energy occ=Occupation factors. ktab
OUTPUT
max_tene=Maximum transition energy. nreh(nsppol)=Number of resonant transitions for each spin.
SIDE EFFECTS
Trans(:,:) input: allocatable array output: Trans(max_nreh,nsppol) stores the correspondence t -> (band,kbz,spin) and the transition energy.
SOURCE
550 subroutine init_transitions(Trans,lomo_spin,humo_spin,ir_cut,uv_cut,nkbz,nbnds,nkibz,nsppol,nspinor,gw_energy,occ,ktab,& 551 minmax_tene,nreh) 552 553 !Arguments ------------------------------------ 554 !scalars 555 integer,intent(in) :: nkbz,nbnds,nkibz,nsppol,nspinor 556 real(dp),intent(in) :: ir_cut,uv_cut 557 real(dp),intent(out) :: minmax_tene(2) 558 type(transition),allocatable,intent(out) :: Trans(:,:) 559 !arrays 560 integer,intent(in) :: lomo_spin(nsppol),humo_spin(nsppol) 561 integer,intent(in) :: ktab(nkbz) 562 integer,intent(out) :: nreh(nsppol) 563 real(dp),intent(in) :: occ(nbnds,nkibz,nsppol) 564 complex(dpc),intent(in) :: gw_energy(nbnds,nkibz,nsppol) 565 566 !Local variables ------------------------------ 567 !scalars 568 integer :: spin,it,ik_bz,ik_ibz,iv,ic,max_occ,sweep,max_nreh,lomo,humo 569 real(dp) :: tene, delta_f,min_tene,max_tene 570 complex(dpc) :: cplx_enet 571 logical :: add_transition 572 573 !************************************************************************ 574 575 ! Find transitions 576 max_occ=2/(nsppol*nspinor) 577 nreh=0 578 min_tene = -one; max_tene = zero 579 ! 580 ! sweep=1 calculats the number of resonants transitions taking into 581 ! account a possible energy cutoff. 582 ! sweep=2 initializes the tables describing the e-h transition. 583 ! 584 do sweep=1,2 585 ! 586 if (sweep==2) then 587 ! Allocate Trans structure. 588 max_nreh = MAXVAL(nreh) 589 ABI_MALLOC(Trans, (max_nreh,nsppol)) 590 end if 591 ! 592 do spin=1,nsppol 593 it=0 594 lomo = lomo_spin(spin) 595 humo = humo_spin(spin) 596 do ik_bz=1,nkbz 597 ik_ibz=ktab(ik_bz) 598 ! 599 do iv=lomo,humo 600 do ic=lomo,humo 601 delta_f = ( occ(ic,ik_ibz,spin)-occ(iv,ik_ibz,spin) ) / max_occ 602 cplx_enet = gw_energy(ic,ik_ibz,spin)-gw_energy(iv,ik_ibz,spin) 603 tene = DBLE(cplx_enet) 604 605 add_transition = & 606 & (tene > tol12) .and. & ! Resonant transition. 607 & ( ABS(delta_f) > tol12) .and. & ! c-v transition. 608 & (tene < uv_cut .and. tene > ir_cut) ! Energy cutoff. 609 610 if (add_transition) then 611 it = it + 1 612 max_tene = MAX(max_tene, tene) 613 min_tene = MAX(min_tene, tene) 614 end if 615 if (add_transition.and.sweep==2) then 616 Trans(it,spin)%k = ik_bz 617 Trans(it,spin)%v = iv 618 Trans(it,spin)%c = ic 619 Trans(it,spin)%en = cplx_enet 620 end if 621 622 end do 623 end do 624 end do ! ik_bz 625 ! Save number of transitions for this spin. 626 if (sweep==1) nreh(spin) = it 627 end do ! spin 628 ! 629 end do ! sweep 630 631 minmax_tene = [min_tene, max_tene] 632 633 end subroutine init_transitions
m_bs_defs/print_bs_files [ Functions ]
[ Top ] [ m_bs_defs ] [ Functions ]
NAME
print_bs_files
FUNCTION
Printout of the content of the excfiles structure.
INPUTS
BS_files<excfiles>=An object of type excfile storing the filenames used in the Bethe-Salpeter code. [unit]=Unit number for output [prtvol]=Verbosity level [mode_paral]=Either "COLL" or "PERS" [header]=String to be printed as header for additional info.
OUTPUT
Only printing.
SOURCE
737 subroutine print_bs_files(BS_files,header,unit,mode_paral,prtvol) 738 739 !Arguments ------------------------------------ 740 class(excfiles),intent(in) :: BS_files 741 integer,optional,intent(in) :: unit,prtvol 742 character(len=4),optional,intent(in) :: mode_paral 743 character(len=*),optional,intent(in) :: header 744 745 !Local variables ------------------------------ 746 !scalars 747 integer :: my_unt,my_prtvol 748 character(len=4) :: my_mode 749 character(len=500) :: msg 750 ! ********************************************************************* 751 752 !@excfiles 753 my_unt =std_out; if (PRESENT(unit )) my_unt =unit 754 my_prtvol=0 ; if (PRESENT(prtvol )) my_prtvol=prtvol 755 my_mode ='COLL' ; if (PRESENT(mode_paral)) my_mode =mode_paral 756 757 msg=' ==== Files used for the Bethe-Salpeter calculation ==== ' 758 if (PRESENT(header)) msg=' ==== '//TRIM(ADJUSTL(header))//' ==== ' 759 call wrtout(my_unt,msg,my_mode) 760 761 if (BS_files%in_hreso /= BSE_NOFILE) then 762 call wrtout(my_unt," Resonant block will be read from: "//trim(BS_files%in_hreso),my_mode) 763 end if 764 765 if (BS_files%in_hcoup /= BSE_NOFILE) then 766 call wrtout(my_unt," Coupling block will be read from: "//trim(BS_files%in_hcoup),my_mode) 767 end if 768 769 if (BS_files%in_eig /= BSE_NOFILE) then 770 call wrtout(my_unt," BS eigenstates will be read from: "//trim(BS_files%in_eig),my_mode) 771 end if 772 773 if (BS_files%in_haydock_basename /= BSE_NOFILE) then 774 call wrtout(my_unt," Haydock restart files have basename: "//trim(BS_files%in_haydock_basename),my_mode) 775 end if 776 777 end subroutine print_bs_files
m_bs_defs/print_bs_parameters [ Functions ]
[ Top ] [ m_bs_defs ] [ Functions ]
NAME
print_bs_parameters
FUNCTION
Printout of the parameters used for the BS calculation.
INPUTS
p<excparam>=Datatype storing the parameters of the Bethe-Salpeter calculation.
OUTPUT
Only printing.
SOURCE
329 subroutine print_bs_parameters(BSp, header, unit, mode_paral, prtvol) 330 331 !Arguments ------------------------------------ 332 !scalars 333 class(excparam),intent(in) :: BSp 334 integer,optional,intent(in) :: unit,prtvol 335 character(len=4),optional,intent(in) :: mode_paral 336 character(len=*),optional,intent(in) :: header 337 338 !Local variables ------------------------------ 339 !scalars 340 integer :: my_unt,my_prtvol,iq,ii,spin 341 character(len=4) :: my_mode 342 character(len=500) :: msg 343 344 ! ********************************************************************* 345 346 my_unt =std_out; if (PRESENT(unit )) my_unt =unit 347 my_prtvol=0 ; if (PRESENT(prtvol )) my_prtvol=prtvol 348 my_mode ='COLL' ; if (PRESENT(mode_paral)) my_mode =mode_paral 349 350 msg=' ==== Parameters of the Bethe-Salpeter run ==== ' 351 if (PRESENT(header)) msg=' ==== '//TRIM(ADJUSTL(header))//' ==== ' 352 call wrtout(my_unt,msg,my_mode) 353 354 select case (Bsp%algorithm) 355 case (BSE_ALGO_NONE) 356 msg = " Algorithm: Build Hamiltonian but skip the calculation of the spectrum." 357 case (BSE_ALGO_DDIAGO) 358 msg = " Algorithm: Direct diagonalization." 359 case (BSE_ALGO_HAYDOCK) 360 msg = " Algorithm: Haydock technique." 361 case (BSE_ALGO_CG) 362 msg = " Algorithm: Conjugate gradient." 363 case default 364 msg = " Algorithm: Unknown!." 365 end select 366 call wrtout(my_unt,msg,my_mode) 367 368 write(msg,'(4(a,i0,a))')& 369 ' Dimension of the v, W matrices, npweps = ',BSp%npweps,ch10,& 370 ' Cutoff for the wavefunctions, npwwfn = ',BSp%npwwfn,ch10,& 371 ' Number of k-points in the IBZ, nkibz = ',BSp%nkibz,ch10,& 372 ' Highest empty band included, nband = ',BSp%nbnds,"" 373 call wrtout(my_unt,msg,my_mode) 374 375 do spin=1,Bsp%nsppol 376 msg = " === Spin UP ==="; if (spin == 2) msg = " === Spin DOWN ===" 377 call wrtout(my_unt,msg,my_mode) 378 write(msg,'(5(a,i0,a))')& 379 ' Number of resonant transitions ',BSp%nreh(spin),ch10,& 380 ' Lowest occupied state ',BSp%lomo_spin(spin),ch10,& 381 ' Highest occupied state ',BSp%homo_spin(spin),ch10,& 382 ' Lowest unoccupied state ',BSp%lumo_spin(spin),ch10,& 383 ' Highest unoccupied state ',BSp%nbnds,"" 384 ! ' Number of valence bands ',BSp%nbndv,ch10,& 385 ! ' Number of conduction bands ',BSp%nbndc,"" 386 call wrtout(my_unt,msg,my_mode) 387 end do 388 389 write(msg,'(3(a,f6.2,a),a,f6.2)')& 390 ' Minimum frequency [eV] Emin = ',BSp%omegai*Ha_eV,ch10,& 391 ' Maximum frequency [eV] Emax = ',BSp%omegae*Ha_eV,ch10,& 392 ' Frequency step [eV] dE = ',BSp%domega*Ha_eV,ch10,& 393 ' Lorentzian broadening [eV] eta = ',BSp%broad*Ha_eV 394 call wrtout(my_unt,msg,my_mode) 395 396 ! Calculation type 397 call bsp_calctype2str(Bsp,msg) 398 call wrtout(my_unt,msg,my_mode) 399 400 if (ABS(Bsp%mbpt_sciss)>tol6) then 401 write(msg,'(a,f5.2)')" Scissors operator energy [eV] = ",Bsp%mbpt_sciss*Ha_eV 402 call wrtout(my_unt,msg,my_mode) 403 end if 404 405 msg=' Local fields effects (v term) excluded' 406 if (BSp%exchange_term>0) msg=' Local fields effects (v term) included' 407 call wrtout(my_unt,msg,my_mode) 408 409 msg=' Excitonic effects (W term) excluded' 410 if (BSp%use_coulomb_term) msg=' Excitonic effects (W term) included' 411 call wrtout(my_unt,msg,my_mode) 412 413 if (BSp%use_coulomb_term) then 414 msg=" Full W_GG' included" 415 if (BSp%use_diagonal_Wgg) msg=' Only diagonal term W_GG included' 416 call wrtout(my_unt,msg,my_mode) 417 if (BSp%wtype==BSE_WTYPE_FROM_SCR) then 418 call wrtout(my_unt," W is read from an external SCR file",my_mode) 419 end if 420 if (BSp%wtype==BSE_WTYPE_FROM_MDL) then 421 call wrtout(my_unt," W is approximated with the model dielectric function",my_mode) 422 end if 423 end if 424 425 msg=' Resonant-only calculation (Hermitian case)' 426 if (BSp%use_coupling>0) msg=' Resonant + Coupling calculation' 427 call wrtout(my_unt,msg,my_mode) 428 429 if(Bsp%use_interp) then 430 call wrtout(my_unt,' Interpolation technique used',my_mode) 431 end if 432 433 if(Bsp%use_interp) then 434 select case (Bsp%interp_mode) 435 case (1) 436 msg = ' Interpolation using WFK on the dense mesh' 437 case (2) 438 msg = ' Interpolation using WFK on the dense mesh + ABC divergence' 439 case (3) 440 msg = ' Interpolation using WFK on the dense mesh + ABC divergence along diagonal' 441 case (4) 442 msg = ' Interpolation using WFK on the dense mesh' 443 case default 444 msg = ' Unknown interpolation technique' 445 end select 446 call wrtout(my_unt,msg,my_mode) 447 448 if(BSp%prep_interp) then 449 call wrtout(my_unt,' Prepare interpolation technique with ABC',my_mode) 450 end if 451 452 if(BSp%interp_method == BSE_INTERP_YG) then 453 write(msg,'(a)') " Use Y. Gillet interpolation with 8 neighbours" 454 call wrtout(my_unt, msg, my_mode) 455 else if(BSP%interp_method == BSE_INTERP_RL2) then 456 write(msg,'(a)') " Use only 2 neighbours to interpolate linearly (Debug mode)" 457 call wrtout(my_unt, msg, my_mode) 458 else if(BSp%interp_method == BSE_INTERP_RL) then 459 write(msg,'(a,i0)') " Use Rohlfing and Louie with nb = ",BSp%rl_nb 460 call wrtout(my_unt, msg, my_mode) 461 end if 462 463 if(BSp%sum_overlaps) then 464 call wrtout(my_unt, " Summing overlaps to 1 in the interpolation",my_mode) 465 end if 466 end if 467 468 write(msg,'(a)')ch10 469 call wrtout(my_unt,msg,my_mode) 470 471 call wrtout(my_unt,' Calculating epsilon_Macro(q-->0,w), along the following directions:',my_mode) 472 do iq=1,BSp%nq 473 write(msg,'(a,3f10.6,2a)')' q = (',(BSp%q(ii,iq),ii=1,3), ') [r.l.u.]' 474 call wrtout(my_unt,msg,my_mode) 475 end do 476 477 !TODO 478 !Add file sizes and size of the buffer used for the matrix. 479 480 end subroutine print_bs_parameters
m_bs_defs/repr_1trans [ Functions ]
[ Top ] [ m_bs_defs ] [ Functions ]
NAME
repr_1trans
FUNCTION
Returns a string with info on the (k,v,c) transition.
INPUTS
Trans<transition>=structure datatype containing indececes and info on the optical transition. [prtvol]=Verbosity level. Defaults to 0.
OUTPUT
str(len=500)=The string representing the transition.
SOURCE
654 pure function repr_1trans(Trans, prtvol) result(str) 655 656 !Arguments ------------------------------------ 657 class(transition),intent(in) :: Trans 658 integer,optional,intent(in) :: prtvol 659 character(len=500) :: str 660 661 !Local variables ------------------------------ 662 !scalars 663 integer :: my_prtvol 664 665 !************************************************************************ 666 667 my_prtvol=0; if (PRESENT(prtvol)) my_prtvol=prtvol 668 669 if (my_prtvol==0) then 670 write(str,'(3(a,i3))')" k= ",Trans%k," v= ",Trans%v," c= ",Trans%c 671 else 672 write(str,'(3(a,i3),a,2f6.2)')" k= ",Trans%k," v= ",Trans%v," c= ",Trans%c," ene= ",Trans%en*Ha_eV 673 end if 674 675 end function repr_1trans
m_bs_defs/repr_2trans [ Functions ]
[ Top ] [ m_bs_defs ] [ Functions ]
NAME
repr_2trans
FUNCTION
Returns a string with info on two transitions
INPUTS
Trans1, Trans2<transition>=structure datatypes containing indececes and info on the optical transitions [prtvol]=Verbosity level. Defaults to 0.
OUTPUT
string(len=500)=The string representing the transition.
SOURCE
696 pure function repr_2trans(Trans1,Trans2,prtvol) result(string) 697 698 !Arguments ------------------------------------ 699 !scalars 700 integer,optional,intent(in) :: prtvol 701 character(len=500) :: string 702 type(transition),intent(in) :: Trans1,Trans2 703 704 !Local variables ------------------------------ 705 !scalars 706 integer :: my_prtvol 707 708 !************************************************************************ 709 710 my_prtvol=0; if (PRESENT(prtvol)) my_prtvol=prtvol 711 string = repr_1trans(Trans1,my_prtvol)//" | "//trim(repr_1trans(Trans2,my_prtvol)) 712 713 end function repr_2trans
m_bs_defs/transition [ Types ]
[ Top ] [ m_bs_defs ] [ Types ]
NAME
transition
FUNCTION
The transition derived data type is used to store the correspondence between the transition index and the set of quantum numbers (ik_bz,v,c) The energy of the transition is stored as well.
SOURCE
79 type,public :: transition 80 integer :: k=0 ! Index of the k-point in the BZ 81 integer :: v=0 ! Valence band index. 82 integer :: c=0 ! Conduction band index. 83 complex(dpc) :: en=huge(one) ! Transition energy 84 end type transition 85 86 public :: init_transitions ! Main creation method. 87 public :: repr_trans ! Returns a string representing the transition or a couple of transitions.