TABLE OF CONTENTS
- ABINIT/m_melemts
- m_melemts/melements_free
- m_melemts/melements_herm
- m_melemts/melements_init
- m_melemts/melements_mpisum
- m_melemts/melements_print
- m_melemts/melements_t
- m_melemts/melements_zero
- m_melemts/melflags_copy
- m_melemts/melflags_reset
- m_melemts/melflags_t
- m_melemts/my_select_melements
- m_sigma/mels_get_exene_core
ABINIT/m_melemts [ Modules ]
NAME
m_melemts
FUNCTION
This module defines an object used as database to store matrix elements of several potentials and operators between two Bloch states. These values are used in the GW part of abinit to evaluate QP energies using the perturbative approach.
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 .
TODO
NOTES
* This module is supposed to be used only in the GW part to facilitate we might decide to use ragged arrays Mels(nkcalc, nsppol*nspinor**2)%data or replaced nkibz with nkcalc with reduce memory
SOURCE
29 #if defined HAVE_CONFIG_H 30 #include "config.h" 31 #endif 32 33 #include "abi_common.h" 34 35 module m_melemts 36 37 use defs_basis 38 use m_errors 39 use m_xmpi 40 use m_abicore 41 42 use m_fstrings, only : tolower, sjoin 43 use m_numeric_tools, only : print_arr 44 45 implicit none 46 47 private
m_melemts/melements_free [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
melements_free
FUNCTION
Free all dynamic memory of the database
SOURCE
281 subroutine melements_free(Mels) 282 283 !Arguments ------------------------------------ 284 class(melements_t),intent(inout) :: Mels 285 286 ! ************************************************************************* 287 288 ! integer arrays 289 ABI_SFREE(Mels%bands_idx) 290 ABI_SFREE(Mels%iscalc) 291 292 ! real arrays 293 ABI_SFREE(Mels%kibz) 294 295 ! complex arrays 296 ABI_SFREE(Mels%kinetic) 297 ABI_SFREE(Mels%hbare) 298 ABI_SFREE(Mels%sxcore) 299 ABI_SFREE(Mels%vhartree) 300 ABI_SFREE(Mels%vlexx) 301 ABI_SFREE(Mels%vu) 302 ABI_SFREE(Mels%vxc) 303 ABI_SFREE(Mels%vxcval) 304 ABI_SFREE(Mels%vxcval_hybrid) 305 306 ! Reset all has_* flags. 307 call Mels%flags%reset() 308 309 end subroutine melements_free
m_melemts/melements_herm [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
melements_herm
FUNCTION
Reconstruc the lower triangle of all calculated arrays. Assuming Hermitian operator. Works both for collinear and non-collinear case.
INPUTS
Mels=The database [aname]=The name of the array to be symmetrized, by default all calculated arrays are filled.
SIDE EFFECTS
All arrays whose flag is 2, are filled assuming an Hermitian operator.
SOURCE
508 subroutine melements_herm(Mels, aname) 509 510 !Arguments ------------------------------------ 511 !scalars 512 class(melements_t),intent(inout) :: Mels 513 character(len=*),optional,intent(in) :: aname 514 515 !Local variables------------------------------- 516 integer :: is,ik,ib,jb,iab,iab_tr,iname 517 integer,pointer :: flag_p 518 character(len=NAMELEN) :: key 519 !arrays 520 integer,parameter :: trsp_idx(2:4) = [2,4,3] 521 complex(dpc),ABI_CONTIGUOUS pointer :: arr_p(:,:,:,:) 522 523 ! ************************************************************************* 524 525 ! === Symmetrize matrix elements === 526 ! * In the collinear case, generate the lower triangle by just doing a complex conjugate. 527 ! * In the noncollinear case do also a transposition since A_{12}^{ab} = A_{21}^{ba}^* 528 ! 2-->2, 3-->4, 4-->3 529 ! 530 do iname=1,NNAMES 531 key = ANAMES(iname) 532 if (PRESENT(aname)) then 533 if (key /= aname) CYCLE 534 end if 535 536 call my_select_melements(Mels, key, flag_p, arr_p) 537 538 if (flag_p>0) then 539 do ik=1,Mels%nkibz 540 do is=1,Mels%nsppol 541 542 do jb=Mels%bmin,Mels%bmax 543 do ib=Mels%bmin,jb ! Upper triangle 544 545 if (ib/=jb) then 546 arr_p(jb,ib,ik,is)=CONJG(arr_p(ib,jb,ik,is)) 547 if (Mels%nspinor==2) then 548 do iab=2,4 549 iab_tr=trsp_idx(iab) 550 arr_p(jb,ib,ik,iab)=CONJG(arr_p(ib,jb,ik,iab_tr)) 551 end do 552 end if 553 else ! For ib==jb force real-valued 554 arr_p(jb,ib,ik,is)=half*(arr_p(jb,ib,ik,is)+CONJG(arr_p(jb,ib,ik,is))) 555 if (Mels%nspinor==2) arr_p(jb,ib,ik,2)=half*(arr_p(ib,jb,ik,2)+CONJG(arr_p(ib,jb,ik,2))) 556 end if 557 558 end do !ib 559 end do !jb 560 561 end do !is 562 end do !ik 563 end if 564 565 end do !inames 566 567 end subroutine melements_herm
m_melemts/melements_init [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
melements_init
FUNCTION
Initialize the database, allocate arrays according to Mflags_in, zeroing the content of the allocated arrays.
INPUTS
nsppol=Number of independent spin polarizations. nspden=Number of spin-density components nspinor=Number of spinor components band_idx=min and Max band index for each ik_ibz and spin
OUTPUT
Mels=The initialized database with dimensions and allocated memory.
SOURCE
399 subroutine melements_init(Mels, Mflags_in, nsppol, nspden, nspinor, nkibz, kibz, bands_idx) 400 401 !Arguments ------------------------------------ 402 !scalars 403 class(melements_t),intent(out) :: Mels 404 integer,intent(in) :: nspinor,nspden,nsppol,nkibz 405 type(melflags_t),intent(in) :: Mflags_in 406 !arrays 407 integer,intent(in) :: bands_idx(2, nkibz, nsppol) 408 real(dp),intent(in) :: kibz(3,nkibz) 409 410 !Local variables------------------------------- 411 integer :: ikibz,isppol,bmin,bmax,b1,b2 412 413 ! ************************************************************************* 414 415 ! Copy flags. 416 call Mflags_in%copy(Mels%flags) 417 418 ! Copy dimensions. 419 Mels%nkibz = nkibz 420 Mels%nsppol = nsppol 421 Mels%nspinor = nspinor 422 Mels%nspden = nspden 423 424 ABI_MALLOC(Mels%bands_idx, (2, nkibz, nsppol)) 425 Mels%bands_idx = bands_idx 426 427 ABI_MALLOC(Mels%iscalc,(nkibz, nsppol)) 428 Mels%iscalc = 0 429 430 bmin = HUGE(1); bmax =-HUGE(1) 431 do isppol=1,Mels%nsppol 432 do ikibz=1,Mels%nkibz 433 if (ANY(Mels%bands_idx(:,ikibz,isppol)/=0)) then 434 b1 = Mels%bands_idx(1,ikibz,isppol) 435 b2 = Mels%bands_idx(2,ikibz,isppol) 436 Mels%iscalc(ikibz, isppol)=1 437 bmin = MIN(bmin,b1) 438 bmax = MAX(bmax,b2) 439 ABI_CHECK(b2 >= b1 .and. b1 > 0, "Wrong b1, b2") 440 end if 441 end do 442 end do 443 444 if (bmin==HUGE(1).or.bmax==-HUGE(1)) then 445 ABI_BUG("Wrong bands_idx") 446 end if 447 448 Mels%bmin = bmin 449 Mels%bmax = bmax 450 451 b1 = Mels%bmin; b2 = Mels%bmax 452 453 ! real arrays 454 ABI_MALLOC(Mels%kibz, (3,nkibz)) 455 Mels%kibz = kibz 456 457 ! complex arrays 458 if (Mels%flags%has_kinetic == 1) then 459 ABI_CALLOC(Mels%kinetic, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 460 end if 461 if (Mels%flags%has_hbare == 1) then 462 ABI_CALLOC(Mels%hbare, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 463 end if 464 if (Mels%flags%has_sxcore == 1) then 465 ABI_CALLOC(Mels%sxcore, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 466 end if 467 if (Mels%flags%has_vhartree == 1) then 468 ABI_CALLOC(Mels%vhartree, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 469 end if 470 if (Mels%flags%has_lexexch == 1) then 471 ABI_CALLOC(Mels%vlexx, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 472 end if 473 if (Mels%flags%has_vu == 1) then 474 ABI_CALLOC(Mels%vu, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 475 end if 476 if (Mels%flags%has_vxc == 1) then 477 ABI_CALLOC(Mels%vxc, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 478 end if 479 if (Mels%flags%has_vxcval == 1) then 480 ABI_CALLOC(Mels%vxcval, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 481 end if 482 if (Mels%flags%has_vxcval_hybrid == 1) then 483 ABI_CALLOC(Mels%vxcval_hybrid, (b1:b2,b1:b2,nkibz,nsppol*nspinor**2)) 484 end if 485 486 end subroutine melements_init
m_melemts/melements_mpisum [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
melements_mpisum
FUNCTION
Perform a collective SUM within the MPI communicator comm of the matrix elements stored in the database.
INPUTS
Mels=The database [aname]=The name of a particular array to be summed, by default all allocated arrays are considered.
SIDE EFFECTS
All arrays whose flag==1 are summed within the MPI communicator comm. In output the corresponding flas is set to 2.
SOURCE
591 subroutine melements_mpisum(Mels, comm, aname) 592 593 !Arguments ------------------------------------ 594 !scalars 595 class(melements_t),intent(inout) :: Mels 596 integer,intent(in) :: comm 597 character(len=*),optional,intent(in) :: aname 598 599 !Local variables------------------------------- 600 integer :: iname,ierr 601 integer,pointer :: flag_p 602 character(len=NAMELEN) :: key 603 !character(len=500) :: msg 604 !arrays 605 complex(dpc),ABI_CONTIGUOUS pointer :: arr_p(:,:,:,:) 606 607 ! ************************************************************************* 608 609 do iname=1,NNAMES 610 key = ANAMES(iname) 611 if (PRESENT(aname)) then 612 if (key /= aname) CYCLE 613 end if 614 615 call my_select_melements(Mels, key, flag_p, arr_p) 616 617 if (flag_p == 1) then 618 call xmpi_sum(arr_p, comm, ierr) 619 flag_p = 2 ! Tag this array as calculated 620 end if 621 end do 622 623 end subroutine melements_mpisum
m_melemts/melements_print [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
melements_print
FUNCTION
Printout of the content of all calculated array. Optionally, it is possible to print the content of a single entry of the database.
INPUTS
Mels=The database [unit]=the unit number for output, defaults to std_out [prtvol]=verbosity level, defaults to 0 [mode_paral]=either "COLL" or "PERS", default to "COLL" [header]=title for info
OUTPUT
Only writing
SOURCE
648 subroutine melements_print(Mels, names_list, header, unit, prtvol, mode_paral) 649 650 !Arguments ------------------------------------ 651 !scalars 652 class(melements_t),intent(in) :: Mels 653 integer,optional,intent(in) :: prtvol,unit 654 character(len=*),optional,intent(in) :: names_list(:) 655 character(len=*),optional,intent(in) :: header 656 character(len=4),optional,intent(in) :: mode_paral 657 658 !Local variables------------------------------- 659 integer :: my_unt,my_prtvol,max_r,max_c,ii 660 integer :: isppol,ikibz,iab,ib,b1,b2,my_nkeys,ikey 661 integer,pointer :: flag_p 662 character(len=4) :: my_mode 663 character(len=NAMELEN) :: key 664 character(len=500) :: msg,str,fmt 665 !arrays 666 integer,allocatable :: tab(:) 667 character(len=NAMELEN),allocatable :: my_keys(:) 668 complex(dpc),allocatable :: mat(:,:) 669 670 type rarr_dpc4 671 complex(dpc),ABI_CONTIGUOUS pointer :: arr_p(:,:,:,:) 672 end type rarr_dpc4 673 type(rarr_dpc4),allocatable :: data_p(:) 674 675 ! ************************************************************************* 676 677 !@melements_t 678 my_unt =std_out; if (PRESENT(unit )) my_unt =unit 679 my_prtvol=0 ; if (PRESENT(prtvol )) my_prtvol=prtvol 680 my_mode ='COLL' ; if (PRESENT(mode_paral)) my_mode =mode_paral 681 682 if (Mels%nspinor == 2) ABI_WARNING("nspinor=2 not coded") 683 684 if (PRESENT(names_list)) then 685 my_nkeys=SIZE(names_list) 686 ABI_MALLOC(my_keys,(my_nkeys)) 687 my_keys = names_list 688 else 689 my_nkeys=NNAMES 690 ABI_MALLOC(my_keys,(NNAMES)) 691 my_keys = ANAMES 692 end if 693 694 ABI_MALLOC(data_p,(my_nkeys)) 695 ABI_MALLOC(tab,(my_nkeys)) 696 tab=0 697 698 my_nkeys=0; str = " ib"; ii=4 699 do ikey=1,SIZE(my_keys) 700 key = my_keys(ikey) 701 call my_select_melements(Mels,key,flag_p,data_p(ikey)%arr_p) 702 if (flag_p==2) then 703 my_nkeys=my_nkeys+1 704 tab(my_nkeys)=ikey 705 str(ii+1:)=" "//TRIM(tolower(key)) 706 ii = ii+MAX(1+LEN_TRIM(key),10) 707 ABI_CHECK(ii<490,"I'm gonna SIGFAULT!") 708 end if 709 end do 710 711 write(msg,'(2a)')ch10,' === Matrix Elements stored in Mels% [eV] === ' 712 if (PRESENT(header)) write(msg,'(4a)')ch10,' === '//TRIM(ADJUSTL(header))//' [eV] === ' 713 call wrtout(my_unt,msg,my_mode) 714 if (Mels%nspinor == 2) then 715 call wrtout(my_unt, "Sum_ab M_ab, M_11, M_22, Re(M_12), IM(Re_12)" ,my_mode) 716 end if 717 718 if (my_nkeys==0) GOTO 10 719 write(fmt,'(a,i4,a)')'(1x,i3,',my_nkeys,'(1x,f9.5))' ! width of 10 chars 720 721 do isppol=1,Mels%nsppol 722 do ikibz=1,Mels%nkibz 723 if (Mels%iscalc(ikibz,isppol)/=1) CYCLE 724 725 write(msg,'(a,3es16.8,a,i2,a)')" kpt= (",Mels%kibz(:,ikibz),") spin=",isppol,":" 726 call wrtout(my_unt,msg,my_mode) 727 728 b1 = Mels%bands_idx(1,ikibz,isppol) 729 b2 = Mels%bands_idx(2,ikibz,isppol) 730 731 if (Mels%flags%only_diago==1.or.my_prtvol==0) then 732 ! Print only the diagonal. 733 write(msg,'(a)')str 734 call wrtout(my_unt,msg,my_mode) 735 do ib=b1,b2 736 if (Mels%nspinor == 1) then 737 write(msg,fmt)ib,(REAL(data_p(tab(ikey))%arr_p(ib,ib,ikibz,1))*Ha_eV, ikey=1,my_nkeys) 738 call wrtout(my_unt,msg,my_mode) 739 else 740 ! Write sum_ab, then diagonal elements, finally Re_12, Im_12 741 write(msg,fmt)ib,(real(sum(data_p(tab(ikey))%arr_p(ib,ib,ikibz,:)))*Ha_eV, ikey=1,my_nkeys) 742 call wrtout(my_unt,msg,my_mode) 743 if (my_prtvol > 0) then 744 write(msg,fmt)ib,(real(data_p(tab(ikey))%arr_p(ib,ib,ikibz,1))*Ha_eV, ikey=1,my_nkeys) 745 call wrtout(my_unt,msg,my_mode) 746 write(msg,fmt)ib,(real(data_p(tab(ikey))%arr_p(ib,ib,ikibz,2))*Ha_eV, ikey=1,my_nkeys) 747 call wrtout(my_unt,msg,my_mode) 748 write(msg,fmt)ib,(real(data_p(tab(ikey))%arr_p(ib,ib,ikibz,3))*Ha_eV, ikey=1,my_nkeys) 749 call wrtout(my_unt,msg,my_mode) 750 write(msg,fmt)ib,(aimag(data_p(tab(ikey))%arr_p(ib,ib,ikibz,3))*Ha_eV, ikey=1,my_nkeys) 751 call wrtout(my_unt,msg,my_mode) 752 end if 753 end if 754 end do 755 756 else 757 ! Print full matrix. 758 max_r = b2-b1+1 759 max_c = MIN(b2-b1+1, 9) 760 ABI_MALLOC(mat,(b1:b2,b1:b2)) 761 do ikey=1,my_nkeys 762 write(msg,'(3a)')" **** Off-diagonal elements of ",TRIM(my_keys(tab(ikey)))," **** " 763 call wrtout(my_unt,msg,my_mode) 764 do iab=1,Mels%nspinor**2 765 mat = data_p(tab(ikey))%arr_p(b1:b2,b1:b2,ikibz,iab)*Ha_eV 766 call print_arr(mat,max_r,max_c,my_unt,my_mode) 767 end do 768 write(msg,'(a)')ch10 769 call wrtout(my_unt,msg,my_mode) 770 end do 771 ABI_FREE(mat) 772 end if 773 774 end do !ikibz 775 end do ! isppol 776 777 10 continue 778 779 ABI_FREE(my_keys) 780 ABI_FREE(data_p) 781 ABI_FREE(tab) 782 783 end subroutine melements_print
m_melemts/melements_t [ Types ]
[ Top ] [ m_melemts ] [ Types ]
NAME
FUNCTION
Structure defining a database to store the matrix elements of operators needed for GW calculations.
SOURCE
95 type,public :: melements_t 96 97 integer :: nkibz 98 ! Number of k-points in the IBZ. 99 100 integer :: nsppol 101 ! Number of independent spin-polarizations. 102 103 integer :: nspinor 104 ! 1 for collinear, 2 for noncollinear. 105 106 integer :: nspden 107 ! Number of independent spin-density components. 108 109 integer :: bmin, bmax 110 ! min and Max band index over k-points and spin. 111 ! Used to dimension the arrays below. 112 113 integer, allocatable :: bands_idx(:,:,:) 114 ! (2, nkibz, nsppol) 115 ! min and Max band index for each k-point and spin. 116 117 integer, allocatable :: iscalc(:,:) 118 ! (nkibz, nsppol) 119 ! 1 if this k-point and spin has been calculated, 0 otherwise. 120 121 real(dp), allocatable :: kibz(:,:) 122 ! (3, nkibz) 123 ! The list of k-points in reduced coordinates. 124 125 complex(dpc), allocatable :: kinetic(:,:,:,:) 126 ! (b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 127 ! Matrix elements of the kinetic energy. 128 129 complex(dpc), allocatable :: hbare(:,:,:,:) 130 ! (b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 131 ! Matrix elements of the bare Hamiltonian. 132 133 complex(dpc), allocatable :: sxcore(:,:,:,:) 134 ! (b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 135 ! Matrix elements of the Fock operator generated by core electrons. 136 137 complex(dpc), allocatable :: vhartree(:,:,:,:) 138 ! (b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 139 ! Matrix elements of the Hartree potential. 140 141 complex(dpc), allocatable :: vlexx(:,:,:,:) 142 ! (b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 143 ! Matrix elements of the local exact exchange potential. 144 145 complex(dpc), allocatable :: vu(:,:,:,:) 146 ! vu(b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 147 ! Matrix elements of the U Hamiltonian. 148 149 complex(dpc), allocatable :: vxc(:,:,:,:) 150 ! (b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 151 ! Matrix elements of XC potential, including model core if present. 152 153 complex(dpc), allocatable :: vxcval(:,:,:,:) 154 ! (b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 155 ! Matrix elements of the XC potential, valence-only contribution. 156 157 complex(dpc), allocatable :: vxcval_hybrid(:,:,:,:) 158 ! (b1:b2, b1:b2, nkibz, nsppol*nspinor**2) 159 ! Matrix elements of the XC potential for hybrid calculations, valence-only contribution. 160 161 type(melflags_t) :: flags 162 163 contains 164 165 procedure :: init => melements_init ! Initialize the object 166 procedure :: free => melements_free ! Free memory 167 procedure :: herm => melements_herm ! Construct the lower triangle from the upper triangle 168 procedure :: mpisum => melements_mpisum ! Perform a collective SUM within the MPI communicator comm 169 procedure :: print => melements_print ! Print matrix elements 170 procedure :: zero => melements_zero ! Set matrix elements connecting states with different irrep to zero. 171 !procedure :: mels_get_exene_core 172 173 end type melements_t
m_melemts/melements_zero [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
melements_zero
FUNCTION
Set matrix elements connecting states with different irreducible representation to zero.
INPUTS
irrep_tab=Array used to select the entries that have to be set to zero. irrep_tab(ib,ik,is)=gives the index of the irreducible representation associated to state (ib,ik,is). [aname]=The name of the array to be symmetrized, by default all calculated arrays are filled.
SIDE EFFECTS
Mels= All arrays elements connecting states belonging to different irreps are set to zero.
SOURCE
806 subroutine melements_zero(Mels, irrep_tab, aname) 807 808 !Arguments ------------------------------------ 809 !scalars 810 class(melements_t),intent(inout) :: Mels 811 character(len=*),optional,intent(in) :: aname 812 !arrays 813 integer,intent(in) :: irrep_tab(Mels%bmin:Mels%bmax,Mels%nkibz,Mels%nsppol) 814 815 !Local variables------------------------------- 816 integer :: is,ik,ib,jb,iname,irrep_j,irrep_i 817 integer,pointer :: flag_p 818 character(len=NAMELEN) :: key 819 !arrays 820 complex(dpc),ABI_CONTIGUOUS pointer :: arr_p(:,:,:,:) 821 822 ! ************************************************************************* 823 824 do iname=1,NNAMES 825 key = ANAMES(iname) 826 if (PRESENT(aname)) then 827 if (key /= aname) CYCLE 828 end if 829 830 call my_select_melements(Mels,key,flag_p,arr_p) 831 832 if (flag_p>0) then 833 do is=1,Mels%nsppol 834 do ik=1,Mels%nkibz 835 836 do jb=Mels%bmin,Mels%bmax 837 irrep_j = irrep_tab(jb,ik,is) 838 do ib=Mels%bmin,Mels%bmax 839 irrep_i = irrep_tab(ib,ik,is) 840 ! 841 ! Set this matrix element to zero if the irreps are known and they differ. 842 if (irrep_i/=irrep_j .and. ALL((/irrep_i,irrep_j/) /=0) ) then 843 !write(std_out,*)"setting to zero ",ib,jb,ik,is 844 if (Mels%nspinor==2) then 845 arr_p(ib,jb,ik,is)=czero 846 else 847 arr_p(ib,jb,ik,:)=czero 848 end if 849 end if 850 851 end do !ib 852 end do !jb 853 854 end do !is 855 end do !ik 856 end if 857 858 end do !inames 859 860 end subroutine melements_zero
m_melemts/melflags_copy [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
melflags_copy
FUNCTION
Copy an object storing the flags.
INPUTS
Mflags_in=The flags to be copied.
OUTPUT
Mflags_out=The new set of flags.
SOURCE
245 subroutine melflags_copy(Mflags_in, Mflags_out) 246 247 !Arguments ------------------------------------ 248 class(melflags_t),intent(in) :: Mflags_in 249 class(melflags_t),intent(inout) :: Mflags_out 250 251 ! ************************************************************************* 252 253 call Mflags_out%reset() 254 255 ! @melflags_t 256 Mflags_out%has_kinetic = Mflags_in%has_kinetic 257 Mflags_out%has_hbare = Mflags_in%has_hbare 258 Mflags_out%has_sxcore = Mflags_in%has_sxcore 259 Mflags_out%has_vhartree = Mflags_in%has_vhartree 260 Mflags_out%has_vu = Mflags_in%has_vu 261 Mflags_out%has_vxc = Mflags_in%has_vxc 262 Mflags_out%has_vxcval = Mflags_in%has_vxcval 263 Mflags_out%has_vxcval_hybrid = Mflags_in%has_vxcval_hybrid 264 Mflags_out%has_lexexch = Mflags_in%has_lexexch 265 Mflags_out%only_diago = Mflags_in%only_diago 266 267 end subroutine melflags_copy
m_melemts/melflags_reset [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
melflags_reset
FUNCTION
Set all flags in melflags_t to 0.
INPUTS
SOURCE
207 subroutine melflags_reset(Mflags) 208 209 !Arguments ------------------------------------ 210 class(melflags_t),intent(inout) :: Mflags 211 212 ! ************************************************************************* 213 214 Mflags%has_kinetic = 0 215 Mflags%has_hbare = 0 216 Mflags%has_sxcore = 0 217 Mflags%has_vhartree = 0 218 Mflags%has_vu = 0 219 Mflags%has_vxc = 0 220 Mflags%has_vxcval = 0 221 Mflags%has_vxcval_hybrid = 0 222 Mflags%has_lexexch = 0 223 Mflags%only_diago = 0 224 225 end subroutine melflags_reset
m_melemts/melflags_t [ Types ]
[ Top ] [ m_melemts ] [ Types ]
NAME
FUNCTION
Container for the flags defining the status of the corresponding pointer defined in the type melements_t. Possible values are: * 0 if the correspondending array is not allocated. * 1 if allocated but not yet calculated. * 2 if allocated and calculated
SOURCE
64 type,public :: melflags_t 65 66 integer :: has_kinetic = 0 67 integer :: has_hbare = 0 68 integer :: has_lexexch = 0 69 integer :: has_sxcore = 0 70 integer :: has_vhartree = 0 71 integer :: has_vu = 0 72 integer :: has_vxc = 0 73 integer :: has_vxcval = 0 74 integer :: has_vxcval_hybrid = 0 75 integer :: only_diago = 0 76 ! 1 if only diagonal elements are calculated 77 78 contains 79 procedure :: reset => melflags_reset ! Reset the value of the flags. 80 procedure :: copy => melflags_copy ! Copy the object 81 end type melflags_t
m_melemts/my_select_melements [ Functions ]
[ Top ] [ m_melemts ] [ Functions ]
NAME
my_select_melements
FUNCTION
Helper function returning a pointer to the array "aname" as well as the status of the array
INPUTS
Mels<melements_t>=The database. aname=String with the name of the array.
OUTPUT
flag_p=Pointer to the integer defining the status of the array, see melflags_t. arr_p=The pointer to the array.
SOURCE
331 subroutine my_select_melements(Mels, aname, flag_p, arr_p) 332 333 !Arguments ------------------------------------ 334 !scalars 335 integer,pointer :: flag_p 336 character(len=*),intent(in) :: aname 337 type(melements_t),target,intent(in) :: Mels 338 !arrays 339 complex(dpc),ABI_CONTIGUOUS pointer :: arr_p(:,:,:,:) 340 341 ! ************************************************************************* 342 343 SELECT CASE (tolower(aname)) 344 CASE ("kinetic") 345 flag_p => Mels%flags%has_kinetic 346 arr_p => Mels%kinetic 347 CASE ("hbare") 348 flag_p => Mels%flags%has_hbare 349 arr_p => Mels%hbare 350 CASE ("sxcore") 351 flag_p => Mels%flags%has_sxcore 352 arr_p => Mels%sxcore 353 CASE ("vhartree") 354 flag_p => Mels%flags%has_vhartree 355 arr_p => Mels%vhartree 356 CASE ("vlexx") 357 flag_p => Mels%flags%has_lexexch 358 arr_p => Mels%vlexx 359 CASE ("vu") 360 flag_p => Mels%flags%has_vu 361 arr_p => Mels%vu 362 CASE ("vxc") 363 flag_p => Mels%flags%has_vxc 364 arr_p => Mels%vxc 365 CASE ("vxcval") 366 flag_p => Mels%flags%has_vxcval 367 arr_p => Mels%vxcval 368 CASE ("vxcval_hybrid") 369 flag_p => Mels%flags%has_vxcval_hybrid 370 arr_p => Mels%vxcval_hybrid 371 CASE DEFAULT 372 ABI_ERROR(sjoin("Wrong aname: ", aname)) 373 END SELECT 374 375 end subroutine my_select_melements
m_sigma/mels_get_exene_core [ Functions ]
[ Top ] [ m_sigma ] [ Functions ]
NAME
mels_get_exene_core
FUNCTION
Compute exchange energy.
INPUTS
mels<melements_t>=Matrix elements. kmesh<kmesh_t>=BZ sampling. bands<band_t>=Bands with occupation factors
SOURCE
879 !pure function mels_get_exene_core(mels,kmesh,bands) result(ex_energy) 880 ! 881 !!Arguments ------------------------------------ 882 !!scalars 883 ! real(dp) :: ex_energy 884 ! type(melements_t),intent(in) :: mels 885 ! type(kmesh_t),intent(in) :: kmesh 886 ! type(ebands_t),intent(in) :: bands 887 ! 888 !!Local variables------------------------------- 889 !!scalars 890 ! integer :: ik,ib,spin 891 ! real(dp) :: wtk,occ_bks 892 ! 893 !! ************************************************************************* 894 ! 895 ! ex_energy = zero 896 ! 897 ! do spin=1,mels%nsppol 898 ! do ik=1,mels%nkibz 899 ! wtk = kmesh%wt(ik) 900 ! do ib=mels%bmin,mels%bmax 901 ! occ_bks = bands%occ(ib,ik,spin) 902 ! if (mels%nspinor==1) then 903 ! ex_energy = ex_energy + half * occ_bks * wtk * mels%sxcore(ib,ib,ik,spin) 904 ! else 905 ! ex_energy = ex_energy + half * occ_bks wtk *SUM(mels%sxcore(ib,ib,ik,:)) 906 ! end if 907 ! end do 908 ! end do 909 ! end do 910 ! 911 !end function mels_get_exene_core 912 !!!*** 913 914 end module m_melemts