The Elk Code
tdbackup.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2018 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
tdbackup
7
use
modmain
8
use
modtddft
9
use
modmpi
10
use
modramdisk
11
implicit none
12
! local variables
13
integer
ik
14
! allocatable arrays
15
complex(8)
,
allocatable
:: evecsvt(:,:)
16
! ensure eigenvectors are written to disk
17
wrtdisk0
=
wrtdisk
18
wrtdisk
=.true.
19
allocate
(evecsvt(
nstsv
,
nstsv
))
20
do
ik=1,
nkpt
21
! distribute among MPI processes
22
if
(mod(ik-1,
np_mpi
) /=
lp_mpi
) cycle
23
! read in time-dependent Kohn-Sham eigenvectors
24
call
getevecsv
(
filext
,ik,
vkl
(:,ik),evecsvt)
25
! write eigenvectors to backup file
26
call
putevecsv
(
'_TD_BACKUP.OUT'
,ik,evecsvt)
27
end do
28
deallocate
(evecsvt)
29
! write the time step backup file
30
if
(
mp_mpi
)
then
31
open
(50,
file
=
'TIMESTEP_BACKUP.OUT'
,form=
'FORMATTED'
)
32
write
(50,
'(I8,G18.10)'
)
itimes
,
times
(
itimes
)
33
close
(50)
34
end if
35
wrtdisk
=
wrtdisk0
36
! synchronise MPI processes
37
call
mpi_barrier(
mpicom
,
ierror
)
38
end subroutine
39
getevecsv
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition:
getevecsv.f90:7
modmain::filext
character(256) filext
Definition:
modmain.f90:1288
modmpi::mp_mpi
logical mp_mpi
Definition:
modmpi.f90:17
modmain::nkpt
integer nkpt
Definition:
modmain.f90:461
modramdisk::wrtdisk
logical wrtdisk
Definition:
modramdisk.f90:15
modramdisk::wrtdisk0
logical wrtdisk0
Definition:
modramdisk.f90:15
modramdisk
Definition:
modramdisk.f90:6
modramdisk::file
type(file_t), dimension(:), allocatable, private file
Definition:
modramdisk.f90:29
modtddft
Definition:
modtddft.f90:6
modmpi::np_mpi
integer np_mpi
Definition:
modmpi.f90:13
modmain
Definition:
modmain.f90:6
modmain::nstsv
integer nstsv
Definition:
modmain.f90:880
modmain::vkl
real(8), dimension(:,:), allocatable vkl
Definition:
modmain.f90:471
tdbackup
subroutine tdbackup
Definition:
tdbackup.f90:7
modmpi
Definition:
modmpi.f90:6
modtddft::times
real(8), dimension(:), allocatable times
Definition:
modtddft.f90:48
modtddft::itimes
integer itimes
Definition:
modtddft.f90:46
modmpi::lp_mpi
integer lp_mpi
Definition:
modmpi.f90:15
putevecsv
subroutine putevecsv(fext, ik, evecsv)
Definition:
putevecsv.f90:7
modmpi::mpicom
integer mpicom
Definition:
modmpi.f90:11
modmpi::ierror
integer ierror
Definition:
modmpi.f90:19
tdbackup.f90
Generated by
1.8.14