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
54
Issues
54
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
1737a5bf
Commit
1737a5bf
authored
Apr 13, 2018
by
Gregor Michalicek
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
Introduce t_mcd type to cdn/cdnval.F90
parent
aea50b04
Changes
5
Hide whitespace changes
Inline
Side-by-side
Showing
5 changed files
with
91 additions
and
88 deletions
+91
-88
cdn/cdnval.F90
cdn/cdnval.F90
+13
-29
cdn/eparas.f90
cdn/eparas.f90
+9
-11
dos/doswrite.f90
dos/doswrite.f90
+13
-32
orbdep/mcd_init.f90
orbdep/mcd_init.f90
+11
-15
types/types_cdnval.f90
types/types_cdnval.f90
+45
-1
No files found.
cdn/cdnval.F90
View file @
1737a5bf
...
...
@@ -124,7 +124,7 @@ CONTAINS
! .. Local Scalars ..
TYPE
(
t_lapw
)::
lapw
INTEGER
::
llpd
REAL
wronk
,
emcd_lo
,
emcd_up
REAL
wronk
INTEGER
i
,
ie
,
iv
,
ivac
,
j
,
k
,
l
,
n
,
ilo
,
isp
,&
nbands
,
noded
,
nodeu
,
noccbd
,
nslibd
,
na
,&
ikpt
,
jsp_start
,
jsp_end
,
ispin
...
...
@@ -132,8 +132,7 @@ CONTAINS
INTEGER
n_size
,
i_rec
,
n_rank
,
ncored
,
n_start
,
n_end
,
noccbd_l
,
nbasfcn
LOGICAL
l_fmpl
,
l_evp
,
l_orbcomprot
,
l_real
! ...Local Arrays ..
INTEGER
n_bands
(
0
:
dimension
%
neigd
),
ncore
(
atoms
%
ntype
)
REAL
e_mcd
(
atoms
%
ntype
,
input
%
jspins
,
dimension
%
nstd
)
INTEGER
n_bands
(
0
:
dimension
%
neigd
)
REAL
eig
(
dimension
%
neigd
)
REAL
vz0
(
2
)
REAL
uuilon
(
atoms
%
nlod
,
atoms
%
ntype
),
duilon
(
atoms
%
nlod
,
atoms
%
ntype
)
...
...
@@ -154,10 +153,10 @@ CONTAINS
REAL
,
ALLOCATABLE
::
sqlo
(:,:,:)
REAL
,
ALLOCATABLE
::
qal
(:,:,:,:),
sqal
(:,:,:),
ener
(:,:,:)
REAL
,
ALLOCATABLE
::
svac
(:,:),
pvac
(:,:)
,
mcd
(:,:,:)
REAL
,
ALLOCATABLE
::
svac
(:,:),
pvac
(:,:)
REAL
,
ALLOCATABLE
::
enerlo
(:,:,:),
qmat
(:,:,:,:)
COMPLEX
,
ALLOCATABLE
::
qstars
(:,:,:,:)
,
m_mcd
(:,:,:,:)
COMPLEX
,
ALLOCATABLE
::
qstars
(:,:,:,:)
TYPE
(
t_orb
)
::
orb
TYPE
(
t_denCoeffs
)
::
denCoeffs
...
...
@@ -165,6 +164,7 @@ CONTAINS
TYPE
(
t_force
)
::
force
TYPE
(
t_slab
)
::
slab
TYPE
(
t_eigVecCoeffs
)
::
eigVecCoeffs
TYPE
(
t_mcd
)
::
mcd
TYPE
(
t_usdus
)
::
usdus
TYPE
(
t_zMat
)
::
zMat
...
...
@@ -219,15 +219,7 @@ CONTAINS
CALL
orb
%
init
(
atoms
,
noco
,
jsp_start
,
jsp_end
)
IF
(
banddos
%
l_mcd
)
THEN
emcd_lo
=
banddos
%
e_mcd_lo
emcd_up
=
banddos
%
e_mcd_up
ALLOCATE
(
m_mcd
(
dimension
%
nstd
,(
3+1
)
**
2
,
3
*
atoms
%
ntype
,
2
)
)
ALLOCATE
(
mcd
(
3
*
atoms
%
ntype
,
dimension
%
nstd
,
dimension
%
neigd
)
)
IF
(
.not.
banddos
%
dos
)
WRITE
(
*
,
*
)
'For mcd-spectra set banddos%dos=T!'
ELSE
ALLOCATE
(
m_mcd
(
1
,
1
,
1
,
1
),
mcd
(
1
,
1
,
1
)
)
ENDIF
CALL
mcd
%
init1
(
banddos
,
dimension
,
input
,
atoms
)
! calculation of core spectra (EELS) initializations -start-
CALL
corespec_init
(
input
,
atoms
,
coreSpecInput
)
...
...
@@ -276,20 +268,16 @@ CONTAINS
END
IF
END
DO
IF
(
banddos
%
l_mcd
)
THEN
CALL
mcd_init
(
atoms
,
input
,
dimension
,&
vTot
%
mt
(:,
0
,:,:),
g
,
f
,
emcd_up
,
emcd_lo
,
n
,
jspin
,&
ncore
,
e_mcd
,
m_mcd
)
ncored
=
max
(
ncore
(
n
),
ncored
)
CALL
mcd_init
(
atoms
,
input
,
dimension
,
vTot
%
mt
(:,
0
,:,:),
g
,
f
,
mcd
,
n
,
jspin
)
ncored
=
max
(
mcd
%
ncore
(
n
),
ncored
)
END
IF
IF
(
l_cs
)
CALL
corespec_rme
(
atoms
,
input
,
n
,
dimension
%
nstd
,&
input
%
jspins
,
jspin
,
results
%
ef
,&
dimension
%
msh
,
vTot
%
mt
(:,
0
,:,:),
f
,
g
)
!
!---> generate the extra wavefunctions for the local orbitals,
!---> if there are any.
!
IF
(
atoms
%
nlo
(
n
)
>
0
)
THEN
DO
ispin
=
jsp_start
,
jsp_end
CALL
radflo
(
atoms
,
n
,
ispin
,
enpara
%
ello0
(
1
,
1
,
ispin
),
vTot
%
mt
(:,
0
,
n
,
ispin
),
f
(
1
,
1
,
0
,
ispin
),&
...
...
@@ -330,9 +318,6 @@ CONTAINS
!--> loop over k-points: each can be a separate task
IF
(
kpts
%
nkpt
<
mpi
%
isize
)
THEN
l_evp
=
.true.
IF
(
banddos
%
l_mcd
)
THEN
mcd
(:,:,:)
=
0.0
ENDIF
ener
(:,:,:)
=
0.0
sqal
(:,:,:)
=
0.0
qal
(:,:,:,:)
=
0.0
...
...
@@ -574,10 +559,10 @@ CONTAINS
!---> atom and angular momentum
IF
(
.not.
sliceplot
%
slice
)
THEN
CALL
eparas
(
ispin
,
atoms
,
noccbd
,
mpi
,
ikpt
,
noccbd
,
we
,
eig
,&
skip_t
,
l_evp
,
eigVecCoeffs
,
usdus
,&
ncore
,
banddos
%
l_mcd
,
m
_mcd
,
enerlo
(
1
,
1
,
ispin
),
sqlo
(
1
,
1
,
ispin
),&
skip_t
,
l_evp
,
eigVecCoeffs
,
usdus
,
mcd
,
&
banddos
%
l
_mcd
,
enerlo
(
1
,
1
,
ispin
),
sqlo
(
1
,
1
,
ispin
),&
ener
(
0
,
1
,
ispin
),
sqal
(
0
,
1
,
ispin
),&
qal
(
0
:,:,:,
ispin
)
,
mcd
)
qal
(
0
:,:,:,
ispin
))
IF
(
noco
%
l_mperp
.AND.
(
ispin
==
jsp_end
))
THEN
CALL
qal_21
(
atoms
,
input
,
noccbd
,
we
,
noco
,
eigVecCoeffs
,
denCoeffsOffdiag
,
qal
,
qmat
)
...
...
@@ -661,7 +646,7 @@ CONTAINS
!--dw now write k-point data to tmp_dos
CALL
write_dos
(
eig_id
,
ikpt
,
jspin
,
qal
(:,:,:,
jspin
),
qvac
(:,:,
ikpt
,
jspin
),
qis
(:,
ikpt
,
jspin
),&
qvlay
(:,:,:,
ikpt
,
jspin
),
qstars
,
ksym
,
jsym
,
mcd
,
slab
%
qintsl
,&
qvlay
(:,:,:,
ikpt
,
jspin
),
qstars
,
ksym
,
jsym
,
mcd
%
mcd
,
slab
%
qintsl
,&
slab
%
qmtsl
(:,:),
qmtp
(:,:),
orbcomp
)
CALL
timestop
(
"cdnval: write_info"
)
...
...
@@ -698,8 +683,7 @@ CONTAINS
CALL
timestart
(
"cdnval: dos"
)
IF
(
mpi
%
irank
==
0
)
THEN
CALL
doswrite
(
eig_id
,
dimension
,
kpts
,
atoms
,
vacuum
,
input
,
banddos
,&
sliceplot
,
noco
,
sym
,
cell
,
banddos
%
l_mcd
,
ncored
,
ncore
,
e_mcd
,&
results
%
ef
,
results
%
bandgap
,
slab
%
nsld
,
oneD
)
sliceplot
,
noco
,
sym
,
cell
,
mcd
,
ncored
,
results
,
slab
%
nsld
,
oneD
)
IF
(
banddos
%
dos
.AND.
(
banddos
%
ndir
.EQ.
-3
))
THEN
CALL
Ek_write_sl
(
eig_id
,
dimension
,
kpts
,
atoms
,
vacuum
,
input
,
jspin
,
sym
,
cell
,
slab
)
END
IF
...
...
cdn/eparas.f90
View file @
1737a5bf
...
...
@@ -24,13 +24,14 @@ MODULE m_eparas
!
CONTAINS
SUBROUTINE
eparas
(
jsp
,
atoms
,
noccbd
,
mpi
,
ikpt
,
ne
,
we
,
eig
,
skip_t
,
l_evp
,
eigVecCoeffs
,&
usdus
,
ncore
,
l_mcd
,
m_mcd
,
enerlo
,
sqlo
,
ener
,
sqal
,
qal
,
mcd
)
usdus
,
mcd
,
l_mcd
,
enerlo
,
sqlo
,
ener
,
sqal
,
qal
)
USE
m_types
IMPLICIT
NONE
TYPE
(
t_usdus
),
INTENT
(
IN
)
::
usdus
TYPE
(
t_mpi
),
INTENT
(
IN
)
::
mpi
TYPE
(
t_atoms
),
INTENT
(
IN
)
::
atoms
TYPE
(
t_eigVecCoeffs
),
INTENT
(
IN
)
::
eigVecCoeffs
TYPE
(
t_mcd
),
INTENT
(
INOUT
)
::
mcd
! ..
! .. Scalar Arguments ..
INTEGER
,
INTENT
(
IN
)
::
noccbd
,
jsp
...
...
@@ -38,14 +39,11 @@ CONTAINS
LOGICAL
,
INTENT
(
IN
)
::
l_mcd
,
l_evp
! ..
! .. Array Arguments ..
INTEGER
,
INTENT
(
IN
)
::
ncore
(
atoms
%
ntype
)
REAL
,
INTENT
(
IN
)
::
eig
(:)
!(dimension%neigd),
REAL
,
INTENT
(
IN
)
::
we
(
noccbd
)
COMPLEX
,
INTENT
(
IN
)
::
m_mcd
(:,:,:,:)
!(dimension%nstd,(3+1)**2,3*ntypd ,2)
REAL
,
INTENT
(
INOUT
)
::
enerlo
(
atoms
%
nlod
,
atoms
%
ntype
),
sqlo
(
atoms
%
nlod
,
atoms
%
ntype
)
REAL
,
INTENT
(
INOUT
)
::
ener
(
0
:
3
,
atoms
%
ntype
),
sqal
(
0
:
3
,
atoms
%
ntype
)
REAL
,
INTENT
(
INOUT
)
::
qal
(
0
:,:,:)
!(0:3,atoms%ntype,dimension%neigd)
REAL
,
INTENT
(
INOUT
)
::
mcd
(:,:,:)
!(3*atoms%ntype,dimension%nstd,dimension%neigd)
! ..
! .. Local Scalars ..
...
...
@@ -64,7 +62,7 @@ CONTAINS
IF
((
ikpt
.LE.
mpi
%
isize
)
.AND..NOT.
l_evp
)
THEN
IF
(
l_mcd
)
THEN
mcd
(:,:,:)
=
0.0
mcd
%
mcd
(:,:,:)
=
0.0
ENDIF
ener
(:,:)
=
0.0
sqal
(:,:)
=
0.0
...
...
@@ -101,14 +99,14 @@ CONTAINS
sumab
=
sumab
+
eigVecCoeffs
%
acof
(
i
,
lm
,
natom
,
jsp
)
*
CONJG
(
eigVecCoeffs
%
bcof
(
i
,
lm
,
natom
,
jsp
))
sumba
=
sumba
+
eigVecCoeffs
%
bcof
(
i
,
lm
,
natom
,
jsp
)
*
CONJG
(
eigVecCoeffs
%
acof
(
i
,
lm
,
natom
,
jsp
))
ENDDO
DO
icore
=
1
,
ncore
(
n
)
DO
icore
=
1
,
mcd
%
ncore
(
n
)
DO
ipol
=
1
,
3
index
=
3
*
(
n
-1
)
+
ipol
mcd
(
index
,
icore
,
i
)
=
mcd
(
index
,
icore
,
i
)
+
fac
*
(&
suma
*
CONJG
(
m
_mcd
(
icore
,
lm
+1
,
index
,
1
))
*
m_mcd
(
icore
,
lm
+1
,
index
,
1
)
+
&
sumb
*
CONJG
(
m
_mcd
(
icore
,
lm
+1
,
index
,
2
))
*
m_mcd
(
icore
,
lm
+1
,
index
,
2
)
+
&
sumab
*
CONJG
(
m
_mcd
(
icore
,
lm
+1
,
index
,
2
))
*
m_mcd
(
icore
,
lm
+1
,
index
,
1
)
+
&
sumba
*
CONJG
(
m
_mcd
(
icore
,
lm
+1
,
index
,
1
))
*
m_mcd
(
icore
,
lm
+1
,
index
,
2
)
)
mcd
%
mcd
(
index
,
icore
,
i
)
=
mcd
%
mcd
(
index
,
icore
,
i
)
+
fac
*
(&
suma
*
CONJG
(
m
cd
%
m_mcd
(
icore
,
lm
+1
,
index
,
1
))
*
mcd
%
m_mcd
(
icore
,
lm
+1
,
index
,
1
)
+
&
sumb
*
CONJG
(
m
cd
%
m_mcd
(
icore
,
lm
+1
,
index
,
2
))
*
mcd
%
m_mcd
(
icore
,
lm
+1
,
index
,
2
)
+
&
sumab
*
CONJG
(
m
cd
%
m_mcd
(
icore
,
lm
+1
,
index
,
2
))
*
mcd
%
m_mcd
(
icore
,
lm
+1
,
index
,
1
)
+
&
sumba
*
CONJG
(
m
cd
%
m_mcd
(
icore
,
lm
+1
,
index
,
1
))
*
mcd
%
m_mcd
(
icore
,
lm
+1
,
index
,
2
)
)
ENDDO
ENDDO
ENDIF
! end MCD
...
...
dos/doswrite.f90
View file @
1737a5bf
...
...
@@ -11,13 +11,8 @@ MODULE m_doswrite
!-- now read data from tmp_dos and write to vacdos&dosinp .. dw
!
CONTAINS
SUBROUTINE
doswrite
(&
&
eig_id
,
DIMENSION
,
kpts
,
atoms
,
vacuum
,&
&
input
,
banddos
,&
&
sliceplot
,
noco
,
sym
,&
&
cell
,&
&
l_mcd
,
ncored
,
ncore
,
e_mcd
,&
&
efermi
,
bandgap
,
nsld
,
oneD
)
SUBROUTINE
doswrite
(
eig_id
,
DIMENSION
,
kpts
,
atoms
,
vacuum
,
input
,
banddos
,&
sliceplot
,
noco
,
sym
,
cell
,
mcd
,
ncored
,
results
,
nsld
,
oneD
)
USE
m_eig66_io
,
ONLY
:
read_dos
,
read_eig
USE
m_evaldos
USE
m_cdninf
...
...
@@ -35,19 +30,13 @@ CONTAINS
TYPE
(
t_cell
),
INTENT
(
IN
)
::
cell
TYPE
(
t_kpts
),
INTENT
(
IN
)
::
kpts
TYPE
(
t_atoms
),
INTENT
(
IN
)
::
atoms
TYPE
(
t_mcd
),
INTENT
(
IN
)
::
mcd
TYPE
(
t_results
),
INTENT
(
IN
)
::
results
! .. Scalar Arguments ..
INTEGER
,
PARAMETER
::
n2max
=
13
INTEGER
,
INTENT
(
IN
)
::
nsld
,
eig_id
INTEGER
,
INTENT
(
IN
)
::
ncored
REAL
,
INTENT
(
IN
)
::
efermi
,
bandgap
LOGICAL
,
INTENT
(
IN
)
::
l_mcd
! ..
! .. Array Arguments ..
INTEGER
,
INTENT
(
IN
)
::
ncore
(
atoms
%
ntype
)
REAL
,
INTENT
(
IN
)
::
e_mcd
(
atoms
%
ntype
,
input
%
jspins
,
ncored
)
!-odim
!+odim
! locals
INTEGER
::
jsym
(
DIMENSION
%
neigd
),
ksym
(
DIMENSION
%
neigd
)
...
...
@@ -105,20 +94,13 @@ CONTAINS
ENDIF
DO
ikpt
=
1
,
kpts
%
nkpt
call
read_eig
(
eig_id
,
ikpt
,
kspin
,&
neig
=
ne
,
eig
=
eig
)
call
read_dos
(
eig_id
,
ikpt
,
kspin
,&
&
qal
(:,:,:,
kspin
),
qvac
(:,:,
ikpt
,
kspin
),&
&
qis
(:,
ikpt
,
kspin
),&
&
qvlay
(:,:,:),
qstars
,
ksym
,
jsym
)
call
read_eig
(
eig_id
,
ikpt
,
kspin
,
neig
=
ne
,
eig
=
eig
)
call
read_dos
(
eig_id
,
ikpt
,
kspin
,
qal
(:,:,:,
kspin
),
qvac
(:,:,
ikpt
,
kspin
),&
qis
(:,
ikpt
,
kspin
),
qvlay
(:,:,:),
qstars
,
ksym
,
jsym
)
CALL
cdninf
(&
&
input
,
sym
,
noco
,
kspin
,
atoms
,&
&
vacuum
,
sliceplot
,
banddos
,
ikpt
,
kpts
%
bk
(:,
ikpt
),&
&
kpts
%
wtkpt
(
ikpt
),
cell
,
kpts
,&
&
ne
,
eig
,
qal
(
0
:,:,:,
kspin
),
qis
,
qvac
,&
&
qvlay
(:,:,:),&
&
qstars
,
ksym
,
jsym
)
CALL
cdninf
(
input
,
sym
,
noco
,
kspin
,
atoms
,
vacuum
,
sliceplot
,
banddos
,
ikpt
,
kpts
%
bk
(:,
ikpt
),&
kpts
%
wtkpt
(
ikpt
),
cell
,
kpts
,
ne
,
eig
,
qal
(
0
:,:,:,
kspin
),
qis
,
qvac
,&
qvlay
(:,:,:),
qstars
,
ksym
,
jsym
)
ENDDO
ENDDO
! end spin loop (kspin = 1,input%jspins)
...
...
@@ -135,10 +117,9 @@ CONTAINS
!
!
IF
(
banddos
%
dos
.AND.
(
banddos
%
ndir
.LT.
0
))
THEN
CALL
evaldos
(&
&
eig_id
,
input
,
banddos
,
vacuum
,
kpts
,
atoms
,
sym
,
noco
,
oneD
,
cell
,&
&
DIMENSION
,
efermi
,
bandgap
,&
&
l_mcd
,
ncored
,
ncore
,
e_mcd
,
nsld
)
CALL
evaldos
(
eig_id
,
input
,
banddos
,
vacuum
,
kpts
,
atoms
,
sym
,
noco
,
oneD
,
cell
,&
DIMENSION
,
results
%
ef
,
results
%
bandgap
,&
banddos
%
l_mcd
,
ncored
,
mcd
%
ncore
,
mcd
%
e_mcd
,
nsld
)
ENDIF
!
! Now write to vacwave if nstm=3
...
...
orbdep/mcd_init.f90
View file @
1737a5bf
MODULE
m_mcdinit
CONTAINS
SUBROUTINE
mcd_init
(
atoms
,
input
,
DIMENSION
,&
vr
,
g
,
f
,
emcd_up
,
emcd_lo
,
itype
,
jspin
,
ncore
,
e_mcd
,
m_mcd
)
SUBROUTINE
mcd_init
(
atoms
,
input
,
DIMENSION
,
vr
,
g
,
f
,
mcd
,
itype
,
jspin
)
!-----------------------------------------------------------------------
!
...
...
@@ -22,6 +21,7 @@ CONTAINS
TYPE
(
t_dimension
),
INTENT
(
IN
)
::
DIMENSION
TYPE
(
t_input
),
INTENT
(
IN
)
::
input
TYPE
(
t_atoms
),
INTENT
(
IN
)
::
atoms
TYPE
(
t_mcd
),
INTENT
(
INOUT
)
::
mcd
INTEGER
,
PARAMETER
::
l_max
=
3
...
...
@@ -29,13 +29,9 @@ CONTAINS
INTEGER
,
INTENT
(
IN
)
::
itype
INTEGER
,
INTENT
(
IN
)
::
jspin
REAL
,
INTENT
(
IN
)
::
emcd_up
,
emcd_lo
REAL
,
INTENT
(
IN
)
::
vr
(
atoms
%
jmtd
,
atoms
%
ntype
,
input
%
jspins
)
REAL
,
INTENT
(
IN
)
::
f
(
atoms
%
jmtd
,
2
,
0
:
atoms
%
lmaxd
,
jspin
:
jspin
)
REAL
,
INTENT
(
IN
)
::
g
(
atoms
%
jmtd
,
2
,
0
:
atoms
%
lmaxd
,
jspin
:
jspin
)
INTEGER
,
INTENT
(
OUT
)
::
ncore
(
atoms
%
ntype
)
REAL
,
INTENT
(
OUT
)
::
e_mcd
(
atoms
%
ntype
,
input
%
jspins
,
DIMENSION
%
nstd
)
COMPLEX
,
INTENT
(
OUT
)
::
m_mcd
(:,:,:,:)
! Locals ...
...
...
@@ -55,7 +51,7 @@ CONTAINS
! core setup
ncore
(
itype
)
=
0
mcd
%
ncore
(
itype
)
=
0
bmu
=
0.0
CALL
setcor
(
itype
,
1
,
atoms
,
input
,
bmu
,
nst
,
kappa
,
nprnc
,
occ
)
...
...
@@ -85,7 +81,7 @@ CONTAINS
CALL
differ
(
fn
,
fl
,
fj
,
c
,
atoms
%
zatom
(
itype
),
atoms
%
dx
(
itype
),
atoms
%
rmsh
(
1
,
itype
),&
rn
,
d
,
DIMENSION
%
msh
,
vrd
,
e
,
a
,
b
,
ierr
)
IF
(
ierr
/
=
0
)
CALL
juDFT_error
(
"error in core-levels"
,
calledby
=
"mcd_init"
)
IF
(
(
e
.LE.
emcd_up
)
.AND.
(
e
.GE.
emcd_lo
)
)
THEN
IF
(
(
e
.LE.
mcd
%
emcd_up
)
.AND.
(
e
.GE.
mcd
%
emcd_lo
)
)
THEN
WRITE
(
*
,
*
)
'good ev = '
,
e
n_core
=
n_core
+
1
j_core
(
n_core
)
=
fj
...
...
@@ -112,7 +108,7 @@ CONTAINS
DO
iri
=
3
*
(
itype
-1
)
+1
,
3
*
(
itype
-1
)
+3
DO
l
=
1
,
(
l_max
+1
)
**
2
DO
icore
=
1
,
DIMENSION
%
nstd
m_mcd
(
icore
,
l
,
iri
,
i
)
=
CMPLX
(
0.0
,
0.0
)
m
cd
%
m
_mcd
(
icore
,
l
,
iri
,
i
)
=
CMPLX
(
0.0
,
0.0
)
ENDDO
ENDDO
ENDDO
...
...
@@ -145,13 +141,13 @@ CONTAINS
! write(*,*) j_core(icore),l_core(icore),l_max,ms
CALL
nabla
(
itype
,
icore
,
atoms
%
jri
(
itype
),
atoms
%
dx
(
itype
),
DIMENSION
%
nstd
,
atoms
%
ntype
,&
j_core
(
icore
),
l_core
(
icore
),
l_max
,
ms
,
atoms
%
rmsh
(:,
itype
),
gc
(:,
icore
,
ispin
),&
gv
(:,
0
:,
ispin
,
i
),
dgv
(:,
0
:,
ispin
,
i
),
m_mcd
(:,:,:,
i
)
)
gv
(:,
0
:,
ispin
,
i
),
dgv
(:,
0
:,
ispin
,
i
),
m
cd
%
m
_mcd
(:,:,:,
i
)
)
ENDDO
DO
i
=
1
,
2
*
icore
*
l_core
(
icore
)
ncore
(
itype
)
=
ncore
(
itype
)
+
1
IF
(
ncore
(
itype
)
>
DIMENSION
%
nstd
)
CALL
juDFT_error
(
"dimension%nstd too small"
,
calledby
=
"mcd_init"
)
e_mcd
(
itype
,
ispin
,
ncore
(
itype
))
=
e_mcd1
(
icore
)
mcd
%
ncore
(
itype
)
=
mcd
%
ncore
(
itype
)
+
1
IF
(
mcd
%
ncore
(
itype
)
>
DIMENSION
%
nstd
)
CALL
juDFT_error
(
"dimension%nstd too small"
,
calledby
=
"mcd_init"
)
mcd
%
e_mcd
(
itype
,
ispin
,
mcd
%
ncore
(
itype
))
=
e_mcd1
(
icore
)
ENDDO
ENDDO
ENDDO
...
...
@@ -164,8 +160,8 @@ CONTAINS
! DO i = 1, 2
! DO iri = 3*(itype-1)+1 , 3*(itype-1)+3
! write (*,*) iri
! DO icore = 1, ncore(itype)
! write (*,'(10f10.5)') (m_mcd(icore,l,iri,i),l=1,9)
! DO icore = 1,
mcd%
ncore(itype)
! write (*,'(10f10.5)') (m
cd%m
_mcd(icore,l,iri,i),l=1,9)
! ENDDO
! ENDDO
! ENDDO
...
...
types/types_cdnval.f90
View file @
1737a5bf
...
...
@@ -142,7 +142,19 @@ PRIVATE
PROCEDURE
,
PASS
::
init
=>
eigVecCoeffs_init
END
TYPE
t_eigVecCoeffs
PUBLIC
t_orb
,
t_denCoeffs
,
t_denCoeffsOffdiag
,
t_force
,
t_slab
,
t_eigVecCoeffs
TYPE
t_mcd
REAL
::
emcd_lo
,
emcd_up
INTEGER
,
ALLOCATABLE
::
ncore
(:)
REAL
,
ALLOCATABLE
::
e_mcd
(:,:,:)
REAL
,
ALLOCATABLE
::
mcd
(:,:,:)
COMPLEX
,
ALLOCATABLE
::
m_mcd
(:,:,:,:)
CONTAINS
PROCEDURE
,
PASS
::
init1
=>
mcd_init1
END
TYPE
t_mcd
PUBLIC
t_orb
,
t_denCoeffs
,
t_denCoeffsOffdiag
,
t_force
,
t_slab
,
t_eigVecCoeffs
,
t_mcd
CONTAINS
...
...
@@ -534,4 +546,36 @@ SUBROUTINE eigVecCoeffs_init(thisEigVecCoeffs,dimension,atoms,noco,jspin,noccbd)
END
SUBROUTINE
eigVecCoeffs_init
SUBROUTINE
mcd_init1
(
thisMCD
,
banddos
,
dimension
,
input
,
atoms
)
USE
m_types_setup
IMPLICIT
NONE
CLASS
(
t_mcd
),
INTENT
(
INOUT
)
::
thisMCD
TYPE
(
t_banddos
),
INTENT
(
IN
)
::
banddos
TYPE
(
t_dimension
),
INTENT
(
IN
)
::
dimension
TYPE
(
t_input
),
INTENT
(
IN
)
::
input
TYPE
(
t_atoms
),
INTENT
(
IN
)
::
atoms
ALLOCATE
(
thisMCD
%
ncore
(
atoms
%
ntype
))
ALLOCATE
(
thisMCD
%
e_mcd
(
atoms
%
ntype
,
input
%
jspins
,
dimension
%
nstd
))
IF
(
banddos
%
l_mcd
)
THEN
thisMCD
%
emcd_lo
=
banddos
%
e_mcd_lo
thisMCD
%
emcd_up
=
banddos
%
e_mcd_up
ALLOCATE
(
thisMCD
%
m_mcd
(
dimension
%
nstd
,(
3+1
)
**
2
,
3
*
atoms
%
ntype
,
2
))
ALLOCATE
(
thisMCD
%
mcd
(
3
*
atoms
%
ntype
,
dimension
%
nstd
,
dimension
%
neigd
)
)
IF
(
.NOT.
banddos
%
dos
)
WRITE
(
*
,
*
)
'For mcd-spectra set banddos%dos=T!'
ELSE
ALLOCATE
(
thisMCD
%
m_mcd
(
1
,
1
,
1
,
1
))
ALLOCATE
(
thisMCD
%
mcd
(
1
,
1
,
1
))
ENDIF
thisMCD
%
e_mcd
=
0.0
thisMCD
%
mcd
=
0.0
thisMCD
%
m_mcd
=
CMPLX
(
0.0
,
0.0
)
END
SUBROUTINE
mcd_init1
END
MODULE
m_types_cdnval
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