Skip to content
GitLab
Menu
Projects
Groups
Snippets
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
west-public
West
Commits
9d3e7f84
Commit
9d3e7f84
authored
Feb 24, 2018
by
Matteo Gerosa
Browse files
Fixed bug in class_bz_grid. Works with non-collinear spin.
parent
fb59a885
Changes
7
Hide whitespace changes
Inline
Side-by-side
DFPT_kernel/apply_sternheimerop_to_m_wfcs.f90
View file @
9d3e7f84
...
...
@@ -63,7 +63,7 @@ SUBROUTINE apply_sternheimerop_to_m_wfcs(nbndval, psi, hpsi, e, alpha, m)
za
=
CMPLX
(
-
e
(
ibnd
),
0._DP
,
KIND
=
DP
)
CALL
ZAXPY
(
npw
,
za
,
psi
(
1
,
ibnd
),
1
,
hpsi
(
1
,
ibnd
),
1
)
ENDDO
IF
(
n
pol
==
2
)
THEN
IF
(
n
oncolin
)
THEN
DO
ibnd
=
1
,
m
za
=
CMPLX
(
-
e
(
ibnd
),
0._DP
,
KIND
=
DP
)
CALL
ZAXPY
(
npw
,
za
,
psi
(
npwx
+1
,
ibnd
),
1
,
hpsi
(
npwx
+1
,
ibnd
),
1
)
...
...
DFPT_kernel/dfpt_module.f90
View file @
9d3e7f84
...
...
@@ -240,7 +240,7 @@ MODULE dfpt_module
!
ENDDO
!
IF
(
npol
==
2
)
THEN
IF
(
noncolin
)
THEN
DO
ibnd
=
1
,
nbndval
!
CALL
single_invfft_k
(
dffts
,
npwkq
,
npwx
,
evckmq
(
npwx
+1
,
ibnd
),
psic
,
'Wave'
,
igk_k
(
1
,
ikqs
))
...
...
@@ -316,7 +316,7 @@ MODULE dfpt_module
!
ENDDO
!
IF
(
npol
==
2
)
THEN
IF
(
noncolin
)
THEN
DO
ibnd
=
1
,
nbndval
!
CALL
single_invfft_k
(
dffts
,
npwkq
,
npwx
,
evckmq
(
npwx
+1
,
ibnd
),
psic
,
'Wave'
,
igk_k
(
1
,
ikqs
))
...
...
DFPT_kernel/linsolve_sternheimer_m_wfcts.f90
View file @
9d3e7f84
...
...
@@ -101,7 +101,7 @@ SUBROUTINE linsolve_sternheimer_m_wfcts ( nbndval, m, b, x, e, eprec, tr2, ierr
DO
ibnd
=
1
,
m
CALL
ZAXPY
(
npw
,(
-1._DP
,
0._DP
),
b
(
1
,
ibnd
),
1
,
g
(
1
,
ibnd
),
1
)
ENDDO
IF
(
n
pol
==
2
)
THEN
IF
(
n
oncolin
)
THEN
DO
ibnd
=
1
,
m
CALL
ZAXPY
(
npw
,(
-1._DP
,
0._DP
),
b
(
npwx
+1
,
ibnd
),
1
,
g
(
npwx
+1
,
ibnd
),
1
)
ENDDO
...
...
Tools/class_bz_grid.f90
View file @
9d3e7f84
...
...
@@ -23,8 +23,8 @@ MODULE class_bz_grid
TYPE
,
PUBLIC
::
bz_grid
!
INTEGER
::
ngrid
(
3
)
=
(/
1
,
1
,
1
/)
! number of points in each direction
INTEGER
::
np
=
1
! total number of points
= ngrid(1) * ngrid(2) * ngrid(3) **** NOOOOOO ******
INTEGER
::
ns
=
1
! total number of spin = nspin
INTEGER
::
np
=
1
! total number of points
INTEGER
::
ns
=
1
! total number of spin = nspin
_lsda
INTEGER
::
nps
=
1
! total number of points and spins = np * ns
INTEGER
,
ALLOCATABLE
::
ip
(:)
! given ips --> ip
INTEGER
,
ALLOCATABLE
::
is
(:)
! given ips --> is
...
...
@@ -51,9 +51,8 @@ MODULE class_bz_grid
USE
cell_base
,
ONLY
:
at
,
bg
USE
klist
,
ONLY
:
xk
,
wk
,
nkstot
USE
start_k
,
ONLY
:
nk1
,
nk2
,
nk3
USE
pwcom
,
ONLY
:
nspin
USE
noncollin_module
,
ONLY
:
nspin
_lsda
USE
control_flags
,
ONLY
:
gamma_only
!USE westcom, ONLY : nq
USE
constants
,
ONLY
:
eps8
USE
westcom
,
ONLY
:
qlist
!
...
...
@@ -76,8 +75,9 @@ MODULE class_bz_grid
! This is a workaround to prevent ngrid(:) to be set to (/ 0, 0, 0 /) in the gamma_only case (espresso default)
IF
(
.NOT.
gamma_only
)
this
%
ngrid
(
1
:
3
)
=
(/
nk1
,
nk2
,
nk3
/)
this
%
np
=
this
%
ngrid
(
1
)
*
this
%
ngrid
(
2
)
*
this
%
ngrid
(
3
)
this
%
ns
=
nspin
this
%
nps
=
nkstot
! = np * ns
this
%
ns
=
nspin_lsda
! = 1 if nspin = 1 (unpolarized) or nspin = 4 (noncollinear)
! = 2 if nspin = 2 (collinear)
this
%
nps
=
nkstot
! = np * ns
!
! generate p-vectors in cart
!
...
...
Wfreq/calc_exx2.f90
View file @
9d3e7f84
...
...
@@ -70,11 +70,11 @@ SUBROUTINE calc_exx2_gamma( sigma_exx, nb1, nb2 )
CALL
io_push_bar
()
!
ALLOCATE
(
pertg
(
ngm
)
)
IF
(
noncolin
)
THEN
ALLOCATE
(
pertr_nc
(
dffts
%
nnr
,
npol
)
)
ELSE
!
IF(noncolin) THEN
!
ALLOCATE( pertr_nc( dffts%nnr, npol ) )
!
ELSE
ALLOCATE
(
pertr
(
dffts
%
nnr
)
)
ENDIF
!
ENDIF
!
CALL
pot3D
%
init
(
'Dense'
,
.FALSE.
,
'gb'
)
!
...
...
@@ -135,40 +135,40 @@ SUBROUTINE calc_exx2_gamma( sigma_exx, nb1, nb2 )
!
DO
ib
=
nb1
,
nb2
!
IF
(
gamma_only
)
THEN
!
IF(gamma_only) THEN
CALL
single_invfft_gamma
(
dffts
,
npw
,
npwx
,
evc
(
1
,
ib
),
psic
,
'Wave'
)
ELSEIF
(
noncolin
)
THEN
CALL
single_invfft_k
(
dffts
,
npw
,
npwx
,
evc
(
1
,
ib
),
psic_nc
(
1
,
1
),
'Wave'
,
igk_k
(
1
,
current_k
))
CALL
single_invfft_k
(
dffts
,
npw
,
npwx
,
evc
(
1
+
npwx
,
ib
),
psic_nc
(
1
,
2
),
'Wave'
,
igk_k
(
1
,
current_k
))
!
ELSE
!
CALL single_invfft_k(dffts,npw,npwx,evc(1,ib),psic,'Wave',igk_k(1,current_k))
ENDIF
!
ELSEIF(noncolin) THEN
!
CALL single_invfft_k(dffts,npw,npwx,evc(1 ,ib),psic_nc(1,1),'Wave',igk_k(1,current_k))
!
CALL single_invfft_k(dffts,npw,npwx,evc(1+npwx,ib),psic_nc(1,2),'Wave',igk_k(1,current_k))
!
ELSE
!
CALL single_invfft_k(dffts,npw,npwx,evc(1,ib),psic,'Wave',igk_k(1,current_k))
!
ENDIF
!
DO
iv
=
1
,
vband
%
nloc
!
iv_glob
=
vband
%
l2g
(
iv
)
!
! Bring it to R-space
IF
(
gamma_only
)
THEN
!
IF(gamma_only) THEN
CALL
single_invfft_gamma
(
dffts
,
npw
,
npwx
,
evc
(
1
,
iv_glob
),
pertr
,
'Wave'
)
DO
ir
=
1
,
dffts
%
nnr
pertr
(
ir
)
=
psic
(
ir
)
*
pertr
(
ir
)
ENDDO
CALL
single_fwfft_gamma
(
dffts
,
ngm
,
ngm
,
pertr
,
pertg
,
'Dense'
)
ELSEIF
(
noncolin
)
THEN
CALL
single_invfft_k
(
dffts
,
npw
,
npwx
,
evc
(
1
,
iv_glob
),
pertr_nc
(
1
,
1
),
'Wave'
,
igk_k
(
1
,
current_k
))
CALL
single_invfft_k
(
dffts
,
npw
,
npwx
,
evc
(
1
+
npwx
,
iv_glob
),
pertr_nc
(
1
,
2
),
'Wave'
,
igk_k
(
1
,
current_k
))
DO
ir
=
1
,
dffts
%
nnr
pertr_nc
(
ir
,
1
)
=
DCONJG
(
psic_nc
(
ir
,
1
))
*
pertr_nc
(
ir
,
1
)
+
DCONJG
(
psic_nc
(
ir
,
2
))
*
pertr_nc
(
ir
,
2
)
ENDDO
CALL
single_fwfft_k
(
dffts
,
ngm
,
ngm
,
pertr_nc
(
1
,
1
),
pertg
,
'Dense'
)
! no igk
!
ELSE
!
CALL single_invfft_k(dffts,npw,npwx,evc(1,iv_glob),pertr,'Wave',igk_k(1,current_k))
!
DO ir=1,dffts%nnr
!
pertr(ir)=DCONJG(psic(ir))*pertr(ir)
!
ENDDO
!
CALL single_fwfft_k(dffts,ngms,ngms,pertr,pertg,'Smooth') ! no igk
ENDIF
!
ELSEIF(noncolin) THEN
!
CALL single_invfft_k(dffts,npw,npwx,evc(1 ,iv_glob),pertr_nc(1,1),'Wave',igk_k(1,current_k))
!
CALL single_invfft_k(dffts,npw,npwx,evc(1+npwx,iv_glob),pertr_nc(1,2),'Wave',igk_k(1,current_k))
!
DO ir=1,dffts%nnr
!
pertr_nc(ir,1)=DCONJG(psic_nc(ir,1))*pertr_nc(ir,1)+DCONJG(psic_nc(ir,2))*pertr_nc(ir,2)
!
ENDDO
!
CALL single_fwfft_k(dffts,ngm,ngm,pertr_nc(1,1),pertg,'Dense') ! no igk
!
ELSE
!
CALL single_invfft_k(dffts,npw,npwx,evc(1,iv_glob),pertr,'Wave',igk_k(1,current_k))
!
DO ir=1,dffts%nnr
!
pertr(ir)=DCONJG(psic(ir))*pertr(ir)
!
ENDDO
!
CALL single_fwfft_k(dffts,ngms,ngms,pertr,pertg,'Smooth') ! no igk
!
ENDIF
!
DO
ig
=
1
,
ngm
pertg
(
ig
)
=
pertg
(
ig
)
*
pot3D
%
sqvc
(
ig
)
...
...
@@ -192,11 +192,11 @@ SUBROUTINE calc_exx2_gamma( sigma_exx, nb1, nb2 )
!
DEALLOCATE
(
pertg
)
! DEALLOCATE( mysqvc )
IF
(
noncolin
)
THEN
DEALLOCATE
(
pertr_nc
)
ELSE
!
IF( noncolin ) THEN
!
DEALLOCATE( pertr_nc )
!
ELSE
DEALLOCATE
(
pertr
)
ENDIF
!
ENDIF
!
END
SUBROUTINE
!
...
...
Wfreq/solve_gfreq.f90
View file @
9d3e7f84
...
...
@@ -673,7 +673,11 @@ SUBROUTINE solve_gfreq_k(l_read_restart)
ENDDO
! KPOINT-SPIN (MATRIX ELEMENT)
!
DEALLOCATE
(
phase
)
DEALLOCATE
(
psick
)
IF
(
noncolin
)
THEN
DEALLOCATE
(
psick_nc
)
ELSE
DEALLOCATE
(
psick
)
ENDIF
DEALLOCATE
(
evck
)
!
CALL
stop_bar_type
(
barra
,
'glanczos'
)
...
...
Wfreq/solve_wfreq.f90
View file @
9d3e7f84
...
...
@@ -1154,6 +1154,11 @@ SUBROUTINE solve_wfreq_k(l_read_restart,l_generate_plot)
CALL
stop_bar_type
(
barra
,
'wlanczos'
)
!
DEALLOCATE
(
evckpq
)
IF
(
noncolin
)
THEN
DEALLOCATE
(
psick_nc
)
ELSE
DEALLOCATE
(
psick
)
ENDIF
DEALLOCATE
(
phase
)
!
! EPS-1 imfreq
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment