The Elk Code
 
Loading...
Searching...
No Matches
rdmwritededn.f90
Go to the documentation of this file.
1
2! Copyright (C) 2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: rdmwritededn
8! !INTERFACE:
9subroutine rdmwritededn(dedn)
10! !USES:
11use modmain
12use modrdm
13! !INPUT/OUTPUT PARAMETERS:
14! dedn : derivative of energy (in,real(nstsv,nkpt))
15! !DESCRIPTION:
16! Writes the derivative of total energy with respect to occupation numbers to
17! file {\tt RDM\_DEDN.OUT}.
18!
19! !REVISION HISTORY:
20! Created 2008 (Sharma)
21!EOP
22!BOC
23implicit none
24! arguments
25real(8), intent(in) :: dedn(nstsv,nkpt)
26! local variables
27integer ik,ist
28open(50,file='RDM_DEDN.OUT',form='FORMATTED')
29write(50,'(I6," : nkpt")') nkpt
30write(50,'(I6," : nstsv")') nstsv
31do ik=1,nkpt
32 write(50,*)
33 write(50,'(I6,3G18.10," : k-point, vkl")') ik,vkl(:,ik)
34 write(50,'(" (state, occupancy and derivative below)")')
35 do ist=1,nstsv
36 write(50,'(I6,4G18.10)') ist,occsv(ist,ik),-dedn(ist,ik)
37 end do
38end do
39close(50)
40end subroutine
41!EOC
42
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
subroutine rdmwritededn(dedn)