File: timestep.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (123 lines) | stat: -rw-r--r-- 3,727 bytes parent folder | download
1
2
3
4
5
6
7
8
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
61
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

! Copyright (C) 2014 K. Krieger, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine timestep
use modmain
use modtddft
use modmpi
use modomp
implicit none
! local variables
integer ik,is,ias
integer i,j,nthd
real(8) dt,t1
complex(8) z1
! allocatable arrays
real(8), allocatable :: vmt(:,:),vir(:),rfmt(:),w(:)
complex(8), allocatable :: evecsv(:,:),evectv(:,:),evecsvt(:,:)
complex(8), allocatable :: a(:,:),b(:,:),c(:,:)
if (itimes.ge.ntimes) then
  write(*,*)
  write(*,'("Error(timestep): itimes >= ntimes : ",2I8)') itimes,ntimes
  write(*,*)
  stop
end if
! convert muffin-tin Kohn-Sham potential to spherical coordinates
allocate(vmt(npcmtmax,natmtot))
call omp_hold(natmtot,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(rfmt,is) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ias=1,natmtot
  allocate(rfmt(npcmtmax))
  is=idxis(ias)
  call rfmtftoc(nrmt(is),nrmti(is),vsmt(:,ias),rfmt)
  call rbsht(nrcmt(is),nrcmti(is),rfmt,vmt(:,ias))
  deallocate(rfmt)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
! multiply interstitial potential by characteristic function
allocate(vir(ngtot))
vir(:)=vsir(:)*cfunir(:)
! backup existing time-dependent Kohn-Sham eigenvectors if required
call tdbackup
! loop over k-points
call omp_hold(nkpt/np_mpi,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(w,evecsv,evectv,evecsvt) &
!$OMP PRIVATE(a,b,c,i,j,dt,t1,z1) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ik=1,nkpt
! distribute among MPI processes
  if (mod(ik-1,np_mpi).ne.lp_mpi) cycle
  allocate(w(nstsv))
  allocate(evecsv(nstsv,nstsv),evectv(nstsv,nstsv),evecsvt(nstsv,nstsv))
  allocate(a(nstsv,nstsv),b(nstsv,nstsv),c(nstsv,nstsv))
! generate the Hamiltonian matrix in the ground-state second-variational basis
  call genhmlt(ik,vmt,vir,evectv)
! diagonalise the Hamiltonian to get third-variational eigenvectors
  if (spinpol.and.(.not.ncmag)) then
! collinear case requires block diagonalisation
    call eveqnz(nstfv,nstsv,evectv,w)
    i=nstfv+1
    call eveqnz(nstfv,nstsv,evectv(i,i),w(i))
    do i=1,nstfv
      do j=1,nstfv
        evectv(i,j+nstfv)=0.d0
        evectv(i+nstfv,j)=0.d0
      end do
    end do
  else
! non-collinear or spin-unpolarised: full diagonalisation
    call eveqnz(nstsv,nstsv,evectv,w)
  end if
! read in ground-state eigenvectors
  call getevecsv('.OUT',ik,vkl(:,ik),evecsv)
! convert third-variational eigenvectors to first-variational basis
  call zgemm('N','N',nstsv,nstsv,nstsv,zone,evecsv,nstsv,evectv,nstsv,zzero,a, &
   nstsv)
  deallocate(evecsv,evectv)
! time propagate instantaneous eigenvectors across one time step
  dt=times(itimes+1)-times(itimes)
  if (tdphi.eq.0.d0) then
! real time evolution
    do i=1,nstsv
      t1=-w(i)*dt
      z1=cmplx(cos(t1),sin(t1),8)
      b(:,i)=z1*a(:,i)
    end do
  else
! complex time evolution
    do i=1,nstsv
      t1=-w(i)*dt
      z1=t1*cmplx(sin(tdphi),cos(tdphi),8)
      z1=exp(z1)
      b(:,i)=z1*a(:,i)
    end do
  end if
! read in time-dependent Kohn-Sham eigenvectors (first-variational basis)
  call getevecsv(filext,ik,vkl(:,ik),evecsvt)
! apply time-evolution operator
  call zgemm('C','N',nstsv,nstsv,nstsv,zone,a,nstsv,evecsvt,nstsv,zzero,c,nstsv)
  call zgemm('N','N',nstsv,nstsv,nstsv,zone,b,nstsv,c,nstsv,zzero,evecsvt,nstsv)
! orthonormalise the eigenvectors if required
  if (tdphi.ne.0.d0) call orthevsv(evecsvt)
! write the new eigenvectors to file
  call putevecsv(filext,ik,evecsvt)
  deallocate(w,a,b,c,evecsvt)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
deallocate(vmt,vir)
! synchronise MPI processes
call mpi_barrier(mpicom,ierror)
return
end subroutine