Source module last modified on Mon, 28 Mar 2005, 21:52;
HTML image of Fortran source automatically generated by
for2html on Mon, 23 May 2005, 21:29.
#//////////////////////////////////////////////////////////////////////
#//// ////
#//// NEutron Scattering Simulation - v.2.0, (c) J.Saroun, 1997 ////
#//// ////
#//////////////////////////////////////////////////////////////////////
#////
#//// Subroutines describing objects - SAMPLES
#////
#////
#//////////////////////////////////////////////////////////////////////
#
logical*4 FUNCTION PWD_GO(pwd,neui,neuf,q)
#
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
INCLUDE 'randvars.inc'
record /SLIT/ pwd
record /NEUTRON/ neui,neuf
logical*4 log1,map(3),SAM_BOARDER
real*8 v(3),k(3),r(3),q,x1(3),x2(3),rt1(3,3),rt2(3,3)
real*8 kf,pp,t1,t2,dz,dz0,stb,ctb,s2tb,c2tb,dal,dfi,sinn,coss
real*4 RAN1
integer*4 i
data map /.true.,.true.,.true./
log1=.true.
pp=1
call SLIT_PRE1(pwd,neui.r,neui.k,v,k)
#/// move neutron to the centre of the SLIT:
neuf.t=neui.t-v(3)/hovm/k(3)
do 10 i=1,2
10 r(i)=v(i)-v(3)/k(3)*k(i)
r(3)=0.
# write(*,*) 'PWD: '
#
if (SAM_BOARDER(pwd,r,k,t1,t2)) then
# DZ0=2*SQRT((PWD.SIZE(1)/2)**2-R(1)**2)
kf=sqrt(k(1)**2+k(2)**2+k(3)**2)
dz0=kf*(t2-t1)
dz=(RAN1()-0.5)*dz0
pp=dz0*0.01 ! suppose sigma=0.1 cm^-1
r(3)=dz
# write(*,*) DZ0,DZ,T1,T2
neuf.t=neuf.t+dz/hovm/k(3)
kf=sqrt(k(1)**2+k(2)**2+k(3)**2)
stb=q/2/kf
if (abs(stb).ge.1) then
neuf.p=0
PWD_GO=.false.
return
endif
ctb=sqrt(1-stb**2)
s2tb=2*stb*ctb
c2tb=sqrt(1-s2tb**2)
if (stb**2.gt.0.5) c2tb=-c2tb
dal=atan(k(1)/k(3))
dfi=atan(k(2)/sqrt(k(1)**2+k(3)**2))
call MK_ROT3(2,dal,rt1)
call MK_ROT3(1,dfi,rt2)
sinn=sin(xrnd(7))
coss=sqrt(1-sinn**2)
if(abs(xrnd(7)).gt.pi/2) coss=-coss
x1(1)=kf*s2tb*coss
x1(2)=kf*s2tb*sinn
x1(3)=kf*c2tb
call M3XV3(-1,map,rt1,x1,x2)
call M3XV3(-1,map,rt2,x2,k)
do i=1,3
neuf.r(i)=r(i)
neuf.k(i)=k(i)
end do
neuf.p=neui.p*pp
neuf.s=neui.s
pwd.count=pwd.count+1
else
log1=.false.
neuf.p=0
end if
PWD_GO=log1
return
end
#
logical*4 FUNCTION SAM_BOARDER(sam,r,k,t1,t2)
# return cross-section times with sample borders in T1, T2
#
implicit none
INCLUDE 'structures.inc'
record /SLIT/ sam
real*8 k(3),r(3),t1,t2,t3,t4,z,a,b,c,dd
# write(*,*) SAM.SIZE(1),SAM.SIZE(2),SAM.SIZE(3)
a=(k(1)*2/sam.size(1))**2+(k(3)*2/sam.size(3))**2
b=2*(k(1)*r(1)*4/sam.size(1)**2+k(3)*r(3)*4/sam.size(3)**2)
c=(r(1)*2/sam.size(1))**2+(r(3)*2/sam.size(3))**2-1.0
dd=b**2-4*a*c
# write(*,*) 'ABCD:',A,B,C,DD
if(dd.gt.0) then
t1=(-b-sqrt(dd))/2/a
t2=(-b+sqrt(dd))/2/a
if(abs(k(2)).gt.1e-10) then
t3=(-sam.size(2)/2-r(2))/k(2)
t4=(sam.size(2)/2-r(2))/k(2)
if(t4.lt.t3) then
z=t3
t3=t4
t4=z
endif
else
t3=-1e31
t4=1e31
endif
if (t1.lt.t3) t1=t3
if (t2.gt.t4) t2=t4
if(t1.lt.t2) then
SAM_BOARDER=.true.
else
SAM_BOARDER=.false.
endif
else
# A=SQRT(R(1)**2+R(3)**2)
# write(*,*) 'SAM_BORDER: det<=0 ',A,R(2)
SAM_BOARDER=.false.
endif
end
#
logical*4 FUNCTION VAN_GO(van,neui,neuf,q)
#
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
INCLUDE 'randvars.inc'
real*8 sscat,sabs
#// scattering and absorption cross-section for k=2.664A-1
parameter (sscat=0.0362,sabs=0.0476)
# PARAMETER (SSCAT=0.0362,SABS=0.)
record /SLIT/ van
record /NEUTRON/ neui,neuf
logical*4 log1,SAM_BOARDER
real*8 v(3),k(3),r(3),q,x1(3)
real*8 kf,pp,p0,t1,t2,delta,stb,ctb,s2tb,c2tb,sinn,coss,ksi
real*4 dt
real*4 RAN1
integer*4 i
log1=.true.
pp=1
call SLIT_PRE1(van,neui.r,neui.k,v,k)
#/// move neutron to the centre of the SLIT:
neuf.t=neui.t-v(3)/hovm/k(3)
do 10 i=1,2
10 r(i)=v(i)-v(3)/k(3)*k(i)
r(3)=0.
if (SAM_BOARDER(van,r,k,t1,t2)) then
#/// move neutron to the sample entry:
neuf.t=neui.t+t1
do i=1,3
r(i)=r(i)+k(i)*t1
enddo
delta=t2-t1
kf=sqrt(k(1)**2+k(2)**2+k(3)**2)
ksi=RAN1()
p0=1.-exp(-sscat*kf*delta) ! scattering probability
dt=-log(1-ksi*p0)/sscat/kf
#/// move to the scattering point:
neuf.t=neuf.t+dt
do i=1,3
r(i)=r(i)+k(i)*dt
enddo
stb=q/2/kf
if (abs(stb).ge.1) then
neuf.p=0
VAN_GO=.false.
return
endif
ctb=sqrt(1-stb**2)
s2tb=sin(xrnd(8)+2*atan(stb/ctb))
c2tb=sqrt(1-s2tb**2)
if (stb**2.gt.0.5) c2tb=-c2tb
sinn=sin(xrnd(7))
coss=sqrt(1-sinn**2)
if(abs(xrnd(7)).gt.pi/2) coss=-coss
x1(1)=kf*s2tb*coss
x1(2)=kf*sinn
x1(3)=kf*c2tb*coss
do i=1,3
neuf.r(i)=r(i)
neuf.k(i)=x1(i)
end do
if (SAM_BOARDER(van,r,x1,t1,t2)) then
pp=exp(-sabs*(dt+t2)*2.664) ! absorption
neuf.p=neui.p*p0*pp/4/pi ! convert Sigma to dSigma/dOmega
neuf.s=neui.s
log1=(neuf.p.gt.0.d0)
if(log1) van.count=van.count+1
endif
else
log1=.false.
neuf.p=0
end if
VAN_GO=log1
return
end
#
logical*4 FUNCTION VAN_TRANS(van,neui,neuf)
#
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
real*8 sscat,sabs
parameter (sscat=0.0362,sabs=0.0476)
# PARAMETER (SSCAT=0.0362,SABS=0.)
record /NEUTRON/ neui,neuf
record /SLIT/ van
logical*4 SAM_BOARDER
real*8 v(3),k(3)
real*8 kf,t1,t2,exp
integer*4 i
call SLIT_PRE1(van,neui.r,neui.k,v,k)
#/// move neutron to the centre of the SLIT:
neuf.t=neui.t-v(3)/hovm/k(3)
do 10 i=1,2
neuf.r(i)=v(i)-v(3)/k(3)*k(i)
10 neuf.k(i)=k(i)
neuf.r(3)=0.
if (SAM_BOARDER(van,neuf.r,k,t1,t2)) then
#/// transmission through the sample:
kf=sqrt(k(1)**2+k(2)**2+k(3)**2)
neuf.p=neui.p*exp(-(sabs+sscat)*kf*(t2-t1))
neuf.s=neui.s
else
neuf.p=neui.p
end if
van.count=van.count+1
VAN_TRANS=.true.
return
end
#
logical*4 FUNCTION SAM_GO(van,neui,neuf)
# Spreads n. to all Q's and energies, to get resolution function R(Q,E)
# except: ISQOM>0 defines the scatterig function
#
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'randvars.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
record /SLIT/ van
record /NEUTRON/ neui,neuf
logical*4 log1,SAM_BOARDER,ESAM_GO
real*8 v(3),k(3),r(3),x1(3)
real*8 kf,pp,t1,t2,delta,s2tb,c2tb,sinn,coss,z
real*4 dt
real*4 RAN1
integer*4 i
real*8 SQE_AMAG,SQE_AMAG1
if (emode.eq.1) then
SAM_GO= ESAM_GO(van,neui,neuf)
return
endif
pp=1
call SLIT_PRE1(van,neui.r,neui.k,v,k)
#/// move neutron to the centre of the SLIT:
neuf.t=neui.t-v(3)/hovm/k(3)
do 10 i=1,2
10 r(i)=v(i)-v(3)/k(3)*k(i)
r(3)=0.
if (SAM_BOARDER(van,r,k,t1,t2)) then
#/// move neutron to the scattering point:
delta=t2-t1
dt=RAN1()*delta
neuf.t=neui.t+t1+dt
do i=1,3
r(i)=r(i)+k(i)*(t1+dt)
enddo
z=xrnd(8)+omega
s2tb=sin(z)
c2tb=sqrt(1-s2tb**2)
if (abs(z).gt.pi/2) c2tb=-c2tb
sinn=sin(xrnd(7))
coss=sqrt(1-sinn**2)
if(abs(xrnd(7)).gt.pi/2) coss=-coss
kf=stp.kf+xrnd(9)
x1(1)=kf*s2tb*coss
x1(2)=kf*sinn
x1(3)=kf*c2tb*coss
do i=1,3
neuf.r(i)=r(i)
neuf.k(i)=x1(i)
end do
# multiply by local thickness * kf/ki
neuf.p=neui.p*abs(delta)*kf
neuf.s=neui.s
if (isqom.eq.1) then ! magnon in antiferomag.
neuf.p=neuf.p*SQE_AMAG(neui.k,neuf.k)
else if (isqom.eq.3) then ! magnon in antiferomag.
neuf.p=neuf.p*SQE_AMAG1(neui.k,neuf.k)
endif
log1=(neuf.p.gt.0.d0)
if(log1) van.count=van.count+1
else
log1=.false.
neuf.p=0
end if
SAM_GO=log1
return
end
#------------------------------------------
logical*4 FUNCTION ESAM_GO(van,neui,neuf)
# Spread n. to all Q's with zero energy transfer, to get elastic resol. function R(Q)
#
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'randvars.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
real*8 sscat,sabs
#// scattering and absorption cross-section for k=2.664A-1
parameter (sscat=0.0362,sabs=0.0476)
# RECORD /RANDFIELD/ RNDLIST
record /SLIT/ van
record /NEUTRON/ neui,neuf
logical*4 log1,SAM_BOARDER
real*8 v(3),k(3),r(3),x1(3)
real*8 kf,ki,pp,t1,t2,delta,s2tb,c2tb,sinn,coss,z
real*4 dt
real*4 RAN1
integer*4 i
#1 format(a10,2x,8(G10.4,2x))
# DBGREF=.TRUE.
pp=1
call SLIT_PRE1(van,neui.r,neui.k,v,k)
#/// move neutron to the centre of the SLIT:
neuf.t=neui.t-v(3)/hovm/k(3)
do 10 i=1,2
10 r(i)=v(i)-v(3)/k(3)*k(i)
r(3)=0.
if (SAM_BOARDER(van,r,k,t1,t2)) then
# if (DBGREF) write(*,1) VAN.NAME,R,K,T1,T2
#/// move neutron to the sample entry:
neuf.t=neui.t+t1
do i=1,3
r(i)=r(i)+k(i)*t1
enddo
delta=t2-t1
dt=RAN1()*delta
#/// move to the scattering point:
neuf.t=neuf.t+dt
do i=1,3
r(i)=r(i)+k(i)*dt
enddo
z=xrnd(8)+omega
s2tb=sin(z)
c2tb=sqrt(1-s2tb**2)
if (abs(z).gt.pi/2) c2tb=-c2tb
sinn=sin(xrnd(7))
coss=sqrt(1-sinn**2)
if(abs(xrnd(7)).gt.pi/2) coss=-coss
ki=sqrt(k(1)**2+k(2)**2+k(3)**2)
#// kf=SQRT(ki**2-STP.E/HSQOV2M)
kf=ki
x1(1)=kf*s2tb*coss
x1(2)=kf*sinn
x1(3)=kf*c2tb*coss
# if (DBGREF) write(*,1) VAN.NAME,XRND(7)/deg,XRND(8)/deg,
# & OMEGA/deg,X1
# if (DBGREF) pause
do i=1,3
neuf.r(i)=r(i)
neuf.k(i)=x1(i)
end do
neuf.p=neui.p*abs(delta)*kf ! multiply by local thickness*kf/ki
if (isqom.eq.2) then ! Vanadium
neuf.p=neuf.p*sscat/(4.d0*pi)
log1=SAM_BOARDER(van,r,neuf.k,t1,t2)
neuf.p=neuf.p*exp(-sabs*(dt+t2)*kf) ! absorption
if (RAN1().le.0.666666) then ! spin flip
neuf.s=-neui.s
else
neuf.s=neui.s
endif
else
neuf.s=neui.s
endif
log1=(neuf.p.gt.0.d0)
if(log1) van.count=van.count+1
else
log1=.false.
neuf.p=0
end if
ESAM_GO=log1
return
end
#--------------------------------------------------------------
logical*4 FUNCTION BRAGG_GO(van,neui,neuf)
# no sample, neutrons just scattered at nominal (Q,E) with P=1
#
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'randvars.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
record /SLIT/ van
record /NEUTRON/ neui,neuf
logical*4 log1,SAM_BOARDER
real*8 v(3),k(3),r(3)
real*8 pp,t1,t2,delta
real*8 vq(3),wq(3)
real*4 dt
real*4 RAN1
integer*4 i,j
pp=1
call SLIT_PRE1(van,neui.r,neui.k,v,k)
#/// move neutron to the centre of the SLIT:
neuf.t=neui.t-v(3)/hovm/k(3)
do 10 i=1,2
10 r(i)=v(i)-v(3)/k(3)*k(i)
r(3)=0.
if (SAM_BOARDER(van,r,k,t1,t2)) then
#/// move neutron to the sample entry:
neuf.t=neui.t+t1
do i=1,3
r(i)=r(i)+k(i)*t1
enddo
delta=t2-t1
dt=RAN1()*delta
#/// move to the scattering point:
neuf.t=neuf.t+dt
do i=1,3
r(i)=r(i)+k(i)*dt
enddo
#* transform from C&N coord.
do i=2,3
wq(i)=0
enddo
wq(1)=stp.q
do i=1,3
vq(i)=0
do j=1,3
vq(i)=vq(i)+mlc(i,j)*wq(j)
enddo
enddo
#// kf=SQRT(ki**2-STP.E/HSQOV2M)
do i=1,3
neuf.r(i)=r(i)
neuf.k(i)=k(i)+vq(i)
end do
# NEUF.P=NEUI.P*ABS(DELTA)*kf ! multiply by local thickness*kf/ki
neuf.p=neui.p
neuf.s=neui.s
log1=(neuf.p.gt.0.d0)
if(log1) van.count=van.count+1
else
log1=.false.
neuf.p=0
end if
BRAGG_GO=log1
return
end
#----------------------------------------------------------
real*8 FUNCTION SQE_AMAG(vki,vkf)
# inelastic scattering cross-section
# spin waves in antiferromagnetics
#----------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'randvars.inc'
real*8 bf
parameter(bf=11.609) !Bose factor conversion kT -> meV
real*8 vki(3),vkf(3)
real*8 kf0,ki0,vq(3),wq(3),tau(3),ve
real*8 r1,nw,eq,z,wqs
integer*4 i,j
real*8 cs,egap,gamma,sig,eps,temp
data cs,egap,gamma,sig,eps,temp /1.d-5,1.65,0.2,1000,0.6,1.6/
# DATA ITMP/0/
kf0=vkf(1)**2+vkf(2)**2+vkf(3)**2
ki0=vki(1)**2+vki(2)**2+vki(3)**2
#* get Q in lab. coord.
do i=1,3
vq(i)=vkf(i)-vki(i)
enddo
#* subtract nominal values from k vectors
vq(3)=vq(3)-stp.kf*comega+stp.ki
vq(1)=vq(1)-stp.kf*somega
#* transform to C&N coord.
do i=1,3
wq(i)=0
do j=1,3
wq(i)=wq(i)+mlc(j,i)*vq(j)
enddo
enddo
#* transform to r.l. coord.
do i=1,3
vq(i)=0
do j=1,3
vq(i)=vq(i)+mrc(i,j)*wq(j)
enddo
enddo
#* get energy transfer
ve=hsqov2m*(ki0-kf0)
#* get TAU
do i=1,3
tau(i)=nint(qhkl(i))
enddo
#* get propagation vector
do i=1,3
vq(i)=vq(i)+qhkl(i)-tau(i)
enddo
#// bose factor
z=exp(-abs(ve)*bf/temp)
nw=z/(1-z)
if (ve.gt.0) nw=nw+1
#// excitation energy
wqs=egap**2*(1+sig*(vq(1)**2+vq(2)**2+eps*vq(3)**2)) ! dispersion law
#// dynamic structure factor
eq=wqs+gamma**2
z=(ve**2-eq)**2+(ve*gamma)**2
r1=ve*gamma/((ve**2-eq)**2+(ve*gamma)**2)
z=cs*nw*r1
SQE_AMAG=z
end
#----------------------------------------------------------
real*8 FUNCTION SQE_AMAG1(vki,vkf)
# inelastic scattering cross-section
# dispersion with finite width and curvature along QH
#----------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'randvars.inc'
real*8 bf
parameter(bf=11.609) !Bose factor conversion kT -> meV
real*8 vki(3),vkf(3)
real*8 kf0,ki0,vq(3),wq(3),tau(3),ve
real*8 r1,nw,eq,z,wqs
integer*4 i,j
real*8 cs,egap,gamma,sig,eps,temp
data cs,egap,gamma,sig,eps,temp /1.d-5,1.65,0.2,1000,0,1.6/
# DATA ITMP/0/
kf0=vkf(1)**2+vkf(2)**2+vkf(3)**2
ki0=vki(1)**2+vki(2)**2+vki(3)**2
#* get Q in lab. coord.
do i=1,3
vq(i)=vkf(i)-vki(i)
enddo
#* subtract nominal values from k vectors
vq(3)=vq(3)-stp.kf*comega+stp.ki
vq(1)=vq(1)-stp.kf*somega
#* transform to C&N coord.
do i=1,3
wq(i)=0
do j=1,3
wq(i)=wq(i)+mlc(j,i)*vq(j)
enddo
enddo
#* transform to r.l. coord.
do i=1,3
vq(i)=0
do j=1,3
vq(i)=vq(i)+mrc(i,j)*wq(j)
enddo
enddo
#* get energy transfer
ve=hsqov2m*(ki0-kf0)
#* get TAU
do i=1,3
tau(i)=nint(qhkl(i))
enddo
#* get propagation vector
do i=1,3
vq(i)=vq(i)+qhkl(i)-tau(i)
enddo
#// bose factor
z=exp(-abs(ve)*bf/temp)
nw=z/(1-z)
if (ve.gt.0) nw=nw+1
#// excitation energy
wqs=egap**2*(1+sig*(vq(1)**2)) ! dispersion law
#// dynamic structure factor
eq=wqs+gamma**2
z=(ve**2-eq)**2+(ve*gamma)**2
r1=ve*gamma/((ve**2-eq)**2+(ve*gamma)**2)
z=cs*nw*r1
SQE_AMAG1=z
end