mapatom.F90 9.19 KB
Newer Older
1 2 3 4 5 6 7
      MODULE m_mapatom
      use m_juDFT
!*******************************************************************
!     determines the group operation which maps the representive
!     atom into its equivalent atoms     c.l.fu
!*******************************************************************
      CONTAINS
8
      SUBROUTINE mapatom(sym,atoms,cell,input,noco)
9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60
!
!     if (l_f) setup multab,invtab,invarop,invarind for force_a12 & 21
!***********************************************************************
! the contribution to the hamiltonian and to the overlap matrix in a
! system with inversion symmetry from one muffin tin is the complex
! conjugate of the contribution from the "invers" muffin tin. this fact
! can be exploited to save cpu-time in hssphn. Thus, it is nessessary to
! know whether an atom can be mapped onto an equivalent atom via 3d
! inversion. where both atoms have to belong to the same unit cell, i.e.
! the are not related to each other by a lattice translation. therefore,
! an array invsatom is set up.
! invsatom(natom) =
! 0 if the atom cannot be mapped onto an eqivalent atom via inversion
! 1 if the atom can be mapped onto an eqivalent atom via inversion, and
!   has a smaller atom index than the related atom
! 2 if the atom can be mapped onto an eqivalent atom via inversion, and
!   has a bigger atom index than the related atom
! p.kurz aug. 1996
!***********************************************************************
!
      USE m_socsym
      USE m_types
      IMPLICIT NONE
      TYPE(t_sym),INTENT(INOUT)   :: sym
      TYPE(t_atoms),INTENT(INOUT) :: atoms
      TYPE(t_cell),INTENT(IN)  :: cell
      TYPE(t_input),INTENT(IN) :: input
      TYPE(t_noco),INTENT(IN)  :: noco

!     .. Local Scalars ..
      REAL s3,norm
      INTEGER i,icount,j,j1,j2,j3,jop,n,na,nat1,nat2,nb,na_r
      INTEGER k,ij,n1,n2,ix,iy,iz,na2
      REAL, PARAMETER :: del = 1.0e-4
!     ..
!     .. Local Arrays ..
      INTEGER mt(3,3),mp(3,3)
      REAL aamat(3,3),sum_tau_lat(3),sum_taual(3)
      REAL gam(3),gaminv(3),gamr(3),sr(3),ttau(3)
      LOGICAL error(sym%nop)
!     ..
!     .. Intrinsic Functions ..
      INTRINSIC real,sqrt
!     ..
      
    !  CALL dotset(&
    ! &            cell,&
    ! &            aamat, )
      aamat=matmul(transpose(cell%amat),cell%amat)
    
      IF (noco%l_soc) THEN  ! check once more here...
        CALL soc_sym(&
61
     &               sym%nop,sym%mrot,noco%theta,noco%phi,cell%amat,&
62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211
     &               error)
      ELSE
        error(:) = .false.
      ENDIF
                               
      WRITE (6,FMT=8000)
 8000 FORMAT (/,/,5x,'group operations on equivalent atoms:')
      nat1 = 1
      DO n = 1,atoms%ntype
         nat2 = nat1 + atoms%neq(n) - 1
         atoms%ngopr(nat1) = 1
!+gu
         na_r = nat1
         DO na = nat1,nat2
            IF (atoms%ntypsy(na).NE.atoms%ntypsy(na_r)) na_r = na
!-gu
            DO i = 1,3
               gam(i) = atoms%taual(i,na)
            END DO
            sym%invarind(na) = 0
            icount = 0
            DO  jop = 1,sym%nop
               DO i = 1,3
                  gamr(i) = 0.
                  DO j = 1,3
                     gamr(i) = gamr(i) + sym%mrot(i,j,jop)*gam(j)
                  END DO
                  gamr(i) = gamr(i) + sym%tau(i,jop)
               END DO
               DO i = 1,3
                  gaminv(i) = gamr(i) - atoms%taual(i,na)
                  gamr(i)   = gamr(i) - atoms%taual(i,nat1) ! cf local_sym
               END DO
               IF (icount.EQ.0) THEN
                  DO j3 = -2,2
                     sr(3) = gamr(3) + real(j3)
                     DO j2 = -2,2
                        sr(2) = gamr(2) + real(j2)
                        DO j1 = -2,2
                           sr(1) = gamr(1) + real(j1)
                           s3 = sqrt(dot_product(matmul(sr,aamat),sr))
                           IF ((s3.LT.del).AND.(.not.error(jop))) THEN
                              icount = icount + 1
                              atoms%ngopr(na) = jop
                           END IF
                        END DO
                     END DO
                  END DO
               END IF
!
! search for operations which leave taual invariant
!
               IF (input%l_f.OR.(atoms%n_u.GT.0)) THEN 
                  DO j3 = -2,2
                     sr(3) = gaminv(3) + real(j3)
                     DO j2 = -2,2
                        sr(2) = gaminv(2) + real(j2)
                        DO j1 = -2,2
                           sr(1) = gaminv(1) + real(j1)
                           s3 = sqrt(dot_product(matmul(sr,aamat),sr))
                           IF (s3.LT.del) THEN
                              sym%invarind(na) = sym%invarind(na) + 1
                              sym%invarop(na,sym%invarind(na)) = jop
                           END IF
                        END DO
                     END DO
                  END DO
               ENDIF
