Source module last modified on Sat, 21 May 2005, 16:26;
HTML image of Fortran source automatically generated by
for2html on Mon, 23 May 2005, 21:29.
#f//////////////////////////////////////////////////////////////////////
#//// ////
#//// NEutron Scattering Simulation - v.1.0, (c) J.Saroun, 1996 ////
#//// ////
#//////////////////////////////////////////////////////////////////////
#////
#//// Subroutines for handling events & simple command interpreter
#////
#//// * SUBROUTINE NESS_LOOP
#//// * SUBROUTINE READCOM(ICOMM,IOE)
#//// * SUBROUTINE NESS(ITASK)
#//// * LOGICAL FUNCTION SAFETY_POOL()
#//// * SUBROUTINE SWPOOL
#//// * SUBROUTINE MAXV_UPD(ITASK)
#//// * SUBROUTINE RANDFILL
#//// * SUBROUTINE RESINT(ICOM,VAL,KI,R,KF)
#//// * SUBROUTINE NESS_RUN(ICOM,NCNT,NEVENT)
#//// * SUBROUTINE VALID(ICOM,NCNT)
#//// * SUBROUTINE ISAMPLE(ICOM,VAL)
#//// * BLOCK DATA
#////
#////
#//////////////////////////////////////////////////////////////////////
#-------------------------------
SUBROUTINE NESSEND
# NESS destructor
#-------------------------------
implicit none
call EVARRAY(-1,0)
call NSTORE_FREE
return
end
#------------------------------------------------------------
SUBROUTINE NESS(itask,dnev)
#
# **** Main procedure for MC ray-tracing ****
#
# ITASK=1 ... inelastic scattering, TAS resolution
# ITASK=2 ... sample -> source
# ITASK=3 ... source -> sample
# ITASK=4 ... sample -> source + sample(powder) -> detector (no analyser)
# ITASK=5 ... source -> sample + sample(powder) -> detector (no analyser)
# ITASK=6 ... sample -> source + sample (Vanad) -> monitor(IMONIT)
# ITASK=7 ... source -> monitor(IMONIT)
#------------------------------------------------------------
implicit none
integer*4 itask,maxev
parameter(maxev=500000)
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'randvars.inc'
record /STATI/ cov_qe
integer*4 nevent,ncnt,nctot
real*8 ntot,nout
common /result/ cov_qe
real*8 e4(4),e16(crnd)
logical*4 verbose
integer*4 nev
common /mcsetting/ verbose,nev
real t1,t2,t3
real*8 volresmc,resmatcn(4,4),resmatrl(4,4)
real*8 bragcn(4),vancn(4),bragrl(4),vanrl(4)
common /resmatpar/ volresmc,resmatcn,resmatrl,bragcn,vancn,
* bragrl,vanrl
record /NEUTRON/ neui,neuf,neui1,neuf1
common /neuif/ neui,neuf,neui1,neuf1
logical*4 timeout,cycle2
integer*4 i,j,it2,imess,iupd,iswpool,mess,itask1,nsum,isum
integer*4 nev1,nev2,nevent1
real*8 z,z1,vol1,vol2,smon1,znev,dnev
real*8 psum(5),psum1(5),sum1,dsum,dsum1
real*8 dcnt,dcntf,dcnt1,dcntf1,scnt1,scntf1
real*8 suma,center,fwhm
real*4 par(3),dpar(3),CHI2SPC,tol
external CHI2SPC
equivalence(e16,e4)
save ncnt
# write(*,*) ' NESS: ',ITASK,DNEV
1 format( ' Incident intensity : ',g10.4, ' +- ',g10.4,
* ', flux : ',g10.4, ' +- ',g12.4,
* ' , <dEi> = ',g10.4,a5)
2 format( ' Incident intensity : ',g10.4, ' +- ',g10.4,
* ' ,<dEi> = ',g10.4,a5)
3 format( ' Wait please ... ')
4 format( ' Time expected: ',f10.1, ' s ')
5 format( ' Timeout reached at ',f10.1, ' s')
6 format( ' Number of events [x 1000] : ',$)
7 format(10(1x,e10.4))
8 format( ' Time spent: ',f10.1, ' s ')
9 format( ' Events: ',i8, ' target, ',i8, ' final, ',f12.0, ' total')
12 format( ' Safety pool hits: ',16(1x,i3))
13 format( ' Monitor',i2, ': Intensity: ',g12.6, ' +- ',g10.4)
14 format( ' suma: ',g10.4, ' +- ',g9.3, ' fwhm: ',g10.4, ' spread: ',
* g10.4, ' center:',g10.4)
15 format( ' Vanad width: ',g10.4, ' shift: ',g10.4, ' ',a5,
* ' Int: ',g10.4, ' cm*sr*meV')
16 format( ' Width: ',g10.4, ' Int: ',g10.4)
17 format( ' Powder peak width: ',g10.4, ' fwhm: ',g10.4,
* ' center:',g10.4, ' Int: ',g10.4)
18 format( ' Gaussian peak fwhm: ',g10.4, ' center: ',g10.4,
* ' Imax: ',g10.4)
19 format( ' suma: ',g10.4, ' fwhm: ',g10.4, ' spread: ',
* g10.4, ' center:',g10.4)
20 format( ' TAS RESOLUTION: '/,
* ' I(mon.)= ',g12.5, ' I(anal.)= ',g12.5,/,
* ' Cooper & Nathans : ',/,
* ' Volume :',g12.5,/,
* ' Bragg: ',4(2x,g12.5),/,
* ' Vanad: ',4(2x,g12.5),/,
* ' Reciprocal lattice: ',/,
* ' Bragg: ',4(2x,g12.5),/,
* ' Vanad: ',4(2x,g12.5))
21 format( ' TAS RESOLUTION: '/,
* ' I(tot.)= ',g12.5,/,
* ' Cooper & Nathans : ',/,
* ' Volume : ',g12.5,/,
* ' Bragg: ',4(2x,g12.5),/,
* ' Vanad: ',4(2x,g12.5),/,
* ' Reciprocal lattice: ',/,
* ' Bragg: ',4(2x,g12.5),/,
* ' Vanad: ',4(2x,g12.5))
22 format( '.',$)
#// Get required number of events (NEV)
if (dnev.le.0) then
if(dnev.lt.0) then
znev=abs(dnev)
else
140 write(*,6)
read(sinp,*,err=140) znev
endif
verbose=.true.
else
znev=dnev
verbose=.false.
endif
nev=int(znev*1000)
nev=max(nev,100)
if(nev.gt.maxev) then
write(*,*) 'Maximum number of events is ',maxev
nev=maxev
endif
nev1=nev
#// reset hit counters
spcn=0
do i=1,rndlist.dim
hit(i)=0
end do
#// adjust report targets
it2=500
imess=1000
iupd=500
iswpool=max(5000,int(nev1/2))
nout=1.d+4*nev1
if(nev1.ge.2*iupd.and.iopt.eq.1) then
mess=0
else
mess=1
endif
#// reset event counters
nevent=0
ncnt=0
ntot=0.d+0
nctot=0
#// no sampling optimization in debug mode
if (idbg.ge.2) then
iupd=nev1+10
iswpool=nev1+10
nout=1.d+7*nev1
endif
#// some tasks require 2-cycle process
cycle2=(.false..or.
& itask.eq.8.or. ! TAS
& itask.eq.9.or. ! PWD
& itask.eq.10) ! PWDS
#// don't switch the pool off for TAS resolution function simulations
if (itask.eq.8.or.itask.eq.9) iswpool=nev1+10
#// Allocate memory for monitor events (monitor counts)
call EVARRAY(0,1,nev1) ! (X,Y,E,time)
call EVARRAY(0,0,nev1) ! (kx..kz,spin)
#// Special storage for some tasks: (r,k,s,p) for initial and final neutrons
if (itask.eq.1.or. ! TAS
& itask.eq.2.or. ! NFLUX
& cycle2) then
call NSTORE_ALLOCATE(nev1)
endif
itask1=itask
if (cycle2) itask1=2 ! cycle 1 accumulates events at the sample
#// Initialize TAS components for cycle 1
call SPEC_INI(0,itask1)
#// reset max. values of random variables
call MAXV_UPD(0)
t1=secnds(0.0) ! timestamp 1
timeout=.false.
#----------------------- Begin cycle 1 ------------------------------
if (verbose) write(*,3)
do while ((ncnt.lt.nev1).and.(.not.timeout))
170 ntot=ntot+1.d+0
#// show progress dots ....
z=ntot/50000.
if(z-int(z).eq.0) write(*,22)
# Z=NTOT/50000.D0
# if (Z-INT(Z).EQ.0) then
# dbgref=(Z-INT(Z).EQ.0)
dbgref=.false.
# if (dbgref) then
# write(*,*) 'NESS ',NEVENT,NTOT,MON.FRAME.COUNT,SOU.COUNT
# CALL GETSTATE(ITASK1,NEVENT)
# WRITE(*,12) (HIT(I),I=1,RNDLIST.DIM)
# pause
# endif
timeout=(ntot.gt.nout.or.((nctot+1)*1d5.lt.ntot))
#// Reset counters if NEVENT=0
if (nevent.eq.0) then
scnt=0.
scntf=0.
scnt1=0.
scntf1=0.
dcnt=0.
dcntf=0.
dcnt1=0.
dcntf1=0.
vol1=1
do i=1,rndlist.dim
vol1=vol1*rndlist.limits(i)
end do
nsum=0
isum=0
do i=1,5
psum(i)=0
psum1(i)=0
enddo
nevent1=0
vol2=0
call SPEC_INI(1,itask1)
endif
#// Copy counter values to the counters for the 1st part of simulation
#// (safety pool ON)
if (ncnt.le.iswpool) then
nev2=nevent
scnt1=scnt
scntf1=scntf
dcnt1=dcnt
dcntf1=dcntf
endif
#// remember old counter values
i=ncnt
j=nevent
z=scnt
z1=scntf
#// Trace neutron through the instrument
call NESS_RUN(ncnt,nevent,itask1)
if (nevent.eq.0.and.(.not.timeout)) goto 170
#// increment counters
isum=isum+ncnt-i
nevent1=nevent1+nevent-j
nctot=nctot+ncnt-i
dcnt=dcnt+scnt-z
dcntf=dcntf+scntf-z1
if(ncnt.eq.it2) t2=secnds(t1) ! timestamp 2
#/// When IUPD events were accumulated, covariance matrix is
#// calculated and simulation restarts with new limits and covariances
if((ncnt.eq.iupd).and.(mess.lt.1)) then
call COV_NEW
call MAXV_UPD(0)
mess=1
nevent=0
ncnt=0
if (verbose) write(*,12) (hit(i),i=1,rndlist.dim)
do i=1,rndlist.dim
hit(i)=0
end do
endif
#// When IMESS events were accumulated, total time is estimated:
if((mess.lt.2).and.(ncnt.ge.imess)) then
mess=2
t3=secnds(t1)
z=t2+(t3-t2)*(nev1-imess*1.)/it2
if (verbose) write(*,4) z
endif
#// After ISWPOOL events, safety pool is switched off.
#// VOL1 and VOL2 are the sampling volumes before and after the switch
if((ncnt.eq.iswpool).and.(mess.lt.3)) then
if (verbose) write(*,*) 'SWP: Initial volume: ',vol1
call MAXV_UPD(2)
call SWPOOL
vol2=1
do i=1,rndlist.dim
vol2=vol2*rndlist.limits(i)
end do
if (verbose) write(*,*) 'SWP: Final volume: ',vol2
mess=3
endif
#// Get partial sums for error estimation
if (isum.eq.int(nev1/5.d0).and.nevent1.gt.0) then
nsum=nsum+1
psum(nsum)=(vol1*dcnt1+vol2*(dcnt-dcnt1))/nevent1
psum1(nsum)=(vol1*dcntf1+vol2*(dcntf-dcntf1))/nevent1
isum=0
dcnt=0.
dcntf=0.
dcnt1=0.
dcntf1=0.
nevent1=0
endif
end do
if (verbose) write(sout,*)
#----------------------- End cycle 1------------------------------
if(ncnt.ge.iswpool) call SWPOOL ! safety pool back ON
#// NSEED index is used to point to an incident neutron taken from KSTACK storage
#// This event is further processed in the cycle 2
nseed=0
t3=secnds(t1) ! timestamp 3
#// timeout reached => print status and exit
if (timeout) then
write(sout,5) t3
write(sout,9) ncnt,nevent,ntot
call GETSTATE(itask1,nevent)
write(*,*)
iinc=0
diinc=0
i3ax=0
di3ax=0
ipwd=0
dipwd=0
return
endif
#// report status
if (verbose) then
write(sout,8) t3
write(sout,9) ncnt,nevent,ntot
write(*,12) (hit(i),i=1,rndlist.dim)
write(*,*)
call GETSTATE(itask1,nevent)
write(*,*)
endif
#// calculate norms
sum=(vol1*scnt1+vol2*(scnt-scnt1))/nevent
sum1=(vol1*scntf1+vol2*(scntf-scntf1))/nevent
sum=sum*1.d14/100 ! mm2-> cm2
sum1=sum1*1.d14/100 ! mm2-> cm2
dsum=0
dsum1=0
if (nsum.gt.1) then
do i=1,nsum
dsum=dsum+(psum(i)*1.d14/100-sum)**2
dsum1=dsum1+(psum1(i)*1.d14/100-sum1)**2
enddo
dsum=sqrt(dsum/nsum/(nsum-1))
dsum1=sqrt(dsum1/nsum/(nsum-1))
endif
#// Report results after cycle 1
#// powder peak parameters (command NFLUX 2)
if (itask1.eq.4.or.itask1.eq.5) then
call GETPEAKPARAM(0,suma,center,fwhm,spread)
ipwd=suma
dipwd=suma*dsum/sum
wpwd=fwhm
spwd=spread
cpwd=center
if (verbose) write(sout,14) suma,suma*dsum/sum,fwhm,spread,
& center
#// beam profile at the sample
else if (itask1.eq.2.or.itask1.eq.3) then
z=0.d0
if (scnt.gt.0) z=r8ln2*sqrt(dei/scnt-(dei0/scnt)**2)
if (verbose) then
write(sout,1) sum,dsum,sum1,dsum1,z,cunit
endif
iinc=sum
diinc=dsum
einc=z
#// intensity at the monitor
else if (itask1.eq.6.or.itask1.eq.7) then
if (verbose) write(sout,13) imonit,sum,dsum
iinc=sum
diinc=dsum
#// beam profile at the detector in double-crystal setup
else if (itask1.eq.11) then
if (verbose) write(sout,13) 99,sum,dsum
iinc=sum
diinc=dsum
i3ax=sum
di3ax=dsum
einc=0
endif
#-------------------------Start cycle 2 ------------------------------
#// 2nd cycle for secondary spectrometer: TAS resolution function
if(itask.eq.8.or.itask.eq.9) then
call SPEC_INI(0,itask)
call MAXV_UPD(0)
ncnt=0
nevent=0
ntot=0.d0
nout=100*nev1
do while ((ncnt.lt.nev1).and.(ntot.lt.nout))
700 if (nevent.eq.0) then
scnt=0.
vol1=1
do i=1,rndlist.dim
vol1=vol1*abs(rndlist.limits(i))
end do
call SPEC_INI(1,itask)
dcnt=0.
nsum=0
isum=0
do i=1,5
psum(i)=0
enddo
nevent1=0
smon1=0
smon=0
endif
ntot=ntot+1.d+0
i=ncnt
j=nevent
z=scnt
z1=smon
call NESS_RUN(ncnt,nevent,itask)
if (nevent.eq.0) goto 700
isum=isum+ncnt-i
nevent1=nevent1+nevent-j
smon1=smon1+smon-z1
dcnt=dcnt+scnt-z
#// Get parcial sums for error estimation
if (isum.eq.nint(nev1/5.d0).and.smon1.gt.0.d0) then
nsum=nsum+1
psum(nsum)=(vol1*dcnt)/smon1
isum=0
dcnt=0.
nevent1=0
smon1=0
endif
enddo
t3=secnds(t1)
# SUM=(VOL1*SCNT)/NEVENT
dsum=0.d0
sum=0.d0
if (smon.gt.0) sum=(vol1*scnt)/smon
if (nsum.gt.1.and.sum.gt.0.d0) then
do i=1,nsum
dsum=dsum+(psum(i)-sum)**2
# write(*,*) PSUM(I),SUM,(PSUM(I)-SUM)/SUM
enddo
dsum=sqrt(dsum/nsum/(nsum-1))
endif
if (verbose) then
write(sout,8) t3
write(sout,9) ncnt,nevent,ntot
write(*,12) (hit(i),i=1,rndlist.dim)
write(*,*)
write(*,*) 'cycle 2 - final volume: ',vol1
call GETSTATE(itask,nevent)
write(*,*)
endif
endif
#// 2nd cycle for PWDS: multidetector scan
if(itask.eq.10) then
spcn=spcm
call SCANMDET(spcx,spcy,spcd,spcm)
call GETPEAKPARAM(3,ipwd,cpwd,fwhm,spwd)
dipwd=diinc/iinc*ipwd
par(1)=ipwd
par(2)=cpwd
par(3)=spwd
do i=1,3
dpar(i)=0.
enddo
tol=0.01
call LMOPT(CHI2SPC,par,3,tol,dpar,2)
t3=secnds(t1)
if (verbose) then
write(sout,8) t3
write(sout,17) spwd,fwhm,cpwd,ipwd
write(sout,18) par(3),par(2),par(1)
else
write(sout,16) par(3),par(1)
endif
endif
#------ End of the 2nd cycle ---------------------------------------
#// ITASK=1,8 do the same thing in a different way:
#// ITASK=1 .. simulate TAS in one cycle
#// ITASK=8 .. split into 2 cycles: sample->source and sample->detector
if (itask1.eq.1.or.itask.eq.8) then
e3ax=0.d0
if (stas.gt.0) e3ax=2.35482d0*sqrt(abs(dei/stas-(dei0/stas)**2))
# I3AX=(VOL1*SCNT)/NEVENT*HSQOV2M*2*STP.KF/10.
if(itask1.eq.1) then
i3ax=sum*hsqov2m*2*stp.kf/10.
di3ax=dsum*hsqov2m*2*stp.kf/10.
endif
if(itask.eq.8) then
i3ax=iinc*sum*hsqov2m*2*stp.kf/10.
di3ax=sqrt((iinc*dsum)**2+(diinc*sum)**2)*hsqov2m*2*stp.kf/10.
endif
if (stas.gt.0) de3ax=dei0/stas
call GETPEAKPARAM(1,suma,center,fwhm,spread)
call RESMAT
write(sout,19) i3ax,fwhm,spread,center
if (verbose) then
if(itask.eq.8) then
write(sout,20) iinc,i3ax,
* volresmc,
* (bragcn(i),i=1,4),(vancn(i),i=1,4),
* (bragrl(i),i=1,4),(vanrl(i),i=1,4)
else
write(sout,21) i3ax,
* volresmc,
* (bragcn(i),i=1,4),(vancn(i),i=1,4),
* (bragrl(i),i=1,4),(vanrl(i),i=1,4)
endif
else
write(sout,16) e3ax,i3ax
endif
endif
#// ITASK=9,10 do similar things:
#// ITASK=9 .. simulate powder curve using ki,kf events
#// ITASK=10 .. simulate powder curve by step-by-step scanning
if (itask.eq.9) then
z=0.d0
ipwd=iinc*(vol1*scnt)/nevent/10. ! mm->cm
dipwd=ipwd*diinc/iinc
call GETPEAKPARAM(2,suma,cpwd,fwhm,spwd)
call THETA_SCAN(spcx,spcy,spcd,65)
par(1)=spcy(32)
par(2)=cpwd
par(3)=spwd
do i=1,3
dpar(i)=0.
enddo
tol=0.01
if (verbose) then
call LMOPT(CHI2SPC,par,3,tol,dpar,2)
write(sout,17) spwd,fwhm,cpwd,ipwd
write(sout,18) par(3),par(2),par(1)
else
write(sout,16) spwd,ipwd
endif
endif
if (verbose) write(sout,*)
#100 FORMAT(' NESS: ',3(2X,E12.5))
return
end
#-------------------------------------------------------------------
logical*4 FUNCTION SAFETY_POOL(lasthit)
# Checks, if the value of any random variable X(I) is found
# in the safety pool. If yes, corresponding limits are relaxed.
# assignment of XRND to physical variables:
# 1-3 .. K(i) vector
# 4-5 .. R(1), R(2)
# 6 .. not used
# 7 .. vertical scattering angle
# 8 .. horizontal scattering angle (for Vanad)
#-------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'randvars.inc'
real*8 z
logical*4 log1
integer*4 i,lasthit
5 format(a10,1x,i2,4(2x,f12.5))
log1=.false.
lasthit=0
do 10 i=1,rndlist.dim
if (rndlist.active(i).gt.0) then
z=abs(2*rndlist.pool(i)*xnorm(i))-rndlist.limits(i)
if (z.gt.0) then
lasthit=i
log1=.true.
# write(*,5) 'HIT: ',I,XNORM(I),XRND(I),RNDLIST.LIMITS(I)
rndlist.limits(i)=rndlist.limits(i)*rndlist.pool(i)
hit(i)=hit(i)+1
#
# write(*,5) 'HIT3: ',i,X(i),XRND(i),RNDLIST.LIMITS(i),TMEAN(i)
#
if(rndlist.limits(7).gt.pi*2) then
rndlist.limits(7)=pi*2
rndlist.active(7)=0
#
endif
endif
endif
10 continue
SAFETY_POOL=log1
return
end
# --------------------------------------------------
SUBROUTINE SWPOOL
# switch safety pool off/on
# --------------------------------------------------
implicit none
INCLUDE 'randvars.inc'
integer*4 mypool(crnd)
integer ipool_off,i
real*8 mylim(crnd)
logical*4 verbose
integer*4 nev
common /mcsetting/ verbose,nev
save ipool_off,mylim,mypool
data ipool_off/0/
# write(*,*) 'ipool=',ipool_off
if(ipool_off.eq.0) then
if (verbose) write(*,*) 'Safety pool OFF'
do 10 i=1,rndlist.dim
mylim(i)=rndlist.limits(i)
mypool(i)=rndlist.active(i)
if (mypool(i).ne.0) then
rndlist.limits(i)=rndlist.limits(i)/rndlist.pool(i)
rndlist.active(i)=0
endif
10 continue
#101 FORMAT('TLIM: ',16(1X,G10.4))
# write(*,101) (RNDLIST.LIMITS(J),J=1,RNDLIST.DIM)
ipool_off=1
else
if (verbose) write(*,*) 'Safety pool ON'
do 20 i=1,rndlist.dim
rndlist.limits(i)=mylim(i)
rndlist.active(i)=mypool(i)
20 continue
ipool_off=0
endif
return
end
#-----------------------------------------------------------------------
SUBROUTINE MAXV_UPD(itask)
# ITASK=0 ... clears MAXV(I) array
# ITASK=1 ... MAXV(I) is compared with X(I) and changed if necessary
# ITASK=2 ... limits are changed according to MAXV(I)
#-----------------------------------------------------------------------
implicit none
INCLUDE 'randvars.inc'
integer*4 i,itask
if(itask.eq.0) then
do i=1,rndlist.dim
maxv(i)=0.d0
enddo
else if(itask.eq.1) then
do i=1,rndlist.dim
if (rndlist.active(i).gt.0) then
if(abs(xnorm(i)).gt.maxv(i)) maxv(i)=abs(xnorm(i))
endif
enddo
else if(itask.eq.2) then
do i=1,rndlist.dim
if (rndlist.active(i).gt.0) then
rndlist.limits(i) = 2.*maxv(i)*rndlist.pool(i)*1.001
endif
enddo
endif
return
end
# -------------------
SUBROUTINE COV_CLR
# -------------------
implicit none
INCLUDE 'randvars.inc'
integer*4 i,j
do 10 i=1,crnd
mval(i)=0
do 20 j=1,i
cov(j,i)=0.d0
cov(i,j)=0.d0
20 continue
10 continue
ncov=0
scov=0.d0
return
end
# -------------------------
SUBROUTINE COV_UPD(x1,p)
# -------------------------
implicit none
INCLUDE 'randvars.inc'
integer*4 i,j
real*8 x1(crnd),p
do 10 i=1,rndlist.dim
if (rndlist.active(i).gt.0) then
mval(i)=mval(i)+x1(i)*1
do 20 j=1,i
if (rndlist.active(j).gt.0) then
cov(i,j)=cov(i,j)+x1(i)*x1(j)*1
cov(j,i)=cov(i,j)
end if
20 continue
endif
10 continue
ncov=ncov+1
scov=scov+1 ! not weighted in this version
return
end
# -------------------------
SUBROUTINE COV_NEW
# -------------------------
implicit none
INCLUDE 'randvars.inc'
integer*4 i,j,nrot
real*8 aux(crnd,crnd),aux1(crnd,crnd),vlim(crnd),z
logical*4 verbose
integer*4 nev
common /mcsetting/ verbose,nev
5 format(a10,5(2x,e11.5))
104 format( ' MAT: ',10(1x,g10.4))
100 format( ' COV: ',10(1x,g10.4))
101 format( ' LIM: ',10(1x,g10.4))
102 format( ' COV - Old volume: ',g12.6)
103 format( ' COV - New volume: ',g12.6)
if (ncov.gt.30) then
do 10 i=1,rndlist.dim
10 if (rndlist.active(i).lt.1) cmat(i,i)=1.d0
do 20 i=1,rndlist.dim
if (rndlist.active(i).gt.0) then
tmean(i)=mval(i)/scov
do 21 j=1,rndlist.dim
if(rndlist.active(j).gt.0) cmat(i,j)=cov(i,j)/scov
21 continue
endif
20 continue
do i=1,rndlist.dim
do j=1,rndlist.dim
if ((rndlist.active(i).gt.0).and.(rndlist.active(j).gt.0)) then
cmat(i,j)=cmat(i,j)-tmean(i)*tmean(j)
endif
enddo
enddo
# filter covariances
do i=1,rndlist.dim
do j=1,rndlist.dim
if (abs(cmat(i,j))/sqrt(cmat(i,i)*cmat(j,j)).lt.0.10)
* cmat(i,j)=0.d0
enddo
enddo
# do i=1,RNDLIST.DIM
# write(*,104) (CMAT(i,j)/SQRT(CMAT(i,i)*CMAT(j,j)),
# * j=1,RNDLIST.DIM)
# enddo
# write(*,*)
# do i=1,5
# write(*,100) (aux(i,j),j=1,5)
# enddo
# write(*,*)
call JACOBI(cmat,aux,rndlist.dim,crnd,vlim,aux1,nrot)
# CALL JACOBI(CMAT,AUX,RNDLIST.DIM,CRND,VLIM,aux1,NROT)
# CALL MXM(1,RNDLIST.DIM,CRND,CMAT,AUX1,AUX)
# CALL MXM(-1,RNDLIST.DIM,CRND,AUX1,AUX,CMAT)
# write(*,*) 'Check JACOBI'
# do i=1,RNDLIST.DIM
# write(*,100) (CMAT(i,j),j=1,RNDLIST.DIM)
# enddo
# write(*,*) 'End'
do 80 i=1,crnd
do 80 j=1,crnd
if (i.ne.j) then
tmat(i,j)=0.d0
else
tmat(i,j)=1.d0
endif
80 continue
z=1
do i=1,rndlist.dim
z=z*rndlist.limits(i)
enddo
if (verbose) write(*,102) z
do i=1,rndlist.dim
if (rndlist.active(i).gt.0) then
rndlist.limits(i)=2*sqrt(6*abs(vlim(i)))
do j=1,rndlist.dim
if (rndlist.active(j).gt.0) then
tmat(i,j)=aux1(j,i)
endif
enddo
endif
enddo
# do i=1,RNDLIST.DIM
# write(*,100) (tmat(i,j),j=1,RNDLIST.DIM)
# enddo
# write(*,101) (RNDLIST.LIMITS(I),I=1,RNDLIST.DIM)
# write(*,*)
z=1
do i=1,rndlist.dim
z=z*rndlist.limits(i)
enddo
if (verbose) write(*,103) z
endif
return
end
#----------------------------------------------------------
SUBROUTINE RANDFILL
#
#----------------------------------------------------------
implicit none
INCLUDE 'randvars.inc'
integer*4 i
real*4 RAN1
do i=1,rndlist.dim
xnorm(i)=rndlist.limits(i)*(RAN1()-0.5)
end do
return
end
#-------------------------------------------------------------------
SUBROUTINE NESS_RUN(ncnt,nevent,itask)
# trace one event, starting at the sample
# ITASK=1 ... inelastic scattering, TAS resolution
# ITASK=2 ... sample -> source
# ITASK=4 ... sample -> source + sample(powder) -> detector (no analyser)
# ITASK=6 ... sample -> source + sample (Vanad) -> monitor(IMONIT)
#-------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'randvars.inc'
INCLUDE 'sim_grf.inc'
record /NEUTRON/ neui,neuf,neui1,neuf1
integer*4 nevent,ncnt,itask,i,ix,iy,i1,i2,ns
real*8 c1,c2,ki
common /neuif/ neui,neuf,neui1,neuf1
logical*4 SPEC_GO,SAFETY_POOL
real*8 V3xV3
real*4 RAN1,dum
# DBGREF=(NCNT.EQ.0)
if(nevent.eq.0) then
#1 format(a,1x,8(G12.6,2x))
# write(*,1) 'NESS_RUN: ',(XRND(i),i=1,RNDLIST.DIM)
scnt=0.d0
scntf=0.d0
dei=0.d0
dei0=0.d0
sum=0.d0
stas=0.d0
smon=0.d0
call COV_CLR
if (itask.eq.4.or.itask.eq.5) then
do ix=1,mimax
do iy=1,mimax
svol(ix,iy)=0.e0
enddo
enddo
endif
endif
nevent=nevent+1
call RANDFILL
call M16XV16(-1,rndlist.dim,tmat,xnorm,xrnd)
do i=1,rndlist.dim
xrnd(i)=xrnd(i)+tmean(i)
enddo
if (itask.eq.8.or.itask.eq.9.or.itask.eq.10) then ! get NEUI from events simulated before
call NSTORE_N(i1,i2,ns)
nseed=nseed+1
if (nseed.gt.i1.or.nseed.le.0) nseed=1
call NSTORE_READ1(nseed,neui)
ki=sqrt(V3XV3(neui.k,neui.k)) ! |k| for incident neutron
if (normmon.eq.1) then
smon=smon+neui.p/ki ! * monitor efficiency
else
smon=smon+neui.p
endif
else
neui.r(1)=xrnd(4)
neui.r(2)=xrnd(5)
neui.r(3)=0.
# DO I=1,3
# NEUI.K(I)=XRND(I)
# END DO
if(abs(xrnd(1)).ge.1.d0.or.abs(xrnd(2)).ge.1.d0.or.
* xrnd(3).le.1.d-2) then
# write(*,*) 'NESS_RUN: ',(XRND(i),i=1,RNDLIST.DIM)
return
endif
c1=sqrt(abs(1-xrnd(1)**2))
c2=sqrt(abs(1-xrnd(2)**2))
neui.k(1)=xrnd(3)*xrnd(1)*c2
neui.k(2)=-xrnd(3)*xrnd(2)
neui.k(3)=xrnd(3)*c1*c2
neui.p=1.d0/c1/c2
neui.t=0
dum=RAN1()
neui.s=2*nint(dum)-1
endif
# if (mod(nevent,1000).EQ.0) then
# CALL GETSTATE(ITASK,NEVENT)
# endif
if(SPEC_GO(itask)) then
call MAXV_UPD(1)
if (SAFETY_POOL(i)) then
call SPEC_INI(1,itask)
nevent=0
ncnt=0
return
endif
ncnt=ncnt+1
call VALID_EVENT(itask,ncnt)
endif
return
end
#------------------------------------------------------------------------
SUBROUTINE VALID_EVENT(itask,ncnt)
# Store event information
# ITASK=1 ... inelastic scattering, TAS resolution
# ITASK=2 ... sample -> source
# ITASK=3 ... source -> sample
# ITASK=4 ... sample -> source + sample(powder) -> detector (no analyser and col4)
# ITASK=5 ... source -> sample + sample(powder) -> detector (no analyser and col4)
# ITASK=6 ... source -> monitor(IMONIT) (using Vanad sample))
# ITASK=7 ...
# ITASK=8 ... sample -> detector, TAS (diffusee inelastic)
# ITASK=9 ... sample -> detector, PD (without anal. and col4), diffuse elastic
# ITASK=10 ... sample -> detector, PD (without anal. and col4), powder smaple
# ITASK=11 ... source -> detector, TAS, bragg scattering
#------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'rescal.inc'
INCLUDE 'randvars.inc'
INCLUDE 'sim_grf.inc'
record /NEUTRON/ neui,neuf,neui1,neuf1,neu
integer*4 i,itask,ix,iy,ncnt
real*8 e(4),kki,z,pp
real*4 V3xV3
common /neuif/ neui,neuf,neui1,neuf1
if (itask.gt.1.and.itask.le.3) neu=neui
if ((itask.ge.4.and.itask.le.11).or.itask.eq.1) neu=neuf1
if (itask.eq.1) then
pp=neui.p*neuf.p
else
pp=neu.p
endif
call COV_UPD(xrnd,pp)
#// intensity and flux at the sample are calculated using SCNT sums
if (itask.eq.8) pp=neui.p*neuf.p
scnt=scnt+pp
if (abs(neu.r(1)).le.5.and.abs(neu.r(2)).le.5) then
scntf=scntf+pp
endif
if ((itask.ge.2.and.itask.le.7).or.(itask.eq.11)) then
# beam profile is accumulated
kki=V3XV3(neu.k,neu.k)
e(1)=neu.r(1)
e(2)=neu.r(2)
e(3)=hsqov2m*(kki-stp.ki**2)
e(4)=neu.t/1000 ! in [ms]
dei=dei+neu.p*e(3)**2
dei0=dei0+e(3)*neu.p
call EVARRAY(1,1,ncnt,e,neu.p)
do i=1,3
e(i)=neu.k(i)
enddo
#
e(4)=neu.s
call EVARRAY(1,0,ncnt,e,neu.p)
# Accumulate gauge volume
if(itask.eq.4.or.itask.eq.5) then
ix=(neuf.r(1)/sam.size(1)+0.5)*mimax+1
iy=(neuf.r(3)/sam.size(3)+0.5)*mimax+1
if(ix.gt.0.and.ix.le.mimax.and.iy.gt.0.and.iy.le.mimax) then
svol(ix,iy)=svol(ix,iy)+neuf.p
endif
endif
endif
if (itask.eq.2) then
call NSTORE_WRITE1(ncnt,neui)
else if (itask.eq.8) then ! get width in E
call NSTORE_WRITE2(ncnt,nseed,neuf)
call NSTORE_GETQE(ncnt,e,pp,0.)
z=e(4)
stas=stas+pp
dei=dei+pp*z**2
dei0=dei0+pp*z
else if (itask.eq.9) then ! get width in |Q|
call NSTORE_WRITE2(ncnt,nseed,neuf)
call NSTORE_GETQE(ncnt,e,pp,0.)
z=sqrt((stp.q+e(1))**2+e(2)**2+e(3)**2)-stp.q
stas=stas+pp
dei=dei+pp*z**2
dei0=dei0+pp*z
else if (itask.eq.1) then ! get width in E
call NSTORE_WRITE1(ncnt,neui)
call NSTORE_WRITE2(ncnt,ncnt,neuf)
call NSTORE_GETQE(ncnt,e,pp,0.)
z=e(4)
stas=stas+pp
dei=dei+pp*z**2
dei0=dei0+pp*z
endif
end
#------------------------------------------------------------------
SUBROUTINE SCANMDET(fx,fy,dfy,n)
#// scan by multidetector (powder)
#// uses incident neutrons previously stored in KSTACK
#------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'rescal.inc'
INCLUDE 'randvars.inc'
integer*4 i,j,k,n,i2,nev1,ialloc,ncnt,nevent,nevent0,j0
real*4 fx(n),fy(n),dfy(n)
real*8 dx,xmax,xmin,omega0,vol1,sump,dth
real*8 ntot,nout
record /NEUTRON/ neui,neuf,neui1,neuf1
common /neuif/ neui,neuf,neui1,neuf1
dth=res_dat(i_da4)
xmax=dth*n/2.
xmin=-xmax
dx=dth
do i=1,n
fx(i)=xmin+(i-0.5d0)*dx
enddo
do i=1,n
fy(i)=0.
dfy(i)=0.
enddo
call NSTORE_N(nev1,i2,ialloc)
if (nev1.gt.0) then
nev1=nev1/5
#// Scan center
ntot=0.d0
nout=1000.*nev1
call MAXV_UPD(0)
call SPEC_INI(0,4)
omega0=sol3.frame.axi
do i=1,6
rndlist.active(i)=0
rndlist.limits(i)=1
enddo
tmat(2,7)=0.
j0=(n+1)/2
sol3.frame.axi=omega0+fx(j0)*minute
call SLIT_INIT(sol3.frame)
ncnt=0
nevent=0
do while ((ncnt.lt.nev1).and.(ntot.lt.nout))
ntot=ntot+1.d0
if (nevent.eq.0) then
vol1=1
do i=1,rndlist.dim
vol1=vol1*rndlist.limits(i)
end do
call SPEC_INI(1,4)
sump=0
endif
call NESS_RUN(ncnt,nevent,10)
sump=sump+neui.p
enddo
if(ntot.ge.nout) then
write(*,*) 'TIMEOUT',nevent,ncnt
return
endif
nevent0=nevent
fy(j0)=iinc/sump*(vol1*scnt)/10. ! mm->cm
call MAXV_UPD(2)
10 format(i3,2(2x,g12.6))
write(*,10) j0,fx(j0),fy(j0)
scnt=1.
do j=j0+1,n
if (scnt.gt.0) then
sol3.frame.axi=omega0+fx(j)*minute
call SLIT_INIT(sol3.frame)
ncnt=0
nevent=0
do while (nevent.lt.nevent0)
if (nevent.eq.0) then
vol1=1
do i=1,rndlist.dim
vol1=vol1*rndlist.limits(i)
end do
call SPEC_INI(1,4)
sump=0
endif
call NESS_RUN(ncnt,nevent,10)
sump=sump+neui.p
enddo
fy(j)=iinc/sump*(vol1*scnt)/10.
dfy(j)=fy(j)/sqrt(1.+ncnt)
write(*,10)j,fx(j),fy(j)
else ! skip rest if last step was FY=0
fy(j)=0.
dfy(j)=0.
endif
enddo
scnt=1.
do k=1,j0-1
j=j0-k
if (scnt.gt.0) then
sol3.frame.axi=omega0+fx(j)*minute
call SLIT_INIT(sol3.frame)
ncnt=0
nevent=0
do while (nevent.lt.nevent0)
if (nevent.eq.0) then
vol1=1
do i=1,rndlist.dim
vol1=vol1*rndlist.limits(i)
end do
call SPEC_INI(1,5)
sump=0
endif
call NESS_RUN(ncnt,nevent,10)
sump=sump+neui.p
enddo
fy(j)=iinc/sump*(vol1*scnt)/10.
dfy(j)=fy(j)/sqrt(1.+ncnt)
write(*,10)j,fx(j),fy(j)
else
fy(j)=0.
dfy(j)=0.
endif
enddo
sol3.frame.axi=omega0
call SLIT_INIT(sol3.frame)
endif
end
#-------------------------------------------------------------------
SUBROUTINE GETPEAKPARAM(ivar,suma,center,fwhm,wspread)
# Get peak parameters
# IVAR=3 .. from anything accumulated in SPCX,Y arrays
# IVAR=2 .. from |Q| distribution (powder curve)
# IVAR=1 .. from energy distribution (Vanad scan)
# otherwise .. spatial profile (along X) at the monitor
#-------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
real*8 e(4),p,s0,s1,s2,suma,center,fwhm,wspread,z,conv
integer*4 i,ncnt,ivar,i1,ialloc
suma=0
center=0
wspread=0
fwhm=0
if (ivar.eq.3.and.spcn.gt.0) then
suma=0.
s0=0
s1=0
s2=0
do i=1,spcn
z=spcx(i)
s0=s0+spcy(i)
s1=s1+spcy(i)*z
s2=s2+spcy(i)*z**2
if(spcy(i).gt.suma) suma=spcy(i)
enddo
if(s0.gt.0) then
center=s1/s0
wspread=r8ln2*sqrt(abs(s2/s0-center**2))
call GETFWHM(spcx,spcy,spcn,fwhm)
endif
else if (ivar.eq.2) then
conv=180*60/pi/sqrt(abs(stp.ki**2-(stp.q/2.)**2))
call NSTORE_N(i1,ncnt,ialloc) ! get number of events NCNT
if (ncnt.gt.0) then
s0=0
s1=0
s2=0
do i=1,ncnt
call NSTORE_GETQE(i,e,p,0.)
z=sqrt((stp.q+e(1))**2+e(2)**2+e(3)**2)-stp.q
z=z*conv
s0=s0+p
s1=s1+p*z
s2=s2+p*z**2
enddo
if(s0.gt.0) then
spcn=spcm
suma=ipwd
center=s1/s0
wspread=r8ln2*sqrt(abs(s2/s0-center**2))
call THETA_SCAN(spcx,spcy,spcd,spcm)
call GETFWHM(spcx,spcy,spcm,fwhm)
endif
endif
else if (ivar.eq.1) then
call NSTORE_N(i1,ncnt,ialloc) ! get number of events NCNT
if (ncnt.gt.0) then
s0=0
s1=0
s2=0
do i=1,ncnt
call NSTORE_GETQE(i,e,p,0.)
s0=s0+p
s1=s1+p*e(4)
s2=s2+p*e(4)**2
enddo
if(s0.gt.0) then
spcn=spcm
suma=i3ax*iinc
center=s1/s0
wspread=r8ln2*sqrt(abs(s2/s0-center**2))
call E_SCAN(spcx,spcy,spcd,spcm)
call GETFWHM(spcx,spcy,spcm,fwhm)
endif
endif
else
call EVARRAY(3,1,ncnt,e,p) ! get number of events NCNT
if (ncnt.gt.0) then
s0=0
s1=0
s2=0
do i=1,ncnt
call EVARRAY(2,1,i,e,p)
s0=s0+p
s1=s1+p*e(1)
s2=s2+p*e(1)**2
enddo
if(s0.gt.0) then
spcn=spcm
suma=sum
center=s1/s0
wspread=r8ln2*sqrt(abs(s2/s0-center**2))
call PSD_ARRAY(spcx,spcy,spcd,spcm)
call GETFWHM(spcx,spcy,spcm,fwhm)
endif
endif
endif
end
#-------------------------------------------------------------------
SUBROUTINE GETFWHM(x,y,n,fwhm)
# Calculate fwhm
#-------------------------------------------------------------------
implicit none
integer*4 n,i,imax,i1,i2
real*4 x(n),y(n),z1,z2
real*8 ymax,x1,x2,fwhm
ymax=y(1)
imax=1
do i=1,n
if(ymax.le.y(i)) then
imax=i
ymax=y(i)
endif
enddo
i1=1
do while ((y(i1).lt.ymax/2.d0).and.(i1.lt.n))
i1=i1+1
enddo
i2=n
do while ((y(i2).lt.ymax/2.d0).and.(i2.gt.0))
i2=i2-1
enddo
if((i1.gt.1).and.(i1.lt.n).and.(i2.gt.1).and.(i2.lt.n).and.
* (i1.le.i2)) then
z1= y(i1)-y(i1-1)
z2=y(i2+1)-y(i2)
if (z1.ne.0.and.z2.ne.0) then
x1=x(i1-1)+(ymax/2.d0-y(i1-1))/(y(i1)-y(i1-1))*(x(i1)-x(i1-1))
x2=x(i2)+(ymax/2.d0-y(i2))/(y(i2+1)-y(i2))*(x(i2+1)-x(i2))
fwhm=x2-x1
else
fwhm=0.
endif
else
fwhm=0
endif
end
#----------------------------------------------------
SUBROUTINE PSD_ARRAY(fx,fy,dfy,n)
# get beam profile integrated along vertical coordinate
#----------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
integer*4 i,j,n,ncnt
real*4 fx(n),fy(n),dfy(n),y0(129)
real*8 xyet(4),p,dx,sum1
real*8 xmax,xmin
call EVARRAY(3,1,ncnt,xyet,p) ! get number of events NCNT
if (ncnt.gt.0) then
xmax=det.frame.size(1)/2.0d0
xmin=-xmax
dx=(xmax-xmin)/n
do i=1,n
fx(i)=xmin+(i-0.5d0)*dx
enddo
do i=1,n
y0(i)=0.
fy(i)=0.
enddo
sum1=0.
do i=1,ncnt
call EVARRAY(2,1,i,xyet,p)
j=int((xyet(1)-xmin)/dx)+1
if (j.ge.1.and.j.le.n) then
fy(j)=fy(j)+p
y0(j)=y0(j)+1.
sum1=sum1+p
endif
enddo
do i=1,n
fy(i)=fy(i)*sum/sum1
if(fy(i).gt.0.) then
dfy(i)=fy(i)/sqrt(y0(i)/2.)
else
dfy(i)=0.
endif
enddo
endif
end
#----------------------------------------------------
SUBROUTINE E_SCAN(fx,fy,dfy,n)
# get Vanad scan profile from Q,E events
#----------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'rescal.inc'
integer*4 i,j,n,ncnt,i1,ialloc
real*4 fx(n),fy(n),dfy(n),y0(129)
real*8 qe(4),p,dx,sum1
real*8 xmax,xmin,de
call NSTORE_N(i1,ncnt,ialloc)
# write(*,*) 'E-SCAN: ',NCNT,STAS,SMON,SCNT
# write(*,*) 'E-SCAN: ',IINC,I3AX,N
if (ncnt.gt.0.and.stas.gt.0.and.scnt.gt.0) then
de=res_dat(i_den)
xmax=de*n/2.
xmin=-xmax
dx=de
do i=1,n
fx(i)=xmin+(i-0.5d0)*dx
enddo
do i=1,n
fy(i)=0.
y0(i)=0.
enddo
sum1=0.
do i=1,ncnt
call NSTORE_GETQE(i,qe,p,0.)
# write(*,*) (QE(J),J=1,4),P
# pause
j=int((qe(4)-xmin)/dx)+1
if (j.ge.1.and.j.le.n) then
fy(j)=fy(j)+p
y0(j)=y0(j)+1.
sum1=sum1+p
endif
enddo
do i=1,n
fy(i)=fy(i)*iinc*i3ax/sum1
if(fy(i).gt.0.) then
dfy(i)=fy(i)/sqrt(y0(i)/2.)
else
dfy(i)=0.
endif
# write(*,*) FX(I),FY(I),DFY(I)
enddo
endif
end
#----------------------------------------------------
SUBROUTINE THETA_SCAN(fx,fy,dfy,n)
# get powder peak profile from Q,E events
#----------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'rescal.inc'
integer*4 i,j,n,ncnt,i1,ialloc
real*4 fx(n),fy(n),dfy(n),y0(129)
real*8 qe(4),p,dx,sum1,conv,qq
real*8 xmax,xmin,dth
call NSTORE_N(i1,ncnt,ialloc)
if (ncnt.gt.0.and.stas.gt.0.and.scnt.gt.0) then
dth=res_dat(i_da4)
xmax=dth*n/2.
xmin=-xmax
dx=dth
do i=1,n
fx(i)=xmin+(i-0.5d0)*dx
enddo
do i=1,n
fy(i)=0.
y0(i)=0.
enddo
sum1=0.
conv=1.d0/sqrt(stp.ki**2-(stp.q/2.)**2)/minute ! conversion factor dQ->2dThetaS
do i=1,ncnt
call NSTORE_GETQE(i,qe,p,0.)
qq=sqrt((qe(1)+stp.q)**2+qe(2)**2+qe(3)**2)-stp.q
qq=qq*conv
# write(*,*) QQ,P
# pause
j=int((qq-xmin)/dx)+1
if (j.ge.1.and.j.le.n) then
fy(j)=fy(j)+p
y0(j)=y0(j)+1.
sum1=sum1+p
endif
enddo
do i=1,n
fy(i)=fy(i)*ipwd/sum1
if(fy(i).gt.0.) then
dfy(i)=fy(i)/sqrt(y0(i)/2.)
else
dfy(i)=0.
endif
enddo
endif
end
#xxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCx
#
# **** Benchmark ****
#
#xxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCxxxxxxxxxCx
#------------------------------------------------------------------------
SUBROUTINE NESS_BENCH(ncnt)
# measures transmission speed for crystal, colimator and powder sample
#------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'randvars.inc'
record /NEUTRON/ neui
integer*4 i,nc,what,ncnt
real*4 secnds
logical*4 BENCH_CR,BENCH_SOL2,BENCH_PWD
real*4 t1(3),t2(3)
character*15 whatch(3)
real*4 RAN1,dum
real*8 c1,c2
data whatch/ 'Crystal: ', 'Colimator 2: ', 'Powder sample: '/
1 format(a15,g12.6, ' events/sec (',g10.4, ' sec)')
do what=1,3
nc=0
call RAN1SEED(100001)
do while (nc.le.ncnt)
call RANDFILL
call M16XV16(-1,rndlist.dim,tmat,xnorm,xrnd)
do i=1,rndlist.dim
xrnd(i)=xrnd(i)+tmean(i)
enddo
neui.r(1)=xrnd(4)
neui.r(2)=xrnd(5)
neui.r(3)=0.
# DO I=1,3
# NEUI.K(I)=XRND(I)
# END DO
if(abs(xrnd(1)).ge.1.d0.or.abs(xrnd(2)).ge.1.d0) then
return
endif
c1=sqrt(abs(1-xrnd(1)**2))
c2=sqrt(abs(1-xrnd(2)**2))
neui.k(1)=xrnd(3)*xrnd(1)*c2
neui.k(2)=-xrnd(3)*xrnd(2)
neui.k(3)=xrnd(3)*c1*c2
neui.p=1.d0/c1/c2
neui.t=0
dum=RAN1()
neui.s=2*nint(dum)-1
if((mon.mag*neui.s).lt.0) neui.p=0
if(what.eq.1) then
if (BENCH_CR(0,neui)) nc=nc+1
else if(what.eq.2) then
if (BENCH_SOL2(0,neui)) nc=nc+1
else if(what.eq.3) then
if (BENCH_PWD(0,neui)) nc=nc+1
endif
end do
t1(what)=secnds(0.0)
nc=0
call RAN1SEED(100001)
do while (nc.le.ncnt)
call RANDFILL
call M16XV16(-1,rndlist.dim,tmat,xnorm,xrnd)
do i=1,rndlist.dim
xrnd(i)=xrnd(i)+tmean(i)
enddo
neui.r(1)=xrnd(4)
neui.r(2)=xrnd(5)
neui.r(3)=0.
# DO I=1,3
# NEUI.K(I)=XRND(I)
# END DO
if(abs(xrnd(1)).ge.1.d0.or.abs(xrnd(2)).ge.1.d0) then
return
endif
c1=sqrt(abs(1-xrnd(1)**2))
c2=sqrt(abs(1-xrnd(2)**2))
neui.k(1)=xrnd(3)*xrnd(1)*c2
neui.k(2)=-xrnd(3)*xrnd(2)
neui.k(3)=xrnd(3)*c1*c2
neui.p=1.d0/c1/c2
neui.t=0
dum=RAN1()
neui.s=2*nint(dum)-1
if((mon.mag*neui.s).lt.0) neui.p=0
if(what.eq.1) then
if (BENCH_CR(1,neui)) nc=nc+1
else if(what.eq.2) then
if (BENCH_SOL2(1,neui)) nc=nc+1
else if(what.eq.3) then
if (BENCH_PWD(1,neui)) nc=nc+1
endif
end do
t2(what)=secnds(0.0)
write(sout,1) whatch(what),
* ncnt/(t2(what)-t1(what)),t2(what)-t1(what)
enddo
return
end
#------------------------------------------------------------------------
SUBROUTINE NESS_ROCK(icr,nc,nth,dth,rth,divh,divv)
# measures rocking curve of the monochromator
#------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
record /NEUTRON/ neui,neu
integer*4 icr,i,i1,i2,imax,nth,nc,istep,iev,nc0
real*8 z1,z2,fwhm,center,sprd,rmax,theta,theta0,sum1,sum2,dth
real*8 rth(nth),gon1,divh,divv
logical*4 CRYST_GO2
real*8 alfh,alfv,sump,pp,kk !,DUM
real*4 RAN1,ksi1,ksi2
record /CRYSTAL/ cr
1 format( 'fwhm: ',e10.3, ' [min] , center: ',
1 e10.3, ' [min] R0: ',e10.3)
2 format( 'fwhm: ',e10.3, ' [steps], center: ',
1 e10.3, ' [steps]')
if (icr.eq.2) then
kk=stp.kf
cr=ana
else
kk=stp.ki
cr=mon
endif
nc0=nc
if(cr.hmos.lt.sec.and.cr.nb.eq.1.and.cr.nh.eq.1) nc0=1
theta0=cr.thb
gon1=cr.frame.gon(1)
do istep=1,nth
rth(istep)=0
theta=(-(nth-1)/2.+istep-1)*dth
cr.frame.gon(1)=gon1+theta
call CRYST_INIT2(cr)
sump=0.
# DUM=0.
do iev=1,nc0
neui.r(1)=2.*cr.dh
neui.r(2)=2.*cr.dv
neui.r(3)=2.*cr.db+cr.frame.dist
pp=1.
alfh=0.
alfv=0.
if (divh.gt.0) then
ksi1=RAN1()
alfh=(2.0*ksi1-1.0)*divh
pp=pp*(1.-abs(alfh/divh))
endif
if (divv.gt.0) then
ksi2=RAN1()
alfv=(2.0*ksi2-1.0)*divv
pp=pp*(1.-abs(alfv/divv))
endif
neui.k(1)=kk*alfh
neui.k(2)=kk*alfv
neui.k(3)=kk
z1=sqrt(neui.k(1)**2+neui.k(2)**2+neui.k(3)**2)
do i=1,3
neui.k(i)=neui.k(i)*kk/z1
enddo
neui.p=pp
neui.t=0
# DUM=DUM+ALFH*180*60/PI
#3 format(6(1x,G11.4))
# write(*,3) ALFH*180*60/PI,ALFV*180*60/PI,PP,DUM/IEV,KSI1,KSI2
# pause
if (CRYST_GO2(cr,neui,neu)) then
rth(istep)=rth(istep)+neu.p
endif
sump=sump+neui.p
end do
rth(istep)=rth(istep)/sump
end do
rmax=0
sum1=0
sum2=0
do i=1,nth
if (rth(i).gt.rmax) then
rmax=rth(i)
imax=i
endif
sum1=sum1+rth(i)
sum2=sum2+rth(i)**2
enddo
center=sum1/nth
sprd=2.35*sqrt(sum2/nth-center**2)
i1=2
i2=nth-1
do i=1,nth
if(i1.eq.2.and.rth(i).gt.rmax/2.) i1=i-1
if(i2.eq.nth-1.and.i.gt.imax.and.rth(i).lt.rmax/2.) i2=i-1
enddo
z1=i1+(rmax/2.-rth(i1))/(rth(i1+1)-rth(i1))
z2=i2+(rmax/2.-rth(i2))/(rth(i2+1)-rth(i2))
fwhm=z2-z1
write(sout,1) fwhm*dth*180*60/pi,center*dth*180*60/pi,rmax
write(sout,2) fwhm,center
cr.thb=theta0
cr.frame.gon(1)=gon1
call CRYST_INIT2(cr)
spcn=nth
if (spcn.gt.65) spcn=65
do i=1,spcn
spcx(i)=(-(nth-1)/2.+i-1)*dth*180*60/pi
spcy(i)=rth(i)
spcd(i)=rth(i)/sqrt(1.d0*(nc0+1))
enddo
return
end
#------------------------------------------------------------------------
SUBROUTINE TEST_SYMMETRY(nc,nth,dth,divh,divv)
# measures rocking curve of the monochromator
#------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
record /NEUTRON/ neui,neu,neu1
integer*4 i,nth,nc,istep,iev,nc0
real*8 z1,theta0,dth
real*8 gon1,divh,divv,dist0
logical*4 CRYST_GO2
real*8 alfh,alfv,sump,pp,tin,tout !,DUM
real*4 RAN1,ksi1,ksi2
nc0=10
if(mon.hmos.lt.sec.and.mon.nb.eq.1.and.mon.nh.eq.1) nc0=1
theta0=mon.thb
dist0=mon.frame.dist
mon.frame.dist=0.
gon1=mon.frame.gon(1)
do istep=1,nth
sump=0.
# DUM=0.
do iev=1,nc0
neui.r(1)=2.*mon.dh
neui.r(2)=2.*mon.dv
neui.r(3)=2.*mon.db
pp=1.
alfh=0.
alfv=0.
if (divh.gt.0) then
ksi1=RAN1()
alfh=(2.0*ksi1-1.0)*divh
pp=pp*(1.-abs(alfh/divh))
endif
if (divv.gt.0) then
ksi2=RAN1()
alfv=(2.0*ksi2-1.0)*divv
pp=pp*(1.-abs(alfv/divv))
endif
neui.k(1)=stp.ki*alfh
neui.k(2)=stp.ki*alfv
neui.k(3)=stp.ki
z1=sqrt(neui.k(1)**2+neui.k(2)**2+neui.k(3)**2)
do i=1,3
neui.k(i)=neui.k(i)*stp.ki/z1
enddo
neui.p=pp
neui.t=0
mon.frame.gon(1)=gon1
call CRYST_INIT2(mon)
if (CRYST_GO2(mon,neui,neu)) then
neu1.p=neui.p
call SLIT_PRE1(mon.frame,neui.r,neui.k,neu1.r,neu1.k)
call CR_BORDER(mon,neu1.r,neu1.k,tin,tout)
do i=1,3
neu1.r(i)=neu1.r(i)+tin*neu1.k(i)
end do
call wrtneu(neu1)
neu1.p=neu.p
call SLIT_PRE1(mon.frame,neu.r,neu.k,neu1.r,neu1.k)
call wrtneu(neu1)
do i=1,3
neu.k(i)=-neu.k(i)
enddo
neu.p=1.
# MON.FRAME.GON(1)=-gon1-2.*MON.CHI
# CALL CRYST_INIT2(MON)
neu1.p=neu.p
call SLIT_PRE1(mon.frame,neu.r,neu.k,neu1.r,neu1.k)
call wrtneu(neu1)
do i=1,10
if (CRYST_GO2(mon,neu,neui)) then
neu1.p=neui.p
call SLIT_PRE1(mon.frame,neui.r,neui.k,neu1.r,neu1.k)
call wrtneu(neu1)
else
write(*,*) 'lost back'
endif
enddo
pause
else
write(*,*) 'lost forth'
endif
end do
end do
mon.thb=theta0
mon.frame.gon(1)=gon1
mon.frame.dist=dist0
call CRYST_INIT2(mon)
end
real*4 FUNCTION OPTMC(par)
#***********************************************************************
# returns value to be minimized when optimizing bending radii by M.C.
# J.S. 5/7/2001
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
INCLUDE 'trax.inc'
INCLUDE 'rescal.inc'
real*8 eps
parameter(eps=1.d-10)
real*4 b(4),par(1),z
integer*4 optpar,optmerit
real*8 optev
common /mcoptim/ optpar,optmerit,optev
integer*4 i
do i=1,4
b(i)=res_dat(i_romh+i-1)
enddo
do i=1,4
if (optpar.eq.i) res_dat(i_romh+i-1)=par(1)
enddo
# CALL BEFORE ! should call this in some cases, not for ROMH..ROAV
call NESS_CONV(0)
#1 format(a,2(2x,I3),2(2x,G12.6),$)
# write(*,1) 'OPTMC: ',OPTPAR,OPTMERIT,OPTEV,PAR(1)
z=0.
if (optmerit.eq.1) then
call NESS(2,optev)
if (iinc.ge.eps) z=1.d8/iinc
else if (optmerit.eq.2) then
call NESS(2,optev)
if (iinc.ge.eps) z=1.d8*einc/iinc
else if (optmerit.eq.3) then
call NESS(2,optev)
if (iinc.ge.eps) z=1.d8*einc**2/iinc
else if (optmerit.eq.4) then
call NESS(9,optev)
if (ipwd.ge.eps) z=1.d8*spwd**2/ipwd
else if (optmerit.eq.5) then
call NESS(4,optev)
if (ipwd.ge.eps) z=1.d8*spwd**2/ipwd
else if (optmerit.eq.6) then
call NESS(8,optev)
if (iinc*i3ax.ge.eps) z=1.d8*e3ax**2/(iinc*i3ax)
endif
OPTMC=z
do i=1,4
res_dat(i_romh+i-1)= b(i)
enddo
# write(*,*) Z
end
#***********************************************************************
real*4 FUNCTION CHI2SPC(par)
# returns value to be minimized when fitting gaussian to the data in SPCX,Y,DY array.
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'structures.inc'
INCLUDE 'ness_common.inc'
real*8 z
real*4 par(3),chi
integer*4 i
chi=0
par(3)=abs(par(3))
if (abs(par(3)).lt.1e-10) par(3)=1e-10
do i=1,spcn
if(spcd(i).ne.0) then
z=par(1)*exp(-0.5*(spcx(i)-par(2))**2/(par(3)/r8ln2)**2)
chi=chi+((spcy(i)-z))**2
endif
enddo
CHI2SPC=chi
end