The Elk Code
tddft.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 tddft
7 use modmain
8 use modtddft
9 use moddftu
10 use modmpi
11 use modomp
12 use modramdisk
13 use modtest
14 implicit none
15 if (tshift) then
16  write(*,*)
17  write(*,'("Error(tddft): use tshift = .false. for the ground-state run")')
18  write(*,*)
19  stop
20 end if
21 ! initialise TDDFT variables
22 call tdinit
23 ! set the stop signal to .false.
24 tstop=.false.
25 !---------------------------------!
26 ! main loop over time steps !
27 !---------------------------------!
28 if (mp_mpi) write(*,*)
29 ! synchronise MPI processes
30 call mpi_barrier(mpicom,ierror)
32  if (mp_mpi) then
33  write(*,'("Info(tddft): time step ",I8," of ",I8,", t = ",G18.10)') &
35  end if
36 ! reset the OpenMP thread variables
37  call omp_reset
38 ! check for STOP file
39  call checkstop
40 ! write all files on last loop
41  if ((itimes == ntimes-1).or.tstop) wrtdsk=.true.
42 ! flag for writing observables at this time step
43  ttswrite=.false.
44  if (ntswrite(1) > 0) then
45  if (mod(itimes-1,ntswrite(1)) == 0) then
46  if ((itimes == 1).or.(itimes >= ntswrite(2))) ttswrite=.true.
47  end if
48  end if
49 ! flag for calculating forces at this time step
50  ttsforce=(mod(itimes-1,ntsforce) == 0)
51 ! evolve the wavefunctions across a single time step
52  call timestep
53 ! generate the density and magnetisation at current time step
54  call rhomag
55 ! compute the gauge-invariant current j(r) if required
56  if (tjr) call genjr
57 ! time step the induced A-field
58  if (tafindt) call afindtstep
59 ! calculate the electric field
60  call genefieldt
61 ! compute the time-dependent Kohn-Sham potential and magnetic field
62  call potkst
63 ! add the fixed spin moment effective field if required
64  call addbfsm
65 ! DFT+U
66  if (dftu /= 0) then
67  call gendmatmt
68  call genvmatmt
69  call vmatmtsc
70  end if
71 ! compute the total energy
72  call energytd
73 ! calculate the atomic forces if required
74  if (tforce.and.ttsforce) call force
75 ! time step the atomic positions for Ehrenfest dynamics using forces calculated
76 ! during the previous TDDFT run
77  if (tatdisp.and.ttsforce) call atptstep(forcet(:,:,itimes))
78 ! write general TDDFT output
79  if (mp_mpi) call writetddft
80 ! write optional TDDFT output
81  if (ttswrite) then
82 ! write time-dependent DOS
83  if (tddos) call writetddos
84 ! write muffin-tin L, S and J if required
85  if (tdlsj) call writetdlsj
86 ! write the k-point dependent total current
87  if (tdjtk) call writetdjtk
88 ! write the k-point dependent excited density and magnetisation
89  if (tdxrmk) call writexrmk
90  end if
91  if (tstop) exit
92 end do
93 filext='.OUT'
94 ! restore original input parameters
97 tjr=tjr0
98 tatdisp=.false.
99 ! write the total current of the last step to test file
100 call writetest(460,'total current of last time step',nv=3,tol=5.d-4,rva=jtot)
101 ! synchronise MPI processes
102 call mpi_barrier(mpicom,ierror)
103 end subroutine
104 
subroutine writetdjtk
Definition: writetdjtk.f90:7
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
logical ttsforce
Definition: modtddft.f90:100
real(8), dimension(:,:,:), allocatable forcet
Definition: modtddft.f90:102
logical tjr
Definition: modmain.f90:620
subroutine energytd
Definition: energytd.f90:7
character(256) filext
Definition: modmain.f90:1301
subroutine timestep
Definition: timestep.f90:7
logical tatdisp
Definition: modmain.f90:59
subroutine writetddos
Definition: writetddos.f90:10
logical mp_mpi
Definition: modmpi.f90:17
logical tjr0
Definition: modmain.f90:620
real(8), dimension(3) jtot
Definition: modmain.f90:748
logical tshift
Definition: modmain.f90:352
logical tdlsj
Definition: modtddft.f90:94
integer ntimes
Definition: modtddft.f90:42
logical tfav0
Definition: modmain.f90:1003
logical tforce0
Definition: modmain.f90:989
Definition: modomp.f90:6
logical tstop
Definition: modmain.f90:1055
subroutine genjr
Definition: genjr.f90:7
subroutine tddft
Definition: tddft.f90:7
logical tafindt
Definition: modtddft.f90:72
subroutine gendmatmt
Definition: gendmatmt.f90:7
subroutine rhomag
Definition: rhomag.f90:7
subroutine checkstop
Definition: checkstop.f90:7
logical tdxrmk
Definition: modtddft.f90:94
subroutine writexrmk
Definition: writexrmk.f90:7
logical tforce
Definition: modmain.f90:989
logical tdjtk
Definition: modtddft.f90:94
logical wrtdsk
Definition: modramdisk.f90:15
subroutine vmatmtsc
Definition: vmatmtsc.f90:7
integer ntsforce
Definition: modtddft.f90:98
subroutine tdinit
Definition: tdinit.f90:7
integer, dimension(2) ntswrite
Definition: modtddft.f90:78
subroutine writetdlsj
Definition: writetdlsj.f90:7
Definition: modmpi.f90:6
logical tfav00
Definition: modmain.f90:1003
integer dftu
Definition: moddftu.f90:32
subroutine force
Definition: force.f90:10
subroutine potkst
Definition: potkst.f90:7
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
integer itimes
Definition: modtddft.f90:46
subroutine omp_reset
Definition: modomp.f90:71
integer itimes0
Definition: modtddft.f90:44
subroutine writetddft
Definition: writetddft.f90:7
logical ttswrite
Definition: modtddft.f90:80
logical tddos
Definition: modtddft.f90:94
subroutine addbfsm
Definition: addbfsm.f90:7
subroutine genefieldt
Definition: genefieldt.f90:7
subroutine atptstep(ft)
Definition: atptstep.f90:7
integer mpicom
Definition: modmpi.f90:11
subroutine genvmatmt
Definition: genvmatmt.f90:10
subroutine afindtstep
Definition: afindtstep.f90:10
integer ierror
Definition: modmpi.f90:19