Source module last modified on Wed, 13 Jul 2005, 16:20;
HTML image of Fortran source automatically generated by
for2html on Mon, 29 May 2006, 15:06.
#//////////////////////////////////////////////////////////////////////
#////
#//// R E S T R A X 4.4
#////
#//// Subroutines from the program TRAX
#////
#//// * SUBROUTINE THRAX
#//// * SUBROUTINE SSTV
#//// * SUBROUTINE TRCAN(VL0,VLI,VLC,S10,S20,NF,XSHI)
#//// * SUBROUTINE TCR
#//// * SUBROUTINE HIRANG(OB,HII,HI)
#//// * SUBROUTINE VANAD(DVN,YVN) J.S. 3/6/97
#//// * FUNCTION OPTV(N,RO) J.S. 3/6/97
#////
#//////////////////////////////////////////////////////////////////////
#********************************
SUBROUTINE THRAX
#********************************
# COMPUTES THE RESOLUTION MATRIX, THE NORMALIZATION FACTORS AND
# DIFFERENT SCAN WIDTHS AND ABSOLUTE INTENSITIES AT DETECTOR
# FOR A THREE-AXIS NEUTRON SPECTROMETER
# implicit REAL*8 (a-h,o-z)
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
# INCLUDE 'rescal.inc'
INCLUDE 'trax.inc'
integer*4 i,j
real*8 right,sim,vkfm,sia,vkim,two,cose,sine,ar,co,si,ttemp,
* am,bet,xaa,cot,com,v1h,v0h,rmh,v0v,v1v,rmv,volph1,erm,
* wdthef,zef,areaef,pathef,detas,fact,flu,ysam,ysam0
real*8 em(4,4),cm1v(5,5),x33(1,1),aux3(3,3),aux4(4,4),au2(2,2),
1 ahp(3,4),aux6(6,6),cnh(6,6),cnv(3,3),au23(2,3),aux2(2,2),
2 ah(4,6),av(3,5),av1(1,3),au3(3,3),ah8(6,8),cki0(3,3),
3 aux5(5,5),cnini(5,5),det1,det2,det3,det4
real*8 aaa(4,4)
common /matrix/ aaa
integer*4 ierr
common /error/ierr
real*8 emmin1(4,4),volres
common /em/emmin1,volres
real*8 dv(5,5),bm,ba,ahs(3,4)
common /ctr/dv,bm,ba,ahs
real*8 cmh(8,8),cmv(5,5),cmhp(4,4),cmvp(3,3)
common /covar/cmh,cmv,cmhp,cmvp
real*8 sh(10,10),sh1(10,10),sv(5,5),sv1(5,5),shp(5,5),svp(3,3)
common /stv/sh,sh1,sv,sv1,shp,svp
real*8 DETERM
# make local variables static in THRAX !!!
save right,sim,vkfm,sia,vkim,two,cose,sine,ar,co,si,ttemp,
* am,bet,xaa,cot,com,v1h,v0h,rmh,v0v,v1v,rmv,volph1,erm,
* wdthef,zef,areaef,pathef,detas,fact,flu,ysam,ysam0
save em,cm1v,x33,aux3,aux4,au2,ahp,aux6,cnh,cnv,au23,aux2,
* ah,av,av1,au3,ah8,cki0, aux5,cnini,det1,det2,det3,det4
data pii,tdr,tmr,r8ln2,hsqovm/3.141592653589793239d0,.01745329,
+2.9089e-4,2.3548,4.144219/
# write(*,*) 'THRAX', NEFIX
if ((im*isc*ia).eq.0) then
do i=1,4
do j=1,4
aaa(i,j)=0.
end do
end do
reltrax=0.
write(*,*) 'IM,ISC,IA: ',im,isc,ia
return
endif
# WRITE(*,*) NGUIDE,GAMACR
right=pii/2.
if(nefix.gt.1)go to 14
vki=sqrt(2./hsqovm*ei0)
slamdi=2.*pii/vki
om=asin(slamdi/2./cryd(1))*im
tetam=om/tdr
sim=sin(om)
vkfm=pii/cryd(2)
vkf=vki*vki-homega*2./hsqovm
if(vkf.le.(vkfm*vkfm))go to 7
vkf=sqrt(vkf)
slamdf=2.*pii/vkf
oa=asin(vkfm/vkf*ia)
tetaa=oa/tdr
ef0=0.5*hsqovm*vkf*vkf
go to 15
14 vkf=sqrt(2./hsqovm*ef0)
slamdf=2.*pii/vkf
oa=asin(slamdf/2./cryd(2))*ia
tetaa=oa/tdr
sia=sin(oa)
vkim=pii/cryd(1)
vki=vkf*vkf+homega*2./hsqovm
if(vki.le.(vkim*vkim))go to 7
vki=sqrt(vki)
slamdi=2.*pii/vki
om=asin(vkim/vki)*im
tetam=om/tdr
ei0=.5*hsqovm*vki*vki
15 continue
# VQ0=Q
#
36 two=.5*(vki/vkf+vkf/vki-vq0*vq0/vki/vkf)
if(abs(two).ge.1.) go to 9
two=acos(two)
os=isc*two/2.
tetas=os/tdr
phi=atan2(-vkf*sin(2.*os),vki-vkf*cos(2.*os))
phim=phi-2.*os
call HIRANG(om,himon,him)
call HIRANG(oa,hiana,hia)
# DIRTAU=DIRTD*TDR
if(nsam.eq.0) go to 24
if(his.le.(right-os)) go to 29
his=his-pii
# DIRTAU=DIRTAU-PII
29 if(his.ge.(-right-os)) go to 24
his=his+pii
# DIRTAU=DIRTAU+PII
24 cose=1.
sine=0.
epsln=0.
#
8 etamef=etm*tmr
etvmef=etm*anrm*tmr
ar=abs(om+him)/(om+him)
co=cos(him)
si=sin(him)
if(rohm.eq.0.) go to 18
ttemp=poiss(1)/(1-poiss(1))
am = 1 - (1+ttemp)*co*co
bm = (1+ttemp)*si*co
#
18 etaaef=eta*tmr
etvaef=eta*anra*tmr
bet=abs(oa+hia)/(oa+hia)
co=cos(hia)
si=sin(hia)
if(roha.eq.0.) go to 21
ttemp=poiss(2)/(1-poiss(2))
xaa = 1 - (1+ttemp)*co*co
ba = (1+ttemp)*si*co
#
21 pm=1. ! peak reflectivity=1 is supposed
pa=1.
#
#
#
call SSTV
#
call TCR
#
call BABT(cmv,dv,5,5,cm1v)
# CM1V IS THE COVARIANCE MATRIX OF THE VERTICAL VARIABLES
# (DELTA0,DELTA1,DELTA2,DELTA3,ZS).
do 1 i=1,6
do 1 j=1,8
1 ah8(i,j)=0.
cot=1./tan(om)
do 30 i=1,3
do 30 j=1,4
30 ahp(i,j)=0.
co=cos(om-him)/vl1*vki
si=vki*sin(om-him)/vl1
ah8(1,1)=cot*(vki*rohm*ar-si)
ah8(1,2)=-cot*co+rohm*ar*vki*(-am+cot*bm)
ah8(2,2)=co
ah8(2,1)=si
co=cos(os+his)
si=sin(os+his)
ah8(3,3)=co
ah8(3,4)=si
ah8(4,3)=-si
ah8(4,4)=co
co=vki*co/vl1
si=vki*si/vl1
ah8(1,3)=-cot*co
ah8(1,4)=-cot*si
ah8(1,7)=vki*cot
ahp(1,3)=-ah8(1,7)/vl1
if(etm.eq.0.) ah8(1,7)=ah8(1,7)*sin(om+him)/abs(sin(om-him))
ahp(1,4)=ah8(1,7)
ahp(1,1)=ah8(1,1)
ahp(1,2)=ah8(1,2)
ahp(2,1)=ah8(2,1)
ahp(2,2)=ah8(2,2)
ahp(3,3)=1.
ahp(2,3)=vki/vl1
ah8(2,3)=co
ah8(2,4)=si
cot=1./tan(oa)
co=vkf*cos(os-his)/vl2
si=vkf*sin(os-his)/vl2
ah8(5,3)=-co*cot
ah8(5,4)=si*cot
ah8(6,3)=-co
ah8(6,4)=si
co=vkf*cos(oa+hia)/vl2
si=vkf*sin(oa+hia)/vl2
ah8(5,5)=cot*(si-roha*vkf*bet)
ah8(5,6)=-cot*co-vkf*roha*bet*(xaa+cot*ba)
ah8(5,8)=-vkf*cot
if(eta.eq.0.) ah8(5,8)=ah8(5,8)*sin(oa-hia)/abs(sin(oa-hia))
ah8(6,5)=si
ah8(6,6)=-co
call BABT(cmh,ah8,8,6,cnh)
# THE COVARIANCE MATRIX OF THE VECTOR (DKI,KI*GAMAI,YS,XS,DKF,
# KF*GAMAF) IS CNH(6,6) AND HAS BEEN COMPUTED
call BABT(cmhp,ahp,4,3,aux3)
do 31 i=1,2
do 31 j=1,2
31 cki0(i,j)=aux3(i,j)
com=DETERM(aux3,3,au3)
# write(*,*) 'GOTO 1 ', COM
if(com.lt.0.) go to 10
v1h=sqrt(com)
do 32 j=1,4
32 ahs(1,j)=ahp(1,j)
call BABT(cmhp,ahs,4,3,aux3)
com=DETERM(aux3,3,au3)
# write(*,*) 'GOTO 2 ', COM
if(com.lt.0.) go to 10
v0h=sqrt(com)
rmh=v0h/v1h
do 33 i=1,2
do 33 j=1,3
33 au23(i,j)=0.
au23(1,1)=-vki/vl0
au23(1,2)=-au23(1,1)
au23(2,1)=1.
call BABT(cmvp,au23,3,2,au2)
com=DETERM(au2,2,aux2)
# write(*,*) 'GOTO 3 ', COM
if(com.lt.0.) go to 10
v0v=sqrt(com)
au23(1,1)=0.
au23(1,2)=-vki/vl1
au23(1,3)=-au23(1,2)
au23(2,1)=0.
au23(2,3)=1.
call BABT(cmvp,au23,3,2,au2)
com=DETERM(au2,2,aux2)
# write(*,*) 'GOTO 4 ', COM
if(com.lt.0.) go to 10
v1v=sqrt(com)
cki0(3,3)=au2(1,1)
rmv=v0v/v1v
volph1=v1h*v1v*2.*pii*sqrt(2.*pii)
# VOLPH1 IS THE PHASE VOLUME AFTER MONOCHROMATOR FOR AN INFINITELY
# EXTENDED SAMPLE
erm=rmh*rmv
do 2 i=1,3
do 2 j=1,5
2 av(i,j)=0.
av(1,2)=vki
av(2,3)=vkf
av(3,5)=1.
call BABT(cm1v,av,5,3,cnv)
# THE COVARIANCE MATRIX OF THE VECTOR (KI*DELTAI,KF*DELTAF,ZS) IS
# CNV(3,3) AND HAS BEEN COMPUTED
do 11 i=1,2
do 11 j=1,2
11 cki(i,j)=cnh(i,j)
do 12 i=1,2
cki(i,3)=0.
ckf(i,3)=0.
ckf(3,i)=0.
cki0(i,3)=0.
cki0(3,i)=0.
12 cki(3,i)=0.
cki(3,3)=cnv(1,1)
# CKI(3,3) IS THE COVARIANCE MATRIX OF THE (DKI,KI*GAMAI,KI*DELTAI)
# VECTOR AND DESCRIBES THE EFFECTIVE MONOCHROMATOR ELLIPSOID
# CKI0(3,3) IS THE SAME THING FOR AN INFINITELY EXTENDED SAMPLE
com=DETERM(cki,3,aux3)
# write(*,*) 'GOTO 5 ', COM
if(com.lt.0.) go to 10
volcki=2.*pii*sqrt(2.*pii*com)
do 34 i=1,2
do 34 j=1,2
ckf(i,j)=cnh(i+4,j+4)
34 continue
ckf(3,3)=cnv(2,2)
#//added/////////////////////
com=DETERM(ckf,3,aux3)
# write(*,*) 'GOTO 6 ', COM
if(com.lt.0.) go to 10
volckf=2.*pii*sqrt(2.*pii*com)
#/////////////////////////////////
wdthef=sqrt(12.*cnh(3,3))
radet=sqrt((sh(4,4)*sh(5,5)-sh(4,5)*sh(5,4))*sv(3,3))
fracv=(cmh(3,3)*cmh(4,4)-cmh(3,4)*cmh(4,3))*cmv(3,3)
# write(*,*) 'GOTO 7 ', FRACV
if(fracv.lt.0.) go to 10
fracv=radet*sqrt(fracv)
if(nsam.gt.0) volsam=wsam*thsam*hsam
if(nsam.eq.0) volsam=.25*pii*diasam*diasam*hsam
volef=fracv*volsam
zef=sqrt(12.*cnv(3,3))
areaef=wdthef*zef
pathef=volef/areaef
do 3 i=1,4
do 3 j=1,6
3 ah(i,j)=0.
com=cos(phim)
sim=sin(phim)
co=cos(phi)
si=sin(phi)
ah(1,1)=co
ah(1,2)=si
ah(1,5)=-com
ah(1,6)=-sim
ah(2,1)=-si
ah(2,2)=co
ah(2,5)=sim
ah(2,6)=-com
ah(4,1)=hsqovm*vki
ah(4,5)=-hsqovm*vkf
call BABT(cnh,ah,6,4,emmin1)
av1(1,1)=1.
av1(1,2)=-1.
av1(1,3)=0.
call BABT(cnv,av1,3,1,x33)
detas=ets*tmr*vq0/r8ln2
detas=detas*detas
emmin1(2,2)=emmin1(2,2)+detas
emmin1(3,3)=x33(1,1)+detas
call INVERT(4,emmin1,4,em,4)
# CALL INV1(EMMIN1,4,AUX4,EM) ! replaced by INVERT, J.S.
# AT THIS MOMENT THE RESOLUTION MATRIX EM(4,4) AND ITS INVERSE, THE
# COVARIANCE MATRIX EMMIN1(4,4) HAVE BEEN COMPUTED.
#
volres=DETERM(emmin1,4,aux4)
# write(*,*) 'GOTO 8 ', VOLRES
if(volres.lt.0.) go to 10
volres=4.*pii*pii*sqrt(volres)
do 17 i=1,3
do 17 j=1,3
17 au3(i,j)=cnh(i,j)
det1=DETERM(au3,3,aux3)
# write(*,*) 'GOTO 9 ',DET1
if(det1.lt.0.) go to 10
det2=cnv(1,1)*cnv(3,3)-cnv(1,3)*cnv(3,1)
fact=erm
co=4.*pii*pii*pm*fact
flu=1.d+13 ! FLX(VKI)
si=sqrt(det1)*sqrt(det2*2.*pii)
ysam=flu*co*si
ysam0=flu*co*volph1/2./pii
det1=DETERM(cnh,6,aux6)
# write(*,*) 'GOTO 10 ',DET1
if(det1.le.0.d0) go to 10
det2=DETERM(cnv,3,aux3)
# write(*,*) 'GOTO 11 ',DET2
if(det2.le.0.d0) go to 10
r0trax=sqrt((2.*pii)**9)*sqrt(det1)*sqrt(det2)
det3=DETERM(cki,3,aux3)
# write(*,*) 'GOTO 12 ',DET3
if(det3.le.0.d0) go to 10
do i=1,5
do j=1,5
if (i.le.3.and.j.le.3) then
cnini(i,j)=cnh(i,j)
else
cnini(i,j)=0
endif
end do
end do
cnini(4,5)=cnv(1,3)
cnini(5,4)=cnv(3,1)
cnini(4,4)=cnv(1,1)
cnini(5,5)=cnv(3,3)
det4=DETERM(cnini,5,aux5)
reltrax=(2*pii)**3*sqrt(det1*det2/det3)
#
# RELTRAX=(2*PII)**2*SQRT(DET1*DET2/DET4)*10 ! cm -> mm
do 100 i=1,4
do 100 j=1,4
100 aaa(i,j)=em(i,j)
return
7 write(smes,98)homega
98 format(5x, 'THE REQUIRED ENERGY TRANSFER,',f8.3,
& ' MILLIEV, IS INACC ESIBLE')
ierr=1
return
9 write(smes,99)vq0
99 format(5x, 'THE REQUIRED WAVE-VECTOR TRANSFER,',f8.3,
& ' 1/A, IS INACCESIBLE')
ierr=1
return
10 write(smes,95)
95 format(5x, 'ERROR IN INVERTING A MATRIX DUE TO A TOO SMALL'/5x,
1 'MOSAIC SPREAD OR COLLIMATOR DIVERGENCE')
ierr=1
return
end
#
#**************************************************************
#
SUBROUTINE SSTV
# THIS SUBROUTINE COMPUTES THE INITIAL PROBABILITY MATRICES OF THE
# HORIZONTAL AND VERTICAL SPATIAL VARIABLES (Y0,LM,GM,LS,GS,LA,GA,YD,
# CSIM,CSIA),AND (Z0,ZM,ZS,ZA,ZD)
implicit real*8 (a-h,o-z)
INCLUDE 'trax.inc'
common/stv/sh(10,10),sh1(10,10),sv(5,5),sv1(5,5),shp(5,5),svp(3,3)
data twel,sixt/12.,16./
# initialization
do 1 i=1,10
do 1 j=1,10
1 sh(i,j)=0.
do 2 i=1,5
do 2 j=1,5
shp(i,j)=0.
2 sv(i,j)=0.
do 3 i=1,3
do 3 j=1,3
3 svp(i,j)=0.
# source dimensions
if(nsou.gt.0)go to 4
sh(1,1)=sixt/diasou/diasou
sv(1,1)=sh(1,1)
go to 5
4 sh(1,1)=twel/wsou/wsou
sv(1,1)=twel/hsou/hsou
# monochromator dimensions
5 sh(2,2)=twel/wmon/wmon
sh(3,3)=twel/thmon/thmon
sv(2,2)=twel/hmon/hmon
do 6 i=1,3
do 6 j=1,3
# copy to matrices for primary part
6 shp(i,j)=sh(i,j)
svp(1,1)=sv(1,1)
svp(2,2)=sv(2,2)
# sample dimensions
if(nsam.gt.0)go to 7
w=sixt/diasam/diasam
sh(4,4)=w
sh(5,5)=w
go to 8
7 sh(4,4)=twel/wsam/wsam
sh(5,5)=twel/thsam/thsam
8 sv(3,3)=twel/hsam/hsam
# analyzer dimensions
sh(6,6)=twel/wana/wana
sh(7,7)=twel/thana/thana
sv(4,4)=twel/hana/hana
# mosaicities
w=r8ln2*r8ln2
sh(9,9)=w/etamef/etamef
shp(5,5)=sh(9,9)
sh(10,10)=w/etaaef/etaaef
# detector
if(ndet) 9,9,10
9 sh(8,8)=sixt/diadet/diadet
sv(5,5)=sh(8,8)
return
10 sh(8,8)=twel/wdet/wdet
sv(5,5)=twel/hdet/hdet
return
end
#
#**************************************************************
#
SUBROUTINE TRCAN(vl0,vli,vlc,s10,s20,nf,xshi)
implicit real*8 (a-h,o-z)
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
# INCLUDE 'trax.inc'
# THIS SUBROUTINE COMPUTES THE TRANSMISSION MATRIX FOR A COARSE COLL
# IMATOR OF CIRCULAR (NF=0) OR RECTANGULAR (NF=1) CROSS SECTION
real*8 xshi(2,2)
data un,con1,con2,con3/1.d0,.866d0,1.3333d0,.6667d0/
if(vlc.eq.0.) vlc=0.1
do 4 i=1,2
do 4 j=1,2
4 xshi(i,j)=0.d0
if(vl0.lt.(vli+vlc))go to 15
alfa=vli/vl0
beta=vlc/vl0
eps=alfa+beta
s1=s10
s2=s20
if(nf)10,10,5
5 s1=con1*s10
s2=con1*s20
10 if(vlc.eq.0.)go to 30
y0k=abs(s1*eps-s2*alfa)/beta
y0j=(eps*s1+s2*alfa)/beta
ai0=2.*(s1/alfa+s2/eps)*(y0j-y0k)-(un/alfa-un/eps)*(y0j**2-y0k**2)
1+4.*y0k*s2/eps
a11=con2*s2*y0k**3/eps+con3*(s2/eps+s1/alfa)*(y0j**3-y0k**3)-
1.5*(un/alfa-un/eps)*(y0j**4-y0k**4)
a22=con3*((s1/alfa)**3+(s2/eps)**3)*(y0j-y0k)-(s1**2*(un-alfa)/
2alfa**3-s2**2*(un-eps)/eps**3)*(y0j**2-y0k**2)+con3*(s1*(un-alfa
3)**2/alfa**3+s2*(un-eps)**2/eps**3)*(y0j**3-y0k**3)-.1667*((un/
4alfa-un)**3-(un/eps-un)**3)*(y0j**4-y0k**4)+con2*y0k*(s2/eps)**
53+con2*s2*(un-eps)**2*(y0k/eps)**3
a12=.5*((s1/alfa)**2-(s2/eps)**2)*(y0j**2-y0k**2)-con3*(s1*(un-
1alfa)/alfa**2+s2*(un-eps)/eps**2)*(y0j**3-y0k**3)+.25*((un/alfa-
2un)**2-(un/eps-un)**2)*(y0j**4-y0k**4)-con2*s2*(un-eps)/eps**2*
3y0k**3
x=a11/ai0
y=a22/ai0
z=a12/ai0
d=x*y-z**2
xshi(1,1)=y/d
xshi(2,2)=x/d
xshi(1,2)=-z/d
xshi(2,1)=xshi(1,2)
go to 25
30 xmo=s1*s1/12.
xshi(1,1)=(un-alfa)**2/xmo
xshi(2,2)=alfa**2/xmo
xshi(1,2)=alfa*(un-alfa)/xmo
xshi(2,1)=xshi(1,2)
go to 25
15 write(smes,20)
20 format(2x, 'WRONG INPUT DATA; THE PRESENCE OF THIS COLLIMATOR',
& ' IS IGNORED')
25 return
end
#
#**************************************************************
SUBROUTINE TCR
# COMPUTES THE PROBABILITY MATRICES OF THE SPATIAL VARIABLES AS
# MODIFIED BY THE BRAGG CONSTRAINTS, NEUTRON GUIDE, SOLLER
# COLLIMATORS AND COARSE COLLIMATORS OR SLITS
implicit real*8 (a-h,o-z)
INCLUDE 'trax.inc'
real*8 xdh(4,10),xgh(10,8),tv(2,5),aco(4,4),xf(2,2),
& dv1(4,5),can(10,10),vcan(5,5),xshi(2,2),ah8(8,8),tv1(1,3),
& canp(5,5),vcanp(3,3),dhp(2,5),acop(2,2),dvp(2,3),ghp(5,4),xf1(1)
common/stv/sh(10,10),sh1(10,10),sv(5,5),sv1(5,5),shp(5,5),svp(3,3)
common /qef/pathef,wdthef,areaef,zef
common /ctr/dv(5,5),bm,ba,ahs(3,4)
common /covar/cmh(8,8),cmv(5,5),cmhp(4,4),cmvp(3,3)
#1 FORMAT(a,10(G10.4,2x))
do 111 i=1,10
do 111 j=1,10
sh1(i,j)=sh(i,j)
111 can(i,j)=0.
do 112 i=1,5
do 112 j=1,5
canp(i,j)=0.
sv1(i,j)=sv(i,j)
112 vcan(i,j)=0.
do 110 i=1,2
do 110 j=1,5
110 dhp(i,j)=0.
do 113 i=1,3
do 113 j=1,3
113 vcanp(i,j)=0.
do 116 i=1,3
do 116 j=1,4
116 ahs(i,j)=0.
if(nfm.lt.0)go to 15
call TRCAN(vl0,vlcanm,vlsm,hdm1,hdm2,nfm,xshi)
co=cos(om+him)
si=sin(om+him)
can(1,1)=xshi(1,1)
can(1,2)=xshi(1,2)*si
can(1,3)=-xshi(1,2)*co
can(2,2)=xshi(2,2)*si*si
can(2,3)=-xshi(2,2)*co*si
can(3,3)=xshi(2,2)*co*co
can(2,1)=can(1,2)
can(3,1)=can(1,3)
can(3,2)=can(2,3)
call TRCAN(vl0,vlcanm,vlsm,vdm1,vdm2,nfm,xshi)
# write(*,1) 'VL0: ',VL0,VLCANM,VLSM,VDM1,VDM2,NFM,XSHI
vcan(1,1)=xshi(1,1)
vcan(2,2)=xshi(2,2)
vcan(1,2)=xshi(1,2)
vcan(2,1)=vcan(1,2)
do 114 i=1,3
do 114 j=1,3
vcanp(i,j)=vcan(i,j)
114 canp(i,j)=can(i,j)
15 if(nfs.lt.0)go to 20
co=cos(om-him)
si=sin(om-him)
a=cos(os+his)
b=sin(os+his)
call TRCAN(vl1,vlcans,vlms,hds1,hds2,nfs,xshi)
can(2,2)=can(2,2)+xshi(1,1)*si*si
can(3,3)=can(3,3)+xshi(1,1)*co*co
can(2,3)=can(2,3)+xshi(1,1)*si*co
can(3,2)=can(2,3)
canp(4,4)=xshi(2,2)
can(4,4)=xshi(2,2)*a*a
can(5,5)=xshi(2,2)*b*b
can(4,5)=xshi(2,2)*a*b
can(5,4)=can(4,5)
can(2,4)=can(2,4)-xshi(1,2)*si*a
can(4,2)=can(2,4)
can(2,5)=can(2,5)-xshi(1,2)*si*b
can(5,2)=can(2,5)
can(3,4)=can(3,4)-xshi(1,2)*co*a
can(4,3)=can(3,4)
can(3,5)=can(3,5)-xshi(1,2)*co*b
can(5,3)=can(3,5)
call TRCAN(vl1,vlcans,vlms,vds1,vds2,nfs,xshi)
# write(*,1) 'VL1: ',VL1,VLCANS,VLMS,VDS1,VDS2,NFS,XSHI
vcan(2,2)=vcan(2,2)+xshi(1,1)
vcan(3,3)=xshi(2,2)
vcan(2,3)=xshi(1,2)
vcan(3,2)=vcan(2,3)
do 115 i=1,3
do 115 j=1,3
canp(i,j)=can(i,j)
115 vcanp(i,j)=vcan(i,j)
20 if(nfa.lt.0)go to 25
co=cos(oa+hia)
si=sin(oa+hia)
a=cos(os-his)
b=sin(os-his)
call TRCAN(vl2,vlcana,vlsa,hda1,hda2,nfa,xshi)
can(4,4)=can(4,4)+xshi(1,1)*a*a
can(5,5)=can(5,5)+xshi(1,1)*b*b
can(4,5)=can(4,5)-xshi(1,1)*a*b
can(5,4)=can(4,5)
can(6,6)=xshi(2,2)*si*si
can(7,7)=xshi(2,2)*co*co
can(6,7)=-xshi(2,2)*si*co
can(7,6)=can(6,7)
can(4,6)=can(4,6)+xshi(1,2)*si*a
can(6,4)=can(4,6)
can(4,7)=can(4,7)-xshi(1,2)*a*co
can(7,4)=can(4,7)
can(5,6)=can(5,6)-xshi(1,2)*si*b
can(6,5)=can(5,6)
can(5,7)=can(5,7)+xshi(1,2)*co*b
can(7,5)=can(5,7)
call TRCAN(vl2,vlcana,vlsa,vda1,vda2,nfa,xshi)
# write(*,1) 'VL2: ',VL2,VLCANA,VLSA,VDA1,VDA2,NFA,XSHI
vcan(3,3)=vcan(3,3)+xshi(1,1)
vcan(4,4)=vcan(4,4)+xshi(2,2)
vcan(3,4)=xshi(1,2)
vcan(4,3)=vcan(3,4)
25 if(nfd.lt.0)go to 30
call TRCAN(vl3,vlcand,vlad,hdd1,hdd2,nfd,xshi)
co=cos(oa-hia)
si=sin(oa-hia)
can(6,6)=can(6,6)+xshi(1,1)*si*si
can(7,7)=can(7,7)+xshi(1,1)*co*co
can(6,7)=can(6,7)+xshi(1,1)*si*co
can(7,6)=can(6,7)
can(8,8)=xshi(2,2)
can(6,8)=-xshi(1,2)*si
can(8,6)=can(6,8)
can(7,8)=-xshi(1,2)*co
can(8,7)=can(7,8)
call TRCAN(vl3,vlcand,vlad,vdd1,vdd2,nfd,xshi)
# write(*,1) 'VL3: ',VL3,VLCAND,VLAD,VDD1,VDD2,NFD,XSHI
vcan(4,4)=vcan(4,4)+xshi(1,1)
vcan(5,5)=xshi(2,2)
vcan(4,5)=xshi(1,2)
vcan(5,4)=vcan(4,5)
30 call SUM(sh,can,10,sh1)
call SUM(sv,vcan,5,sv1)
call SUM(shp,canp,5,shp)
call SUM(svp,vcanp,3,svp)
# THE PRESENCE OF COARSE COLLIMATORS OR SLITS WAS ACCOUNTED FOR
do 1001 i=1,4
do 1001 j=1,10
1001 xdh(i,j)=0.
do 2 i=1,5
do 2 j=1,5
2 dv(i,j)=0.
do 3 i=1,10
do 3 j=1,8
3 xgh(i,j)=0.
xgh(2,1)=1.
xgh(3,2)=1.
xgh(4,3)=1.
xgh(5,4)=1.
xgh(6,5)=1.
xgh(7,6)=1.
xgh(9,7)=1.
xgh(10,8)=1.
co=cos(om+him)
si=sin(om+him)
sss=abs(si)/si
a=1./vl0
xdh(1,1)=-a
xdh(1,3)=-a*co
xdh(1,2)=a*si
dv(1,1)=-a
dv(1,2)=a
a=1./vl1
xgh(1,1)=si-2.*vl0*rohm*sss
xgh(1,2)=-co-2.*vl0*bm*rohm*sss
co=cos(om-him)
si=sin(om-him)
xdh(2,3)=a*co
xdh(2,2)=a*si
dv(2,2)=-a
dv(2,3)=a
xgh(1,1)=xgh(1,1)+vl0*a*si
xgh(1,2)=xgh(1,2)+vl0*a*co
co=cos(os+his)
si=sin(os+his)
xdh(2,5)=a*si
xdh(2,4)=a*co
dhp(2,4)=a
a=vl0*a
xgh(1,3)=a*co
xgh(1,4)=a*si
co=cos(os-his)
si=sin(os-his)
a=1./vl2
xdh(3,5)=a*si
xdh(3,4)=-a*co
dv(3,3)=-a
dv(3,4)=a
xgh(8,3)=a*vl3*co
xgh(8,4)=-a*vl3*si
co=cos(oa+hia)
si=sin(oa+hia)
sss=abs(si)/si
xdh(3,7)=-a*co
xdh(3,6)=a*si
xgh(8,5)=-vl3*a*si+2.*vl3*roha*sss
xgh(8,6)=vl3*a*co+2.*vl3*ba*roha*sss
a=1./vl3
co=cos(oa-hia)
si=sin(oa-hia)
xdh(4,7)=a*co
xdh(4,6)=a*si
xdh(4,8)=a
dv(4,4)=-a
dv(4,5)=a
dv(5,3)=1.
xgh(8,5)=xgh(8,5)-si
xgh(8,6)=xgh(8,6)-co
xgh(1,7)=-2.*vl0
if(etm.eq.0.) xgh(1,7)=xgh(1,7)*sin(him)*cos(om)/abs(sin(om-him))
xgh(8,8)=2.*vl3
if(eta.eq.0.) xgh(8,8)=xgh(8,8)*sin(hia)*cos(oa)/abs(sin(oa-hia))
do 16 i=1,5
do 16 j=1,4
16 ghp(i,j)=xgh(i,j)
ghp(1,3)=vl0/vl1
ghp(1,4)=xgh(1,7)
a=r8ln2*r8ln2
if(etm.ne.0.)xf(1,1)=a/etvmef/etvmef
if(etm.eq.0.)xf(1,1)=a/etamef/etamef*abs(sin(om-him)/sin(om+him))
if(eta.ne.0.)xf(2,2)=a/etvaef/etvaef
if(eta.eq.0.)xf(2,2)=a/etaaef/etaaef*abs(sin(oa-hia)/sin(oa+hia))
xf(1,2)=0.
xf(2,1)=0.
do 4 j=1,5
tv(1,j)=.5*(dv(1,j)-dv(2,j))/sin(om)
4 tv(2,j)=.5*(dv(3,j)-dv(4,j))/sin(oa)
tv(1,2)=tv(1,2)-rovm*abs(om)/om*abs(cos(him))
tv(2,4)=tv(2,4)-rova*abs(oa)/oa*abs(cos(hia))
call BTAB(xf,tv,2,5,vcan)
call SUM(sv1,vcan,5,sv1)
do 18 j=1,3
18 tv1(1,j)=tv(1,j)
xf1(1)=xf(1,1)
call BTAB(xf1,tv1,1,3,vcanp)
call SUM(svp,vcanp,3,svp)
do 5 i=1,4
do 5 j=1,4
5 aco(i,j)=0.
co=0.
if(nguide.ne.0) si=3./(gamacr*slamdi)**2
do 6 i=1,4
co=co+alpha(i)
6 if(alpha(i).ne.0.) aco(i,i)=a/alpha(i)/alpha(i)/tmr/tmr
if(nguide.ne.0) aco(1,1)=aco(1,1)+si
co=co+nguide
if(co.eq.0.) go to 7
call BTAB(aco,xdh,4,10,can)
call SUM(sh1,can,10,sh1)
co=0.
acop(1,2)=0.
acop(2,1)=0.
do 11 i=1,2
co=co+alpha(i)
11 acop(i,i)=aco(i,i)
co=co+nguide
if(co.eq.0.) go to 7
do 13 i=1,2
do 13 j=1,3
13 dhp(i,j)=xdh(i,j)
call BTAB(acop,dhp,2,5,vcan)
call SUM(shp,vcan,5,shp)
7 do 9 i=1,4
do 9 j=1,4
9 aco(i,j)=0.
co=0.
do 10 i=1,4
co=co+beta(i)
10 if(beta(i).ne.0.) aco(i,i)=a/beta(i)/beta(i)/tmr/tmr
if(nguide.ne.0) aco(1,1)=aco(1,1)+si
co=co+nguide
if(co.eq.0.)go to 8
do 12 i=1,4
do 12 j=1,5
12 dv1(i,j)=dv(i,j)
call BTAB(aco,dv1,4,5,vcan)
call SUM(sv1,vcan,5,sv1)
co=0.
do 14 i=1,2
co=co+beta(i)
14 acop(i,i)=aco(i,i)
co=co+nguide
if(co.eq.0.) go to 8
do 17 i=1,2
do 17 j=1,3
17 dvp(i,j)=dv(i,j)
call BTAB(acop,dvp,2,3,vcanp)
call SUM(svp,vcanp,3,svp)
8 continue
do 117 j=1,4
ahs(3,j)=ghp(1,j)
117 ahs(2,j)=-vki/vl0*ghp(1,j)
ahs(2,1)=ahs(2,1)+vki/vl0*sin(om+him)
ahs(2,2)=ahs(2,2)-vki/vl0*cos(om+him)
call BTAB(shp,ghp,5,4,aco)
call INVERT(4,aco,4,cmhp,4)
# CALL INV1(ACO,4,AU4,CMHP) ! replaced by INVERT, J.S.
# SOLLER COLLIMATORS AND NEUTRON GUIDE HAVE BEEN ACCOUNTED FOR
call BTAB(sh1,xgh,10,8,ah8)
call INVERT(8,ah8,8,cmh,8)
# CALL INV1(AH8,8,AU8,CMH) ! ! replaced by INVERT, J.S.
call INVERT(5,sv1,5,cmv,5)
# CALL INV1(SV1,5,VCAN,CMV) ! ! replaced by INVERT, J.S.
call INVERT(3,svp,3,cmvp,3)
# CALL INV1(SVP,3,VCANP,CMVP) ! ! replaced by INVERT, J.S.
# write(*,1) 'SVP: ',SVP
# write(*,1) 'VCANP: ',VCANP
# write(*,1) 'CMVP: ',CMVP
# BRAGG CONSTRAINTS HAVE BEEN ACCOUNTED FOR
# CMH AND CMV ARE THE MODIFIED COVARIANCE MATRICES OF THE
# INDEPENDENT VARIABLES (HORIZONTAL AND VERTICAL, RESPECTIVELY)
return
end
#**************************************************************
#
SUBROUTINE HIRANG(ob,hii,hi)
# MAKES SURE THAT THE CRYSTAL INCLINATION ANGLE IS IN THE CORRECT
# RANGE (MINUS TETA BRAGG TO PI MINUS TETA BRAGG FOR POSITIVE TETA,
# MINUS PI MINUS TETA BRAGG TO MINUS TETA BRAGG FOR NEGATIVE TETA)
implicit real*8 (a-h,o-z)
INCLUDE 'trax.inc'
hi=tdr*hii
if(ob.lt.0.) go to 1
if(hi.lt.(-ob)) hi=hi+pii
if(hi.gt.(pii-ob)) hi=hi-pii
go to 2
1 if(hi.lt.(-pii-ob)) hi=hi+pii
if(hi.gt.(-ob)) hi=hi-pii
2 hii=hi/tdr
return
end
#
SUBROUTINE VANAD(a,dvn,yvn)
#***********************************************************************
# returns Vanad scan width and intensity (in rel. units)
# added by J.S. 3/6/1997
#***********************************************************************
implicit none
INCLUDE 'trax.inc'
real*8 a(4,4),emmin1(4,4),volres,aux(4,4)
real*8 dvn,yvn
common /em/emmin1,volres
call INVERT(4,a,4,aux,4)
dvn=r8ln2*sqrt(aux(4,4))
yvn=r0trax/sqrt(emmin1(4,4)) ! intensity in rel. units
end
#
real*4 FUNCTION OPTV(ro)
#***********************************************************************
# returns value to be minimized when optimizing bending radii
# added by J.S. 3/6/1997
#***********************************************************************
implicit real*8 (a-h,o-z)
INCLUDE 'trax.inc'
dimension a(4,4)
common /matrix/a
parameter(eps=1.d-30)
real*4 ro(4),b(4)
b(1)=rohm
b(2)=rovm
b(3)=roha
b(4)=rova
rohm=ro(1)/100
rovm=ro(2)/100
roha=ro(3)/100
rova=ro(4)/100
call THRAX
call VANAD(a,dvn,yvn)
if (abs(yvn).ge.eps) then
OPTV=dvn*dvn/yvn
else
OPTV=0.
endif
# write(*,*) DVN*DVN/YVN,DVN,YVN
rohm=b(1)
rovm=b(2)
roha=b(3)
rova=b(4)
end