File: potclulr.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 (91 lines) | stat: -rw-r--r-- 2,537 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

! Copyright (C) 2018 T. Mueller, 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 potclulr
use modmain
use modulr
use modomp
implicit none
! local variables
integer iq,ifq,is,ias
integer nr,nri,ir
integer npc,i,nthd
! allocatable arrays
complex(8), allocatable :: zrhomt(:,:),zvclmt(:,:),zfmt(:)
call omp_hold(nqpt,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(zrhomt,zvclmt,zfmt) &
!$OMP PRIVATE(ifq,ias,is,nr,nri) &
!$OMP PRIVATE(i,ir,npc) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do iq=1,nqpt
  allocate(zrhomt(npmtmax,natmtot),zvclmt(npmtmax,natmtot))
  allocate(zfmt(npcmtmax))
  ifq=iqfft(iq)
! convert the complex density from a coarse to fine radial mesh
  do ias=1,natmtot
    is=idxis(ias)
    call zcopy(npcmt(is),rhoqmt(:,ias,ifq),1,zrhomt(:,ias),1)
  end do
  call zfmtctof(zrhomt)
! add the core density for Q=0
  if (iq.eq.1) then
    do ias=1,natmtot
      is=idxis(ias)
      nr=nrmt(is)
      nri=nrmti(is)
      i=1
      do ir=1,nri
        zrhomt(i,ias)=zrhomt(i,ias)+rhocr(ir,ias,1)
        i=i+lmmaxi
      end do
      do ir=nri+1,nr
        zrhomt(i,ias)=zrhomt(i,ias)+rhocr(ir,ias,1)
        i=i+lmmaxo
      end do
    end do
  end if
! solve the complex Poisson's equation in the muffin-tins
  call genzvclmt(nrmt,nrmti,nrspmax,rsp,npmtmax,zrhomt,zvclmt)
! add the nuclear monopole potentials for Q=0
  if (iq.eq.1) then
    do ias=1,natmtot
      is=idxis(ias)
      nr=nrmt(is)
      nri=nrmti(is)
      i=1
      do ir=1,nri
        zvclmt(i,ias)=zvclmt(i,ias)+vcln(ir,is)
        i=i+lmmaxi
      end do
      do ir=nri+1,nr
        zvclmt(i,ias)=zvclmt(i,ias)+vcln(ir,is)
        i=i+lmmaxo
      end do
    end do
  end if
! solve Poisson's equation in the entire unit cell
  call zpotcoul(nrmt,nrmti,npmt,npmti,nrspmax,rsp,ngridg,igfft,ngvec, &
   gqc(:,iq),gclgq(:,iq),ngvec,jlgqrmt(:,:,:,iq),ylmgq(:,:,iq),sfacgq(:,:,iq), &
   rhoqir(:,ifq),npmtmax,zvclmt,vsqir(:,ifq))
  do ias=1,natmtot
    is=idxis(ias)
    npc=npcmt(is)
! convert from fine to coarse radial mesh
    call zfmtftoc(nrmt(is),nrmti(is),zvclmt(:,ias),zfmt)
! convert to spherical coordinates
    call zbsht(nrcmt(is),nrcmti(is),zfmt,vsqmt(:,ias,ifq))
! multiply by the phase factor function exp(-iQ.r)
    vsqmt(1:npc,ias,ifq)=vsqmt(1:npc,ias,ifq)*conjg(expqmt(1:npc,ias,iq))
  end do
  deallocate(zrhomt,zvclmt,zfmt)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
return
end subroutine