Line data Source code
1 : MODULE m_hsoham
2 : !
3 : !*********************************************************************
4 : ! set up spin-orbit contribution to hamiltonian
5 : !
6 : ! OpenMP parrallelization restored (U.Alekseeva 12.2019)
7 : !*********************************************************************
8 : !
9 :
10 : #ifdef CPP_MPI
11 : use mpi
12 : #endif
13 : CONTAINS
14 548 : SUBROUTINE hsoham(&
15 1096 : atoms,noco,input,nsz,neigd,chelp,rsoc,ahelp,bhelp,&
16 : nat_start,nat_stop,n_rank,n_size,SUB_COMM,&
17 548 : hsomtx)
18 :
19 : USE m_types
20 : IMPLICIT NONE
21 : #ifdef CPP_MPI
22 : INTEGER ierr(3)
23 : #endif
24 : TYPE(t_input),INTENT(IN) :: input
25 : TYPE(t_noco),INTENT(IN) :: noco
26 : TYPE(t_atoms),INTENT(IN) :: atoms
27 : TYPE(t_rsoc),INTENT(IN) :: rsoc
28 : ! ..
29 : ! .. Scalar Arguments ..
30 : ! ..
31 : INTEGER, INTENT (IN) :: nat_start,nat_stop,n_rank,n_size,SUB_COMM,neigd
32 : ! .. Array Arguments ..
33 : INTEGER, INTENT (IN) :: nsz(:)!(input%jspins)
34 : COMPLEX, INTENT (IN) :: ahelp((atoms%lmaxd+2)*atoms%lmaxd,nat_stop-nat_start+1,neigd,input%jspins)
35 : COMPLEX, INTENT (IN) :: bhelp((atoms%lmaxd+2)*atoms%lmaxd,nat_stop-nat_start+1,neigd,input%jspins)
36 : COMPLEX, INTENT (IN) :: chelp(-atoms%llod:atoms%llod,neigd,atoms%nlod,nat_stop-nat_start+1,input%jspins)
37 : COMPLEX, INTENT (OUT):: hsomtx(neigd,neigd,2,2)
38 : ! ..
39 : ! .. Local Scalars ..
40 : COMPLEX c_1,c_2,c_3,c_4,c_5
41 : INTEGER i,j,jsp,jsp1,l,lwn,m1,n,na,nn,i1,j1,ilo,ilop,m,nat_l,na_g,lm,ll1,lm1
42 : ! ..
43 : ! .. Local Arrays ..
44 548 : COMPLEX, ALLOCATABLE :: c_b(:,:),c_a(:,:),c_c(:,:,:),c_buf(:)
45 : ! ..
46 : !
47 : !---------------------------------------------------------------------
48 : ! ss' _
49 : ! H = \ (xhelp(s,i,na,l,m) conjg(yhelp(s',j,na,l,m')*rsoxy(na,l,s,s')
50 : ! *<slm|L*S|s'lm'>
51 : ! ij /_
52 : ! na,l,m,m'
53 : ! x,y = a,b
54 : !---------------------------------------------------------------------
55 : !
56 548 : nat_l = nat_stop - nat_start + 1 ! atoms processed by this pe
57 : !
58 : !---> update hamiltonian matrices: upper triangle
59 : !
60 548 : if (nat_l>0)THEN
61 1587 : DO i1 = 1,2
62 1058 : jsp = i1
63 1058 : IF (input%jspins.EQ.1) jsp = 1
64 3703 : DO j1 = 1,2
65 2116 : jsp1 = j1
66 2116 : IF (input%jspins.EQ.1) jsp1 = 1
67 : !$OMP PARALLEL DEFAULT(none)&
68 : !$OMP& PRIVATE(na,na_g,n,nn,l,ll1,m,lm,m1,lm1,ilo,i,lwn,ilop)&
69 : !$OMP& PRIVATE(c_a,c_b,c_c,c_1,c_2,c_3,c_4,c_5) &
70 : !$OMP& SHARED(hsomtx,i1,jsp,j1,jsp1,nsz,atoms)&
71 : !$OMP& SHARED(ahelp,bhelp,chelp,noco,nat_start,nat_stop,nat_l)&
72 3174 : !$OMP& SHARED(rsoc)
73 : ALLOCATE ( c_b((atoms%lmaxd+2)*atoms%lmaxd,nat_l),&
74 : c_a((atoms%lmaxd+2)*atoms%lmaxd,nat_l),&
75 : c_c(-atoms%llod:atoms%llod, atoms%nlod, nat_l) )
76 : !$OMP DO
77 : DO j = 1,nsz(jsp1)
78 : !
79 : ! prepare \sum_m' conjg( xhelp(m',l,na,j,jsp1) ) * soangl(l,m,i1,l,m',j1)
80 : !
81 : na = 0 ; na_g = 0
82 : DO n = 1,atoms%ntype
83 : DO nn = 1, atoms%neq(n)
84 : na_g = na_g + 1
85 : IF ((na_g.GE.nat_start).AND.(na_g.LE.nat_stop)) THEN
86 : na = na + 1
87 : !--> regular part
88 : DO l = 1,atoms%lmax(n)
89 : ll1 = l*(l+1)
90 : DO m = -l,l
91 : lm = ll1 + m
92 : c_a(lm,na) = CMPLX(0.,0.)
93 : c_b(lm,na) = CMPLX(0.,0.)
94 : DO m1 = -l,l
95 : lm1 = ll1 + m1
96 : c_a(lm,na) = c_a(lm,na) + rsoc%soangl(l,m,i1,l,m1,j1)&
97 : * CONJG(ahelp(lm1,na,j,jsp1))
98 : c_b(lm,na) = c_b(lm,na) + rsoc%soangl(l,m,i1,l,m1,j1)&
99 : * CONJG(bhelp(lm1,na,j,jsp1))
100 : ENDDO
101 : ENDDO
102 : ENDDO
103 : !--> LO contribution
104 : DO ilo = 1,atoms%nlo(n)
105 : l = atoms%llo(ilo,n)
106 : IF (l.GT.0) THEN
107 : DO m = -l,l
108 : c_c(m,ilo,na) = CMPLX(0.,0.)
109 : DO m1 = -l,l
110 : c_c(m,ilo,na) = c_c(m,ilo,na) + CONJG(&
111 : chelp(m1,j,ilo,na,jsp1))*rsoc%soangl(l,m,i1,l,m1,j1)
112 : ENDDO
113 : ENDDO
114 : ENDIF
115 : ENDDO
116 : ! end lo's
117 : ENDIF
118 : ENDDO ! nn
119 : ENDDO ! n
120 : !
121 : ! continue loop structure
122 : !
123 : DO i = 1,nsz(jsp)
124 : hsomtx(i,j,i1,j1) = CMPLX(0.,0.)
125 : na = 0 ; na_g = 0
126 : !
127 : !---> loop over each atom type
128 : !
129 : DO n = 1,atoms%ntype
130 : lwn = atoms%lmax(n)
131 : !
132 : !---> loop over equivalent atoms
133 : !
134 : DO nn = 1,atoms%neq(n)
135 : na_g = na_g + 1
136 : IF ((na_g.GE.nat_start).AND.(na_g.LE.nat_stop)) THEN
137 : na = na + 1
138 : DO l = 1,lwn
139 : ll1 = l*(l+1)
140 : DO m = -l,l
141 : lm = ll1 + m
142 : c_1 = rsoc%rsopp(n,l,i1,j1) * ahelp(lm,na,i,jsp) +&
143 : rsoc%rsopdp(n,l,i1,j1) * bhelp(lm,na,i,jsp)
144 : c_2 = rsoc%rsoppd(n,l,i1,j1) * ahelp(lm,na,i,jsp) +&
145 : rsoc%rsopdpd(n,l,i1,j1) * bhelp(lm,na,i,jsp)
146 : hsomtx(i,j,i1,j1) = hsomtx(i,j,i1,j1) +&
147 : c_1*c_a(lm,na) + c_2*c_b(lm,na)
148 : ENDDO
149 : !
150 : ENDDO
151 : !--> LO contribution
152 : DO ilo = 1,atoms%nlo(n)
153 : l = atoms%llo(ilo,n)
154 : ll1 = l*(l+1)
155 : IF (l.GT.0) THEN
156 : DO m = -l,l
157 : lm = ll1 + m
158 : c_3 = rsoc%rsopplo(n,ilo,i1,j1) *ahelp(lm,na,i,jsp) +&
159 : rsoc%rsopdplo(n,ilo,i1,j1) *bhelp(lm,na,i,jsp)
160 : c_4 = rsoc%rsoplop(n,ilo,i1,j1) *chelp(m,i,ilo,na,jsp)
161 : c_5 =rsoc%rsoplopd(n,ilo,i1,j1) *chelp(m,i,ilo,na,jsp)
162 : hsomtx(i,j,i1,j1) = hsomtx(i,j,i1,j1) + &
163 : c_4*c_a(lm,na) + c_5*c_b(lm,na) +&
164 : c_3*c_c(m,ilo,na)
165 : ENDDO
166 : DO ilop = 1,atoms%nlo(n)
167 : IF (atoms%llo(ilop,n).EQ.l) THEN
168 : DO m = -l,l
169 : hsomtx(i,j,i1,j1) = hsomtx(i,j,i1,j1) + &
170 : rsoc%rsoploplop(n,ilop,ilo,i1,j1) * &
171 : chelp(m,i,ilop,na,jsp) * c_c(m,ilo,na)
172 : ENDDO
173 : ENDIF
174 : ENDDO
175 : ENDIF
176 : ENDDO
177 : ! end lo's
178 : ENDIF
179 : ENDDO
180 : ENDDO ! atoms
181 : ENDDO
182 : !!i
183 : ENDDO
184 : !!j
185 : !$OMP END DO
186 : DEALLOCATE ( c_b,c_a,c_c)
187 : !$OMP END PARALLEL
188 : ENDDO
189 : !!jsp1
190 : ENDDO
191 : !!jsp
192 : !
193 : !---> update hamiltonian matrices: lower triangle
194 : !
195 27001 : DO i = 1,nsz(1)
196 1405817 : DO j = 1,nsz(input%jspins)
197 1405288 : hsomtx(j,i,2,1) = CONJG(hsomtx(i,j,1,2))
198 : ENDDO
199 : ENDDO
200 : else
201 398533 : hsomtx=0.0 !nothing to be done here
202 : endif
203 : #ifdef CPP_MPI
204 548 : CALL MPI_BARRIER(SUB_COMM,ierr(1))
205 548 : n = 4*nsz(1)*nsz(input%jspins)
206 1644 : ALLOCATE(c_buf(n))
207 548 : CALL MPI_REDUCE(hsomtx,c_buf,n,MPI_DOUBLE_COMPLEX,MPI_SUM,0,SUB_COMM,ierr(2))
208 548 : IF (n_rank.EQ.0) THEN
209 276 : CALL zcopy(n,c_buf,1,hsomtx,1)
210 : ENDIF
211 548 : DEALLOCATE(c_buf)
212 : #endif
213 : !
214 548 : RETURN
215 548 : END SUBROUTINE hsoham
216 : END MODULE m_hsoham
|