Source module last modified on Mon, 8 May 2006, 23:39;
HTML image of Fortran source automatically generated by
for2html on Mon, 29 May 2006, 15:06.
#//////////////////////////////////////////////////////////////////////
#//// $Id: res_hist.f,v 1.9 2006/05/08 21:39:57 saroun Exp $
#////
#//// R E S T R A X 4.73
#////
#//// Subroutines for handling histograms + encapsulating fitting routine
#////
#//////////////////////////////////////////////////////////////////////
# in filling histograms, WHATHIS flag is set to signal:
# bit 1 .. histogram ready
# bit 2 .. using ray-tracing
# bit 4 .. using EXCI module
#---------------------------------------------------------------------------
SUBROUTINE HISTINIT
#///***** Initialize histograms *****
#/// this is the default partitioning of histogram arrays:
#/// space is made for current data (mf_cur) and all other active datasets
#-----------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'restrax.inc'
integer*4 i,j,nh0
do i=1,mdat
if ((mf_active(i).and.i.le.mf_max).or.(i.eq.mf_cur)) then
nhist(i)=nhist(i-1) + nhi
nh0=(nhi+1)/2
do j=nhist(i-1)+1,nhist(i)
xhist(j)=(j-nh0-nhist(i-1))
rhist(j)=0
ihist(j)=i
end do
shist(i)=0
else
nhist(i)=nhist(i-1)
shist(i)=0
endif
end do
whathis=iand(whathis,254) ! set bit1=0 => no RHIST ready
# write(*,*) 'HISTINIT ',WHATHIS
end
#---------------------------------------------------------------------
real*8 FUNCTION GETSQOM(imin,imax)
# Create the vector of all the (Q,E) events
# Arguments:
# IMIN,IMAX ... min. and max. dataset number
# Return: CHKQOM
#
# arrays created:
# QOM(4,I) ... Q(3),E of the I-th event in r.l.u. coordinates
# PQOM(I) ... weight of the I-th event
# NQOM(item) ... NQOM(item)-NQOM(item-1) = number of events for item-th data set
# IQOM(I) ... index of the data set the I-th event belongs to
# IMPORTANT ! QOM contains full values of Qhkl,E, while the events
# returned by GETQE are only relative to the scan centre, Qhkl(0) E(0)
#---------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'exciimp.inc'
# RECORD /QOMEGA/ rq
real*8 one,eps
parameter (one=1.,eps=1e-20)
integer*4 i,j,k,nph,item,imin,imax
real*8 x(4),xx(4),var(4),aux(4),aux4(4,4),p,chkbase
real*8 amat(4,4),arc(4,4),ada(4,4),b(4,4)
real*4 xr(4,nxr)
common /rnum/ xr
real*4 GASDEV
#1 format('GETSQOM, ',a,6(1x,G12.6))
#------------------------------------------------------------------------
# HNORM is calculated integrated intensity
# Integrated intensity is normalized to 10^6 counts per 1cm^2 of the monitor
# supposing unit monitor efficiency at ki=1A^-1
# and unit value of S(Q,E) = 1 /cm/meV/ster
# RELTR=(2*PI)**3*SQRT(DET(ki,r,kf))/DET(ki)) ~ V*d(Omega)*d(E)
# RELMC=the same, but calculated from M.C. events
# ZNORM = 10^6*HSQOVM*SQRT(24/PI) is caclulated in RT_CONVRT (factor kf/ki included)
# SQRT(24/PI) is a correction for the gaussian approximation of the sample volume.
# SUMAMC is the sum of events describing R(Q,omega)
# SHIST is the sum of the histogram (events are weighted by S(Q,omega))
# from EXCI and by 1/dE when accumulating RHIST!!!)
#------------------------------------------------------------------------
# write(*,*) 'GETSQOM START: '
#/// for R(Q,E) from TRAX:
#/// Initialize random numbers corresponding to exp(-.5 x^2)
if(swraytr.eq.0) then
nph=nxr
if(xr(1,1).eq.0) then ! generate random numbers only once
do i=1,nxr
xr(1,i)=GASDEV()
xr(2,i)=GASDEV()
xr(3,i)=GASDEV()
xr(4,i)=GASDEV()
enddo
endif
whathis=iand(whathis,253) ! set bit2=0 => resol. is from TRAX
else
whathis=ior(whathis,2) ! set bit2=1 => resol. is from ray-tracing
endif
#// set initial pointers NQOM(I)=0 ! NOTE: NQOM(0)=0 always
rq.nqom(0)=0
do item=1,imin-1
rq.nqom(item)=rq.nqom(item-1)
end do
#/// Create the vector of all the (Q,w) events
#/// QOM(4,I) is Q(3),E of an event in r.l.u. coordinates with weight PQOM(I)
#///-------------------------------------------------------------------------
rq.chkqom=0.d0
chkbase=exp(1.d0)
#// Start cycle through ITEM ... No. of a data set
do 30 item=imin,imax
# fill arrays for all data IF(.NOT.mf_active(ITEM)) goto 30 ! only active data sets
sumamc(item)=0
#// Fill HNORM array for all data sets
if(swraytr.eq.0) then
#// Elipsoid main axes are calculated and stored in ADA(i,i) and B(i,j) ***
#// !!! For each ITEM, take the resol. matrix from mf_A(J,K,ITEM) array !!!
do j=1,4
do k=1,4
amat(j,k)=mf_a(j,k,item)
arc(j,k)=mf_mrc(j,k,item)
end do
end do
call CN2RLU_MF(amat,aux4,item)
call DIAG(aux4,ada,b)
do j=1,4
var(j)= 1./sqrt(abs(ada(j,j)))
end do
hnorm(item)=reltr(item)*znorm
else
call KSTACK_N(nph,item) ! Get number of events from ray-tracing
hnorm(item)=relmc(item)*znorm/1000.
do j=1,4
do k=1,4
arc(j,k)=mf_mrc(j,k,item) ! get conversion matrix C&N->RLU
end do
end do
endif
do j=1,4
rq.qom0(j,item)=mf_par(i_qh-1+j,item) ! store Qhkl,E to share with EXCI
end do
rq.nqom(item)=nph+rq.nqom(item-1) ! Increment NQOM by the number of events
#// Start cycle through all events of ITEM-th resol. function
do 29 i=rq.nqom(item-1)+1,rq.nqom(item)
#// !!! Take Q,E from mf_par(k,item) array, not from QH,QK,... !!!
do k=1,4
rq.qom(k,i) = mf_par(i_qh+k-1,item)
end do
if (swraytr.eq.0) then ! events generated by MC/Trax
do k=1,4
xx(k) = xr(k,i-rq.nqom(item-1))*var(k) ! Transform point in unit sphere to res. elipsoid
end do
do k=1,4
do j=1,4
rq.qom(k,i) = rq.qom(k,i)+ b(k,j)*xx(j) !!!!! Attention: B(I,J) is transposed !!!!!!!!!!!!
end do
end do
rq.pqom(i)=1.d0
# write(*,1) 'QOM: ',(rq.QOM(K,I),K=1,4)
# write(*,1) 'VAR : ',(XX(K),K=1,4)
# write(*,1) 'XX : ',(XX(K),K=1,4)
# pause
else ! events are generated by NESS
call GETQE(i-rq.nqom(item-1),item,x,p)
rq.pqom(i)=p
call M4XV4_3(arc,x,aux) ! transform C&N to r.l.u.
do k=1,4
rq.qom(k,i) = rq.qom(k,i)+aux(k)
end do
endif
rq.chkqom=rq.chkqom+chkbase*i*(rq.qom(1,i)+rq.qom(2,i)**2+rq.qom(3,i)**3)
sumamc(item)=sumamc(item)+rq.pqom(i)
rq.iqom(i)=item
29 continue
30 continue
#// set other pointers NQOM
do item=imax+1,mdat
rq.nqom(item)=rq.nqom(item-1)
end do
# pass QOMEGA data to EXCI
call setqomega(rq)
# write(*,*) 'GETSQOM END: ',rq.NQOM(MDAT),rq.CHKQOM
# return rq.CHKQOM as the result
GETSQOM=rq.chkqom
# qom0(1)=2.5
# qom0(3)=3.5
# call getqomarrays(qom0)
# write(*,1) 'qom0: ',qom0(1:4,1)
# write(*,1) 'rq.qom0: ',rq.qom0(1:4,1)
end
#---------------------------------------------------------------------
SUBROUTINE FILLQOMARRAY(port,indx,nimx,nimy,tm)
# Fill QOM array with QHKL,E values from all branches of EXCI model
# QHKL,E values are defined by the 2-dim array according to the PORT
# Fill PQOM with adequate S(Q,E) values returned by EXCI.
# (QOM,PQOM are defined in exci.inc)
# PORT ... defines input QHKL array
# INDX ... dataset index
# NIMX,NIMY ... number of pixels
# TM ... The transformation matrix TM must convert Q from the PORT
# coordinates to rec. lat. units
# CALLS: EXCI, INITEXCI
# CALLED BY: FILLSQ
#---------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'res_grf.inc'
INCLUDE 'exciimp.inc'
# RECORD /QOMEGA/ rq
record /MODEL/ rm
integer*4 nimx,nimy,ix,iy,j,k,i,l,m
integer*4 indx,idata
real*8 tm(4,4),qe(4),qe1(4),dimx,dimy,v6(6),w6(6)
record /VIEWSET/ port
real*8 chkbase
call getmodel(rm)
if (nimx*nimy*rm.nbr.gt.mqom) then
write(smes,*) 'FILLQOMARRAY: dimension of QOM array exceeded!'
stop
endif
do j=1,4
rq.qom0(j,indx)=mf_par(i_qh-1+j,indx) ! store Qhkl,E to share with EXCI
end do
rq.chkqom=0.d0
chkbase=exp(1.d0)
1 format( 'filling array ',$)
2 format( '.',$)
idata=indx
if(idata.le.0.or.idata.gt.mf_max) idata=mf_cur
rq.ndatqom=idata ! inform EXCI about current data index
dimx=(port.wx2-port.wx1)/nimx
dimy=(port.wy2-port.wy1)/nimy
ix=port.ix
iy=port.iy
do j=1,3
qe(j)=0.
enddo
qe(4)=mf_par(i_en,idata)
#! Partitioning of QOM array
rq.nqom(1)=nimx*nimy
do i=2,mdat
rq.nqom(i)=rq.nqom(i-1)
end do
#// fill Q-values in QOM array first
do j=1,nimx
do k=1,nimy
qe(ix)=(j-0.5)*dimx+port.wx1
qe(iy)=(k-0.5)*dimy+port.wy1
call M4XV4_3(tm,qe,qe1)
l=k+(j-1)*nimy
do m=1,3
rq.qom(m,l)=qe1(m)
enddo
#! don't forget to update CHKQOM - check sum for EXCI
rq.chkqom=rq.chkqom+chkbase*l*(rq.qom(1,l)+rq.qom(2,l)**2+rq.qom(3,l)**3)
rq.iqom(l)=idata
# write(*,*) 'FILLQOMARRAY: ',IX,IY,QE(IX),QE(IY),rq.QOM(1:3,L)
# pause
enddo
enddo
#// initialize EXCI with new QOM values
call setqomega(rq)
call INITEXCI(0,0) ! do not call GETSQOM inside INITEXCI (arg2=0)
#// then fill OMEXC and SQOM values into SQOM(4) and PQOM
#// use array positions above i=NIMX*NIMY to store higher branches
write(sout,1)
do j=1,nimx
if (mod(j,8).eq.0) write(sout,2)
do k=1,nimy
l=k+(j-1)*nimy
do m=1,3
qe1(m)=rq.qom(m,l)
enddo
call EXCI(l,qe1,v6,w6)
do i=1,rm.nbr ! cycle through branches
m=l+(i-1)*nimx*nimy
rq.qom(4,m)=v6(i)
rq.pqom(m)=w6(i)
end do
enddo
enddo
call setqomega(rq)
write(sout,*) ' done.'
end
#------------------------------------------------------------
real*4 FUNCTION HIST_LIN(x,y,nx,nmax,item)
#// creates histogram using planar disp. surface
#// works with data set indexed as ITEM
#// - called by PHON command
#C added in ver. 4.77: use event sweeping technique for a3 scans
#------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'exciimp.inc'
# RECORD /QOMEGA/ rq
integer*4 i,j,k,nh,nx(0:mdat),item,ib,ibs,np,nps,nmax
real*8 dum,gnorm,omexc,den,grad(3),qstep(3),enn
real*4 x(nmax),y(nmax),ampl,backg
real*8 omd(3),qq0(3),qom3(3),enstep,zmean,z,da3,aux(3)
real*8 QxQ
real*8 gama(3,0:mhis),eps(0:mhis),edist,edist0,flag,flag0
real*8 eps0
real*4 omd4(3)
call getqomega(rq)
#/// ATTENTION !
#/// if X,Y<>XHIST,RHIST => SHIST, DHIST will not contain up-to-date values
shist(item)=0
do i=nx(item-1)+1,nx(item)
y(i) = 0
do j=1,6
dhist(j,i) = 0
end do
end do
zmean=0.
#//// copy gradient, QE-step and QHKLE for ITEM-th channel
do i=1,3
grad(i)=mf_par(i_gh-1+i,item)
qstep(i)=mf_par(i_dqh-1+i,item)
qq0(i)=mf_par(i_qh-1+i,item)
enddo
den=mf_par(i_den,item)
enn=mf_par(i_en,item)
da3=mf_par(i_da3,item) ! get step in a3
#//// Norm of G(3) in r.l.u.
call QNORM(grad,gnorm,dum)
if(gnorm.eq.0) goto 99 ! no direction => no scan ....
# Dispersion sheet coefficients OMD are derived from G(3) and
# GMOD of RESCAL, they are given in Energy/r.l.u.
do i=1,3
omd(i)=grad(i)*mf_par(i_gmod,item)/gnorm
end do
enstep=-QxQ(omd,qstep)+den ! energy step relative to the disp. sheet
ib=nx(item-1)+1 ! IB is the base index for ITEM-th histogram
# *** Do nothing if the scan is parallel to dispersion surface
# write(*,*) 'DA3 = ',DA3,' OMD = ',(OMD(I),I=1,3)
# pause
if (da3.ne.0) then
if(mf_par(i_gmod,item).eq.0) goto 99
else if (enstep.eq.0) then
goto 99
endif
if(ib.gt.nx(item)) goto 99 ! No space allocated for ITEM-th histogram
# *** prepare fields for A3 scan:
# gama(i,k) .. gradient of disp. after rotation by (X(k)+0.5D0)*DA3
# EPS(k) .. scalar product gama(i,k)*QHKL(k)-(X(k)+0.5D0)*DEN
# dispersion plane (GHKL) is rotated at each step, thus we avoid rotations for each event
if(da3.ne.0) then
eps0=QxQ(omd,qq0) ! scalar product GRAD*QHKL, used for a3 scans
do i=1,3
omd4(i)=omd(i) ! copy to REAL*4 array
enddo
# get OMD vectors rotated to each of a3 values corresponding to bin edges
call ROTA3(omd4(1),-(x(ib)-0.5d0)*da3,gama(1,0))
eps(0)=eps0-enn
do i=1,3
eps(0)=eps(0)+(x(j)-0.5d0)*den
enddo
do j=ib,nx(item)
k=j-ib+1
call ROTA3(omd4(1),-(x(j)+0.5d0)*da3,gama(1,k))
eps(k)=eps0-enn
do i=1,3
eps(k)=eps(k)+(x(j)+0.5d0)*den
enddo
# write(*,1) k,-(X(k+IB-1)-0.5D0)*DA3,GAMA(1,k),GAMA(3,k),EPS(k)
enddo
endif
# *** Accumulate histogram in differences along the SCAN direction ***
do i=rq.nqom(item-1)+1,rq.nqom(item)
# *** A3 scan ***
# ignore steps in QHKL and do the scan in A3 (sample rotation) ...
# => scan is non-linear in QHKL, sweeping algorithm is applied
if (da3.ne.0) then
#1 format(I,6(2x,G12.6))
#2 format(6(2x,G12.6))
# write(*,2) (QOM(j,i),j=1,4)
# write(*,*) 'k DA3 GH GL ENDIST'
k=0 ! go to the left edge of the bins range
do j=1,3
aux(j)=rq.qom(j,i) ! have a copy in REAL*8 array
enddo
edist0=QxQ(gama(1,k),aux(1))-eps(k)-rq.qom(4,i)
flag0 = sign(1.d0,edist0)
# write(*,1) k,-(X(IB)-0.5D0)*DA3,GAMA(1,k),GAMA(3,k),EDIST0
do k=1,nx(item)-ib+1 ! sweep through histogram bins
edist=QxQ(gama(1,k),aux(1))-eps(k)-rq.qom(4,i)
flag = sign(1.d0,edist)
if(flag.ne.flag0) then ! crossed disp. branch
enstep=abs(edist-edist0)
flag0 = flag
nh=k+ib-1
y(nh) = y(nh)+rq.pqom(i)/enstep
shist(item)=shist(item)+rq.pqom(i)/enstep
# write(*,*) 'BIN ',NH,PQOM(i)/ENSTEP,X(NH),Y(NH)
endif
edist0=edist
# write(*,1) k,-(X(k+IB-1)+0.5D0)*DA3,GAMA(1,k),GAMA(3,k),EDIST
enddo
# pause
else
# QHKL,E scan
# or take linear scan in QHKLE ...
do k=1,3
qom3(k)=rq.qom(k,i)-mf_par(i_qh-1+k,item)
end do
omexc=QxQ(omd,qom3)+enn
z=rq.qom(4,i)-omexc
nh = (-z/enstep+0.5-x(ib))+ib
if(nh.ge.ib.and.nh.le.nx(item)) then
y(nh) = y(nh)+rq.pqom(i)/abs(enstep)
shist(item)=shist(item)+rq.pqom(i)/abs(enstep)
endif
endif
end do
do j=ib,nx(item)
y(j) =y(j)*hnorm(item)/sumamc(item)
end do
# *** Linear fit is used to determine scale factor and background automatically ***
# *** only for current spectrum mf_cur
if(npt(item).gt.npt(item-1)) then
ibs=npt(item-1)+1 ! this is base index from which data are taken
nps=npt(item)-npt(item-1)
np=nx(item)-nx(item-1) ! number of points for ITEM-th histogram
call LINFIT(x(ib),y(ib),np,spx(ibs),spy(ibs),spz(ibs),nps,
* ampl,backg,chisqr)
fpar(1)=ampl
fpar(2)=backg
nfpar=2
# parname(1)='Scale'
# parname(2)='Background'
else
fpar(1)=1
fpar(2)=0
nfpar=0
endif
do j=ib,nx(item)
y(j) = fpar(1)*y(j)+fpar(2)
end do
whathis=iand(whathis,251) ! set bit3=0 => planar disp. was used to produce RHIST
whathis=ior(whathis,1) ! set bit1=1 => RHIST is updated
# write(*,*) 'HISTLIN ',WHATHIS
HIST_LIN=1.
return
99 HIST_LIN=-1.
return
end
#------------------------------------------------------------------------------------
real*4 FUNCTION HIST(x,y,nx,nmax,fitpar,np)
# *** Creates simulated histograms ****
# X,Y ... created histograms
# NX ... pointers to last point number for each data set in X,Y arrays
# i.e. the partitioning of the X,Y arrays.
# NMAX ... dimension of X,Y
# FITPAR ... model parameters
# NP ... dimension of fitpar
# X(I) are supposed to be dimensionless steps in units of real scan steps
# i.e. X(I)=-5,-4,...+4,+5 for equidistant data centered at nominal QHKL,E setting
# Added in ver. 4.77: makes scan in A3 if DA3<>0 (step in QHKLE is ignored)
#-------------------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'exciimp.inc'
# RECORD /QOMEGA/ rq
record /MODEL/ rm
integer*4 nx(0:mdat),np,nmax
real*4 x(nmax),y(nmax),fitpar(mpar)
real*8 omexc(6),omexc1(6),sqom(6),qomsc(4),flag0(6),flag(6)
real*8 edist(6),edist0(6)
integer*4 i,j,k,icom,is,ib,item,item0
real*8 dee,dd(3),domexc,da3
10 format( '.',$)
910 format( 'wait ...',$)
# write(*,*) 'HIST: before getqomega'
# pause
call getqomega(rq)
call getmodel(rm)
# write(*,*) 'HIST: after getqomega'
if (rq.nqom(mdat).gt.20000) write(sout,910)
HIST=0.e0
# *** Pass parameters to EXCI
do j=1,np
rm.param(j) = fitpar(j)
# write(*,*) 'HIST, rm.PARAM: ',J,rm.PARAM(J)
enddo
call setmodel(rm)
# *** Clear SCAN histogram ***
do i=1,mdat
shist(i)=0.d0
enddo
do i=1,nx(mdat)
y(i) = 0
do j=1,rm.nbr
dhist(j,i) = 0.e0
end do
end do
icom = 1
item0=0
i=0
# write (*,*) 'START NQOM: ',IQOM(1),NX(1),mf_active(1)
# *** Accumulate histogram in differences along the SCAN direction ***
# START MAIN CYCLE through all events
#------------------------------------------------------------------------------
# write (*,*) 'NQOM(MDAT): ',rq.NQOM(MDAT),rm.nbr
do 39 while (i.lt.rq.nqom(mdat))
i=i+1
icom = i ! passes event number into exci
item=rq.iqom(i) ! ITEM is index to histogram = index of the data set
if(item.le.0) goto 39 ! Exit if there are no data in QOM arrays
#* step in QHKL,E must be determined for each new histogram:
#------------------------------------------------------------------------------
if(item0.ne.item) then ! next histogram begins
if (item.gt.1) write(sout,10)
ib=nx(item-1)+1 ! IB is the base index for ITEM-th histogram
# write (*,*) 'NEW ITEM: ',I,ITEM,IB,NX(ITEM),mf_active(ITEM)
if(ib.gt.nx(item).or.
& (.not.(mf_active(item)).and.item.ne.mf_cur)) then
#* No space has been allocated for the ITEM-th histogram or the data
#* are not active ==> go to the next spectrum
i=rq.nqom(item)
# write (*,*) 'NO SPACE: ',IB,NX(ITEM),mf_active(ITEM)
goto 39
endif
do j=1,3
dd(j)=mf_par(i_dqh-1+j,item)
end do
dee=mf_par(i_den,item)
da3=mf_par(i_da3,item)
item0=item
endif
#1 format(a,I5,4(2x,G12.5))
# if (mod(icom,1000).EQ.0)
# if (icom.gt.0) write (*,1) 'DE: ',ICOM,(DD(j),j=1,3),DEE
#* Start filling the ITEM-th histograms
#--------------------------------------------------------------------------------
if (da3.ne.0) then ! ignore steps in QHKL and do the scan in A3 (sample rotation)
call ROTA3(rq.qom(1,i),(x(ib)-0.5d0)*da3,qomsc(1))
qomsc(4) = rq.qom(4,i)+(x(ib)-.5d0)*dee
else
do k=1,3
qomsc(k) = rq.qom(k,i)+(x(ib)-.5)*dd(k)
enddo
qomsc(4) = rq.qom(4,i)+(x(ib)-.5)*dee
endif
call exci(icom,qomsc,omexc,sqom)
#1 format('HIST, ',a,I5,6(2x,G12.5))
# write(*,1) 'Q: ',NX(ITEM),qomsc,qomsc(4)-omexc(1),SQOM(1)
#// Set flag0=+1/-1 if the event lies ABOVE/BELOW the dispersion surface
do 36 j=1,rm.nbr
edist0(j)=qomsc(4)-omexc(j)
36 flag0(j) = sign(1.d0,edist0(j))
do is=ib,nx(item) ! cycle through histogram bins
if (da3.ne.0) then ! ignore steps in QHKL and do the scan in A3 (sample rotation)
call ROTA3(rq.qom(1,i),(x(is)+0.5d0)*da3,qomsc(1))
qomsc(4) = rq.qom(4,i)+(x(is)+.5d0)*dee
else
do k=1,3
qomsc(k) = rq.qom(k,i)+(x(is)+.5)*dd(k)
enddo
qomsc(4) = rq.qom(4,i)+(x(is)+.5)*dee
endif
# write (*,1) 'QOM: ',IS,(rq.qom(j,i),j=1,4)
call exci(icom,qomsc,omexc1,sqom)
# write (*,1) 'Q: ',IS,(qomsc(j),j=1,4),qomsc(4)-omexc1(1),SQOM(1)
#11 format(a,6(2x,G12.6))
# write(*,11) 'qomsc,EDIST,SQ: ',
# * (qomsc(k),k=1,4),OMEXC1(1)-qomsc(4),SQOM(1)
do j=1,rm.nbr ! cycle through branches
edist(j)=qomsc(4)-omexc1(j)
flag(j) = sign(1.d0,edist(j))
if((flag(j).ne.flag0(j)).and.(rm.wen(j).le.0)) then ! crossed disp. branch
domexc=abs(edist(j)-edist0(j))
flag0(j) = flag(j)
dhist(j,is) = dhist(j,is)+sqom(j)*rq.pqom(i)/domexc
y(is) = y(is)+sqom(j)*rq.pqom(i)/domexc
shist(item)=shist(item)+rq.pqom(i)/domexc
else if (rm.wen(j).gt.0) then ! finite disp. width
dhist(j,is) = dhist(j,is)+sqom(j)*rq.pqom(i)
y(is) = y(is)+sqom(j)*rq.pqom(i)
shist(item)=shist(item)+rq.pqom(i)
endif
edist0(j)=edist(j)
end do
# write (*,1) 'Y: ',IS,Y(is),rq.PQOM(I),rm.wen(1)
end do
# pause
39 continue
# write (*,*) 'HIST loop passed ',SUMAMC(1)
do i=1,mdat
if(sumamc(i).gt.0.and.nx(i).gt.nx(i-1)) then
# write (*,1) 'HNORM,SUMAMC: ',I,HNORM(I)/SUMAMC(I),FITPAR(1)
do j=nx(i-1)+1,nx(i)
# write (*,1) 'Y: ',J,Y(J)
y(j) = fitpar(1)*y(j)*hnorm(i)/sumamc(i)+fitpar(2)
ihist(j)=i
end do
endif
end do
# pause
whathis=ior(whathis,4) ! set bit3=1 => EXCI module was used to produce RHIST
whathis=ior(whathis,1) ! set bit1=1 => RHIST is updated
# write(*,*) 'HIST ',WHATHIS
HIST = 1. ! return 1 if the histogram has been filled
if (rq.nqom(mdat).gt.20000) write(sout,*)
# write (*,*) 'HIST END, WHATHIS= ',WHATHIS
return
end
#
#*********************************************************************************
SUBROUTINE RESFIT(ncmax)
# Perform data fitting with EXCI model and all loaded (and active) data sets
# Result is stored in XHIST,RHIST,NHIST arrays
# Arguments:
# NCMAX ... max. number of iteration steps
# (for NCMAX=0, calculate only current model curve and exit)
#
# Upgraded RESFIT version (28/4/2002)
# Internal dialogs have been removed, the routine is called from
# command handler FIT_CMD, based on LINP library (structured command-line interpreter).
#*********************************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'res_grf.inc'
INCLUDE 'restrax.inc'
real*4 deltaa(mpar),delmin(mpar)
integer*4 ncmax,nh0
integer*4 i,j,nfree,nc,maxitem,irun
real*4 eval,flamda,chiold,coef
real*4 HIST,FCHISQ
# EXTERNAL HIST
800 format( 'CHISQR = ',g12.5)
801 format( 'CHISQR(I) : ',10(1x,g11.4))
906 format( 'NC=',i4, ' CHISQR=',g12.5, ' LAMBDA=',g12.5)
908 format( 'A(I): ',10(1x,g10.4))
909 format( 'Fit not finished. Reached ',i3, ' iterations.')
# write(*,*) 'call RESFIT ',NCMAX
#* initialize EXCI
call INITEXCI(0,1) ! call GETSQOM inside INITEXCI (arg2=1)
# write(*,*) 'call RESFIT after INITEXCI'
# ***** CONSTANTS ******
jfit=1 ! fit is running
mode=1 ! Weighting by 1/SPZ**2
#* No data loaded: only fill in the model histogram
if (npt(mdat).le.0) then
if (ncmax.gt.0) write(sout,*) 'No data to fit ...'
call HISTINIT
eval = hist(xhist,rhist,nhist,nhist(mdat),fpar,nfpar)
jfit=0
return
endif
#* partitioning of histograms has to be compatible with data
do i=1,mdat
nhist(i)=npt(i)
nh0=(nhist(i)-nhist(i-1)+1)/2
do j=nhist(i-1)+1,nhist(i)
xhist(j)=spx(j)
rhist(j)=0
ihist(j)=i
enddo
shist(i)=0
end do
whathis=iand(whathis,254) ! set bit1=0 => RHIST not ready
# write(*,*) 'RESFIT ',WHATHIS
#* calculate number of fitted points and number of degrees of freedom
nfree=0
maxitem=0
do i=1,npt(mdat)
if (ipt(i).gt.0) nfree=nfree+1
if (maxitem.lt.ipt(i)) maxitem=ipt(i)
enddo
nfree=nfree-nfpar
# If NCMAX=0, only calculate histogram and/or CHISQ and exit
if (ncmax.le.0) then
eval = hist(xhist,rhist,nhist,nhist(mdat),fpar,nfpar)
if (nfree.gt.0) then
chisqr=FCHISQ(spy,spz,ipt,npt(mdat),nfree,mode,rhist,dchisq)
else
chisqr=0.d0
endif
jfit=2
return
endif
#* Check the number of free parameters
if(nfree.le.0) then
write(smes,*) 'Not enough data points !'
return
else
write(sout,*) 'Starting to fit ',nfree+nfpar, ' data points...'
endif
#* calculate initial histograms
eval = hist(spx,rhist,npt,npt(mdat),fpar,nfpar)
#* calculate initial CHISQRs
chisqr = FCHISQ(spy,spz,ipt,npt(mdat),nfree,mode,rhist,dchisq)
write(sout,800) chisqr
write(sout,801) (dchisq(i),i=1,maxitem)
#* set initial parameter values
do j=1,nfpar
fpari(j) = fpar(j)
enddo
#* set initial values of FLAMDA, etc..
nc = 1
flamda = fitlam0
chiold = 1.e12
irun=1
# ** get minimum increments
coef=0.01
do i=1,nfpar
delmin(i) = abs(fpar(i))*coef
enddo
# write(*,*) 'RESFIT: start fitting cycle'
# ************ FITTING CYCLE ************
grfarg(0)=5 ! use PLOT_MDAT on plot refresh
do while(irun.gt.0.and.nc.le.ncmax)
#* Adjust increments for the calculation of derivatives
coef = min(flamda,0.01)
do i=1,nfpar
deltaa(i) = abs(fpar(i))*coef
if(deltaa(i).lt.delmin(i)) deltaa(i)=delmin(i) ! avoid floating overflows
deltaa(i) =deltaa(i)*jfixed(i)
enddo
#* make one iteration step
call CURF3T(deltaa,flamda,nfree)
#* inform about CHISQ and FPAR values
write(sout,906) nc,chisqr,flamda
write(sout,801) (dchisq(i),i=1,maxitem)
write(sout,908) (fpar(i),i=1,nfpar)
write(sout,*)
#* show progress, except for flat-cone mode
if (cfgmode.eq.0) call PLOTOUT
#* stop on too small change
if(abs(chiold-chisqr)/abs(chiold).gt.fittol) then
chiold = chisqr
if (abs(chiold).lt.1.e-10) irun=0 ! just for sure, it should not happen
else
irun=0
endif
nc = nc+1
#* stop on iteration limit
if(nc.gt.ncmax.and.irun.ne.0) then
write(sout,909) nc
irun=0
endif
enddo
jfit=2 ! fit finished
# ************ END OF FITTING CYCLE ************
#* get final histogram
call HISTINIT
eval = hist(xhist,rhist,nhist,nhist(mdat),fpar,nfpar)
end
#*****************************************************************************
SUBROUTINE LISTFITPAR
# list fitting parameters in the format:
# number name value [fix]
#*****************************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'exciimp.inc'
INCLUDE 'restrax.inc'
record /MODEL/ rm
integer*4 i
character*3 fix
character*4 ch
character*5 s,CONCAT
1 format(a, ' ',a10, ' ',g13.5,1x,a3,1x,g13.5)
2 format(i4)
call getmodel(rm)
do i=1,rm.nterm
if (rm.fixparam(i).eq.0) then
fix= 'fix'
else
fix= ' '
endif
s= ' '
write(ch,2) i
s=CONCAT( 'a',ch)
write(sout,1) s,rm.parname(i),rm.param(i),fix !,FPAR(I)
enddo
end
#*********************************************************************************
SUBROUTINE SIMDATA(suma,nev,icom)
# Simulate spectra using EXCI, as if they were read from data file(s)
# SUMA ... total counts, if SUMA=0 then use fpar(1) as amplitude
# NEV ... number of events for simulation
# ICOM=0 ... only generate new histogram and calculate chi^2
#*********************************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'res_grf.inc'
integer*4 icom,i,j,item,ib,nev
real*4 eval
real*4 HIST,bcg,FCHISQ
real*8 suma,z,ampl
if (icom.eq.0) goto 10
do item=1,mf_max
do i=1,4
qe0(i,item)=mf_par(i_qh+i-1,item)
dqe0(i,item)=mf_par(i_dqh+i-1,item)
enddo
do i=5,6
dqe0(i,item)=mf_par(i_dqh+i-1,item)
enddo
enddo
#* run MC
call RUNMC(1,nev)
#* fill QOM array
# Z=GETSQOM(1,mf_max,1)
#* initialize EXCI
call INITEXCI(0,1) ! call GETSQOM inside INITEXCI (arg2=1)
#* partitioning of histograms has to be compatible with data
do i=1,mdat
nhist(i)=npt(i)
# NH0=(NHIST(I)-NHIST(I-1)+1)/2
do j=nhist(i-1)+1,nhist(i)
xhist(j)=spx(j)
rhist(j)=0.
ihist(j)=i
enddo
shist(i)=0
end do
#* create histogram for scale=1 and background=0
ampl=fpar(1)
fpar(1)=1. ! scale=1
bcg=fpar(2)
fpar(2)=0. ! background=0
eval = hist(spx,rhist,npt,npt(mdat),fpar,nfpar)
fpar(2)=bcg
# write(*,*) 'suma=', suma
do item=1,mf_max
ib=npt(item-1)+1
if (item.eq.1) then
z=0.
do i=ib,npt(item)
z=z+rhist(i)
enddo
endif
if (z.gt.0) then
if (suma.ne.0.d0) then ! set FPAR(1)=AMPL or normalize on SUM
fpar(1)=suma/z
else
fpar(1)=ampl
endif
do i=ib,npt(item)
spy(i)=rhist(i)*fpar(1)+fpar(2)
spz(i)=sqrt(abs(spy(i))+1.)
spy(i)=spy(i)+erhist(i)*spz(i)
enddo
else
do i=ib,npt(item)
spy(i)=0.
spz(i)=1.
enddo
endif
end do
jfit=2
#* get Chi^2
10 eval = hist(spx,rhist,npt,npt(mdat),fpar,nfpar)
# GRFARG(0)=4
# CALL PLOTOUT
chisqr = FCHISQ(spy,spz,ipt,npt(mdat),npt(mdat),1,rhist,dchisq)
end
#***********************************************************************
real*4 FUNCTION OPTTAS(p)
# returns value to be minimized when optimizing TAS configuration
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'restrax.inc'
real*8 eps
parameter(eps=1.d-30)
integer*4 i
real*4 p(1)
real*8 b,c0,c1,c2,der2,z
real*8 suma
data suma/1.d4/
if (optdpar.le.eps) then
OPTTAS=0.
return
endif
1 format(/,a,6(2x,g10.3))
# write(*,1) 'OPTTAS: ',P(1),RES_DAT(OPTPAR)
# A=RES_DAT(OPTPAR)
call RAN1SEED(10001)
res_dat(optpar)=p(1)
call BEFORE
do i=1,mf_max
mf_par(optpar,i)=p(1)
enddo
call SIMDATA(suma,optev,1)
c0=chisqr
b=fpar(optfpar)
fpar(optfpar)=b+optdpar
call RESFIT(0)
# CALL SIMDATA(SUMA,OPTEV,0)
c2=chisqr
fpar(optfpar)=b-optdpar
call RESFIT(0)
# CALL SIMDATA(SUMA,OPTEV,0)
c1=chisqr
fpar(optfpar)=b
call RESFIT(0)
# CALL SIMDATA(SUMA,OPTEV,0)
# C0=CHISQR
# RES_DAT(OPTPAR)=A
#* 2nd derivative from chi^2 along FPAR(OPTFPAR)
der2=abs(c1+c2-2*c0)/(optdpar**2)
if (optmerit.eq.1) then ! w^2
z=1./sqrt(der2)
else if (optmerit.eq.2) then ! w^2/Intensity
z=fpar(1)/sqrt(der2)/vkiness
endif
# write(*,1) 'OPTTAS: ',P(1),Z,DER2,C0,C1,C2
# write(*,1) 'SWITCHES: ',WHATHIS,SWRAYTR,STP.KF
OPTTAS=z
end