mod_overlap Module

Calculates the overlap integral of test function PHI with regular or irregular wavefunction (PZ, QZ: spherical part of the large component, PQNS: nonspherical part of the large component) for LDA+U. The overlap is given out in array RESULT

In the inner region the non-spherical wavefunctions are approximated as:

  • The regular one (ir < irmin = irws-irns) : where is the regular wavefunction of the spherically symmetric part of the potential and the alpha matrix (see sub. regns)

  • The irregular one (ir < irmin) : where is the regular and is the irregular wavefunction of the spherically symmetric part of the potential and , the matrices calculated at the point (see sub. irwns)

The integrand is convoluted with the shape functions: Finally, the result should then be transformed to complex spherical harmonics basis!

Note

Here only large component is used, but was PHI normalised according to large + small component? (qldau)

Warning

last index of pns,qns is not the spin index but sra index.


Uses

    • mod_datatypes

Subroutines

public subroutine overlap(result, phi, pz, qz, pqns, acr, dr, lirreg, ipan, ircut, drdi, irmin, lphi, ipand, lmaxd, lmmaxd, mmaxd, lmpotd, irmind, irmd)

Author
Ph. Mavropoulos
License
Creative Commons License
Category
lda+u, KKRhost

C a l c u l a t e s

t h e

o v e r l a p

i n t e g r a l

o f

t e s t

f u n c t i o n

P H I

w i t h

r e g u l a r

o r

i r r e g u l a r

w a v e f u n c t i o n

Read more…

Arguments

Type IntentOptional Attributes Name
complex(kind=dp), intent(out), dimension(mmaxd, mmaxd) :: result

Overlap matrix

complex(kind=dp), intent(in), dimension(irmd) :: phi
complex(kind=dp), intent(in), dimension(irmd, 0:lmaxd) :: pz
complex(kind=dp), intent(in), dimension(irmd, 0:lmaxd) :: qz
complex(kind=dp), intent(in), dimension(lmmaxd, lmmaxd, irmind:irmd, 2) :: pqns
complex(kind=dp), intent(in), dimension(lmmaxd, lmmaxd) :: acr
complex(kind=dp), intent(in), dimension(lmmaxd, lmmaxd) :: dr
logical, intent(in) :: lirreg

Is true if the irrregular wavefunction is to be used

integer, intent(in) :: ipan

Number of panels in non-MT-region

integer, intent(in), dimension(0:ipand) :: ircut

r points of panel borders

real(kind=dp), intent(in), dimension(irmd) :: drdi

Derivative dr/di

integer, intent(in) :: irmin

Max R for spherical treatment

integer, intent(in) :: lphi

Is the angular momentum of PHI. l-value for LDA+U

integer, intent(in) :: ipand

Number of panels in non-spherical part

integer, intent(in) :: lmaxd

Maximum l component in wave function expansion

integer, intent(in) :: lmmaxd

(KREL+KORBIT+1)(LMAX+1)*2

integer, intent(in) :: mmaxd
integer, intent(in) :: lmpotd

(lpot+1)**2

integer, intent(in) :: irmind

irmd - irnsd

integer, intent(in) :: irmd

Maximum number of radial points