The Elk Code
moldyn.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2023 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 moldyn
7 use modmain
8 use modtddft
9 use modmpi
10 use modomp
11 use moddelf
12 implicit none
13 ! local variables
14 integer is,ia
15 ! initialise universal variables
16 call init0
17 ! store original parameters
18 atposl0(:,:,:)=atposl(:,:,:)
19 atposc0(:,:,:)=atposc(:,:,:)
23 ! no shifting of atomic basis allowed
24 tshift=.false.
25 ! calculate atomic forces
26 tforce=.true.
27 ! average force can be non-zero (allow for translation of atomic basis)
28 tfav0=.false.
29 ! generate the time step grid
30 call gentimes
31 ! flag for starting at t=0 or a restart
32 tdt0=(task == 420)
33 if (tdt0) then
34 ! start from t=0
35  itimes0=1
36 else
37 ! restart if required
38  call readtimes
40  trdatdv=.true.
41 end if
42 if (trdatdv) then
43 ! read the atomic displacements and velocities
44  call readatdvc
45 else
46 ! set the displacements and velocities to zero
47  atdvc(:,:,:,:)=0.d0
48 end if
49 trdstate=.false.
50 if (tdt0.and.mp_mpi) then
51  call delfile('TOTENERGY_TD.OUT')
52  if (spinpol) then
53  call delfile('MOMENT_TD.OUT')
54  call delfile('MOMENTM_TD.OUT')
55  call delfile('MOMENTMT_TD.OUT')
56  call delfile('MOMENTIR_TD.OUT')
57  end if
58  call delfile('FORCETOT_TD.OUT')
59  call delfile('FORCEMAX_TD.OUT')
60  call delfile('ATDISPL_TD.OUT')
61  call delfile('ATDISPC_TD.OUT')
62 end if
64  if (mp_mpi) then
65  write(*,'("Info(moldyn): time step ",I8," of ",I8,", t = ",G18.10)') &
67  end if
68 ! reset the OpenMP thread variables
69  call omp_reset
70 ! add the displacements to the atomic positions
71  do is=1,nspecies
72  do ia=1,natoms(is)
73  atposc(:,ia,is)=atposc0(:,ia,is)+atdvc(:,0,ia,is)
74  call r3mv(ainv,atposc(:,ia,is),atposl(:,ia,is))
75  end do
76  end do
77 ! calculate the ground-state and atomic forces
78  call gndstate
79 ! subsequent calculations will read in the potential from STATE.OUT
80  trdstate=.true.
81 ! time step the atomic positions within the adiabatic approximation
82  call atptstep(forcetot)
83  if (mp_mpi) then
84 ! write the time step to file
85  call writetimes
86 ! write time-dependent total energy
87  call writetdengy
88 ! write spin moments if required
89  if (spinpol) call writemomtd
90 ! write time-dependent atomic forces
91  call writetdforces
92 ! write the time-dependent atomic displacements
93  call writeatdisp
94 ! write the XCrysden animation file crystal.axsf
95  call writeaxsf
96  end if
97  if (tstop) exit
98 end do
99 ! restore original input parameters
100 atposl(:,:,:)=atposl0(:,:,:)
104 ! synchronise MPI processes
105 call mpi_barrier(mpicom,ierror)
106 end subroutine
107 
subroutine gndstate
Definition: gndstate.f90:10
real(8), dimension(3, maxatoms, maxspecies) atposl0
Definition: modmain.f90:51
subroutine writetimes
Definition: writetimes.f90:7
integer task
Definition: modmain.f90:1299
logical mp_mpi
Definition: modmpi.f90:17
subroutine writetdengy
Definition: writetdengy.f90:7
logical tshift0
Definition: modmain.f90:352
subroutine writetdforces
logical spinpol
Definition: modmain.f90:228
logical tshift
Definition: modmain.f90:352
integer ntimes
Definition: modtddft.f90:42
logical tfav0
Definition: modmain.f90:1003
logical tforce0
Definition: modmain.f90:989
real(8), dimension(3, 3) ainv
Definition: modmain.f90:14
real(8), dimension(3, 0:1, maxatoms, maxspecies) atdvc
Definition: modmain.f90:64
Definition: modomp.f90:6
logical tstop
Definition: modmain.f90:1055
subroutine moldyn
Definition: moldyn.f90:7
subroutine readtimes
Definition: readtimes.f90:7
real(8), dimension(:,:), allocatable forcetot
Definition: modmain.f90:993
subroutine writemomtd
Definition: writemomtd.f90:7
logical tforce
Definition: modmain.f90:989
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
integer ntsforce
Definition: modtddft.f90:98
subroutine gentimes
Definition: gentimes.f90:7
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
subroutine readatdvc
Definition: readatdvc.f90:7
subroutine delfile(fname)
Definition: moddelf.f90:15
Definition: modmpi.f90:6
logical tfav00
Definition: modmain.f90:1003
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
integer itimes
Definition: modtddft.f90:46
subroutine omp_reset
Definition: modomp.f90:71
integer nspecies
Definition: modmain.f90:34
integer itimes0
Definition: modtddft.f90:44
logical trdstate
Definition: modmain.f90:682
subroutine init0
Definition: init0.f90:10
logical tdt0
Definition: modtddft.f90:50
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
subroutine writeaxsf
Definition: writeaxsf.f90:7
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54
subroutine atptstep(ft)
Definition: atptstep.f90:7
real(8), dimension(3, maxatoms, maxspecies) atposc0
Definition: modmain.f90:54
integer mpicom
Definition: modmpi.f90:11
subroutine writeatdisp
Definition: writeatdisp.f90:7
logical trdatdv
Definition: modmain.f90:62
integer ierror
Definition: modmpi.f90:19