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_abcrot2
8 : PRIVATE
9 : PUBLIC :: abcrot2
10 : CONTAINS
11 0 : SUBROUTINE abcrot2(itype,na,atoms,banddos,eigVecCoeffs,jsp,acof,bcof,ccof)
12 : USE m_dwigner
13 : USE m_types
14 : IMPLICIT NONE
15 : INTEGER,INTENT(IN) :: itype,na
16 : TYPE(t_atoms),INTENT(IN) :: atoms
17 : TYPE(t_banddos),INTENT(IN) :: banddos
18 : TYPE(t_eigVecCoeffs),INTENT(IN) :: eigVecCoeffs
19 : COMPLEX, ALLOCATABLE,INTENT(INOUT) :: acof(:,:)
20 : COMPLEX, ALLOCATABLE,INTENT(INOUT) :: bcof(:,:)
21 : COMPLEX, ALLOCATABLE,INTENT(INOUT) :: ccof(:,:,:)
22 : ! ..
23 : ! .. Scalar Arguments ..
24 : INTEGER, INTENT (IN) :: jsp
25 : ! ..
26 : ! .. Local Scalars ..
27 : INTEGER ilo,i,l
28 : REAL amx(3,3,1),imx(3,3)
29 0 : COMPLEX d_wgn(-atoms%lmaxd:atoms%lmaxd,-atoms%lmaxd:atoms%lmaxd,1:atoms%lmaxd,1)
30 :
31 :
32 0 : CALL euler(banddos%alpha(na),banddos%beta(na),banddos%gamma(na), amx)
33 :
34 0 : imx(:,:) = 0. ; imx(1,1) = 1. ; imx(2,2) = 1. ; imx(3,3) = 1.
35 :
36 0 : CALL d_wigner(1,amx,imx,atoms%lmaxd, d_wgn)
37 :
38 0 : DO l = 1, atoms%lmax(itype)
39 0 : DO i = 1, size(acof,1)
40 : acof(i,l**2:l*(l+2)) = MATMUL(CONJG(d_wgn(-l:l,-l:l,l,1)),&
41 0 : eigVecCoeffs%abcof(i,l**2:l*(l+2),0,na,jsp))
42 : bcof(i,l**2:l*(l+2)) = MATMUL(CONJG(d_wgn(-l:l,-l:l,l,1)),&
43 0 : eigVecCoeffs%abcof(i,l**2:l*(l+2),1,na,jsp))
44 : ENDDO
45 : ENDDO
46 0 : DO ilo = 1, atoms%nlo(itype)
47 0 : l = atoms%llo(ilo,itype)
48 0 : IF (l.GT.0) THEN
49 0 : DO i = 1 ,size(acof,1)
50 : ccof(-l:l,i,ilo) = MATMUL(CONJG(d_wgn(-l:l,-l:l,l,1)),&
51 0 : eigVecCoeffs%ccof(-l:l,i,ilo,na,jsp))
52 : ENDDO
53 : ENDIF
54 : ENDDO
55 0 : END SUBROUTINE abcrot2
56 :
57 : !********************************************************************
58 : !********************************************************************
59 0 : SUBROUTINE euler(alpha,beta,gamma,amx)
60 : IMPLICIT NONE
61 :
62 : REAL, INTENT (IN) :: alpha,beta,gamma
63 : REAL, INTENT (OUT) :: amx(3,3,1)
64 :
65 : REAL alph,bet,gamm
66 : REAL bmx(3,3),cmx(3,3),dmx(3,3),hmx(3,3)
67 : INTEGER nwf,i,j,ii
68 :
69 : !..define the D,C,B-matrices
70 0 : amx(:,:,:)=0.
71 :
72 0 : alph = alpha ; bet = beta ; gamm = gamma
73 :
74 0 : dmx(1,1) = COS(alph) ; dmx(1,2) = SIN(alph) ; dmx(1,3) = 0.
75 0 : dmx(2,1) =-SIN(alph) ; dmx(2,2) = COS(alph) ; dmx(2,3) = 0.
76 0 : dmx(3,1) = 0. ; dmx(3,2) = 0. ; dmx(3,3) = 1.
77 :
78 0 : cmx(1,1) = 1. ; cmx(1,2) = 0. ; cmx(1,3) = 0.
79 0 : cmx(2,1) = 0. ; cmx(2,2) = COS(bet) ; cmx(2,3) = SIN(bet)
80 0 : cmx(3,1) = 0. ; cmx(3,2) =-SIN(bet) ; cmx(3,3) = COS(bet)
81 :
82 0 : bmx(1,1) = COS(gamm) ; bmx(1,2) = SIN(gamm) ; bmx(1,3) = 0.
83 0 : bmx(2,1) =-SIN(gamm) ; bmx(2,2) = COS(gamm) ; bmx(2,3) = 0.
84 0 : bmx(3,1) = 0. ; bmx(3,2) = 0. ; bmx(3,3) = 1.
85 :
86 0 : hmx(:,:) = 0.
87 0 : DO i = 1,3
88 0 : DO j = 1,3
89 0 : DO ii = 1,3
90 0 : hmx(i,j) = hmx(i,j) + cmx(i,ii)*dmx(ii,j)
91 : ENDDO
92 : ENDDO
93 : ENDDO
94 :
95 0 : DO i = 1,3
96 0 : DO j = 1,3
97 0 : DO ii = 1,3
98 0 : amx(i,j,1) = amx(i,j,1) + bmx(i,ii)*hmx(ii,j)
99 : ENDDO
100 : ENDDO
101 : ENDDO
102 :
103 0 : END SUBROUTINE euler
104 :
105 : END MODULE m_abcrot2
|