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.
#//////////////////////////////////////////////////////////////////////
#//// ////
#//// NEutron Scattering Simulation - v.1.2, (c) J.Saroun, 1997 ////
#//// update October 1998 ////
#//////////////////////////////////////////////////////////////////////
#////
#//// Subroutines describing objects - SLIT, SOLLER,CRYSTAL
#////
#//// * LOGICAL*4 FUNCTION INSIDE(SLIT,R)
#//// * LOGICAL*4 FUNCTION CR_INSIDE(CRYST,R)
#//// * SUBROUTINE SLIT_INIT(SLIT)
#//// * LOGICAL*4 FUNCTION SLIT_GO(SLIT,NEUI,NEUF)
#//// * LOGICAL*4 FUNCTION SOLLER_GO(SOLLER,NEUI,NEUF)
#//// * LOGICAL*4 FUNCTION BENDER_GO(BENDER,NEUI,NEUF)
#//// * SUBROUTINE CRYST_INIT(OBJECT)
#//// * LOGICAL*4 FUNCTION CRYST_GO(CRYST,NEUI,NEUF,DKK)
#//// * SUBROUTINE SLIT_PRE(SLIT,R0,K0,R,K)
#//// * SUBROUTINE SLIT_POST(SLIT,R0,K0,R,K)
#//// * SUBROUTINE SLIT_PRE1(SLIT,R0,K0,R,K)
#//// * SUBROUTINE SLIT_POST1(SLIT,R0,K0,R,K)
#////
#////
#//////////////////////////////////////////////////////////////////////
#
logical*4 FUNCTION CR_INSIDE(cryst,r)
# INSIDE function for CRSYTAL object ... takes into
#
#
implicit none
INCLUDE 'structures.inc'
record /CRYSTAL/ cryst
real*8 r(3)
logical*4 INSIDE
#
#
#
CR_INSIDE=INSIDE(cryst.frame,r)
return
end
# ////////////////// End of definition - SOLLER ///////////////////
#---------------------------------------
SUBROUTINE CRYST_INIT(cr)
#---------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
integer*4 i,j
record /CRYSTAL/ cr
#/// THB and LAMBDA must be specified before !
call SLIT_INIT(cr.frame)
cr.gtot=2*pi/abs(cr.dhkl)
cr.g(1)=cr.gtot*sin(cr.chi)
cr.g(2)=0
cr.g(3)=cr.gtot*cos(cr.chi)
cr.stmch=sin(cr.thb-cr.chi)
cr.ctmch=cos(cr.thb-cr.chi)
# G-gradient for elastically bent crystal
do i=1,3
cr.mapg(i)=.true.
do j=1,3
cr.dg_dr(i,j)=0
enddo
enddo
if(cr.hmos.lt.sec.or.cr.nh.eq.1) then
cr.dg_dr(1,1)=-cos(cr.chi)*cr.gtot*cr.rh
cr.dg_dr(1,3)=sin(cr.chi)*cr.gtot*cr.rh
cr.dg_dr(3,1)=sin(cr.chi)*cr.gtot*cr.rh
cr.dg_dr(2,2)=0 ! no vertical bending
cr.dg_dr(3,3)=-cr.poi*cos(cr.chi)*cr.gtot*cr.rh
cr.mapg(1)=.true.
cr.mapg(3)=.true.
endif
#* vertical bending considered only in the case of single segment
#* and zero mosaicity
if(cr.hmos.le.sec.and.cr.nv.eq.1) then
cr.dg_dr(2,2)=-cos(cr.chi)*cr.gtot*cr.rv
cr.mapg(2)=.true.
endif
cr.qhkl=(cr.fhkl/cr.vol*cr.lambda)**2*
1 abs(cr.dhkl)*1e-5
cr.dext=cr.vol*cos(pi/2-cr.thb-cr.chi)/
1 cr.lambda/cr.fhkl*0.1
cr.ref=1 ! CRYST_REF(CR)
cr.deta=3.*cr.hmos
cr.gama(1)=cos(cr.chi)
cr.gama(3)=-sin(cr.chi)
cr.gama(2)=0
cr.nb=1
cr.rb=0.d0
return
end
#
logical*4 FUNCTION CRYST_GO(cryst,neui,neuf,dkk)
#
implicit none
INCLUDE 'const.inc'
INCLUDE 'ness_common.inc'
record /CRYSTAL/ cryst
record /NEUTRON/ neui,neuf
logical*4 INSIDE,emod,CRYST_SIMPLE
real*8 v(3),k(3),r(3),v1(3)
real*8 g0(3),g(3),dgr(3),tnorm
real*8 dkk,c,eta1,eta2,dt,gg,alfv,alfh,dah,dav,kk,tin,tout
integer*4 i,ih,iv
real*4 GASDEV1,rn,RAN1,z,z1
real*8 v3xv3
common /mode/ emod
5 format(a10,5(g12.5))
#/// CRYST.DNRND must be specified elsewhere. X(..) is a random number
#/// in the interval (-0.5,0.5)
rn=rndx(cryst.dnrnd+1)
#/// if thb=0 => special device:
neuf=neui
if (cryst.thb.eq.0) then
call SLIT_PRE(cryst.frame,neui.r,neui.k,v,k)
neuf.t=neui.t-v(3)/hovm/k(3)
do i=1,2
r(i)=v(i)-v(3)/k(3)*k(i)
end do
r(3)=0.
call SLIT_POST(cryst.frame,r,k,neuf.r,neuf.k)
cryst.frame.count=cryst.frame.count+1
CRYST_GO=.true.
#// primitive velocity selector (triangular distr.)
if (cryst.nv.eq.0) then
z1=RAN1()-0.5
z=(z1+rn)*cryst.poi ! triangular distribution, fwhm=CR.POI
kk=sqrt(v3xv3(neui.k,neui.k))
dkk=-1.+2*pi/cryst.lambda/kk*(1+z)
do i=1,3
neuf.k(i)=neuf.k(i)*(1.+dkk)
end do
# NEUF.P=NEUF.P*(1.+DKK)**2
neuf.t=neuf.t/(1.+dkk)
return
#// just go through
else if(emod) then
dkk=0.
return
else
#// crystalline filter
c=2*pi/sqrt(v3xv3(neui.k,neui.k))
dkk=-1.+(rn+0.501)*c/cryst.dhkl ! -.5<RN<.5
do i=1,3
neuf.k(i)=neuf.k(i)*(1.+dkk)
end do
# NEUF.P=NEUF.P*(1.+DKK)**2
neuf.t=neuf.t/(1.+dkk)
return
endif
endif
#/// no reflection if spin does not match magnetization
if (cryst.mag*neui.s.lt.0) goto 30
#* for analyzer in elastic mode, use another version
if (cryst.typ.eq.1) then
CRYST_GO=CRYST_SIMPLE(cryst,neui,neuf)
dkk=0.d0
return
endif
#* V&H mosaic block deviation is taken in random
#* gaussian distribution of mosaic blocs is limitted to +-3*sigma
eta1=cryst.hmos*GASDEV1(0.,3.)
eta2=cryst.vmos*GASDEV1(0.,3.)
#* G vector is corrected for the mosaic block orientation
g0(1)=cryst.g(1)+cryst.g(3)*eta1
g0(3)=cryst.g(3)-cryst.g(1)*eta1
g0(2)=cryst.g(2)+cryst.gtot*eta2
#* get nominal flight length through the crystal
do i=1,3
k(i)=0.
r(i)=0.
enddo
k(3)=1.
call SLIT_PRE(cryst.frame,r,k,v,v1)
call CR_BORDER(cryst,v,v1,tin,tout)
tnorm=(tout-tin)
#* transform neutron coordinates to the local system
call SLIT_PRE(cryst.frame,neui.r,neui.k,v,k)
#* The depth in the crystal when the reflection takes place is chosen in random
call CR_BORDER(cryst,v,k,tin,tout)
if(tin.ge.tout) goto 30 ! no intersection with the crystal
dt=tin+(rn+5.e-1)*(tout-tin)
#* move neutron to the point of reflection
neuf.t=neui.t+dt/hovm
do i=1,3
r(i)=v(i)+dt*k(i)
enddo
#* get the position and orientation of reflecting segment
ih=1
if(cryst.hmos.gt.sec) then ! mosaict crystal
ih=int((0.5+r(1)/cryst.frame.size(1))*cryst.nh)+1
dah=cryst.frame.size(1)*cryst.rh/cryst.nh
alfh=(-(cryst.nh-1)/2+ih-1)*dah
else ! elastically bent crystal
alfh=0
endif
iv=int((0.5+r(2)/cryst.frame.size(2))*cryst.nv)+1
dav=cryst.frame.size(2)*cryst.rv/cryst.nv
alfv=(-(cryst.nv-1)/2+iv-1)*dav
#* go out if the neutron doesn't hit the crystal
if((iv.lt.1).or.(iv.gt.cryst.nv).
* or.(ih.lt.1).or.(ih.gt.cryst.nh)) goto 30
#* get new orientation of G vector including segment orientation
g(1)=g0(1)-g0(3)*alfh
g(3)=g0(3)+g0(1)*alfh
if (abs(g0(3)).lt.0.01*cryst.gtot) then
g(2)=g0(2)-g0(1)*alfv ! only in Laue case
else
g(2)=g0(2)-g0(3)*alfv
endif
c=sqrt(V3xV3(g,g))
gg=cryst.gtot/c
do i=1,3
g(i)=g(i)*gg ! renormalize
enddo
gg=cryst.gtot**2
# if elastically bent crystal, include deformation
if(cryst.hmos.le.sec) then
call M3xV3_M(1,cryst.mapg,cryst.dg_dr,r,dgr)
call V3AV3(1,g,dgr,g) ! G now includes elastic deformation
gg=V3XV3(g,g)
endif
# calculate scalar product 2*(K*G)
c=2.*V3XV3(k,g)
#* wavelength and time-of-flight are corrected in order to meet
#* the Bragg condition:
dkk=-gg/c-1
do 15 i=1,3
15 k(i)=k(i)*(1+dkk)
neuf.t=neuf.t/(1+dkk)
neuf.phi=neuf.phi/(1+dkk)
# if (CRYST.FRAME.NAME(1:1).EQ.'m') then
# write(*,5) 'M: ',RN,TOUT-TIN,TNORM,DKK
# write(*,5) 'R: ',(R(I),I=1,3)
# if (mod(cryst.frame.count,10).eq.0) pause
# endif
if (INSIDE(cryst.frame,r)) then
call v3av3(1,k,g,k)
neuf.p=neui.p*(tout-tin)/tnorm
# NEUF.P=NEUI.P
if (neuf.p.ne.0) then
call SLIT_POST(cryst.frame,r,k,neuf.r,neuf.k)
cryst.frame.count=cryst.frame.count+1
CRYST_GO=.true.
else
CRYST_GO=.false.
endif
else
# write(*,*) CRYST.FRAME.NAME,' POZOR!'
# pause
neuf.p=0
CRYST_GO=.false.
endif
return
30 neuf.p=0
CRYST_GO=.false.
return
end
#
logical*4 FUNCTION CRYST_SIMPLE(cryst,neui,neuf)
#
implicit none
INCLUDE 'const.inc'
INCLUDE 'ness_common.inc'
record /CRYSTAL/ cryst
record /NEUTRON/ neui,neuf
real*8 v(3),k(3),r(3)
real*8 g0(3),dgk(3),kag(3)
real*8 b,c,eta1,eta2,dt,pp,g2,tin,tout,tau
integer*4 i
real*4 GASDEV1,rn
real*8 v3xv3,HMOS_DIST
# LOGICAL*4 FLAG
#5 format(a10,5(2x,G12.5))
#/// CRYST.DNRND must be specified elsewhere. X(..) is a random number
#/// in the interval (-0.5,0.5)
rn=rndx(cryst.dnrnd+1)
#/// if thb=0 => special device:
if(neuf.p.eq.0) goto 30
# FLAG=(CRYST.FRAME.COUNT.GT.8000)
# IF (FLAG) write(*,5) 'RND: ',CRYST.DNRND,RN
# IF (FLAG) write(*,5) 'RND: ',(RNDX(I),I=1,5)
# IF (FLAG) write(*,5) 'RND: ',(RNDX(I),I=6,10)
#* transform neutron coordinates to the local system
call SLIT_PRE(cryst.frame,neui.r,neui.k,v,k)
#* V mosaic block deviation is taken in random
#* gaussian distribution of mosaic blocs is limitted to +-3*sigma
eta2=cryst.vmos*GASDEV1(0.,3.)
#* Move neutron inside the crystal, depth is chosen in random
call CR_BORDER(cryst,v,k,tin,tout)
# IF (FLAG) THEN
# write(*,5) 'V',(V(I),I=1,3)
# write(*,5) 'K',(K(I),I=1,3)
# write(*,5) 'T1,T2',TIN,TOUT
# ENDIF
if(tin.ge.tout) goto 30 ! no intersection with the crystal
dt=tin+rn*(tout-tin)
do i=1,3
r(i)=v(i)+k(i)*dt
enddo
# IF (FLAG) write(*,5) 'R',(R(I),I=1,3),DT,RN
tout=tout-tin-rn*(tout-tin)
tin=-rn*(tout-tin)
#* get local G vector, assume no horizontal mosaic spread
call LOCALG(cryst,r,0.d0,eta2,g0)
g2=V3xV3(g0,g0)
c=g2+2*V3xV3(k,g0)
# elastically bent crystal
if(cryst.hmos.le.sec) then
call M3xV3_M(1,cryst.mapg,cryst.dg_dr,k,dgk)
call V3AV3(1,k,g0,kag)
b=2*V3xV3(kag,dgk)
if(abs(b).lt.1d-30) b=1d-30
tau=-c/b
if(tau.gt.tin.and.tau.lt.tout) then
do i=1,3
r(i)=r(i)+k(i)*tau
g0(i)=g0(i)+tau*dgk(i)
enddo
dt=dt+tau
endif
pp=1.d0
else
# mosaic crystal
eta1=-c/g2
pp=HMOS_DIST(cryst,eta1)
if(pp.le.1.d-5) goto 30
call LOCALG(cryst,r,eta1,eta2,g0)
# IF (FLAG) THEN
# write(*,5) 'ETA',ETA1,ETA2
# write(*,5) 'R',(R(I),I=1,3)
# write(*,5) 'K',(K(I),I=1,3)
# write(*,5) 'G',(G0(I),I=1,3)
# CALL V3AV3(1,K,G0,KaG)
# write(*,5) 'K+G',(KaG(I),I=1,3)
# ENDIF
endif
call V3AV3(1,k,g0,kag)
neuf.p=neui.p*pp
neuf.t=neui.t+dt/hovm
call SLIT_POST(cryst.frame,r,kag,neuf.r,neuf.k)
# IF (FLAG) THEN
# write(*,5) 'RF',(NEUF.R(I),I=1,3)
# write(*,5) 'KF',(NEUF.K(I),I=1,3)
# pause
# ENDIF
cryst.frame.count=cryst.frame.count+1
CRYST_SIMPLE=.true.
# write(*,*) 'CR OK '
# pause
return
30 neuf.p=0
CRYST_SIMPLE=.false.
# write(*,*) 'CR false '
# pause
return
end
# -------------------------------------------------------------
SUBROUTINE CR_BORDER(cr,r,k,tin,tout)
# Returns times of intersection with the crystal assembly borders,
# started at current position R and measured along K.
# All in crystal local coordinate.
# !! Time is in units [sec*h/m] i.e. length=time*K
#--------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
record /CRYSTAL/ cr
real*8 r(3),k(3),tin,tout,t1(3),t2(3),dum
integer*4 i
do i=1,3
if (abs(k(i)).gt.1.0d-8) then
t2(i)=(cr.frame.size(i)/2.d0 - r(i))/k(i)
t1(i)=(-cr.frame.size(i)/2.d0 - r(i))/k(i)
if (t1(i).gt.t2(i)) then
dum=t1(i)
t1(i)=t2(i)
t2(i)=dum
endif
else
t2(i)=1.0d30
t1(i)=-1.0d30
endif
end do
tin=max(t1(1),t1(2),t1(3))
tout=min(t2(1),t2(2),t2(3))
if (tin.gt.tout) then
tin=1d30
tout=1d30
endif
end
# -------------------------------------------------------------
SUBROUTINE LOCALG(cr,r,eta,phi,g)
# Calculate local G-vector at R
# I0(3) segment coordinates
# R0(3) segment center physical coordinates
# ETA .. horizontal tilt angle of mosaic domain
# PHI .. vertical tilt angle of mosaic domain
# -------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
record /CRYSTAL/ cr
integer*4 i,i1,i0(3)
real*8 r(3),r0(3),g(3),w(3),at(3),g0(3)
real*8 z,eta,phi,gabs
call SEGCOORD(cr,r,r0,i0)
# Calculate local G-vector (only G-gradient)
do i=1,3
g0(i)=cr.g(i)
if (cr.mapg(i)) then
w(i)=0.d0
do i1=1,3
w(i)=w(i)+cr.dg_dr(i,i1)*(r(i1)-r0(i1))
enddo
g0(i)=g0(i)+w(i)
endif
end do
gabs=sqrt(g0(1)**2+g0(2)**2+g0(3)**2)
# Add segment tilt angle and vertical mosaic spread
call SEGTILT(cr,i0,at)
g(1)=g0(1)-g0(3)*(at(1)+at(3))
g(3)=g0(3)+g0(1)*(at(1)+at(3))
g(2)=g0(2)-g0(3)*at(2)+gabs*phi
# Add the angle of the mosaic block
do i=1,3
g(i)=g(i)+gabs*cr.gama(i)*eta
end do
# Renormalize
z=gabs/sqrt(g(1)**2+g(2)**2+g(3)**2)
do i=1,3
g(i)=g(i)*z
enddo
end
# ---------------------------------------------------------------------------
SUBROUTINE SEGCOORD(cr,r,r0,i0)
# return coordinates of the segment R0, in which the particle at R resides
# ---------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
record /CRYSTAL/ cr
real*8 r(3),r0(3),half
integer*4 i0(3)
parameter(half=0.5d0)
# ih,iv,ib are the integre-coordinates of the closest segment
i0(1)=nint((r(1)/cr.frame.size(1)+half)*cr.nh+half)
i0(2)=nint((r(2)/cr.frame.size(2)+half)*cr.nv+half)
i0(3)=nint((r(3)/cr.frame.size(3)+half)*cr.nb+half)
# get physical coordinates of the segment center
r0(1)=cr.frame.size(1)*(1.d0*(i0(1)-half)/cr.nh-half)
r0(2)=cr.frame.size(2)*(1.d0*(i0(2)-half)/cr.nv-half)
r0(3)=cr.frame.size(3)*(1.d0*(i0(3)-half)/cr.nb-half)
end
# ----------------------------------------------------------
SUBROUTINE SEGTILT(cr,i0,at)
# return tilt angles
# ----------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
record /CRYSTAL/ cr
real*8 at(3),da
integer*4 i0(3)
if(i0(1).gt.0.and.i0(1).le.cr.nh) then
da=cr.frame.size(1)*cr.rh/cr.nh
at(1)=(-(cr.nh-1)/2.d0+i0(1)-1.d0)*da
else
at(1)=0.d0
endif
if(i0(2).gt.0.and.i0(2).le.cr.nv) then
da=cr.frame.size(2)*cr.rv/cr.nv
at(2)=(-(cr.nv-1)/2.d0+i0(2)-1.d0)*da
else
at(2)=0.d0
endif
if(i0(3).gt.0.and.i0(3).le.cr.nb) then
da=cr.frame.size(3)*cr.rb/cr.nb
at(3)=(-(cr.nb-1)/2.d0+i0(3)-1.d0)*da
else
at(3)=0.d0
endif
end
#---------------------------------------------------------------
real*8 FUNCTION HMOS_DIST(cr,x)
# Distribution of horizontal angular deviation of mosaic blocks
# --------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
record /CRYSTAL/ cr
real*8 x
if (cr.hmos.gt.sec) then
HMOS_DIST=exp(-0.5*(x/cr.hmos)**2)
else
if(abs(x).le.sec) then
HMOS_DIST=1
else
HMOS_DIST=0
endif
endif
end