File: dwavefmt.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 (102 lines) | stat: -rw-r--r-- 2,774 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
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102

! Copyright (C) 2013 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.

subroutine dwavefmt(lrstp,ias,ngp,ngpq,apwalmq,dapwalm,evecfv,devecfv,dwfmt)
use modmain
use modphonon
implicit none
! arguments
integer, intent(in) :: lrstp,ias,ngp,ngpq
complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw,natmtot)
!**** remove natmtot

complex(8), intent(in) :: dapwalm(ngkmax,apwordmax,lmmaxapw)
complex(8), intent(in) :: evecfv(nmatmax),devecfv(nmatmax)
real(8), intent(out) :: dwfmt(2,*)
! local variables
integer is,ldi,ldo,io,ilo
integer nrc,nrci,nrco,iro
integer l,m,lm,npc,npci,i
complex(8) z1
! external functions
complex(8) zdotu
external zdotu
is=idxis(ias)
ldi=2*lmmaxi
ldo=2*lmmaxo
iro=nrmti(is)+lrstp
if (lrstp.eq.1) then
  nrc=nrmt(is)
  nrci=nrmti(is)
  npc=npmt(is)
  npci=npmti(is)
else if (lrstp.eq.lradstp) then
  nrc=nrcmt(is)
  nrci=nrcmti(is)
  npc=npcmt(is)
  npci=npcmti(is)
else
  write(*,*)
  write(*,'("Error(dwavefmt): invalid lrstp : ",I8)') lrstp
  write(*,*)
  stop
end if
nrco=nrc-nrci
! zero the wavefunction derivative
dwfmt(:,1:npc)=0.d0
!-----------------------!
!     APW functions     !
!-----------------------!
lm=0
do l=0,lmaxo
  do m=-l,l
    lm=lm+1
    i=npci+lm
    do io=1,apword(l,is)
      z1=zdotu(ngpq,devecfv,1,apwalmq(:,io,lm,ias),1)
      if (ias.eq.iasph) then
        z1=z1+zdotu(ngp,evecfv,1,dapwalm(:,io,lm),1)
      end if
      if (abs(dble(z1)).gt.1.d-14) then
        if (l.le.lmaxi) then
          call daxpy(nrci,dble(z1),apwfr(1,1,io,l,ias),lrstp,dwfmt(1,lm),ldi)
        end if
        call daxpy(nrco,dble(z1),apwfr(iro,1,io,l,ias),lrstp,dwfmt(1,i),ldo)
      end if
      if (abs(aimag(z1)).gt.1.d-14) then
        if (l.le.lmaxi) then
          call daxpy(nrci,aimag(z1),apwfr(1,1,io,l,ias),lrstp,dwfmt(2,lm),ldi)
        end if
        call daxpy(nrco,aimag(z1),apwfr(iro,1,io,l,ias),lrstp,dwfmt(2,i),ldo)
      end if
    end do
  end do
end do
!---------------------------------!
!     local-orbital functions     !
!---------------------------------!
do ilo=1,nlorb(is)
  l=lorbl(ilo,is)
  do m=-l,l
    lm=idxlm(l,m)
    i=npci+lm
    z1=devecfv(ngpq+idxlo(lm,ilo,ias))
    if (abs(dble(z1)).gt.1.d-14) then
      if (l.le.lmaxi) then
        call daxpy(nrci,dble(z1),lofr(1,1,ilo,ias),lrstp,dwfmt(1,lm),ldi)
      end if
      call daxpy(nrco,dble(z1),lofr(iro,1,ilo,ias),lrstp,dwfmt(1,i),ldo)
    end if
    if (abs(aimag(z1)).gt.1.d-14) then
      if (l.le.lmaxi) then
        call daxpy(nrci,aimag(z1),lofr(1,1,ilo,ias),lrstp,dwfmt(2,lm),ldi)
      end if
      call daxpy(nrco,aimag(z1),lofr(iro,1,ilo,ias),lrstp,dwfmt(2,i),ldo)
    end if
  end do
end do
return
end subroutine