Line data Source code
1 : !--------------------------------------------------------------------------------
2 : ! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
3 : ! This file is part of FLEUR and available as free software under the conditions
4 : ! of the MIT license as expressed in the LICENSE file in more detail.
5 : !--------------------------------------------------------------------------------
6 :
7 : module m_wann_socmat_vec
8 : use m_juDFT
9 : contains
10 0 : subroutine wann_socmat_vec(
11 : > fmpi,enpara,input,noco,nococonv,atoms,
12 : > lmaxd,natd,neigd,
13 : > llod,jspd,
14 : > jspins,irank,
15 0 : > vr,
16 0 : > acof,bcof,chelp,
17 0 : < hsomtx_vec)
18 : c***********************************************************************
19 : c Calculate matrix elements of the spin-orbit interaction for those
20 : c Bloch functions out of which the Wannier functions are
21 : c constructed. From these matrix elements the spin-orbit Hamiltonian
22 : c in the basis set of Wannier functions may be constructed.
23 : c
24 : c Frank Freimuth
25 : c***********************************************************************
26 : USE m_types
27 : USE m_spnorb
28 : USE m_hsoham
29 : implicit none
30 :
31 : TYPE(t_mpi),INTENT(IN) :: fmpi
32 :
33 : TYPE(t_enpara),INTENT(IN) :: enpara
34 : TYPE(t_input),INTENT(IN) :: input
35 : TYPE(t_noco),INTENT(IN) :: noco
36 : TYPE(t_nococonv),INTENT(IN) :: nococonv
37 : TYPE(t_atoms),INTENT(IN) :: atoms
38 :
39 : integer, intent(in) :: lmaxd
40 : integer, intent(in) :: natd
41 : integer, intent(in) :: neigd
42 :
43 : integer, intent(in) :: llod
44 : integer, intent(in) :: jspd
45 :
46 : integer, intent(in) :: jspins
47 : integer, intent(in) :: irank
48 :
49 : real, intent(in) :: vr(:,0:,:,:) ! jmtd,0:nlhd,ntypd,jspd
50 :
51 : complex, intent(in) :: acof(:,0:,:,:) !acof(noccbd,0:lmd,natd,jspd)
52 : complex, intent(in) :: bcof(:,0:,:,:) !bcof(noccbd,0:lmd,natd,jspd)
53 : complex, intent(in) :: chelp(-llod:,:,:,:,:) !chelp(-llod:llod,neigd,nlod,natd,jspd)
54 :
55 : complex, intent(out):: hsomtx_vec(:,:,:,:,:) !(2,2,neigd,neigd,3)
56 :
57 : integer :: n,l,nwdd,nw,ispin,ie,na,ll1,m,lm,i,j,i1,i2,nsz(2),compo
58 : real :: s(3),r2
59 : logical :: l_all
60 : CHARACTER*3 chntype
61 : real :: theta,phi,pi
62 :
63 0 : COMPLEX,ALLOCATABLE :: ahelp(:,:,:,:)
64 0 : complex,allocatable :: bhelp(:,:,:,:)
65 0 : COMPLEX,ALLOCATABLE :: hsomtx(:,:,:,:)
66 :
67 0 : TYPE(t_rsoc)::rsoc
68 0 : TYPE(t_usdus):: usdus
69 :
70 0 : call timestart("wann_socmat_vec")
71 0 : nwdd=1
72 0 : nw=1
73 0 : nsz=neigd
74 :
75 0 : theta= 0.0
76 0 : phi= 0.0
77 :
78 0 : ALLOCATE ( ahelp(lmaxd*(lmaxd+2),natd,neigd,input%jspins) )
79 0 : ALLOCATE ( bhelp(lmaxd*(lmaxd+2),natd,neigd,input%jspins) )
80 0 : allocate ( hsomtx(neigd,neigd,2,2) )
81 :
82 0 : do ispin=1,input%jspins
83 0 : DO ie = 1, neigd
84 0 : DO na = 1, natd
85 0 : DO l = 1, lmaxd
86 0 : ll1 = l*(l+1)
87 0 : DO m = -l,l
88 0 : lm = ll1 + m
89 0 : ahelp(lm,na,ie,ispin) = (acof(ie,lm,na,ispin))
90 0 : bhelp(lm,na,ie,ispin) = (bcof(ie,lm,na,ispin))
91 : ENDDO !m
92 : ENDDO !l
93 : ENDDO !na
94 : ENDDO !ie
95 : enddo !ispin
96 :
97 : ALLOCATE(usdus%us(0:atoms%lmaxd,atoms%ntype,jspd),
98 : + usdus%dus(0:atoms%lmaxd,atoms%ntype,jspd),
99 : + usdus%uds(0:atoms%lmaxd,atoms%ntype,jspd),
100 : + usdus%duds(0:atoms%lmaxd,atoms%ntype,jspd),
101 : + usdus%ddn(0:atoms%lmaxd,atoms%ntype,jspd),
102 : + usdus%ulos(atoms%nlod,atoms%ntype,jspd),
103 : + usdus%dulos(atoms%nlod,atoms%ntype,jspd),
104 : + usdus%uulon(atoms%nlod,atoms%ntype,jspd),
105 0 : + usdus%dulon(atoms%nlod,atoms%ntype,jspd))
106 :
107 :
108 : CALL spnorb(
109 : > atoms,noco,nococonv,input,fmpi, enpara,
110 0 : > vr,usdus,rsoc,.true.)
111 :
112 0 : do compo=1,3
113 0 : rsoc%soangl= cmplx(0.0,0.0)
114 :
115 :
116 : CALL spnorb_angles(atoms,fmpi,nococonv%theta,
117 0 : > nococonv%phi,rsoc%soangl,compo)
118 :
119 : CALL hsoham(atoms,noco,input,nsz,neigd,chelp,rsoc,ahelp,
120 : > bhelp,1,natd,fmpi%n_rank,fmpi%n_size,fmpi%SUB_COMM,
121 0 : < hsomtx)
122 :
123 0 : do i1 = 1,2
124 0 : do i2 = 1,2
125 0 : do i = 1,neigd
126 0 : do j = 1,neigd
127 0 : hsomtx_vec(i1,i2,i,j,compo) = hsomtx(i,j,i1,i2)
128 : enddo
129 : enddo
130 : enddo
131 : enddo
132 :
133 : enddo !compo
134 0 : call timestart("wann_socmat_vec")
135 0 : end subroutine wann_socmat_vec
136 : end module m_wann_socmat_vec
|