Source module last modified on Tue, 9 May 2006, 20:05;
HTML image of Fortran source automatically generated by
for2html on Mon, 29 May 2006, 15:06.
#//////////////////////////////////////////////////////////////////////
#//// $Id: res_cal.f,v 1.5 2006/05/09 18:05:06 saroun Exp $
#////
#//// R E S T R A X 4.4
#////
#//// Transformations between coordinate systems etc.
#////
#////
#//////////////////////////////////////////////////////////////////////
#--------------------------------------------------------------------
SUBROUTINE RECLAT
# Calculate transformation matrices for reciprocal lattice
#
# Input:
# AQ(3) ... unit cell base vectors
# ALFA(3) ... unit cell angles
# A1(3), A2(3) ... perpendicular vectors in scattering plane
#
# Output:
# SMAT(i,j) ... transforms Q(hkl) from r.l.u. to
# the orthogonal coordinates given by base vectors A1(i), A2(i), A3(i)
# COSB(3) ... direction cosines for unit cell base vectors
# revised 22/2/2005, J.S.
#--------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'lattice.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'exciimp.inc'
record /RECLATTICE/ rl
real*8 eps,zd,rd,xcc,c1,c2,c3
integer*4 i,j,k,l,m,ier
parameter(eps=1.0d-8)
real*8 aq(3),alfa(3),a1(3),a2(3)
real*8 a(3),xbb(3,3),v1(3),v2(3),v3(3),u(3,3),cosa(3),sina(3),sinb(3)
# REAL*8 W1(3),W2(3),W3(3),XBBI(3,3),
real*8 a3(3),aux(3,3),zn(4),dum
real*8 s(3,3),sd(3,3),b(3)
common /error/ier
# REAL*8 QxQ
equivalence (aq(1),res_dat(i_as))
equivalence (alfa(1),res_dat(i_aa))
equivalence (a1(1),res_dat(i_ax))
equivalence (a2(1),res_dat(i_bx))
# write(*,*) 'RECLAT BEFORE:'
10 format( 'res_cal.f, RECLAT, ',a,1x,10(1x,g10.4))
#1000 format(a15,3(2x,E10.3))
ier=0
zd=2.d0*pi
rd=pi/180.d0
xcc=0.d0
do i=1,3
a(i)=aq(i)/zd
if(abs(a(i)).lt..00001) then ! check lattice spacing
ier=1
write(sout,30)
return
endif
cosa(i)=cos(alfa(i)*rd)
sina(i)=sin(alfa(i)*rd)
xcc=xcc+cosa(i)**2
enddo
xcc=1.+2.*cosa(1)*cosa(2)*cosa(3)-xcc
if(xcc.le.0.) then ! check lattice angles
ier=2
write(sout,31)
return
endif
xcc=sqrt(xcc)
cellvol=xcc*aq(1)*aq(2)*aq(3) ! this is unit cell volume
j=2
k=3
do i=1,3
b(i)=sina(i)/(a(i)*xcc) ! length of the r.l. axes an AngŻ1
cosb(i)=(cosa(j)*cosa(k)-cosa(i))/(sina(j)*sina(k))
sinb(i)=sqrt(1.-cosb(i)*cosb(i)) ! angles btw. r.l. axes
j=k
k=i
enddo
# check that A,B are perpendicular
# NOTE: QxQ uses only COSB, not SMAT ...
# IF (ABS(QxQ(A1(1),A2(1))).GT.EPS) THEN
# write(sout,33)
# ENDIF
# SD(i,j) is a matrix defined so that
# Qhkl = SQRT((hkl)*SD*(hkl)) (in A^-1)
do i=1,3
sd(i,i)=b(i)**2
end do
sd(1,2)=(cosa(1)*cosa(2)-cosa(3))/aq(1)/aq(2)*(2*pi/xcc)**2
sd(1,3)=(cosa(1)*cosa(3)-cosa(2))/aq(1)/aq(3)*(2*pi/xcc)**2
sd(2,3)=(cosa(2)*cosa(3)-cosa(1))/aq(2)/aq(3)*(2*pi/xcc)**2
sd(2,1)=sd(1,2)
sd(3,1)=sd(1,3)
sd(3,2)=sd(2,3)
# XBB(i,j) are projections of rec. lat. base vectors on the orthonormal base:
# Let (a,b,c) and (a*,b*,c*) are direct and reciprocal lattice base vectors
# assume a parallel to a*
# XBB(i,1) ... parallel to a*
# XBB(i,2) ... parallel to c x a*
# XBB(i,3) ... parallel to c
# i.e. the columns are a*,b*,c* in the new orthonormal base
xbb(1,1)=b(1)
xbb(2,1)=0.d0
xbb(3,1)=0.d0
xbb(1,2)=b(2)*cosb(3)
xbb(2,2)=b(2)*sinb(3)
xbb(3,2)=0.d0
xbb(1,3)=b(3)*cosb(2)
xbb(2,3)=-b(3)*sinb(2)*cosa(1)
xbb(3,3)=1/a(3)
# CALL INVERT(3,XBB,3,XBBI,3)
# convert A1,A2 to the new orthonormal system:
do i=1,3
v1(i)=0.d0
v2(i)=0.d0
do j=1,3
v1(i)=v1(i)+xbb(i,j)*a1(j)
v2(i)=v2(i)+xbb(i,j)*a2(j)
enddo
enddo
# get V3 perpendicular to V1,V2
v3(1)=v1(2)*v2(3)-v1(3)*v2(2)
v3(2)=v1(3)*v2(1)-v1(1)*v2(3)
v3(3)=v1(1)*v2(2)-v1(2)*v2(1)
# get V2 perpendicular to V3,V1
v2(1)=v3(2)*v1(3)-v3(3)*v1(2)
v2(2)=v3(3)*v1(1)-v3(1)*v1(3)
v2(3)=v3(1)*v1(2)-v3(2)*v1(1)
# get norms of V1,V2,V3
c1=v1(1)**2+v1(2)**2+v1(3)**2
c2=v2(1)**2+v2(2)**2+v2(3)**2
c3=v3(1)**2+v3(2)**2+v3(3)**2
c1=sqrt(c1)
c2=sqrt(c2)
c3=sqrt(c3)
# convert V1,V2 back to r.l.:
# DO I=1,3
# W1(I)=0.D0
# W2(I)=0.D0
# W3(I)=0.D0
# DO J=1,3
# W1(I)=W1(I)+XBBI(I,J)*V1(J)
# W2(I)=W2(I)+XBBI(I,J)*V2(J)
# W3(I)=W3(I)+XBBI(I,J)*V3(J)
# ENDDO
# ENDDO
if (c1*c2*c3.lt.eps) then ! check scattering plane
ier=3
write(sout,32)
return
endif
# U(i,j) is the orthonormal system made from V1,V2,V3 (also called AB system in RESTRAX)
do i=1,3
u(1,i)=v1(i)/c1
u(2,i)=v2(i)/c2
u(3,i)=v3(i)/c3
enddo
do k=1,3
do m=1,3
s(k,m)=0.d0
do l=1,3
s(k,m)=s(k,m)+u(k,l)*xbb(l,m)
enddo
enddo
enddo
# SMAT converts Q from r.l.u. to AB coordinates (in A^-1)
# SINV is inverse to SMAT
# the matrices are 4-dimensional so that they can operate in (h,k,l,energy) space
call INVERT(3,s,3,aux,3)
do i=1,4
do j=1,4
if(i.le.3.and.j.le.3) then
smat(i,j)=s(i,j)
sinv(i,j)=aux(i,j)
else if (i.eq.4.and.j.eq.4) then
smat(i,j)=1.d0
sinv(i,j)=1.d0
else
smat(i,j)=0.d0
sinv(i,j)=0.d0
endif
enddo
enddo
# write(*,10) 'SMAT: ',SMAT(1:3,1)
# write(*,10) 'SMAT: ',SMAT(1:3,2)
# write(*,10) 'SMAT: ',SMAT(1:3,3)
# write(*,10) 'SINV: ',SINV(1:3,1)
# write(*,10) 'SINV: ',SINV(1:3,2)
# write(*,10) 'SINV: ',SINV(1:3,3)
# pause
#// MABR converts Q from r.l.u. to A,B coordinates, taking A,B as base vectors
#// MRAB is inverse to MABR
do i=1,3
a3(i)=sinv(i,3)*c1
enddo
call QNORM(a1,dum,zn(1)) ! get norms of vectors A,B,C in Ang^-1
call QNORM(a2,dum,zn(2))
call QNORM(a3,dum,zn(3))
zn(4)=1.d0
do i=1,4
do j=1,4
mrab(i,j)=sinv(i,j)*zn(j)
mabr(i,j)=smat(i,j)/zn(i)
enddo
enddo
# pass rec. lattice parameters to EXCI
do i=1,3
rl.cell(i)=aq(i)
rl.cell(i+3)=alfa(i)
rl.veca(i)=a1(i)
rl.vecb(i)=a2(i)
rl.cosb(i)=cosb(i)
enddo
do i=1,4
do j=1,4
rl.smat(i,j)=smat(i,j)
rl.sinv(i,j)=sinv(i,j)
rl.mrab(i,j)=mrab(i,j)
rl.mabr(i,j)=mabr(i,j)
enddo
enddo
rl.cellvol=cellvol
call setreclat(rl)
return
30 format( ' RECLAT: Check Lattice Spacings (AS,BS,CS)')
31 format( ' RECLAT: Check Cell Angles (AA,BB,XCC)')
32 format( ' RECLAT: Check Scattering Plane (AX....BZ)')
33 format( ' RECLAT: Vectors A,B are not perpendicular !')
end
#--------------------------------------------------------------------
SUBROUTINE RLU2AB(a,ab)
# Converts A matrix from r.l.u. to coordinates given by
# the base vectors AX..BZ in the scat. plane (assumed orthogonal)
#--------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'lattice.inc'
real*8 a(4,4),ab(4,4)
call BTAB4(a,mrab,ab)
end
#***********************************************************************
SUBROUTINE CN2RLU(a,ar)
# Converts A from Cooper&Nathans to recip. lattice coordinates
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 a(4,4),ar(4,4)
call BTAB4(a,mcr,ar)
end
#***********************************************************************
SUBROUTINE CN2RLU_MF(a,ar,item)
# Version of CN2RLU for ITEM-th data set
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'restrax.inc'
real*8 a(4,4),ar(4,4)
integer*4 item
call BTAB4(a,mf_mcr(1,1,item),ar)
end
#--------------------------------------------------------------------
SUBROUTINE GET_SCANGLES(ki,kf,q,ss,om,psi,ier)
# Calculates scattering angle (OM) and angle between Q and KI (PSI)
# Input is KI,KF,Q,SS, which determine the scattering triangle
# IER=2 ... error in getting angle, cannot close triangle
#--------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
real*8 ki,kf,q,ss,om,psi
integer*4 ier
real*8 co
ier=2
if (abs(q).lt.1e-6) then
om=0.
psi=pi/2.
ier=0
return
endif
#/// calculates scattering angle omega=2*thetaS
co=(ki**2+kf**2-q**2)/(2*ki*kf)
if(abs(co).gt.1) goto 99
om=sign(1.d0,ss)*abs(acos(co))
#/// calculates scattering angle omega=2*thetaS
co=(kf**2-ki**2-q**2)/(2*ki*q)
if(abs(co).gt.1) goto 99
psi=sign(1.d0,ss)*abs(acos(co))
ier=0
return
99 return
end
#-----------------------------------------------------------------------------
SUBROUTINE ANGSCAN(da3,da4)
# set a scan in QHKL,E equivalent to angular scan (small range) by da3 or da4
#-----------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 a3,a4,da3,da4,vq(4),wq(4),si,co
integer*4 i
1 format( ' (DH,DK,DL,EN) = ',4(2x,f7.4))
2 format( ' VQ = ',4(2x,f7.4))
# write(*,*) 'AGSCAN: ',DA3,DA4
if(da3.ne.0.or.da4.ne.0) then
a4=da4*pi/180
si=sin(a4)
co=cos(a4)
# VQ is the change in QHKL due to rotation by A4, Lab. coordinates (z || ki, y vertical)
vq(1)=kf0*(comega*si+somega*(co-1.d0))
vq(3)=kf0*(-somega*si+comega*(co-1.d0))
vq(2)=0.
vq(4)=0.
# write(*,2) (VQ(I), I=1,4)
call MXV(-1,4,4,mlc,vq,wq) ! from Lab to CN coord.
# write(*,2) (VQ(I), I=1,4)
# add an increment due to sample rotation (dA3)
a3=da3*pi/180
si=sin(a3)
co=cos(a3)
wq(1)=wq(1)+q0*(co-1.d0)
wq(2)=wq(2)-q0*si
# write(*,2) (VQ(I), I=1,4)
# convert from CN to r.l.u.
call MXV(1,4,4,mrc,wq,vq)
# write(*,2) (VQ(I), I=1,4)
do i=1,4
res_dat(i_dqh+i-1)=vq(i)
enddo
endif
end
#*********************************************************************************
SUBROUTINE ROTA3(qi,a3,qf)
# transform QI(3) to QF(3) through the rotation by axis A3
# Q are given in r.l.u.
# WARNING! QI is REAL*4 !!
#*********************************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'lattice.inc'
integer*4 i,j
real*8 a3,qf(3)
real*4 qi(3)
real*8 da3 ,vq(3),wq(2),si,co
# transform from r.l.u to AX..BZ
do j=1,3
vq(j)=0.d0
enddo
do i=1,3
do j=1,3
vq(j)=vq(j)+smat(j,i)*qi(i)
enddo
enddo
# rotate by A3 around vertical axis (i=3)
da3=a3*deg
si=sin(da3)
co=cos(da3)
wq(1)=vq(1)*co+vq(2)*si
wq(2)=vq(2)*co-vq(1)*si
# transform back from AX..BZ to r.l.u
do i=1,3
# QF(I)=SINV(1,I)*WQ(1)+SINV(2,I)*WQ(2)
qf(i)=sinv(i,1)*wq(1)+sinv(i,2)*wq(2)+sinv(i,3)*vq(3)
enddo
end
#------------------------------------------------------------------------------------
SUBROUTINE ROTA4(a4,q)
# Increment scattering angle by A4 [rad]
# Only changes QHKL, DOESN'T UPDATE DEPENDENT FIELDS !!
#------------------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 a4,vq(4),wq(4),q(3),si,co
integer*4 i
if(a4.ne.0.) then
si=sin(a4)
co=cos(a4)
#// VQ is the change of Q in Lab. coordinates
vq(1)=kf0*(comega*si+somega*(co-1.0))
vq(3)=kf0*(comega*(co-1.0)-somega*si)
vq(2)=0.
vq(4)=0.
call MXV(-1,4,4,mlc,vq,wq) ! from Lab to CN coord.
call MXV(1,4,4,mrc,wq,vq) ! from CN to r.l.u.
#// Increment QHKL
do i=1,3
res_dat(i_qh+i-1)=res_dat(i_qh+i-1)+vq(i)
enddo
# swith SS to -SS if needed
if ((omega+a4)*omega.le.0) res_dat(i_ss)=-res_dat(i_ss)
# CALL BEFORE ! don't call BEFORE
endif
end
#--------------------------------------------------------------------
SUBROUTINE GET_A3A4(id,q,a3,a4,ier)
# returns A3,A4 angles for given Q(4).
# ID-th data set is taken as reference (where A3=0)
# if IR>0, error has occured (cannot construct the triangle)
#--------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
integer*4 id,ier
real*8 q(4),a3,a4
real*8 vki,vkf,qq
real*8 kfix,kom,dum,psi,psiref,a4ref,qref
ier=1
# Get correct KI,KF,Q values
kom=q(4)/hsqov2m
call QNORM(q(1),dum,qq)
kfix=mf_par(i_kfix,id)
if (mf_par(i_fx,id).eq.1) then
vki=kfix
vkf=kfix**2-kom
if (vkf.le.0) goto 99
vkf=sqrt(vkf)
else
vkf=kfix
vki=kfix**2+kom
if (vki.le.0) goto 99
vki=sqrt(vki)
endif
#/// calculates scattering angle (A4) and angle between Ki and Q (PSI)
call QNORM(mf_par(i_qh,id),dum,qref)
call GET_SCANGLES(vki,vkf,qref,mf_par(i_ss,id),a4ref,psiref,ier)
if (ier.ne.0) goto 99
call GET_SCANGLES(vki,vkf,qq,mf_par(i_ss,id),a4,psi,ier)
if (ier.ne.0) goto 99
#11 format(a,6(1x,G12.6))
# write(*,*) 'Qref, Q: ',Qref,Q
# write(*,*) 'REF A4, PSI: ',A4ref/deg,PSIref/deg
# write(*,*) 'CUR A4, PSI: ',A4/deg,PSI/deg
#/// calculates sample rotation with respect to
#/// the reference channel (=ID)
ier=3
#// angle between Q and reference (=QHKL from IDth channel)
call GET_ANGLE(mf_par(i_qh,id),q(1),a3)
# write(*,*) 'Q vs. Qref: ',A3/deg
a3=psi-psiref-a3
# call M4xV4_3(mf_MCR(1,1,ID),Q,VQ) ! convert QHKL to C&N
# CA3=VQ(1)/QQ
# SA3=VQ(2)/QQ
# A3=SIGN(1.D0,SA3)*ACOS(CA3)
ier=0
99 return
end
#--------------------------------------------------------------------
SUBROUTINE SCATTRIANGLE
# - define scatteing triangle and associated fields
# - create transformation matrices between C&N and rec. lattice (MCR,MRC)
# - create transformation matrix Lab. coord. -> C&N
# Lab. coordinates are defined with Z//Ki, Y vertical
# W(er.lat.) = MRC*V(C&N)
# MCR = MRC^-1
# RECLAT must be called before !!!!
#--------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'lattice.inc'
real*8 si,co,kfix,kom,dum,xq,yq,qq,TRANS
integer*4 ier,i,j
character*40 mater(3)
data mater /
1 ' Check scaterring triangle ',
2 ' Check monochromator dhkl ',
3 ' Check analyzer dhkl '/
600 format(a20,4(2x,e12.6))
601 format(1x, 'KI: ',e12.6,2x, 'KF: ',e12.6,2x, 'Q: ',e12.6)
602 format(1x, 'DM: ',e12.6,2x, 'KI: ',e12.6,2x)
603 format(1x, 'DA: ',e12.6,2x, 'KF: ',e12.6,2x)
501 format(a)
# Get correct KI,KF,Q values
ier=1
# Get correct KI,KF,Q values
kom=res_dat(i_en)/hsqov2m
call QNORM(res_dat(i_qh),dum,q0)
kfix=res_dat(i_kfix)
if (res_dat(i_fx).eq.1) then
ki0=kfix
kf0=kfix**2-kom
kf0=sign(1.d0,kf0)*sqrt(abs(kf0))
else
kf0=kfix
ki0=kfix**2+kom
ki0=sign(1.d0,ki0)*sqrt(abs(ki0))
endif
if ((res_dat(i_sm).ne.0).and.(ki0.le.pi/res_dat(i_dm))) then
ier=2
goto 99
endif
if ((res_dat(i_sa).ne.0).and.(kf0.le.pi/res_dat(i_da))) then
ier=3
goto 99
endif
#/// calculates scattering angle omega=2*thetaS
comega=(ki0**2+kf0**2-q0**2)/(2*ki0*kf0)
if(abs(comega).gt.1) goto 99
if(res_dat(i_ss).gt.0) then
somega=sqrt(1-comega**2)
omega=acos(comega)
else
somega=-sqrt(1-comega**2)
omega=-acos(comega)
endif
#/// trans. matrix C&N --> Lab
do i=1,4
do j=1,4
mlc(i,j)=0.d0
enddo
enddo
mlc(4,4)=1.d0
# angle between Ki and Q
co=(kf0**2-ki0**2-q0**2)/(2*ki0*q0)
if(abs(co).gt.1) goto 99 ! should never happen !
if(res_dat(i_ss).gt.0) then
si=sqrt(1.d0-co**2)
else
si=-sqrt(1.d0-co**2)
endif
mlc(1,1)=si
mlc(1,2)=co
mlc(2,3)=1.d0
mlc(3,1)=co
mlc(3,2)=-si
#/// trans. matrix r.l.u. --> C&N
xq=TRANS(qhkl,1)
yq=TRANS(qhkl,2)
qq=sqrt(xq**2+yq**2)
co=xq/qq
si=yq/qq
do 11 i=1,3
mcr(1,i)=smat(1,i)*co+smat(2,i)*si
mcr(2,i)=smat(2,i)*co-smat(1,i)*si
mcr(3,i)=smat(3,i)
mcr(4,i)=0.
11 mcr(i,4)=0.
mcr(4,4)=1. ! (MCR): r.l.u. --> C&N
call INVERT(4,mcr,4,mrc,4) ! (MRC): C&N --> r.l.u.
ier=0
return
99 continue
write(sout,501) mater(ier)
if (ier.eq.1) write(sout,601) ki0,kf0,q0
if (ier.eq.2) write(sout,602) res_dat(i_dm),ki0
if (ier.eq.3) write(sout,603) res_dat(i_da),kf0
return
end
#--------------------------------------------------------------------
SUBROUTINE TRANSMAT
# creates transformation matrices for the four coordinate systems:
# R ... reciprocal lattice
# C ... Cooper & Nathans coordinates (C&N)
# G ... C&N rotated so that X(1) // grad E(Q)
# D ... G rotated so that X(4) // normal to the disp. surface
#--------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 mgd(4,4),phi,zra,gnr,gna
integer*4 i,j
call GETMAT(grd,mcg)
# (MGD): grad <-- disp
do i=1,4
do j=1,4
mgd(i,j)=0
enddo
mgd(i,i)=1
enddo
call QNORM(grd,gnr,gna)
zra=gnr/gna
phi=atan(res_dat(i_gmod)*zra) ! units of GMOD are Energy/r.l.u.!
mgd(1,1)=cos(phi)
mgd(4,4)=cos(phi)
mgd(1,4)=-sin(phi)
mgd(4,1)=sin(phi)
# (MCD): C&N <-- disp
call MXM(1,4,4,mcg,mgd,mcd)
# (MDR): r.l.u. --> disp
call MXM(-1,4,4,mcd,mcr,mdr)
end
#--------------------------------------------------------------------
SUBROUTINE GETMAT(v,mat)
# creates rotation matrix converting C&N coordinates to
# coordinates with X//V
#--------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 mat(4,4),v(3)
real*8 xv,yv,zv,xq,yq,vv,vv1,qq,si1,co1, si2,co2,TRANS
integer*4 i,j
# V to AB coordinates
xv=TRANS(v,1)
yv=TRANS(v,2)
zv=TRANS(v,3)
# QHKL to AB coordinates
xq=TRANS(qhkl,1)
yq=TRANS(qhkl,2)
vv=sqrt(xv**2+yv**2+zv**2)
vv1=sqrt(xv**2+yv**2)
qq=sqrt(xq**2+yq**2)
#/// Transformation from (V//x) to C&N
# special case: GRD vertical:
# define Gy//Qy:
if (abs(vv1).le.1e-7) then
do i=1,3
do j=1,3
mat(i,j)=0.d0
enddo
enddo
mat(3,1)=-1.d0
mat(2,2)=1.d0
mat(1,3)=1.d0
else
# sagital angle of V
si1=zv/vv
co1=sqrt(1.d0-si1**2)
# angle between Q and V(component in the scatt. plane)
# sign is positive from QHKL to V
si2=(xq*yv-yq*xv)/qq/vv1
co2=(xv*xq+yv*yq)/qq/vv1
mat(1,1)= co1*co2
mat(2,1)= co1*si2
mat(3,1)= +si1
mat(1,2)= -si2
mat(2,2)= co2
mat(3,2)=0
mat(1,3)=-si1*co2
mat(2,3)=-si1*si2
mat(3,3)=co1
endif
do i=1,3
mat(i,4)=0
mat(4,i)=0
enddo
mat(4,4)=1
end
#
SUBROUTINE GETNORMS(icom,vi,vf,vres,r0phon,r0brag)
#***********************************************************************
# returns volumes of VI (ki), VF (kf) and VRES (Q,E)
# returns normalizing factors for phonons and Bragg scans
# R0PHON=VI*VF/VRES, R0BRAG=R0PHON*(Bragg width)
#
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'trax.inc'
INCLUDE 'restrax.inc'
integer*4 icom
real*8 vi,vf,vres,r0phon,r0brag
real*8 wb,wv,wv1,dtm
real*8 aux1(4,4)
real*8 DETERM
if (icom.eq.1) then ! TRAX
vi=volcki
vf=volckf
# CALL INVERT(4,ATRAX,4,AUX,4)
dtm=DETERM(atrax,4,aux1)
call GETRESSIZE(atrax,qhkl,0,wb,wv,wv1)
else ! Monte Carlo
vi=vkiness
vf=vkfness
# CALL INVERT(4,ANESS,4,AUX,4)
dtm=DETERM(aness,4,aux1)
call GETRESSIZE(aness,qhkl,0,wb,wv,wv1)
endif
wb=wb*sqrt(2*pi)/sqrt8ln2
if (dtm.gt.0) then
vres=(2.*pii)**2/sqrt(dtm)
r0phon=vi*vf/vres
r0brag=r0phon*wb
else
vres=0
r0phon=0
r0brag=0
endif
end
#---------------------------------------------------------
SUBROUTINE GETRESSIZE(a,dir,uni,BRAG,evan,van)
# Returns Bragg and EVANadium (dE=0) fwhm in direction DIR
# Units:
# A-1 (UNI=0) or
# r.l.u. (UNI=1) or
# steps (UNI=2) step is defined by DIR in r.l.u.
#---------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 c8ln2
parameter (c8ln2=5.545177444)
real*8 mat(4,4),a(4,4),aux(4,4),a1(4,4),dir(3)
real*8 xrl,xcn,BRAG,evan,van,GETFWHM
integer*4 uni
call QNORM(dir,xrl,xcn)
if (xrl.gt.0) then
call GETMAT(dir,mat)
call MXM(1,4,4,a,mat,aux)
call MXM(-1,4,4,mat,aux,a1)
call INVERT(4,a1,4,aux,4)
van=sqrt(c8ln2*aux(1,1))
evan=GETFWHM(a1,1)
BRAG=sqrt(c8ln2/a1(1,1))
if(uni.eq.1) then
call QNORM(dir,xrl,xcn)
van=van*xrl/xcn
evan=evan*xrl/xcn
BRAG=BRAG*xrl/xcn
else if (uni.eq.2) then
call QNORM(dir,xrl,xcn)
van=evan/xcn
evan=evan/xcn
BRAG=BRAG/xcn
endif
else
evan=0.
BRAG=0.
van=0.
endif
end
#----------------------------------------------------
SUBROUTINE GETPHONWIDTH(ares,wa,wr,ws)
#// returns phonon fwhm in [Energy]x[A-1](WA) , [r.l.u] (WR)
#// and [steps] (WS)
#----------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 a(4,4),aux(4,4),ares(4,4),scd(4),wa,wr,ws,dnr,dna
# scan direction to coordinates with axis 4 normal to the disp. surface (D-coordinates)
call MXV(1,4,4,mdr,delq,scd)
call QNORM(delq,dnr,dna)
dnr=sqrt(dnr**2+res_dat(i_den)**2)
dna=sqrt(dna**2+res_dat(i_den)**2)
if(abs(scd(4)).lt.1.e-8) goto 999
call MXM(1,4,4,ares,mcd,aux)
call MXM(-1,4,4,mcd,aux,a) ! Res. matrix to D-coordinates
call INVERT(4,a,4,aux,4)
ws=sqrt(aux(4,4))*sqrt8ln2/abs(scd(4)) ! fwhm in steps
wa=ws*dna
wr=ws*dnr
return
999 wa=1.d10
wr=1.d10
ws=1.d10
end
#--------------------------------------------------------------------
SUBROUTINE GetProj(a,qe,a1,qe1,mcn)
#
# Transforms the resolution matrix A and vector QE to A1 and QE1,
# in the coordinates with X(1) parallel to grad(E) and scaled in R.L.U.
# MCN transforms vectors(Q,E) to the same coordinates
# (!! axes 2,3 are not transformed !!)
#--------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 a(4,4),a1(4,4),qe(4),qe1(4)
real*8 aux(4,4),mcn(4,4)
real*8 zra,gnr,gna
integer*4 i,j
call QNORM(grd,gnr,gna)
zra=gnr/gna
call MXM(1,4,4,a,mcg,aux)
call MXM(-1,4,4,mcg,aux,a1)
do i=1,4
a1(1,i)=a1(1,i)/zra
a1(i,1)=a1(1,i)
end do
a1(1,1)=a1(1,1)/zra
do 10 i=1,4
do 10 j=1,4
mcn(i,j)=mcg(i,j)
if(j.eq.1) mcn(i,j)=mcn(i,j)*zra
10 continue
call MXV(-1,4,4,mcn,qe,qe1)
end
#***********************************************************************
SUBROUTINE FWHM(icom)
# writes fwhm of section and projection in arbitrary direction
# (i.e. Bragg width and 'vanadium width' at dE=0)
# *** J.S. February 1999
# ICOM=1 .. TRAX
# ICOM=2 .. Monte Carlo (added by J.S. Sept. 2002)
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
real*8 brlu,bang,vrlu,vang,v(3),wa,wr
integer*4 i,icom
900 format( 'Analytical result (TRAX):')
901 format( 'Monte Carlo result (NESS):')
1 format(1x, 'Bragg width: ',g12.5, ' [A-1]',
* 3x,g12.5, ' [r.l.u]')
2 format(1x, 'Elastic ' 'Vanad' ' width (dE=0): ',g12.5, ' [A-1]',
* 3x,g12.5, ' [r.l.u]')
3 format(1x, 'Inelastic ' 'Vanad' ' width: ',g12.5, ' [A-1]',
* 3x,g12.5, ' [r.l.u]')
if (nos.ge.3) then
if (abs(ret(1))+abs(ret(2))+abs(ret(3)).le.1.d-30) then
write(sout,*) 'Direction vector has zero length !'
return
endif
do i=1,3
v(i)=ret(i)
enddo
if (icom.eq.1) then
call GETRESSIZE(atrax,v,0,bang,vang,wa)
call GETRESSIZE(atrax,v,1,brlu,vrlu,wr)
write(sout,900)
else
call GETRESSIZE(aness,v,0,bang,vang,wa)
call GETRESSIZE(aness,v,1,brlu,vrlu,wr)
write(sout,901)
endif
write(sout,1) bang,brlu
write(sout,2) vang,vrlu
write(sout,3) wa,wr
return
else
write(sout,*) 'Specify direction in reciprocal lattice.'
endif
end
SUBROUTINE BRAG(icom)
#***********************************************************************
# if ICOM=0 ... TRAX matrix
# if ICOM=1 ... M.C. matrix
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
real*8 eps1
parameter(eps1=1.d-6)
integer*4 i,j,icom
real*8 a(4,4),huitlog2
logical*4 lzer
real*8 zv1,zv2,zv3,dqx,dqy,dqz,dee,dvn,dum
real*8 GETFWHM
if (icom.eq.1) then
do i=1,4
do j=1,4
a(i,j)=aness(i,j)
enddo
enddo
else
do i=1,4
do j=1,4
a(i,j)=atrax(i,j)
enddo
enddo
endif
huitlog2=8.d0*log(2.d0)
lzer=.false.
do 11 i=1,4
if (abs(a(i,i)).lt.eps1) lzer=.true.
11 continue
if (lzer) then
write(smes,*) ' BRAG: problem with matrix -> zeros on diagonal'
return
endif
#---------------------------------------------------------------------------
# BRAGG WIDTHS IN RECIPROCRAL ANGSTROMS
dqx=sqrt(huitlog2/a(1,1))
dqy=sqrt(huitlog2/a(2,2))
dqz=sqrt(huitlog2/a(3,3))
# ENERGY BRAGG AND VANADIUM WIDTHS.
dee=sqrt(huitlog2/a(4,4))
call VANAD(a,dvn,dum)
# VANADIUM WIDTHS IN RECIPROCRAL ANGSTROMS
zv1=GETFWHM(a,1)
zv2=GETFWHM(a,2)
zv3=GETFWHM(a,3)
write(sout,5) dqx,dqy,dqz,
* zv1,zv2,zv3,
* cunit,dvn,dee
5 format( ' Bragg widths (radial,tangential,vertical) [A-1]'/
1 ' DQR=',f9.5, ' DQT=',f9.5, ' DQV=',f9.5,/,
4 ' '/
2 ' ' 'Vanad' ' widths (from section at dE=0) [A-1]'/
3 ' DQR=',f9.5, ' DQT=',f9.5, ' DQV=',f9.5,/,
4 ' '/
5 ' Energy widths (Vanad, Bragg) ',a,/,
6 ' DVN=',f9.5, ' DEE=',f9.5)
return
end
#
#***********************************************************************
SUBROUTINE RESOL(icom,iarg)
# Calculate resolutuion matrices and related parameters
# ICOM=1 ... analytically (TRAX)
# ICOM=2 ... Monte Carlo
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
real*8 ada(4,4),b(4,4),az(4,4),a(4,4)
real*8 tmp,vi,vf,vres,r0phon,r0brag
integer*4 icom,iarg,k,i,j
character*2 iax(4),iaxp(4)
character*4 FWHM
data iax/ 'X', 'Y', 'Z', 'W'/,FWHM/ 'FWHM'/
data iaxp/ 'X' '', 'Y' '', 'Z' '', 'W' ''/
900 format( 'Analytical result (TRAX):')
901 format( 'Monte Carlo result (NESS):')
910 format( 'Vol(ki)= ',g12.4, ' Vol(kf)=',g12.4,
* ' Vol(QE)= ',g12.4)
911 format( 'Norm factors: R0(C&N)= ',g12.4)
920 format( 'Resolution Matrix, Cooper & Nathans coordinates [A^-1]')
921 format(2x,4(10x,a))
922 format(2x,a,4e12.4)
930 format( 'Resolution Matrix, rec. lat. coordinates')
940 format( 'Diagonalised Resolution Matrix, [r.l.u.]')
941 format(2x,5(10x,a))
942 format(2x,a,5e12.4)
945 format( 'Direction Cosines (w.r.t. reciprocal lattice)')
if (icom.eq.1) then
do i=1,4
do j=1,4
a(i,j)=atrax(i,j)
enddo
enddo
write(sout,900)
else
do i=1,4
do j=1,4
a(i,j)=aness(i,j)
enddo
enddo
write(sout,901)
endif
if (iarg.le.0.or.iarg.gt.4) then
goto 11
else
go to (11,12,13,14),iarg
endif
return
# --------------------------------------------------------------------
# RESOL 1, get volumes VI VF etc. (J.S. 6/6/97)
11 call GETNORMS(icom,vi,vf,vres,r0phon,r0brag)
write(sout,910) vi,vf,vres
write(sout,911) r0phon ! ,R0BRAG
return
# --------------------------------------------------------------------
# RESOL 2, resolution matrix in C&N coord.
12 write(sout,920)
write(sout,921) iax
write(sout,922) (iax(i),(a(i,j),j=1,4),i=1,4)
return
# --------------------------------------------------------------------
# RESOL 3, resolution matrix in rec. lat. coord.
13 call CN2RLU(a,az)
write(sout,930)
write(sout,921) iax
write(sout,922) (iax(i),(az(j,i),j=1,4),i=1,4)
return
# --------------------------------------------------------------------
# RESOL 4, Diagonalized matrix
14 call CN2RLU(a,az)
call DIAG(az,ada,b)
write(sout,940)
write(sout,921) iaxp
do i=1,4
do k=1,4
if (ada(i,k).lt.1.d-20) ada(i,k)=0.d0
enddo
write(sout,922) iaxp(i),(ada(i,k),k=1,4)
enddo
# *** At Cooper&Nathans the res. function is exp(-.5*MijXiXj) !!!!!!
tmp=8.d0*log(2.d0)
do i=1,4
if (ada(i,i).gt.0.d0) ada(i,i)= sqrt(tmp/abs(ada(i,i)))
enddo
write(sout,941) iax,FWHM
do i=1,4
write(sout,942) iaxp(i),(b(j,i),j=1,4),ada(i,i)
enddo
return
end
#***********************************************************************
#//
#// D E P O S I T - UNITS NOT USED IN CURRENT VERSION
#//
#***********************************************************************
#------------------------------------------------------------------------------------
SUBROUTINE A3_TO_Q(a3,q)
# Get QHKL equivalent to the sample rotation by A3 (nominal QHKL corresponds to A3=0)
# Uses current transformation matrix, MRC !
#------------------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 a3,da3,vq(4),wq(4),q(3),si,co
integer*4 i
if(a3.ne.0.) then
da3=a3*deg
si=sin(da3)
co=cos(da3)
wq(1)=q0*(co-1.d0)
wq(2)=-q0*si
wq(3)=0.
wq(4)=0.
call MXV(1,4,4,mrc,wq,vq) ! from CN to r.l.u.
do i=1,3
q(i)=res_dat(i_qh+i-1)+vq(i)
enddo
else
do i=1,3
q(i)=res_dat(i_qh+i-1)
enddo
endif
end
#------------------------------------------------------------------
SUBROUTINE PSDSCAN(qxmin,qxmax,qymin,qymax,emin,emax)
# calculates direction and step size of the scan performed
# by a linear position-sensitive detector
#------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'trax.inc'
INCLUDE 'rescal.inc'
integer*4 i
real*8 qxmin,qxmax,qymin,qymax,emin,emax
real*8 tbr,coa,sigs,siga,sia, sin1,sin2,ro,xdmax,vkx,vkz,hir
real*8 dqe(4),rqe(4)
1 format(2x, 'PSD scan range (dQ,dE) : ',4(1x,f7.3), ' [meV]')
2 format(2x, 'scale: ',f8.3, ' mm/range')
tbr=tetaa*tdr
call HIRANG(tbr,hiana,hir)
coa=(q0**2-2*homega/hsqovm)/(2*kf0*q0)
sigs=sign(1.d0,tetas)
siga=sign(1.d0,tetaa)
sia=sqrt(1-coa**2)
sin1=sin(abs(tbr-hir))
sin2=sin(abs(tbr+hir))
if(sin1.lt.1.d-6) sin1=1.d-6
if(sin2.lt.1.d-6) sin2=1.d-6
ro=roha
xdmax=wana*(sin1 - 2.*ro*vl3 + vl3/vl2*sin2)/2.
vkx=-vkf*sin2/vl2*wana/2.
vkz=siga*vkf*(ro-sin2/vl2)/abs(tan(tbr))*wana/2.
qxmax= vkz*coa + vkx*sia*sigs
qymax=-vkz*sia*sigs + vkx*coa
emax=-2*vkz/kf0*ef0
qxmin=-qxmax
qymin=-qymax
emin=-emax
dqe(1)=(qxmax-qxmin)
dqe(2)=(qymax-qymin)
dqe(3)=0.
dqe(4)=(emax-emin)
call MXV(1,4,4,mrc,dqe,rqe)
write(sout,*)
write(sout,1) (rqe(i),i=1,4)
write(sout,2) xdmax*20.
write(sout,*)
return
end