The Elk Code
atpstep.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: atpstep
8 ! !INTERFACE:
9 subroutine atpstep
10 ! !USES:
11 use modmain
12 use modmpi
13 ! !DESCRIPTION:
14 ! Makes a geometry optimisation step and updates the current atomic positions
15 ! according to the force on each atom. If ${\bf r}_{ij}^m$ is the position and
16 ! ${\bf F}_{ij}^m$ is the force acting on it for atom $j$ of species $i$ and
17 ! after time step $m$, then the new position is calculated by
18 ! $$ {\bf r}_{ij}^{m+1}={\bf r}_{ij}^m+\tau_{ij}^m\left({\bf F}_{ij}^m
19 ! +{\bf F}_{ij}^{m-1}\right), $$
20 ! where $\tau_{ij}^m$ is a parameter governing the size of the displacement.
21 ! If ${\bf F}_{ij}^m\cdot{\bf F}_{ij}^{m-1}>0$ then $\tau_{ij}^m$ is
22 ! increased, otherwise it is decreased.
23 !
24 ! !REVISION HISTORY:
25 ! Created June 2003 (JKD)
26 !EOP
27 !BOC
28 implicit none
29 ! local variables
30 integer is,ia,ias,n
31 real(8) t1
32 do is=1,nspecies
33  do ia=1,natoms(is)
34  ias=idxas(ia,is)
35 ! compute the dot-product between the current and previous total force
36  t1=dot_product(forcetot(1:3,ias),forcetotp(1:3,ias))
37 ! if the force is in the same direction then increase step size parameter
38  if (t1 > 0.d0) then
39  tauatp(ias)=tauatp(ias)+tau0atp
40  else
41  tauatp(ias)=tau0atp
42  end if
43 ! make atomic position step
44  atposc(1:3,ia,is)=atposc(1:3,ia,is)+tauatp(ias)*(forcetot(1:3,ias) &
45  +forcetotp(1:3,ias))
46  end do
47 end do
48 ! each MPI process should have identical atomic positions
50 call mpi_bcast(atposc,n,mpi_double_precision,0,mpicom,ierror)
51 do is=1,nspecies
52  do ia=1,natoms(is)
53 ! compute the lattice coordinates of the atomic positions
54  call r3mv(ainv,atposc(:,ia,is),atposl(:,ia,is))
55  end do
56 end do
57 end subroutine
58 !EOC
59 
integer, parameter maxspecies
Definition: modmain.f90:30
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
real(8), dimension(3, 3) ainv
Definition: modmain.f90:14
real(8), dimension(:,:), allocatable forcetot
Definition: modmain.f90:993
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
real(8), dimension(:), allocatable tauatp
Definition: modmain.f90:1013
integer, parameter maxatoms
Definition: modmain.f90:32
real(8) tau0atp
Definition: modmain.f90:1011
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
subroutine atpstep
Definition: atpstep.f90:10
real(8), dimension(:,:), allocatable forcetotp
Definition: modmain.f90:995
Definition: modmpi.f90:6
integer nspecies
Definition: modmain.f90:34
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19