Skip to content
GitLab
Projects
Groups
Snippets
Help
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
fleur
Project overview
Project overview
Details
Activity
Releases
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Issues
52
Issues
52
List
Boards
Labels
Service Desk
Milestones
Operations
Operations
Incidents
Packages & Registries
Packages & Registries
Container Registry
Analytics
Analytics
Repository
Value Stream
Wiki
Wiki
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Create a new issue
Commits
Issue Boards
Open sidebar
fleur
fleur
Commits
ac51edc4
Commit
ac51edc4
authored
Oct 20, 2017
by
Gregor Michalicek
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
Move reading of density matrix mixing parameter out of SCF loop
parent
6e605703
Changes
8
Hide whitespace changes
Inline
Side-by-side
Showing
8 changed files
with
91 additions
and
78 deletions
+91
-78
global/types.F90
global/types.F90
+3
-0
io/cdn_io.F90
io/cdn_io.F90
+34
-10
main/cdngen.F90
main/cdngen.F90
+1
-1
main/fleur.F90
main/fleur.F90
+9
-1
main/fleur_init.F90
main/fleur_init.F90
+3
-0
main/mix.F90
main/mix.F90
+8
-19
mix/u_mix.f90
mix/u_mix.f90
+27
-42
mpi/mpi_bc_all.F90
mpi/mpi_bc_all.F90
+6
-5
No files found.
global/types.F90
View file @
ac51edc4
...
...
@@ -582,6 +582,9 @@ MODULE m_types
TYPE
(
t_efield
)::
efield
LOGICAL
::
l_core_confpot
LOGICAL
::
l_useapw
LOGICAL
::
ldauLinMix
REAL
::
ldauMixParam
REAL
::
ldauSpinf
END
TYPE
t_input
TYPE
t_sliceplot
...
...
io/cdn_io.F90
View file @
ac51edc4
...
...
@@ -1159,16 +1159,17 @@ MODULE m_cdn_io
SUBROUTINE
readDensityMatrix
(
input
,
atoms
,
n_mmp
,
l_error
)
TYPE
(
t_input
),
INTENT
(
IN
)
::
input
TYPE
(
t_atoms
),
INTENT
(
IN
)
::
atoms
COMPLEX
,
INTENT
(
OUT
)
::
n_mmp
(
-
lmaxU_const
:
lmaxU_const
,
-
lmaxU_const
:
lmaxU_const
,&
-
lmaxU_const
:
lmaxU_const
,
-
lmaxU_const
:
lmaxU_const
,&
atoms
%
n_u
,
input
%
jspins
)
LOGICAL
,
INTENT
(
OUT
)
::
l_error
INTEGER
::
mode
,
ioStatus
LOGICAL
::
l_exist
CHARACTER
(
LEN
=
30
)
::
filename
TYPE
(
t_input
),
INTENT
(
INOUT
)
::
input
TYPE
(
t_atoms
),
INTENT
(
IN
)
::
atoms
COMPLEX
,
INTENT
(
OUT
)
::
n_mmp
(
-
lmaxU_const
:
lmaxU_const
,
-
lmaxU_const
:
lmaxU_const
,&
-
lmaxU_const
:
lmaxU_const
,
-
lmaxU_const
:
lmaxU_const
,&
atoms
%
n_u
,
input
%
jspins
)
LOGICAL
,
INTENT
(
OUT
)
::
l_error
INTEGER
::
mode
,
ioStatus
,
numLines
,
iofl
,
i
LOGICAL
::
l_exist
CHARACTER
(
LEN
=
30
)
::
filename
REAL
::
alpha
,
spinf
l_error
=
.FALSE.
ioStatus
=
0
...
...
@@ -1201,7 +1202,30 @@ MODULE m_cdn_io
OPEN
(
69
,
file
=
TRIM
(
ADJUSTL
(
filename
)),
status
=
'unknown'
,
form
=
'formatted'
)
READ
(
69
,
'(7f20.13)'
,
IOSTAT
=
ioStatus
)
n_mmp
REWIND
(
69
)
numLines
=
0
DO
READ
(
69
,
*
,
iostat
=
iofl
)
IF
(
iofl
<
0
)
EXIT
numLines
=
numLines
+
1
END
DO
IF
(
MOD
(
numLines
,
14
*
input
%
jspins
)
==
1
)
THEN
!read in alpha, spinf in this case
REWIND
(
69
)
DO
i
=
1
,
numLines
-1
READ
(
69
,
*
,
iostat
=
iofl
)
END
DO
READ
(
69
,
'(2(6x,f5.3))'
,
IOSTAT
=
ioStatus
)
alpha
,
spinf
IF
(
input
%
ldauLinMix
.AND.
(
input
%
ldauMixParam
.LT.
0.0
))
THEN
input
%
ldauMixParam
=
alpha
input
%
ldauSpinf
=
spinf
END
IF
ELSE
IF
(
MOD
(
numLines
,
14
*
input
%
jspins
)
.NE.
0
)
THEN
CALL
juDFT_error
(
"strange n_mmp_mat-file..."
,
calledby
=
"readDensityMatrix"
)
ELSE
IF
(
input
%
ldauMixParam
.LT.
0.0
)
THEN
input
%
ldauLinMix
=
.FALSE.
END
IF
CLOSE
(
69
)
! IF (ioStatus.NE.0) THEN
! l_error = .TRUE.
! RETURN
...
...
main/cdngen.F90
View file @
ac51edc4
...
...
@@ -162,7 +162,7 @@ SUBROUTINE cdngen(eig_id,mpi,input,banddos,sliceplot,vacuum,&
END
DO
! lda+u
IF
((
atoms
%
n_u
.GT.
0
)
.and.
(
mpi
%
irank
.EQ.
0
))
CALL
u_mix
(
atoms
,
input
%
jspin
s
,
inDen
%
mmpMat
,
outDen
%
mmpMat
)
IF
((
atoms
%
n_u
.GT.
0
)
.and.
(
mpi
%
irank
.EQ.
0
))
CALL
u_mix
(
input
,
atom
s
,
inDen
%
mmpMat
,
outDen
%
mmpMat
)
!+t3e
IF
(
mpi
%
irank
.EQ.
0
)
THEN
...
...
main/fleur.F90
View file @
ac51edc4
...
...
@@ -196,10 +196,18 @@ CONTAINS
IF
(
isDensityMatrixPresent
()
.AND.
atoms
%
n_u
>
0
)
THEN
CALL
readDensityMatrix
(
input
,
atoms
,
inDen
%
mmpMat
,
l_error
)
IF
(
l_error
)
CALL
juDFT_error
(
'Error in reading density matrix!'
,
calledby
=
'fleur'
)
ELSE
ELSE
IF
(
atoms
%
n_u
>
0
)
THEN
IF
(
input
%
ldauLinMix
)
THEN
input
%
ldauMixParam
=
0.05
input
%
ldauSpinf
=
1.0
END
IF
inDen
%
mmpMat
=
CMPLX
(
0.0
,
0.0
)
END
IF
END
IF
#ifdef CPP_MPI
CALL
MPI_BCAST
(
input
%
ldauMixParam
,
1
,
MPI_DOUBLE_PRECISION
,
0
,
mpi
%
mpi_comm
,
ierr
)
CALL
MPI_BCAST
(
input
%
ldauSpinf
,
1
,
MPI_DOUBLE_PRECISION
,
0
,
mpi
%
mpi_comm
,
ierr
)
#endif
! Initialize and load inDen density (end)
! Initialize potentials (start)
...
...
main/fleur_init.F90
View file @
ac51edc4
...
...
@@ -143,6 +143,9 @@
CALL
initWannierDefaults
(
wann
)
input
%
minDistance
=
0.0
input
%
ldauLinMix
=
.TRUE.
input
%
ldauMixParam
=
-1.0
input
%
ldauSpinf
=
1.0
kpts
%
ntet
=
1
kpts
%
numSpecialPoints
=
1
IF
(
input
%
l_inpXML
)
THEN
...
...
main/mix.F90
View file @
ac51edc4
...
...
@@ -101,32 +101,21 @@ SUBROUTINE mix(stars,atoms,sphhar,vacuum,input,sym,cell,noco,oneD,&
2
*
vacuum
%
nmzd
*
vacuum
%
nvac
END
IF
! LDA+U (start)
n_mmpTemp
=
inDen
%
mmpMat
n_u_keep
=
atoms
%
n_u
INQUIRE
(
file
=
'n_mmp_mat'
,
exist
=
l_ldaU
)
IF
(
l_ldaU
)
THEN
IF
(
isDensityMatrixPresent
())
THEN
!In an LDA+U caclulation, also the density matrix is included in the
!supervectors (sm,fsm) if no mixing factors are in the n_mmp_mat-file
OPEN
(
69
,
file
=
'n_mmp_mat'
,
status
=
'old'
,
form
=
'formatted'
)
i
=
0
DO
READ
(
69
,
*
,
iostat
=
iofl
)
IF
(
iofl
<
0
)
EXIT
i
=
i
+
1
END
DO
CLOSE
(
69
)
IF
(
MOD
(
i
,
14
*
input
%
jspins
)
==
1
)
THEN
! was already mixed in u_mix
!supervectors (sm,fsm) if no linear mixing is performed on it.
IF
(
input
%
ldauLinMix
)
THEN
atoms
%
n_u
=
0
ELSE
IF
(
MOD
(
i
,
28
*
input
%
jspins
)
==
0
)
THEN
! mix here
mmap
=
mmap
+
7
*
7
*
2
*
atoms
%
n_u
*
input
%
jspins
! add 7*7 complex numbers per atoms%n_u and spin
ELSE
CALL
juDFT_error
(
"strange n_mmp_mat-file..."
,
calledby
=
"mix"
)
mmap
=
mmap
+
7
*
7
*
2
*
atoms
%
n_u
*
input
%
jspins
! add 7*7 complex numbers per atoms%n_u and spin
END
IF
ELSE
atoms
%
n_u
=
0
END
IF
! l_ldaU
atoms
%
n_u
=
0
END
IF
! LDA+U (end)
ALLOCATE
(
sm
(
mmap
),
fsm
(
mmap
))
...
...
mix/u_mix.f90
View file @
ac51edc4
...
...
@@ -12,19 +12,20 @@ MODULE m_umix
! --------------------------------------------------------
! Extension to multiple U per atom type by G.M. 2017
CONTAINS
SUBROUTINE
u_mix
(
atoms
,
jspin
s
,
n_mmp_in
,
n_mmp_out
)
SUBROUTINE
u_mix
(
input
,
atom
s
,
n_mmp_in
,
n_mmp_out
)
USE
m_types
USE
m_cdn_io
USE
m_nmat_rot
USE
m_xmlOutput
! ... Arguments
IMPLICIT
NONE
TYPE
(
t_input
),
INTENT
(
IN
)
::
input
TYPE
(
t_atoms
),
INTENT
(
IN
)
::
atoms
INTEGER
,
INTENT
(
IN
)
::
jspins
COMPLEX
,
INTENT
(
INOUT
)
::
n_mmp_out
(
-3
:
3
,
-3
:
3
,
atoms
%
n_u
,
jspins
)
COMPLEX
,
INTENT
(
INOUT
)
::
n_mmp_in
(
-3
:
3
,
-3
:
3
,
atoms
%
n_u
,
jspins
)
COMPLEX
,
INTENT
(
INOUT
)
::
n_mmp_out
(
-3
:
3
,
-3
:
3
,
atoms
%
n_u
,
input
%
jspins
)
COMPLEX
,
INTENT
(
INOUT
)
::
n_mmp_in
(
-3
:
3
,
-3
:
3
,
atoms
%
n_u
,
input
%
jspins
)
!
! ... Locals ...
INTEGER
j
,
k
,
iofl
,
l
,
itype
,
ios
,
i_u
,
jsp
,
lty
(
atoms
%
n_u
)
...
...
@@ -52,13 +53,13 @@ CONTAINS
END
DO
CLOSE
(
68
)
zero
=
0.0
CALL
nmat_rot
(
zero
,
-
theta
,
-
phi
,
3
,
atoms
%
n_u
,
jspins
,
lty
,
n_mmp_out
)
CALL
nmat_rot
(
zero
,
-
theta
,
-
phi
,
3
,
atoms
%
n_u
,
input
%
jspins
,
lty
,
n_mmp_out
)
END
IF
! Write out n_mmp_out to out.xml file
CALL
openXMLElementNoAttributes
(
'ldaUDensityMatrix'
)
DO
jsp
=
1
,
jspins
DO
jsp
=
1
,
input
%
jspins
DO
i_u
=
1
,
atoms
%
n_u
l
=
atoms
%
lda_u
(
i_u
)
%
l
itype
=
atoms
%
lda_u
(
i_u
)
%
atomType
...
...
@@ -78,26 +79,22 @@ CONTAINS
END
DO
CALL
closeXMLElement
(
'ldaUDensityMatrix'
)
!
! check for LDA+U and open density-matrix - file
!
INQUIRE
(
file
=
'n_mmp_mat'
,
exist
=
n_exist
)
! Check whether density matrix is present and open density-matrix - file
n_exist
=
isDensityMatrixPresent
()
OPEN
(
69
,
file
=
'n_mmp_mat'
,
status
=
'unknown'
,
form
=
'formatted'
)
IF
(
n_exist
)
THEN
ALLOCATE
(
n_mmp
(
-3
:
3
,
-3
:
3
,
atoms
%
n_u
,
jspins
))
READ
(
69
,
9000
)
n_mmp
(:,:,:,:)
n_mmp
=
CMPLX
(
0.0
,
0.0
)
READ
(
69
,
'(2(6x,f5.3))'
,
IOSTAT
=
iofl
)
alpha
,
spinf
IF
(
iofl
==
0
)
THEN
!
! mix here straight with given mixing factors
!
REWIND
(
69
)
IF
(
input
%
ldauLinMix
)
THEN
! mix here straight with given mixing factors
ALLOCATE
(
n_mmp
(
-3
:
3
,
-3
:
3
,
atoms
%
n_u
,
input
%
jspins
))
alpha
=
input
%
ldauMixParam
spinf
=
input
%
ldauSpinf
sum1
=
0.0
IF
(
jspins
.EQ.
1
)
THEN
IF
(
input
%
jspins
.EQ.
1
)
THEN
DO
i_u
=
1
,
atoms
%
n_u
DO
j
=
-3
,
3
DO
k
=
-3
,
3
...
...
@@ -135,18 +132,11 @@ CONTAINS
WRITE
(
69
,
9000
)
n_mmp
WRITE
(
69
,
'(2(a6,f5.3))'
)
'alpha='
,
alpha
,
'spinf='
,
spinf
n_mmp_in
=
n_mmp
DEALLOCATE
(
n_mmp
)
ELSE
! input%ldauLinMix
ELSEIF
(
iofl
>
0
)
THEN
!
! read error ; stop
!
WRITE
(
6
,
*
)
'ERROR READING mixing factors in n_mmp_mat'
WRITE
(
6
,
'(2(a6,f5.3))'
)
'alpha='
,
alpha
,
'spinf='
,
spinf
CALL
juDFT_error
(
"ERROR READING n_mmp_mat"
,
calledby
=
"u_mix"
)
ELSE
!
! calculate distance and write new n_mmp to mix in broyden.F
!
sum1
=
0.0
DO
i_u
=
1
,
atoms
%
n_u
DO
j
=
-3
,
3
...
...
@@ -155,7 +145,7 @@ CONTAINS
END
DO
END
DO
END
DO
IF
(
jspins
.EQ.
1
)
THEN
IF
(
input
%
jspins
.EQ.
1
)
THEN
WRITE
(
6
,
'(a16,f12.6)'
)
'n_mmp distance ='
,
sum1
ELSE
sum2
=
0.0
...
...
@@ -175,20 +165,15 @@ CONTAINS
WRITE
(
6
,
'(14f12.6)'
)
(
n_mmp_out
(
k
,
j
,
1
,
2
),
k
=
-3
,
3
)
END
DO
END
IF
REWIND
(
69
)
WRITE
(
69
,
9000
)
n_mmp_in
WRITE
(
69
,
9000
)
n_mmp_out
END
IF
! iofl == 0
DEALLOCATE
(
n_mmp
)
ELSE
!
END
IF
! input%ldauLinMix
ELSE
! isDensityMatrixPresent()
! first time with lda+u; write new n_mmp
!
WRITE
(
69
,
9000
)
n_mmp_out
WRITE
(
69
,
'(2(a6,f5.3))'
)
'alpha='
,
0.05
,
'spinf='
,
1.0
WRITE
(
69
,
'(2(a6,f5.3))'
)
'alpha='
,
input
%
ldauMixParam
,
'spinf='
,
input
%
ldauSpinf
n_mmp_in
=
n_mmp_out
END
IF
END
IF
! isDensityMatrixPresent()
9000
FORMAT
(
7f20.13
)
...
...
mpi/mpi_bc_all.F90
View file @
ac51edc4
...
...
@@ -40,8 +40,8 @@ CONTAINS
REAL
rdum
! .. Local Arrays ..
INTEGER
i
(
39
),
ierr
(
3
)
REAL
r
(
3
0
)
LOGICAL
l
(
4
4
)
REAL
r
(
3
2
)
LOGICAL
l
(
4
5
)
! ..
! .. External Subroutines..
#ifdef CPP_MPI
...
...
@@ -65,7 +65,7 @@ CONTAINS
r
(
19
)
=
cell
%
volint
;
r
(
20
)
=
hybrid
%
gcutm1
;
r
(
21
)
=
hybrid
%
tolerance1
;
r
(
22
)
=
0.0
r
(
23
)
=
0.0
;
r
(
24
)
=
input
%
delgau
;
r
(
25
)
=
input
%
tkb
;
r
(
26
)
=
input
%
efield
%
vslope
r
(
27
)
=
0.0
;
r
(
28
)
=
0.0
!r(27)=aMix_VHSE() ; r(28)=omega_VHSE()
r
(
29
)
=
input
%
minDistance
;
r
(
30
)
=
obsolete
%
chng
r
(
29
)
=
input
%
minDistance
;
r
(
30
)
=
obsolete
%
chng
;
r
(
31
)
=
input
%
ldauMixParam
;
r
(
32
)
=
input
%
ldauSpinf
l
(
1
)
=
input
%
eonly
;
l
(
2
)
=
input
%
l_useapw
;
l
(
3
)
=
input
%
secvar
;
l
(
4
)
=
sym
%
zrfs
;
l
(
5
)
=
input
%
film
l
(
6
)
=
sym
%
invs
;
l
(
7
)
=
sym
%
invs2
;
l
(
8
)
=
input
%
l_bmt
;
l
(
9
)
=
input
%
l_f
;
l
(
10
)
=
input
%
cdinf
...
...
@@ -77,7 +77,7 @@ CONTAINS
l
(
34
)
=
banddos
%
l_mcd
;
l
(
35
)
=
input
%
sso_opt
(
1
)
l
(
36
)
=
input
%
sso_opt
(
2
)
;
l
(
37
)
=
obsolete
%
pot8
;
l
(
38
)
=
input
%
efield
%
l_segmented
l
(
39
)
=
sym
%
symor
;
l
(
40
)
=
input
%
frcor
;
l
(
41
)
=
input
%
tria
;
l
(
42
)
=
input
%
efield
%
dirichlet
l
(
43
)
=
input
%
efield
%
l_dirichlet_coeff
;
l
(
44
)
=
input
%
l_coreSpec
l
(
43
)
=
input
%
efield
%
l_dirichlet_coeff
;
l
(
44
)
=
input
%
l_coreSpec
;
l
(
45
)
=
input
%
ldauLinMix
ENDIF
!
CALL
MPI_BCAST
(
i
,
SIZE
(
i
),
MPI_INTEGER
,
0
,
mpi
%
mpi_comm
,
ierr
)
...
...
@@ -99,6 +99,7 @@ CONTAINS
vacuum
%
locx
(
1
)
=
r
(
11
);
vacuum
%
locx
(
2
)
=
r
(
12
);
vacuum
%
locy
(
1
)
=
r
(
13
);
vacuum
%
locy
(
2
)
=
r
(
14
)
sliceplot
%
e1s
=
r
(
6
)
;
sliceplot
%
e2s
=
r
(
7
)
;
noco
%
theta
=
r
(
8
)
;
noco
%
phi
=
r
(
9
)
;
vacuum
%
tworkf
=
r
(
10
)
cell
%
omtil
=
r
(
1
)
;
cell
%
area
=
r
(
2
)
;
vacuum
%
delz
=
r
(
3
)
;
cell
%
z1
=
r
(
4
)
;
input
%
alpha
=
r
(
5
)
input
%
ldauMixParam
=
r
(
31
)
;
input
%
ldauSpinf
=
r
(
32
)
!
CALL
MPI_BCAST
(
l
,
SIZE
(
l
),
MPI_LOGICAL
,
0
,
mpi
%
mpi_comm
,
ierr
)
input
%
efield
%
l_dirichlet_coeff
=
l
(
43
)
;
input
%
l_useapw
=
l
(
2
)
...
...
@@ -114,7 +115,7 @@ CONTAINS
sym
%
invs
=
l
(
6
)
;
sym
%
invs2
=
l
(
7
)
;
input
%
l_bmt
=
l
(
8
)
;
input
%
l_f
=
l
(
9
)
;
input
%
cdinf
=
l
(
10
)
input
%
eonly
=
l
(
1
)
;
input
%
secvar
=
l
(
3
)
;
sym
%
zrfs
=
l
(
4
)
;
input
%
film
=
l
(
5
)
input
%
efield
%
l_segmented
=
l
(
38
)
;
sym
%
symor
=
l
(
39
);
input
%
efield
%
dirichlet
=
l
(
40
)
input
%
efield
%
l_dirichlet_coeff
=
l
(
41
)
;
input
%
l_coreSpec
=
l
(
44
)
input
%
efield
%
l_dirichlet_coeff
=
l
(
41
)
;
input
%
l_coreSpec
=
l
(
44
)
;
input
%
ldauLinMix
=
l
(
45
)
!
! -> Broadcast the arrays:
IF
(
input
%
efield
%
l_segmented
)
THEN
...
...
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