!
! end of operations
          ENDDO
            IF (icount.LE.0) THEN
             write(6,*) "Mapping failed for atom:",nat1
             write(6,*) "No of symmetries tested:",sym%nop
             CALL juDFT_error("mapatom",calledby="mapatom")
           ENDIF
            WRITE (6,FMT=8010) nat1,na,atoms%ngopr(na)
 8010       FORMAT (5x,'atom',i3,' can be mapped into atom',i3,&
     &             ' through group  operation',i4)
!
! end of equivalent atoms
       ENDDO
!
         nat1 = nat1 + atoms%neq(n)
!
! end of different types of atoms
    ENDDO

!------------------------- FORCE PART -------------------------------
!+gu this is the remainder of spgset necessary for force calculations
!
      IF (input%l_f.OR.(atoms%n_u.GT.0)) THEN

      WRITE (6,FMT=&
     &  '(//,"list of operations which leave taual invariant",/)')
      DO na = 1,nat2
         WRITE (6,FMT='("atom nr.",i3,3x,(t14,"ops are:",24i3))') na,&
     &     (sym%invarop(na,nb),nb=1,sym%invarind(na))
      END DO

      ENDIF
!------------------------- FORCE PART ENDS --------------------------
!
!     check closure  ; note that:  {R|t} tau = R^{-1} tau -  R^{-1} t
!
!--->    loop over all operations
!
      WRITE (6,FMT=8040)
 8040 FORMAT (/,/,' multiplication table',/,/)
      sym%multab = 0
      DO j=1,sym%nop

!--->    multiply {R_j|t_j}{R_i|t_i}
         DO i=1,sym%nop
            mp = matmul( sym%mrot(:,:,j) , sym%mrot(:,:,i) )
            ttau = sym%tau(:,j) + matmul( sym%mrot(:,:,j) , sym%tau(:,i) )
            ttau = ttau - anint( ttau - 1.e-7 )

!--->    determine which operation this is
            DO k=1,sym%nop
              IF( all( mp(:,:) == sym%mrot(:,:,k) ) .AND.&
     &            ALL( abs( ttau(:)-sym%tau(:,k) ) < 1.e-7 ) ) THEN
                 IF (sym%multab(j,i) .EQ. 0 ) THEN
                    sym%multab(j,i) = k
                    IF (k .EQ. 1) sym%invtab(j)=i
                 ELSE
                    WRITE(6,'(" Symmetry error: multiple ops")')
                     CALL juDFT_error("Multiple ops",calledby="mapatom")
                 ENDIF
              ENDIF
            ENDDO

            IF (sym%multab(j,i).EQ.0) THEN
               WRITE (6,'(" Group not closed")')
               WRITE (6,'("  j , i =",2i4)') j,i
               CALL juDFT_error("mapatom: group not closed",calledby&
     &              ="mapatom")
            ENDIF
         ENDDO
      ENDDO

      DO n1 = 1,sym%nop
         WRITE (6,FMT=8060) (sym%multab(n1,n2),n2=1,sym%nop)
      END DO
 8060 FORMAT (1x,48i3)
      WRITE (6,FMT='(//," inverse operations",//)')
      DO n1 = 1,sym%nop
         WRITE (6,FMT=8060) n1,sym%invtab(n1)
      END DO

Daniel Wortmann's avatar
Daniel Wortmann committed
212
      DO na = 1,atoms%nat
213 214 215
         atoms%invsat(na) = 0
         sym%invsatnr(na) = 0
      END DO
216

217
      IF (.not.(noco%l_soc.and.atoms%n_u>0)) THEN
218 219 220 221 222 223
      IF (sym%invs) THEN
         WRITE (6,FMT=*)
         nat1 = 1
         DO n = 1,atoms%ntype
            nat2 = nat1 + atoms%neq(n) - 1
            DO na = nat1,nat2 - 1
224
               IF (atoms%invsat(na).EQ.0.AND..NOT.noco%l_noco) THEN
225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254
                  naloop:DO na2 = na + 1,nat2
                     DO i = 1,3
                        sum_taual(i) = atoms%taual(i,na) + atoms%taual(i,na2)
                     END DO
                     DO ix = -2,2
                       sum_tau_lat(1) = sum_taual(1) + real(ix)
                       DO iy = -2,2
                         sum_tau_lat(2) = sum_taual(2) + real(iy)
                         DO iz = -2,2
                           sum_tau_lat(3) = sum_taual(3) + real(iz)
                           norm = sqrt(dot_product(matmul(sum_tau_lat,aamat),sum_tau_lat))
                           IF (norm.LT.del) THEN
                              atoms%invsat(na) = 1
                              atoms%invsat(na2) = 2
                              sym%invsatnr(na)  = na2
                              sym%invsatnr(na2) = na
                              WRITE (6,FMT=9000) n,na,na2
                              cycle naloop
                           END IF
                        END DO
                      END DO
                    END DO
               END DO naloop
               END IF
            END DO
            nat1 = nat1 + atoms%neq(n)
         END DO
      WRITE (6,FMT=*) atoms%invsat
 9000 FORMAT ('atom type',i3,': atom',i3,' can be mapped into atom',i3,&
     &       ' via 3d inversion')
255
      END IF
256
      END IF
257
 
258 259
      END  SUBROUTINE mapatom
      END  MODULE m_mapatom