The Elk Code
timestep.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 K. Krieger, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine timestep
7 use modmain
8 use modtddft
9 use modmpi
10 use modomp
11 implicit none
12 ! local variables
13 integer ik,i,nthd
14 real(8) ca,dt,t1
15 complex(8) z1,z2
16 ! automatic arrays
17 real(8) w(nstsv)
18 ! allocatable arrays
19 real(8), allocatable :: vmt(:,:),vir(:),bmt(:,:,:),bir(:,:)
20 complex(8), allocatable :: evecsv(:,:),evectv(:,:),evecsvt(:,:)
21 complex(8), allocatable :: kmat(:,:),pmat(:,:,:),a(:,:),b(:,:),c(:,:)
22 if (itimes >= ntimes) then
23  write(*,*)
24  write(*,'("Error(timestep): itimes >= ntimes :",2(X,I0))') itimes,ntimes
25  write(*,*)
26  stop
27 end if
28 allocate(vmt(npcmtmax,natmtot),vir(ngtc))
29 if (spinpol) then
30  allocate(bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag))
31 else
32  allocate(bmt(1,1,1),bir(1,1))
33 end if
34 ! generate the Kohn-Sham potential and magnetic field in spherical coordinates
35 ! and multiply by the radial integration weights; also multiply the interstitial
36 ! potential with the characteristic function
37 call vblocal(vmt,vir,bmt,bir)
38 ! time step length
39 dt=times(itimes+1)-times(itimes)
40 ! zero the kinetic energy
41 engykn=0.d0
42 ! zero the total current
43 jtot(1:3)=0.d0
44 ! loop over k-points
45 call holdthd(nkpt/np_mpi,nthd)
46 !$OMP PARALLEL DEFAULT(SHARED) &
47 !$OMP PRIVATE(evecsv,evectv,evecsvt) &
48 !$OMP PRIVATE(kmat,pmat,w,a,b,c) &
49 !$OMP PRIVATE(i,t1,z1,z2) &
50 !$OMP NUM_THREADS(nthd)
51 allocate(evecsv(nstsv,nstsv),evectv(nstsv,nstsv),evecsvt(nstsv,nstsv))
52 allocate(kmat(nstsv,nstsv),pmat(nstsv,nstsv,3))
53 allocate(a(nstsv,nstsv),b(nstsv,nstsv),c(nstsv,nstsv))
54 !$OMP DO SCHEDULE(DYNAMIC)
55 do ik=1,nkpt
56 ! distribute among MPI processes
57  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
58 ! get the kinetic matrix elements from file
59  call getkmat(ik,kmat)
60 ! get the momentum matrix elements from file
61  call getpmat(vkl(:,ik),pmat)
62 ! generate the Hamiltonian matrix in the ground-state second-variational basis
63  call genhmlt(ik,vmt,vir,bmt,bir,kmat,pmat,evectv)
64 ! diagonalise the Hamiltonian to get third-variational eigenvectors
65  if (spinpol.and.(.not.ncmag)) then
66 ! collinear case requires block diagonalisation
67  call zheevdi(nstfv,nstsv,evectv,w)
68  evectv(nstfv+1:nstsv,1:nstfv)=0.d0
69  evectv(1:nstfv,nstfv+1:nstsv)=0.d0
70  i=nstfv+1
71  call zheevdi(nstfv,nstsv,evectv(i,i),w(i))
72  else
73 ! non-collinear or spin-unpolarised: full diagonalisation
74  call zheevdi(nstsv,nstsv,evectv,w)
75  end if
76 ! read in ground-state eigenvectors
77  call getevecsv('.OUT',ik,vkl(:,ik),evecsv)
78 ! convert third-variational eigenvectors to first-variational basis
79  call zgemm('N','N',nstsv,nstsv,nstsv,zone,evecsv,nstsv,evectv,nstsv,zzero,a, &
80  nstsv)
81 ! time propagate instantaneous eigenvectors across one time step
82  if (tdphi == 0.d0) then
83 ! real time evolution
84  do i=1,nstsv
85  t1=-w(i)*dt
86  z1=cmplx(cos(t1),sin(t1),8)
87  b(1:nstsv,i)=z1*a(1:nstsv,i)
88  end do
89  else
90 ! complex time evolution
91  z2=cmplx(sin(tdphi),cos(tdphi),8)
92  do i=1,nstsv
93  t1=-w(i)*dt
94  z1=exp(t1*z2)
95  b(1:nstsv,i)=z1*a(1:nstsv,i)
96  end do
97  end if
98 ! read in time-dependent Kohn-Sham eigenvectors (first-variational basis)
99  call getevecsv(filext,ik,vkl(:,ik),evecsvt)
100 ! apply time evolution operator
101  call zgemm('C','N',nstsv,nstsv,nstsv,zone,a,nstsv,evecsvt,nstsv,zzero,c,nstsv)
102  call zgemm('N','N',nstsv,nstsv,nstsv,zone,b,nstsv,c,nstsv,zzero,evecsvt,nstsv)
103 ! orthonormalise the eigenvectors if required
104  if (ntsorth > 0) then
105  if (mod(itimes-1,ntsorth) == 0) call unitary(nstsv,evecsvt)
106  end if
107 ! add to the kinetic energy
108  call engyknk(ik,kmat,evecsv,evecsvt)
109 ! add to the total current
110  call jtotk(ik,pmat,evecsv,evecsvt)
111 ! write the new eigenvectors to file
112  call putevecsv(filext,ik,evecsvt)
113 end do
114 !$OMP END DO
115 deallocate(evecsv,evectv,evecsvt,kmat,pmat,a,b,c)
116 !$OMP END PARALLEL
117 call freethd(nthd)
118 deallocate(vmt,vir,bmt,bir)
119 ! add the kinetic energy and total current from each process and redistribute
120 if (np_mpi > 1) then
121  call mpi_allreduce(mpi_in_place,engykn,1,mpi_double_precision,mpi_sum,mpicom,&
122  ierror)
123  call mpi_allreduce(mpi_in_place,jtot,3,mpi_double_precision,mpi_sum,mpicom, &
124  ierror)
125 end if
126 ! add the core kinetic energy
128 ! coupling constant of the external A-field (-1/c)
129 ca=-1.d0/solsc
130 ! add the diamagnetic current to total
131 jtot(1:3)=jtot(1:3)+ca*afieldt(1:3,itimes)*(chgtot-chgstot(1:3))
132 ! symmetrise the vector
133 call symvec(jtot)
134 ! total current magnitude
135 jtotm=norm2(jtot(1:3))
136 ! write the time step to file
137 if (mp_mpi) call writetimes
138 ! backup existing time-dependent Kohn-Sham eigenvectors if required
139 if (ntsbackup > 0) then
140  if (mod(itimes-1,ntsbackup) == 0) call tdbackup
141 end if
142 ! synchronise MPI processes
143 call mpi_barrier(mpicom,ierror)
144 end subroutine
145 
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
integer npcmtmax
Definition: modmain.f90:216
character(256) filext
Definition: modmain.f90:1288
subroutine timestep
Definition: timestep.f90:7
subroutine writetimes
Definition: writetimes.f90:7
logical mp_mpi
Definition: modmpi.f90:17
integer ntsorth
Definition: modtddft.f90:105
integer ngtc
Definition: modmain.f90:392
real(8), dimension(3) jtot
Definition: modmain.f90:746
logical spinpol
Definition: modmain.f90:228
integer nkpt
Definition: modmain.f90:461
integer ntimes
Definition: modtddft.f90:42
integer ndmag
Definition: modmain.f90:238
real(8) engykn
Definition: modmain.f90:942
subroutine getpmat(vpl, pmat)
Definition: getpmat.f90:7
Definition: modomp.f90:6
complex(8), parameter zone
Definition: modmain.f90:1227
subroutine getkmat(ik, kmat)
Definition: getkmat.f90:7
integer np_mpi
Definition: modmpi.f90:13
subroutine jtotk(ik, pmat, evecsv, evecsvt)
Definition: jtotk.f90:7
subroutine zheevdi(n, ld, a, w)
Definition: zheevdi.f90:7
subroutine engyknk(ik, kmat, evecsv, evecsvt)
Definition: engyknk.f90:7
real(8) jtotm
Definition: modmain.f90:746
subroutine symvec(vc)
Definition: symvec.f90:7
real(8) engykncr
Definition: modmain.f90:944
real(8) solsc
Definition: modmain.f90:1240
complex(8), parameter zzero
Definition: modmain.f90:1227
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
subroutine tdbackup
Definition: tdbackup.f90:7
Definition: modmpi.f90:6
real(8), dimension(:,:), allocatable afieldt
Definition: modtddft.f90:58
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
integer itimes
Definition: modtddft.f90:46
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
real(8) chgtot
Definition: modmain.f90:724
subroutine putevecsv(fext, ik, evecsv)
Definition: putevecsv.f90:7
subroutine genhmlt(ik, vmt, vir, bmt, bir, kmat, pmat, h)
Definition: genhmlt.f90:7
real(8), dimension(3) chgstot
Definition: modtddft.f90:84
integer natmtot
Definition: modmain.f90:40
integer ntsbackup
Definition: modtddft.f90:54
subroutine unitary(n, a)
Definition: unitary.f90:10
logical ncmag
Definition: modmain.f90:240
real(8) tdphi
Definition: modtddft.f90:52
subroutine vblocal(vmt, vir, bmt, bir)
Definition: vblocal.f90:7
integer nstfv
Definition: modmain.f90:878
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19