2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
9 subroutine plot2d(fnum
,nf
,lmax
,ld
,rfmt
,rfir
)
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! fnum : plot file number (in,integer)
14 ! nf : number of functions (in,integer)
15 ! lmax : maximum angular momentum (in,integer)
16 ! ld : leading dimension (in,integer)
17 ! rfmt : real muffin-tin function (in,real(ld,nrmtmax,natmtot,nf))
18 ! rfir : real intersitial function (in,real(ngrtot,nf))
20 ! Produces a 2D plot of the real functions contained in arrays {\tt rfmt} and
21 ! {\tt rfir} on the parallelogram defined by the corner vertices in the global
22 ! array {\tt vclp2d}. See routine {\tt rfarray}.
25 ! Created June 2003 (JKD)
30 integer, intent(in
) :: fnum
31 integer, intent(in
) :: nf
32 integer, intent(in
) :: lmax
33 integer, intent(in
) :: ld
34 real(8), intent(in
) :: rfmt(ld
,nrmtmax
,natmtot
,nf
)
35 real(8), intent(in
) :: rfir(ngrtot
,nf
)
38 real(8) vl1(3),vl2(3),vc1(3),vc2(3)
39 real(8) d1
,d2
,d12
,t1
,t2
,t3
,t4
41 real(8), allocatable
:: vpl(:,:)
42 real(8), allocatable
:: fp(:,:)
43 if ((nf
.lt
.1).or
.(nf
.gt
.4)) then
45 write(*,'("Error(plot2d): invalid number of functions : ",I8)') nf
49 ! allocate local arrays
50 allocate(vpl(3,np2d(1)*np2d(2)))
51 allocate(fp(np2d(1)*np2d(2),nf
))
53 vl1(:)=vclp2d(:,2)-vclp2d(:,1)
54 vl2(:)=vclp2d(:,3)-vclp2d(:,1)
55 vc1(:)=vl1(1)*avec(:,1)+vl1(2)*avec(:,2)+vl1(3)*avec(:,3)
56 vc2(:)=vl2(1)*avec(:,1)+vl2(2)*avec(:,2)+vl2(3)*avec(:,3)
57 d1
=sqrt(vc1(1)**2+vc1(2)**2+vc1(3)**2)
58 d2
=sqrt(vc2(1)**2+vc2(2)**2+vc2(3)**2)
59 if ((d1
.lt
.epslat
).or
.(d2
.lt
.epslat
)) then
61 write(*,'("Error(plot2d): zero length plotting vectors")')
65 d12
=(vc1(1)*vc2(1)+vc1(2)*vc2(2)+vc1(3)*vc2(3))/(d1
*d2
)
70 t1
=dble(ip1
)/dble(np2d(1))
71 t2
=dble(ip2
)/dble(np2d(2))
72 vpl(:,ip
)=t1
*vl1(:)+t2
*vl2(:)+vclp2d(:,1)
75 ! evaluate the functions at the grid points
77 call rfarray(lmax
,ld
,rfmt(:,:,:,i
),rfir(:,i
),ip
,vpl
,fp(:,i
))
79 ! write the functions to file
80 write(fnum
,'(2I6," : grid size")') np2d(:)
85 t1
=dble(ip1
)/dble(np2d(1))
86 t2
=dble(ip2
)/dble(np2d(2))
88 t4
=t2
*d2
*sqrt(abs(1.d0
-d12
**2))
89 write(fnum
,'(6G18.10)') t3
,t4
,(fp(ip
,i
),i
=1,nf
)