The Elk Code
 
Loading...
Searching...
No Matches
tdinit.f90
Go to the documentation of this file.
1
2! Copyright (C) 2022 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
6subroutine tdinit
7use modmain
8use modtddft
9use moddftu
10use modmpi
11use modomp
12use modramdisk
13implicit none
14! local variables
15integer ik
16! allocatable arrays
17real(8), allocatable :: evalfv(:,:)
18complex(8), allocatable :: evecfv(:,:,:),evecsv(:,:)
19! flag for starting at t=0 or a restart
20tdt0=((task == 460).or.(task == 462))
21! determine the static density and charge and write to file
22if (tdt0) call rhostatic
23! store original parameters
27! currents should be calculated with forces
28if (tforce) tjr=.true.
29! average force can be non-zero (allow for translation of atomic basis)
30tfav0=.false.
31! ensure eigenvectors are written to disk during initialisation
33wrtdsk=.true.
34! initialise global variables
35call init0
36call init1
37! read the charge density and potentials from file
38call readstate
39call genvsig
40call gencore
41call energykncr
42call linengy
43call genapwlofr
44call gensocfr
45if (tdt0) then
46! generate eigenvalues and eigenvectors only at t=0 (not for the restart) for
47! the k-point set reduced with the symmetries which leave A(t) invariant for all
48! time steps
49 call genevfsv
50else
51! only read in the second-variational eigenvalues for restarts
52 call readevalsv
53! read in the static density and charge
54 call readrhos
55end if
56! compute the occupation numbers
57call occupy
58! DFT+U
59if (dftu /= 0) then
60 call gendmatmt
61 call genvmatmt
62 call vmatmtsc
63 if (tmwrite) call genwkpr0
64end if
65! generate the kinetic matrix elements in the second-variational basis
66call genkmat(.false.,.false.)
67! write the momentum matrix elements in the second-variational basis
68call genpmat
69! copy EVALFV.OUT, EVECFV.OUT, OCCSV.OUT and EVECSV.OUT to _TD.OUT extension
70if (tdt0) then
71 allocate(evalfv(nstfv,nspnfv),evecfv(nmatmax,nstfv,nspnfv))
72 allocate(evecsv(nstsv,nstsv))
73 do ik=1,nkpt
74! distribute among MPI processes
75 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
76 call getevalfv('.OUT',ik,vkl(:,ik),evalfv)
77 call putevalfv('_TD.OUT',ik,evalfv)
78 call getevecfv('.OUT',ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
79 call putevecfv('_TD.OUT',ik,evecfv)
80 call putoccsv('_TD.OUT',ik,occsv(:,ik))
81 call getevecsv('.OUT',ik,vkl(:,ik),evecsv)
82! randomise eigenvectors at t=0 if required
83 call rndevsv(rndevt0,evecsv)
84 call putevecsv('_TD.OUT',ik,evecsv)
85 end do
86 deallocate(evalfv,evecfv,evecsv)
87end if
88! set global file extension
89filext='_TD.OUT'
90! output the new k-point set to file
91if (mp_mpi) call writekpts
92! synchronise MPI processes
93call mpi_barrier(mpicom,ierror)
94! Ehrenfest dynamics
95if ((task == 462).or.(task == 463)) then
96! forces should not be calculated
97 tforce=.false.
98! enable small amplitude displacements
99 tatdisp=.true.
100! zero the displacements and velocities
101 atdvc(:,:,:,:)=0.d0
102! generate the gradient of the nucleus and static density Coulomb potential
103 call gengvnsmt
104end if
105if (tdt0) then
106! start from t=0
107 itimes0=1
108else
109! restart if required
110 call tdrestart
111end if
112! deallocate the static density if not required
113if (.not.tjr) deallocate(rhosmt,rhosir)
114! read the forces calculated during the previous TDDFT run
115if (tatdisp) call readforcet
116! read the atomic displacements and velocities
117if (trdatdv) call readatdvc
119! make the time-dependent Kohn-Sham orbitals strictly orthogonal after each
120! time step for complex time evolution
121if (tdphi /= 0.d0) ntsorth=1
122
123end subroutine
124
subroutine energykncr
Definition energykncr.f90:7
subroutine genapwlofr
Definition genapwlofr.f90:7
subroutine gencore
Definition gencore.f90:10
subroutine gendmatmt
Definition gendmatmt.f90:7
subroutine genevfsv
Definition genevfsv.f90:7
subroutine gengvnsmt
Definition gengvnsmt.f90:7
subroutine genkmat(tfv, tvclcr)
Definition genkmat.f90:10
subroutine genpmat
Definition genpmat.f90:7
subroutine gensocfr
Definition gensocfr.f90:7
subroutine genvmatmt
Definition genvmatmt.f90:10
subroutine genvsig
Definition genvsig.f90:10
subroutine genwkpr0
Definition genwkpr0.f90:7
subroutine getevalfv(fext, ikp, vpl, evalfv)
Definition getevalfv.f90:7
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition getevecfv.f90:10
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine linengy
Definition linengy.f90:10
logical tmwrite
Definition moddftu.f90:66
integer dftu
Definition moddftu.f90:32
logical tjr0
Definition modmain.f90:620
logical tfav00
Definition modmain.f90:1000
character(256) filext
Definition modmain.f90:1300
logical tjr
Definition modmain.f90:620
logical trdatdv
Definition modmain.f90:62
integer nkpt
Definition modmain.f90:461
integer nspnfv
Definition modmain.f90:289
logical tfav0
Definition modmain.f90:1000
real(8), dimension(:,:,:,:), allocatable vgkl
Definition modmain.f90:503
integer nmatmax
Definition modmain.f90:848
integer task
Definition modmain.f90:1298
integer nstsv
Definition modmain.f90:886
logical tatdisp
Definition modmain.f90:59
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
logical tforce
Definition modmain.f90:986
integer nstfv
Definition modmain.f90:884
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
real(8), dimension(3, 0:1, maxatoms, maxspecies) atdvc
Definition modmain.f90:64
logical tforce0
Definition modmain.f90:986
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
logical wrtdsk0
logical wrtdsk
real(8) tdphi
Definition modtddft.f90:52
real(8), dimension(:,:,:), allocatable rhosmt
Definition modtddft.f90:82
integer ntsorth
Definition modtddft.f90:105
logical tdt0
Definition modtddft.f90:50
real(8), dimension(:,:), allocatable rhosir
Definition modtddft.f90:82
integer itimes0
Definition modtddft.f90:44
real(8) rndevt0
Definition modtddft.f90:96
subroutine occupy
Definition occupy.f90:10
subroutine putevalfv(fext, ik, evalfv)
Definition putevalfv.f90:7
subroutine putevecfv(fext, ik, evecfv)
Definition putevecfv.f90:7
subroutine putevecsv(fext, ik, evecsv)
Definition putevecsv.f90:7
subroutine putoccsv(fext, ik, occsvp)
Definition putoccsv.f90:7
subroutine readatdvc
Definition readatdvc.f90:7
subroutine readevalsv
Definition readevalsv.f90:7
subroutine readforcet
Definition readforcet.f90:7
subroutine readrhos
Definition readrhos.f90:7
subroutine readstate
Definition readstate.f90:10
subroutine rhostatic
Definition rhostatic.f90:7
subroutine rndevsv(rndm, evecsv)
Definition rndevsv.f90:7
subroutine tdinit
Definition tdinit.f90:7
subroutine tdrestart
Definition tdrestart.f90:7
subroutine vmatmtsc
Definition vmatmtsc.f90:7
subroutine writekpts
Definition writekpts.f90:10