The Elk Code
Loading...
Searching...
No Matches
writelambda.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2010 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
subroutine
writelambda
(gq)
7
use
modmain
8
use
modphonon
9
implicit none
10
! arguments
11
real
(8),
intent(in)
:: gq(nbph,nqpt)
12
! local variables
13
integer
iq,i
14
real
(8) t1,t2
15
open
(50,file=
'LAMBDAQ.OUT'
,form=
'FORMATTED'
)
16
write
(50,*)
17
write
(50,
'(I4," : total number of atoms")'
)
natmtot
18
write
(50,
'(I6," : number of q-points")'
) nqpt
19
write
(50,*)
20
do
iq=1,nqpt
21
write
(50,
'(I6," : q-point")'
) iq
22
write
(50,
'(3G18.10," : q-vector (lattice coordinates)")'
)
vql
(:,iq)
23
write
(50,
'(3G18.10," : q-vector (Cartesian coordinates)")'
)
vqc
(:,iq)
24
do
i=1,nbph
25
t1=
pi
*
fermidos
*
wphq
(i,iq)**2
26
if
(t1 > 1.d-8)
then
27
t2=gq(i,iq)/t1
28
else
29
t2=0.d0
30
end if
31
write
(50,
'(I4,G18.10)'
) i,t2
32
end do
33
write
(50,*)
34
end do
35
close
(50)
36
write
(*,*)
37
write
(*,
'("Info(writelambda):")'
)
38
write
(*,
'(" wrote electron-phonon coupling constants for all q-points to &
39
&LAMBDAQ.OUT")'
)
40
end subroutine
41
modmain
Definition
modmain.f90:6
modmain::pi
real(8), parameter pi
Definition
modmain.f90:1229
modmain::vqc
real(8), dimension(:,:), allocatable vqc
Definition
modmain.f90:547
modmain::natmtot
integer natmtot
Definition
modmain.f90:40
modmain::vql
real(8), dimension(:,:), allocatable vql
Definition
modmain.f90:545
modmain::fermidos
real(8) fermidos
Definition
modmain.f90:910
modphonon
Definition
modphonon.f90:6
modphonon::wphq
real(8), dimension(:,:), allocatable wphq
Definition
modphonon.f90:31
writelambda
subroutine writelambda(gq)
Definition
writelambda.f90:7
writelambda.f90
Generated by
1.9.8