The Elk Code
 
Loading...
Searching...
No Matches
writexrmk.f90
Go to the documentation of this file.
1
2! Copyright (C) 2024 Peter Elliott, J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine writexrmk
7use modmain
8use modtddft
9use modmpi
10use modomp
11use moddelf
12implicit none
13! local variables
14integer ik,ist,jst,i
15integer lp,nthd
16real(8) xchg,t1,t2
17complex(8) z1
18! automatic arrays
19real(8) xrk(nkpt),xmk(ndmag,nkpt)
20complex(8) zv(nstsv)
21character(32) fext
22! allocatable arrays
23complex(8), allocatable :: evecsv(:,:),evecsvt(:,:)
24! external functions
25complex(8), external :: zdotc
26call holdthd(nkpt/np_mpi,nthd)
27!$OMP PARALLEL DEFAULT(SHARED) &
28!$OMP PRIVATE(evecsv,evecsvt,zv) &
29!$OMP PRIVATE(ist,jst,z1,i,t1,t2) &
30!$OMP NUM_THREADS(nthd)
31allocate(evecsv(nstsv,nstsv),evecsvt(nstsv,nstsv))
32!$OMP DO SCHEDULE(DYNAMIC)
33do ik=1,nkpt
34! distribute among MPI processes
35 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
36! read in ground-state eigenvectors
37 call getevecsv('.OUT',ik,vkl(:,ik),evecsv)
38! read in time-dependent Kohn-Sham eigenvectors (first-variational basis)
39 call getevecsv(filext,ik,vkl(:,ik),evecsvt)
40! zero the excited density and magnetisation
41 xrk(ik)=0.d0
42 if (spinpol) xmk(1:ndmag,ik)=0.d0
43 do ist=1,nstsv
44 if (evalsv(ist,ik) > efermi) cycle
45! project time-dependent Kohn-Sham state onto ground-state conduction bands
46 zv(1:nstsv)=0.d0
47 do jst=1,nstsv
48 if (evalsv(jst,ik) < efermi) cycle
49 z1=zdotc(nstsv,evecsv(:,jst),1,evecsvt(:,ist),1)
50 zv(1:nstsv)=zv(1:nstsv)+z1*evecsv(1:nstsv,jst)
51 end do
52! calculate the excited density and magnetisation
53 if (spinpol) then
54 i=nstfv+1
55 t1=dble(zdotc(nstfv,zv,1,zv,1))
56 t2=dble(zdotc(nstfv,zv(i),1,zv(i),1))
57 xrk(ik)=xrk(ik)+t1+t2
58 xmk(ndmag,ik)=xmk(ndmag,ik)+t1-t2
59 if (ncmag) then
60 z1=zdotc(nstfv,zv,1,zv(i),1)
61 xmk(1,ik)=xmk(1,ik)+2.d0*z1%re
62 xmk(2,ik)=xmk(2,ik)+2.d0*z1%im
63 end if
64 else
65 xrk(ik)=xrk(ik)+occmax*dble(zdotc(nstsv,zv,1,zv,1))
66 end if
67 end do
68end do
69!$OMP END DO
70deallocate(evecsv,evecsvt)
71!$OMP END PARALLEL
72call freethd(nthd)
73! broadcast arrays to every MPI process
74if (np_mpi > 1) then
75 do ik=1,nkpt
76 lp=mod(ik-1,np_mpi)
77 call mpi_bcast(xrk(ik),1,mpi_double_precision,lp,mpicom,ierror)
78 call mpi_bcast(xmk(:,ik),3,mpi_double_precision,lp,mpicom,ierror)
79 end do
80end if
81! compute the excited charge
82xchg=sum(wkpt(1:nkpt)*xrk(1:nkpt))
83if (mp_mpi) then
84! file extension
85 write(fext,'("_TS",I8.8,".OUT")') itimes
86! write k-point dependent excited density to file
87 open(50,file='XRHOK'//trim(fext),form='FORMATTED',action='WRITE')
88 do ik=1,nkpt
89 write(50,'(I6,4G18.10)') ik,vkl(:,ik),xrk(ik)
90 end do
91 close(50)
92! write k-point dependent excited magnetisation to file
93 if (spinpol) then
94 open(50,file='XMAGK'//trim(fext),form='FORMATTED',action='WRITE')
95 do ik=1,nkpt
96 write(50,'(I6,6G18.10)') ik,vkl(:,ik),xmk(:,ik)
97 end do
98 close(50)
99 end if
100! write the excited charge to file
101 if (itimes <= 1) call delfile('XCHARGE_TD.OUT')
102 open(50,file='XCHARGE_TD.OUT',form='FORMATTED',position='APPEND')
103 write(50,'(2G18.10)') times(itimes),xchg
104 close(50)
105end if
106! synchronise MPI processes
107call mpi_barrier(mpicom,ierror)
108end subroutine
109
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine delfile(fname)
Definition moddelf.f90:15
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
logical ncmag
Definition modmain.f90:240
real(8) efermi
Definition modmain.f90:904
character(256) filext
Definition modmain.f90:1300
logical spinpol
Definition modmain.f90:228
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8) occmax
Definition modmain.f90:898
integer nstfv
Definition modmain.f90:884
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
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
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
real(8), dimension(:), allocatable times
Definition modtddft.f90:48
integer itimes
Definition modtddft.f90:46
subroutine writexrmk
Definition writexrmk.f90:7