File: rdmvaryc.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (60 lines) | stat: -rw-r--r-- 1,448 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60

! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: rdmvaryc
! !INTERFACE:
subroutine rdmvaryc
! !USES:
use modmain
use modrdm
! !DESCRIPTION:
!   Calculates new {\tt evecsv} from old by using the derivatives of the total
!   energy w.r.t. {\tt evecsv}. A single step of steepest-descent is made.
!
! !REVISION HISTORY:
!   Created 2009 (Sharma)
!EOP
!BOC
implicit none
! local variables
integer ik,ist,jst
real(8) t1
complex(8) z1
! allocatable arrays
complex(8), allocatable :: dedc(:,:,:),evecsv(:,:),x(:)
! external functions
real(8) dznrm2
complex(8) zdotc
external dznrm2,zdotc
! compute the derivative w.r.t. evecsv
allocate(dedc(nstsv,nstsv,nkpt))
call rdmdedc(dedc)
allocate(evecsv(nstsv,nstsv),x(nstsv))
do ik=1,nkpt
! get the eigenvectors from file
  call getevecsv(filext,ik,vkl(:,ik),evecsv)
! calculate new evecsv
  evecsv(:,:)=evecsv(:,:)-taurdmc*dedc(:,:,ik)
! othogonalise evecsv (Gram-Schmidt)
  do ist=1,nstsv
    x(:)=evecsv(:,ist)
    do jst=1,ist-1
      z1=zdotc(nstsv,evecsv(:,jst),1,evecsv(:,ist),1)
      x(:)=x(:)-z1*evecsv(:,jst)
    end do
    t1=dznrm2(nstsv,x,1)
    t1=1.d0/t1
    evecsv(:,ist)=t1*x(:)
  end do
! write new evecsv to file
  call putevecsv(filext,ik,evecsv)
! end loop over k-points
end do
deallocate(dedc,evecsv,x)
return
end subroutine
!EOC