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_dipole_takehome
8 : use m_judft
9 : contains
10 0 : subroutine wann_dipole_takehome(
11 : > jspins,l_absolute,
12 : > amat,bmat,omtil,
13 : < electronic_moment)
14 : c***********************************************************
15 : c-----Check if the polarization may be reduced by adding
16 : c-----contributions due to electrons shifted by primitive
17 : c-----lattice translations. This is the solution to the
18 : c-----problem of Wannier functions being pulled away from
19 : c-----their host atom to a mirror host atom during the
20 : c-----Wannierization process.
21 : c Frank Freimuth
22 : c***********************************************************
23 : use m_constants, only: pimach
24 : implicit none
25 :
26 : integer, intent(in) :: jspins
27 : logical, intent(in) :: l_absolute
28 : real, intent(in) :: amat(3,3),bmat(3,3),omtil
29 : real, intent(inout) :: electronic_moment(3,2)
30 :
31 : real,parameter :: elemchargmu=1.60217646e-13
32 : real,parameter :: bohrtocm=0.529177e-8
33 : real :: electronic_polarization(3,2)
34 : integer :: jspin,k
35 : real :: int_electronic_moment(3,2)
36 : real :: tpi
37 :
38 0 : call timestart("wann_dipole_takehome")
39 0 : tpi=2.0*pimach()
40 :
41 0 : write(666,*)"*********************************"
42 0 : write(666,*)" Try to minimize electronic term "
43 0 : write(666,*)"*********************************"
44 0 : write(*,*) "*********************************"
45 0 : write(*,*) " Try to minimize electronic term "
46 0 : write(*,*) "*********************************"
47 :
48 0 : do jspin=1,jspins
49 : int_electronic_moment(:,jspin)=
50 0 : = matmul( bmat,electronic_moment(:,jspin) )/tpi
51 0 : do k=1,3
52 : int_electronic_moment(k,jspin)=
53 : = int_electronic_moment(k,jspin) -
54 0 : - nint( int_electronic_moment(k,jspin) )
55 : enddo
56 : electronic_moment(:,jspin) =
57 0 : = matmul(amat,int_electronic_moment(:,jspin))
58 0 : write(*, fmt=555)jspin,electronic_moment(:,jspin)
59 0 : write(666,fmt=555)jspin,electronic_moment(:,jspin)
60 0 : if(.not.l_absolute)then
61 : electronic_polarization=
62 0 : & electronic_moment/omtil*elemchargmu/((bohrtocm)**2)
63 0 : write(*, fmt=777)jspin,electronic_polarization(:,jspin)
64 0 : write(666,fmt=777)jspin,electronic_polarization(:,jspin)
65 : endif
66 : enddo
67 :
68 : 555 format("spin ",i1.1," electronic moment = (",3f12.6,") a.u.")
69 : 777 format("spin ",i1.1,
70 : & " electronic polarization = (",3f12.6,") uC/cm**2")
71 :
72 0 : call timestop("wann_dipole_takehome")
73 0 : end subroutine wann_dipole_takehome
74 : end module m_wann_dipole_takehome
75 :
|