The Elk Code
 
Loading...
Searching...
No Matches
writetddos.f90
Go to the documentation of this file.
1
2! Copyright (C) 2015 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!BOP
7! !ROUTINE: writetddos
8! !INTERFACE:
9subroutine writetddos
10! !USES:
11use modmain
12use modtddft
13use modmpi
14! !DESCRIPTION:
15! Calculates the time-dependent density of states (DOS). This is defined as
16! $$ {\rm DOS}(\omega,t)=\frac{\Omega}{(2\pi)^3}\int d^3k \sum_i
17! \delta(\varepsilon_{i{\bf k}}-\omega) F_{i{\bf k}}(t), $$
18! where
19! $$ F_{i{\bf k}}(t)=\sum_j f_{j{\bf k}}|\langle\varphi_{i{\bf k}}|
20! \phi_{j{\bf k}}(t)\rangle|^2, $$
21! with occupation numbers $f_{j{\bf k}}$, ground-state orbitals
22! $\varphi_{i{\bf k}}$ and time-dependent orbitals $\phi_{j{\bf k}}(t)$.
23!
24! !REVISION HISTORY:
25! Created April 2015 (JKD)
26!EOP
27!BOC
28implicit none
29! local variables
30integer ik,ist,jst,lp
31real(8) sm,t1
32complex(8) z1
33character(256) fext
34! allocatable arrays
35real(8), allocatable :: occsvp(:,:)
36complex(8), allocatable :: evecsv(:,:),evecsvt(:,:)
37! external functions
38complex(8), external :: zdotc
39! file extension
40write(fext,'("_TS",I8.8,".OUT")') itimes
41allocate(occsvp(nstsv,nkpt))
42allocate(evecsv(nstsv,nstsv),evecsvt(nstsv,nstsv))
43do ik=1,nkpt
44! distribute among MPI processes
45 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
46! read in ground-state eigenvectors
47 call getevecsv('.OUT',ik,vkl(:,ik),evecsv)
48! read in the time evolving eigenvectors
49 call getevecsv('_TD.OUT',ik,vkl(:,ik),evecsvt)
50! determine the time-dependent projected occupation numbers
51 do ist=1,nstsv
52 sm=0.d0
53 do jst=1,nstsv
54 t1=occsv(jst,ik)
55 if (abs(t1) < epsocc) cycle
56 z1=zdotc(nstsv,evecsv(:,ist),1,evecsvt(:,jst),1)
57 sm=sm+t1*(z1%re**2+z1%im**2)
58 end do
59 occsvp(ist,ik)=sm
60 end do
61! write projected occupation numbers to file
62 call putoccsv('P'//trim(fext),ik,occsvp(:,ik))
63end do
64deallocate(evecsv,evecsvt)
65! broadcast projected occupation numbers to every MPI process
66if (np_mpi > 1) then
67 do ik=1,nkpt
68 lp=mod(ik-1,np_mpi)
69 call mpi_bcast(occsvp(:,ik),nstsv,mpi_double_precision,lp,mpicom,ierror)
70 end do
71end if
72if (mp_mpi) then
73! compute the effective electronic temperature and write to file
74 call tdtemp(occsvp)
75! write the DOS to file
76 call dos(fext,.true.,occsvp)
77end if
78deallocate(occsvp)
79! synchronise MPI processes
80call mpi_barrier(mpicom,ierror)
81end subroutine
82!EOC
83
subroutine dos(fext, tocc, occsvp)
Definition dos.f90:10
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
integer nkpt
Definition modmain.f90:461
integer nstsv
Definition modmain.f90:886
real(8) epsocc
Definition modmain.f90:900
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
logical mp_mpi
Definition modmpi.f90:17
integer itimes
Definition modtddft.f90:46
subroutine putoccsv(fext, ik, occsvp)
Definition putoccsv.f90:7
subroutine tdtemp(occsvp)
Definition tdtemp.f90:7
subroutine writetddos