The Elk Code
Loading...
Searching...
No Matches
atptstep.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2020 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
6
subroutine
atptstep
(ft)
7
use
modmain
8
use
modtddft
9
use
modmpi
10
implicit none
11
! arguments
12
real
(8),
intent(in)
:: ft(3,natmtot)
13
! local variables
14
integer
itn,is,ia,ias
15
real
(8) dt,t1,t2
16
! next time step when the forces will be calculated
17
itn=
itimes
+
ntsforce
18
if
(itn >
ntimes
)
return
19
! time step length
20
dt=
times
(itn)-
times
(
itimes
)
21
do
is=1,
nspecies
22
t1=1.d0-
atdfc
*dt
23
if
(t1 < 0.d0) t1=0.d0
24
t2=dt/
spmass
(is)
25
do
ia=1,
natoms
(is)
26
ias=
idxas
(ia,is)
27
! add to the atomic velocities
28
atdvc
(1:3,1,ia,is)=t1*
atdvc
(1:3,1,ia,is)+t2*ft(1:3,ias)
29
! add to the atomic displacements
30
atdvc
(1:3,0,ia,is)=
atdvc
(1:3,0,ia,is)+
atdvc
(1:3,1,ia,is)*dt
31
end do
32
end do
33
! write the atomic displacements and velocities to file
34
if
(
mp_mpi
)
call
writeatdvc
35
end subroutine
36
atptstep
subroutine atptstep(ft)
Definition
atptstep.f90:7
modmain
Definition
modmain.f90:6
modmain::atdfc
real(8) atdfc
Definition
modmain.f90:66
modmain::natoms
integer, dimension(maxspecies) natoms
Definition
modmain.f90:36
modmain::idxas
integer, dimension(maxatoms, maxspecies) idxas
Definition
modmain.f90:42
modmain::nspecies
integer nspecies
Definition
modmain.f90:34
modmain::atdvc
real(8), dimension(3, 0:1, maxatoms, maxspecies) atdvc
Definition
modmain.f90:64
modmain::spmass
real(8), dimension(maxspecies) spmass
Definition
modmain.f90:101
modmpi
Definition
modmpi.f90:6
modmpi::mp_mpi
logical mp_mpi
Definition
modmpi.f90:17
modtddft
Definition
modtddft.f90:6
modtddft::times
real(8), dimension(:), allocatable times
Definition
modtddft.f90:48
modtddft::ntsforce
integer ntsforce
Definition
modtddft.f90:98
modtddft::itimes
integer itimes
Definition
modtddft.f90:46
modtddft::ntimes
integer ntimes
Definition
modtddft.f90:42
writeatdvc
subroutine writeatdvc
Definition
writeatdvc.f90:7
atptstep.f90
Generated by
1.9.8