Source module last modified on Wed, 10 May 2006, 16:36;
HTML image of Fortran source automatically generated by
for2html on Mon, 29 May 2006, 15:06.
#////////////////////////////////////////////////////////////////////////
# $Id: restrax_main.f,v 1.13 2006/05/10 14:36:26 saroun Exp $
#
# *** R E S T R A X ***
# written by J.Saroun and J.Kulda, Institute Laue-Langevin, Grenoble,
# and Nuclear Physics Institute, Rez near Prague, November 1995-2005
#
# The program package for simulating the neutron optics of three-axis
# spectrometers (TAS), optimization of their resolution and luminosity
# and evaluation of experimental data collected with them.
# The program code includes both a high-speed analytical (Gaussian) convolution
# algorithm and a Monte Carlo ray-tracing method providing enhanced accuracy
# in description of most of the spectrometer components.
#
# Using in part:
#
# 1) RESCAL - provides all except of the graphics subroutines and
# the subroutine AFILL, which calculates the resolution matrix.
#
# - written by M. Hergreave & P. Hullah, P.C.L. London, July 1979
# - updated by P. Frings, ILL Grenoble, October 1986
#
# 2) TRAX - provides a kernel replacing the subroutine AFILL and taking
# into account real dimensions of the spectrometer components and
# focusing by curved perfect and mosaic crystals.
#
# - written by (*)M.Popovici, A.D.Stoica and I.Ionita, Institute for
# Nuclear Power Reactors, Pitesti, 1984-1986
# (*) University of Missouri
# [1] M.Popovici, A.D.Stoica and I.Ionita, J.Appl.Cryst. 20 (1987),90.
# [2] M.Popovici et al., Nucl.Instrum. Methods A338 (1994), 99.
#
#
#****************************************************************************
# *** For all additional information contact the authors: ***
# *** ***
# *** kulda@ill.fr saroun@ujf.cas.cz ***
# *** ***
#****************************************************************************
#-------------------------------------------------------------
SUBROUTINE BEFORE
# Call whenever the setup may have changed
# Updates calculated auxilliary fields and runs TRAX
#-------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
logical*4 log1
integer*4 ier
common /error/ier
# ATTENTION : order of following routines is important !!!
ier=0
# write(*,*) 'BEFORE 1'
call RECLAT ! compute reciprocal lattice parameters and matrices
# write(*,*) 'BEFORE 2'
call SCATTRIANGLE ! compute and check KI,KF,Q and tras. matrix Lab -> CN
# write(*,*) 'BEFORE 3'
call ANGSCAN(res_dat(i_da3),0.d0) ! scan in DA3 => adjust DH,DK,DL and set DE=0
# write(*,*) 'BEFORE 4'
call TRANSMAT ! create transformation matrices for coordinate systems
# write(*,*) 'BEFORE 6'
call TAS_TO_NESS ! convert RESCAL+CFG to NESS
# write(*,*) 'BEFORE 7'
call TAS_TO_TRAX ! convert RESCAL+NESS to TRAX
# write(*,*) 'BEFORE 8'
call THRAX ! calculate TRAX res. matrix
# write(*,*) 'BEFORE 9'
call GetXSpec ! prepare arrays for x-axis of current measured-spectrum
# write(*,*) 'BEFORE 10'
# check whether setup has been modified
call SPEC_UNCHANGED(log1)
ischanged=.not.log1
call SPEC_GETCHK(checksum)
call MFIT_GET(mf_cur)
needbefore=.false.
# write(*,*) 'BEFORE set ',mf_cur,' done: ',mf_done(mf_cur),
# * ' changed: ',isCHANGED
end
#-----------------------------
SUBROUTINE LOGO
#-----------------------------
implicit none
INCLUDE 'config.inc'
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
1 format(2x, '-----------------------------------------------',/,
& 2x, 'R E S T R A X - Monte Carlo simulation of TAS ',/,
& 2x, 'Version: ',a40,/,
& 2x, 'Build: ',a40,/,
& 2x, '-----------------------------------------------',/,
& 2x, '(C) J.Saroun & J.Kulda',/,
& 2x, 'ILL Grenoble, NPI Rez near Prague',/,
& 2x, ' ',/,
& 2x, 'Using (in part):',/,
& 2x, 'RESCAL by M.Hargreave, P.Hullah & P. Frings',/,
& 2x, 'TRAX by M.Popovici, A.D.Stoica & I.Ionita',/,
& 2x, '-----------------------------------------------',/,
& 2x, 'type ? for command list',/,/)
write(sout,1) package_version,package_date
return
end
#-----------------------------------------------------------------------
SUBROUTINE SETPATH(sarg)
# select search path for data files
#-----------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
character*(*) sarg
character*128 prompt,answer
prompt= ' Path to data files'
answer= ' Data will be searched in'
call DLG_SETPATH(sarg,prompt,answer,0,datpath)
end
#-----------------------------------------------------------------------
SUBROUTINE SETRESPATH(sarg)
# select search path for configuration files
#-----------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
character*(*) sarg
character*128 prompt,answer
prompt= ' Additional search path for configuration files'
answer= ' Configurations will be searched in'
call DLG_SETPATH(sarg,prompt,answer,0,respath)
end
#--------------------------------
SUBROUTINE DOSHELL(comm)
# execute shell command
#--------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
character*(*) comm
character*256 comm1
integer*4 is,l
1 format( ' Command : ',$)
2 format(a)
comm1= ' '
if((comm(1:1).eq. ' ').or.(comm(2:2).eq.char(0))) then
write(sout,1)
read(sinp,2) comm1
else
l=len(comm)
if(l.gt.256) l=256
comm1=comm(1:l)
endif
call BOUNDS(comm1,is,l)
if (l.gt.0) then
call system(comm1(is:is+l-1))
endif
end
#-----------------------------------------------------
SUBROUTINE LISTCFG
# list available configuration files.
# creates appropriate shell script and executes
#-----------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
integer*4 is,il
1 format( '#!/bin/sh',/,
* 'for F in `ls ',a, '*.cfg` ; do',/,
* 'echo $F " ... " `sed -n -e' '2p' ' $F`',/,
* 'done',/)
open(file= '~tmprestrax',unit=22,status= 'UNKNOWN',err=99)
call BOUNDS(respath,is,il)
write(22,1,err=98) ' '
if (il.gt.0) then
write(22,1,err=98) respath(is:is+il-1)
endif
close(22)
call DOSHELL( 'chmod 755 ~tmprestrax')
call DOSHELL( './~tmprestrax')
call DOSHELL( 'rm -f ./~tmprestrax')
return
98 close(22)
99 return
end
#-------------------------------------------------------------
SUBROUTINE SETCFG(sarg,iread)
# Read configuration file
# IF IREAD>0, prepare also all calculated fields and run TRAX!
#-------------------------------------------------------------
implicit none
character*(*) sarg
integer*4 iread
if (iread.gt.0) then
call TAS_READCFG(sarg) ! read CFG
call BEFORE ! prepare all and run TRAX
else
call TAS_READCFG(sarg) ! only read CFG
endif
end
#-----------------------------------------------------------------------
SUBROUTINE LIST
#-----------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'restrax_cmd.inc'
integer*4 i1,i2,i
19 format(3( ' ',a5, ' = ',g14.7,1x))
51 format( ' ',2(a4, ' = ',f10.5,1x),/, ! DM,DA
2 ' ',3(a4, ' = ',f10.2,1x),/, ! ETAM,ETAA,ETAS
3 ' ',3(a4, ' = ',f10.0,1x),/, ! SM,SA,SS
4 ' ',a4, ' = ',f10.5,1x, a4, ' = ',f10.0,1x,/, ! KFIX,FX
5 ' ',4(a4, ' = ',f10.2,1x),/, ! ALF1..4
6 ' ',4(a4, ' = ',f10.2,1x),/, ! BET1..4
7 4( ' ',3(a4, ' = ',f10.4,1x)/), ! AS,AA,AX,BX
1 ' ',4(a4, ' = ',f10.4,1x),/, ! QH..EN
2 ' ',6(a4, ' = ',f10.4,1x),/, ! DQH..DE,DA3,DA4
3 ' ',4(a4, ' = ',f10.4,1x),/, ! GH..GL,GMOD
4 ' ',4(a4, ' = ',f10.4,1x),/, ! ROMH..ROAV
5 ' ',2(a4, ' = ',f10.2,1x)) ! SDI,SHI
if (nos.ge.2) then
i1=nint(ret(1))
i2=nint(ret(2))
write(sout,19)(res_nam(i),res_dat(i),i=i1,i2)
else
write(sout,51) (res_nam(i),res_dat(i),i=1,res_nvar)
nos=0
call SET_3AX(1)
call SET_3AX(2)
call SET_3AX(3)
call SET_3AX(4)
call SET_3AX(5)
call SET_3AX(6)
endif
end
#***********************************************************************
SUBROUTINE UNITS(sarg)
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
character*(*) sarg
2 format( '[THz]')
5 format( '[meV]')
4 format( ' Units are ',a)
if(sarg(1:1).eq. 'T') then
euni=0.24181
write(cunit,2)
else
euni=1.
write(cunit,5)
endif
write(sout,4) cunit
end
#-----------------------------------------------------------
SUBROUTINE GETROANAL(ro)
# generates "optimal" monochromator and analyzer curvatures
# using analytical expression:
# - focus parallel beam at the sample by monochromator
# - set monochromatic focusing on analyzer
#-----------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'trax.inc'
real*4 ro(4)
integer*4 ierr
real*8 tha,chia,thm,chim
real*8 r(res_nvar)
common /error/ierr
equivalence(r(1),res_dat(1))
if (r(i_sm).ne.0) then
thm=pi/sqrt(ei0/hsqov2m)/r(i_dm)
thm=abs(asin(thm))
chim=-himon*tdr
ro(1)=sin(thm+chim)/vl1*100
ro(2)=1./vl1/(2.*sin(thm)*cos(chim))*100
endif
if (r(i_sa).ne.0) then
tha=pi/sqrt(ef0/hsqov2m)/r(i_da)
tha=abs(asin(tha))
chia=-hiana*tdr
# RO(3)=(VL2*SIN(THA+CHIA) + VL3*SIN(THA-CHIA))/2./VL2/VL3*100
ro(3)=sin(tha-chia)/vl2*100
ro(4)=(1./vl2+1./vl3)/(2.*sin(tha)*cos(chia))*100
endif
end
#-----------------------------------------------------------
SUBROUTINE GETRO(icom)
# generates "optimal" monochromator and analyzer curvatures
# ICOM=0 ... focus parallel beam at the sample and from the
# sample to the detector
# *** modified J.S. 3/6/1997:
# ICOM=1 ... optimize Vanad resolution and intensity
# *** modified J.S. 10/12/2002:
# ICOM=0 ... optimize for Rowland focusing on analyzer
#-----------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
character*10 remark,vname(4)
real*4 ro(4),romi(4),roma(4)
real*4 OPTV,tol
integer*4 icom,i,ierr
real*8 r(res_nvar)
common /error/ierr
equivalence(r(1),res_dat(1))
external OPTV
data vname / ' ROMH = ', ' ROMV = ', ' ROAH = ', ' ROAV = '/
data romi /-10.0,-10.0,-10.0,-10.0/
data roma /10.0,10.0,10.0,10.0/
1 format(a,f8.4,a20)
5 format( ' Numerical optimization failed ! ')
call GETROANAL(ro)
tol=0.001
if (r(i_sa)*r(i_sm).ne.0.and.icom.eq.1) then
call LMOPT(OPTV,ro,romi,roma,4,tol,1)
write(sout,*)
endif
if (ierr.lt.0) write(smes,5)
do i=1,4
if((nos.eq.0).or.((nos.ge.i).and.(ret(i).eq.1.))) then
res_dat(i+i_romh-1)=ro(i)
remark= ' changed'
else
remark= ' '
endif
write(sout,1) vname(i),ro(i), ' [m-1] '//remark
enddo
do i=i_romh,i_roav
mf_par(i,mf_cur)=res_dat(i)
enddo
return
end
#--------------------------------------------------------------
SUBROUTINE OPTINSTR
# Optimise TAS with respect to the IPAR-th parameter
# IFPAR ... parameter from EXCI model used as a figure of merit
#
#--------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'res_grf.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
INCLUDE 'exciimp.inc'
record /MODEL/ rm
character*128 line
real*4 a,a0,ami,ama
real*8 lim(2)
real*4 OPTTAS,tol
integer*4 ipar,ifpar,GETICOM,ierr
common /error/ierr
external OPTTAS
1 format(a)
2 format(a, ' = ',f10.4)
call getmodel(rm)
call GENDT
20 call DLG_STRING( 'TAS parameter to be optimised (name)',line,0)
ipar=GETICOM(line)
if (ipar.le.0.or.ipar.gt.res_nvar) goto 20
write(sout,1) res_nam(ipar)// ' will be optimised'
30 call DLG_INPUT( 'lower limit:upper limit',lim,0)
if (lim(1).ge.lim(2)) goto 30
ami=lim(1)
ama=lim(2)
call DLG_INTEGER( 'EXCI parameter as figure of merit',ifpar,0,3,rm.nterm)
write(sout,1) rm.parname(ifpar)// ' is the figure of merit'
silent=2
optmerit=1
optev=10000
tol=0.1
if (nos.gt.0) then
optev=nint(1000*ret(1)) ! number of events
endif
#// parameter 2 is used by GENDT as the integral for simulated data
if (nos.gt.2) then
optmerit=nint(ret(3)) ! figure of merit (see OPTTAS)
endif
if (nos.gt.3) then
tol=ret(4) ! TOL
endif
write(sout,*) 'FM=',optmerit
optfpar=ifpar
optpar=ipar
optdpar=fpar(ifpar)*0.01
a=res_dat(optpar)
a0=a
call LMOPT(OPTTAS,a,ami,ama,1,tol,0)
select case (ierr)
case (1)
write(sout,*) 'OK, exact match'
case (2)
write(sout,*) 'OK, exit on TOL limit'
case (-1)
write(sout,*) 'not finished, stalled on lambda'
case (-2)
write(sout,*) 'not finished, iteration number limit'
end select
write(sout,2) res_nam(optpar),a
silent=dsilent ! return to default
grfarg(0)=4
call PLOTOUT ! plot the result
end
#**********************************************************************
SUBROUTINE REINP(sarg)
# redirection of input
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
character*(*) sarg
character*128 inpfile
integer*4 iuini,iuext,ires
data iuini,iuext/0,10/
if(iuini.eq.0) iuini=sinp
if (sinp.ne.iuini) close(sinp)
if((sarg(1:1).ne. ' ').and.(sarg(1:1).ne.char(0))) then
call OPENRESFILE(sarg, 'inp',iuext,0,0,inpfile,ires)
if (ires.ne.0) goto 2002
sinp=iuext
else
sinp=iuini
endif
return
2002 write(smes,*) 'Cannot open input file '//sarg
sinp=iuini
end
#**********************************************************************
SUBROUTINE REOUT(sarg)
# redirection of input
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
character*(*) sarg
integer*4 iuini,iuext
data iuini,iuext/0,11/
if(iuini.eq.0) iuini=sout
if (sout.ne.iuini) close(sout)
if((sarg(1:1).ne. ' ').and.(sarg(1:1).ne.char(0))) then
open (unit=iuext,err=2002,file=sarg, status= 'UNKNOWN')
sout=iuext
else
sout=iuini
endif
# WRITE(*,*) SOUT
return
2002 write(smes,*) 'Cannot open output file '//sarg
sout=iuini
end
#**********************************************************************
SUBROUTINE EMODE
# Switch on/off elastic mode (for SA=0 only)
# !! J.S. 5/4/2002, allow EMOD also for SA<>0
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax_cmd.inc'
logical*4 emod
common /mode/ emod
data emod /.false./
if(nos.eq.0) then
if(emod) then
write(sout,*) 'Elastic scattering'
else
write(sout,*) 'Inelastic scattering'
endif
else if (ret(1).eq.1.) then
# IF(SA.EQ.0) THEN
emod=.true.
write(sout,*) 'Elastic scattering mode is On'
# ELSE
# EMOD=.FALSE.
# write(sout,*) 'Elastic scattering mode is possible only'//
# * ' with SA=0 !'
# ENDIF
else if (ret(1).eq.0.) then
emod=.false.
write(sout,*) 'Elastic scattering mode is Off'
endif
end
#**********************************************************************
logical*4 FUNCTION CMDFILTER(icmd)
# filter for commands for special modes. Set .false. if the command is forbidden
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
integer*4 icmd,sm,sa
logical*4 log,emod
common /mode/ emod
sm=res_dat(i_sm)
sa=res_dat(i_sa)
log=.true.
if (emod.or.(sa*sm.eq.0)) then ! forbid for elastic and SA,SM=0 modes
#// forbid all commands using TRAX
if (res_nam(icmd).eq. 'BRAG') log=.false. ! BRAG
if (res_nam(icmd).eq. 'PHON') log=.false. ! PHON
if (res_nam(icmd).eq. 'RES') log=.false. ! RES
if (res_nam(icmd).eq. 'FIT') log=.false. ! FIT
if (sa*sm.eq.0.and.(.not.log)) then
write(sout,*) 'SM or SA = 0 => Only Monte Carlo is accepted'
endif
if (emod.and.(.not.log)) then
write(sout,*) 'Elastic mode => Only Monte Carlo is accepted'
endif
else
CMDFILTER=.true.
endif
end
#**********************************************************************
SUBROUTINE MAKEMC(sarg)
# calls Monte Carlo if required
# also set the SWRAYTR switch
#***********************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'res_grf.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
integer*4 n,nmin,is,il
logical*4 log1,emod,lastmode
character*(*) sarg
character*128 s
common /mode/ emod
data lastmode /.false./
save lastmode
is=1
call FINDPAR(sarg,1,is,il)
s=sarg(is:is+il-1)
log1=.false.
#// Set LOG=.TRUE. if the command requires MC simulation
#// commands with 0 arguments
nmin=0
if (s(1:il).eq. 'MBRAG') log1=.true.
if (s(1:il).eq. 'MPHON') log1=.true.
if (s(1:il).eq. 'MFIT') log1=.true.
if (s(1:il).eq. 'OPTAS') log1=.true.
if (log1) goto 10
#// commands with 1 arguments
nmin=1
if (s(1:il).eq. 'PLOT') then
log1=(
1 (grfarg(0).eq.-5).or.
1 (grfarg(0).eq.-4).or.
1 (grfarg(0).eq.-3).or.
1 (grfarg(0).eq.2).or.
2 (grfarg(0).eq.3).or.
3 (grfarg(0).eq.4.and.iand(whathis,2).eq.2).or.
4 (grfarg(0).eq.5.and.iand(whathis,2).eq.2).or.
5 (grfarg(0).eq.9).or.
5 (grfarg(0).eq.15).or.
5 (grfarg(0).eq.16).or.
& .false.)
endif
if (s(1:il).eq. 'MRES') log1=.true.
if (log1) goto 10
#// commands with 3 arguments
nmin=3
if (s(1:il).eq. 'MFWHM') log1=.true.
#// execute ray-tracing
10 if (log1) then
swraytr=1 ! switch ray-tracing on
if (nos.gt.nmin) then
n=nint(ret(nos)*1000)
else
n=lastnev
endif
il=0
if(lastmode.xor.emod) il=1 ! mode has changed, run MonteCarlo anyway
call RUNMC(il,n) ! call MonteCarlo
lastmode=emod
# write(*,*) 'MAKEMC: ',S(1:IL),SWRAYTR
else
swraytr=0 ! switch ray-tracing off
# write(*,*) 'MAKEMC: ',S(1:IL),SWRAYTR
endif
end
#-----------------------------------------------------------
SUBROUTINE MCPHON
# Wrapper for MCPROFIL: scans through planar dispersion (GH..GMOD)
# Returns phonon widths in gaussian approximation
# Plots PAGE 2
# ICOM=1 .. TRAX
# ICOM=2 .. Ray tracing
#-----------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'res_grf.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
real*8 z,GETSQOM
real*4 eval,HIST_LIN
real*8 phon_po,phon_por,phon_pos,phon_mc,phon_mcr,phon_mcs
common /phon/ phon_po,phon_por,phon_pos,phon_mc,phon_mcr,phon_mcs
# Calculate and show scan profile
z=GETSQOM(1,mf_max) ! fill QOM array with events
call HISTINIT ! initialize histogram
eval = HIST_LIN(xhist,rhist,nhist,nhist(mdat),mf_cur) ! generate histogram using planar dispersion
# Print phonon peak parameters
if(swraytr.eq.1) then
call GETPHONWIDTH(aness,phon_mc,phon_mcr,phon_mcs)
write(sout,*) 'Phonon FWHM, ray-tracing '//
1 'resolution function:'
write(sout,906) phon_mc,cunit,phon_mcr,cunit,phon_mcs
write(sout,907) hnorm(mf_cur) ! norm from M.C.
else
call GETPHONWIDTH(atrax,phon_po,phon_por,phon_pos)
write(sout,*) 'Phonon FWHM, analytical (gaussian) '//
1 'resolution function:'
write(sout,906) phon_po,cunit,phon_por,cunit,phon_pos
write(sout,907) hnorm(mf_cur) ! norm from TRAX
endif
906 format(f8.4, ' [A-1]',a5, ' ',f8.4, ' [r.l.u.]'
1 ,a5, ' ',f8.4, ' [steps]')
907 format(1x, 'Integrated intensity:',g11.5)
grfarg(0)=4
call PLOTOUT
end
#
#-----------------------------------------------------------
SUBROUTINE GENDT
# simulation of scan data incl. counting errors
#-----------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'res_grf.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
integer*4 ne
real*8 suma
data suma/1.d4/ ! set default and remember
swraytr=1 ! ray-tracing on
ne=10000 ! default number of events
if (nos.ge.1) ne=int(ret(1)*1000) ! number of events from the 1st argument
if (nos.ge.2) suma=ret(2) ! set optionally the sum of counts as the 2nd argument
call mfit_set(1)
call DELDATA(1,mf_max) ! clear channels
call ADDDATA( 'channel',1,mf_cur,2) ! add a dataset
call SIMDATA(suma,ne,1) ! simulate data incl. errors
grfarg(0)=4
call PLOTOUT ! plot the result
# CALL WriteHist(' ') ! save results
end
#-----------------------------------------------------------------------------
SUBROUTINE EXPORT_RES(arg)
# Exports resolution function to an ASCII file or std. output (SOUT)
# J. Saroun, 19/5/2000
#-----------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'restrax.inc'
character*(*) arg
integer*4 iu,n,i,j,l
real*8 qe(4),p
1 format( 'Qx',12x, 'Qy',12x, 'Qz',12x, 'dE',12x, 'P')
2 format(5(g12.5,2x))
3 format(i8,3x,3(g12.6,2x), ' ... events, (KI,Q,E)' )
call KSTACK_N(n,mf_cur)
if(n.gt.0) then
if((arg(1:1).ne. ' ').and.(arg(1:1).ne.char(0))) then
open (unit=24,err=2002,file=arg, status= 'UNKNOWN')
iu=24
else
iu=sout
endif
write(iu,3) n,stp.ki,stp.q,stp.e
write(iu,1)
do i=1,n
call GETQE(i,mf_cur,qe,p)
write(iu,2) (qe(j),j=1,4),p
enddo
if(iu.eq.24) close(iu)
endif
return
2002 l=len(arg)
write(smes,*) 'Cannot open resol. function, file ',arg(1:l)
return
end
#-----------------------------------------------------------------------------
SUBROUTINE IMPORT_RES(arg)
# Imports resolution function from an ASCII file or std. input (SINP)
# J. Saroun, 19/5/2000
#-----------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'restrax.inc'
character*(*) arg
integer*4 iu,n,i,j,l
real*8 qe(4),p,ki,qq,e0,eps
parameter (eps=1e-5)
1 format( 'WARNING: R(Q,E) for diferent scattering triangle !')
if((arg(1:1).ne. ' ').and.(arg(1:1).ne.char(0))) then
open (unit=24,err=2002,file=arg, status= 'OLD')
iu=24
else
iu=sinp
endif
read(iu,*,err=2001) n, ki,qq,e0
if((abs(ki-stp.ki).gt.eps*stp.ki).or.
* (abs(qq-stp.q).gt.eps*stp.q).or.
* (abs(e0-stp.e).gt.eps*abs(stp.e))) write(smes,1)
read(iu,*,err=2001)
call KSTACK_ALLOCATE(n,mf_cur)
do i=1,n
read(iu,*,err=2001)(qe(j),j=1,4),p
call SETQE(i,mf_cur,qe,p)
enddo
if(iu.eq.24) close(iu)
call GETCOV_QE(n)
call BEFORE
call mfit_get(mf_cur)
return
2001 write(smes,*) 'Format error while loading resol. function '
return
2002 l=len(arg)
write(smes,*) 'Cannot open resol. function, file ',arg(1:l)
return
end
#----------------------------------------------------------------------
SUBROUTINE PRINTOUT
# Prints text with resolution parameters etc.
#----------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'config.inc'
INCLUDE 'inout.inc'
INCLUDE 'res_grf.inc'
INCLUDE 'restrax.inc'
integer*4 stmp,iuser,inow,i
character*20 user,now
500 format(1x,70( '-'))
1 format( ' *** R E S T R A X ',a, ' ***')
2 format( ' configuration: ',a60)
3 format( ' user: ',a60)
4 format( ' results have been saved in restrax.txt')
call pgqinf( 'USER',user,iuser)
call pgqinf( 'NOW',now,inow)
open(unit=20,file= 'restrax.txt',status= 'UNKNOWN')
write(20,1) trim(package_version)
write(20,*)
write(20,2) cfgname
write(20,*) 'user: '//user(1:iuser)// ' '//now(1:inow)
write(20,500)
stmp=sout
sout=20
i=swplot
swplot=0
call BRAG(1)
write(20,500)
call BRAG(2)
write(20,500)
call MCPHON(1)
write(20,500)
call MCPHON(2)
write(20,500)
call RESOL(1,1)
write(20,500)
call RESOL(1,2)
write(20,500)
call RESOL(1,3)
write(20,500)
call RESOL(1,4)
write(20,500)
call RESOL(2,4)
write(20,500)
close(20)
sout=stmp
call PRINTFILE( 'restrax.txt')
swplot=i
write(sout,4)
end
#----------------------------------------------------------------------
SUBROUTINE PRINTFILE(sarg)
# Print file named SARG
# try commands in following order:
# 1) $PGPLOT_ILL_PRINT_CMD sarg
# 2) $PRINTER sarg
# 3) lpr sarg
#----------------------------------------------------------------------
implicit none
character*(*) sarg
character*40 prn_command
integer*4 is,il
call BOUNDS(sarg,is,il)
if (il.gt.0) then
prn_command= ' '
call getenv( 'PGPLOT_ILL_PRINT_CMD',prn_command)
if (prn_command.eq. ' ') call getenv( 'PRINTER',prn_command)
if (prn_command.eq. ' ') prn_command= 'lpr'
call system(prn_command// ' '//sarg(is:is+il-1))
endif
end
#*****************************************************************************
SUBROUTINE REPORTOMEXC
# Print excitation energies and S(Q,E) from the EXCI module.
# Input [h k l] is either command option or QH .. QL for current dataset
#*****************************************************************************
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
INCLUDE 'exciimp.inc'
record /MODEL/ rm
integer*4 i,initdone
real*8 dum4(4),dum6(6),dum61(6)
data initdone/0/
6 format( 'OMEXC at Q=(',3(1x,g9.3), ') E=',g9.3, ':')
7 format( 'E[meV] : ',6(1x,g11.4))
8 format( 'S(Q,E) : ',6(1x,g11.4))
# IF (INITDONE.EQ.0) THEN
call INITEXCI(0,0)
# INITDONE=1
# ENDIF
if (nos.ge.3) then
do i=1,3
dum4(i)=ret(i)
enddo
else
do i=1,3
dum4(i)=mf_par(i_qh+i-1,mf_cur)
enddo
endif
if (nos.ge.4) then
dum4(4)=ret(4)
else
dum4(4)=mf_par(i_en,mf_cur)
endif
call EXCI(-2,dum4,dum6,dum61)
call getmodel(rm)
write(sout,6) (dum4(i),i=1,4)
write(sout,7) (dum6(i),i=1,rm.nbr)
write(sout,8) (dum61(i),i=1,rm.nbr)
end
#-------------------------------------------
SUBROUTINE FCONE_INI
# Set the flat-cone analyzer mode on/off
# depending on command argument RET(1)=1/0
# Creates analyzer channels as multiple datasets
# Dialog arguments:
# DLGARG(1) ... number of channels
#-------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'rescal.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
character*2 chn
integer*4 i,np
data np/32/ ! default number of channels
2 format(i2)
3 format( 'DA4=0, zero step between channels !',//,
& 'Set DA4<>0 and repeat the command')
7 format( ' Analyzer part is in normal position')
71 format( ' Analyzer part is turned up')
72 format( ' Analyzer part is turned down')
73 format( ' Analyzer part is turned up/down')
if(nos.ge.1) then
if (ret(1).eq.1) then
cfgmode=1
else
cfgmode=0
endif
mf_cfgmode(mf_cur)=cfgmode
endif
if (cfgmode.eq.0) write(sout, 7)
if (cfgmode.eq.1.and.stp.sa.gt.0) write(sout,71)
if (cfgmode.eq.1.and.stp.sa.lt.0) write(sout,72)
if (cfgmode.eq.1.and.stp.sa.eq.0) write(sout,73)
if (nos.le.0) return ! only report state
if (cfgmode.eq.1) then
# get number of channels
if (cmdmode.eq.1) then ! interactive mode
call DLG_INTEGER( 'number of channels',np,1,1,mdat)
else ! take NP from the argument array
np=nint(dlgarg(1))
if (np.gt.mdat) np=mdat
if (np.lt.1) np=1
endif
# check DA4<>0, otherwise switch to normal mode and exit
if (np.gt.1.and.res_dat(i_da4).eq.0.d0) then
write(smes,3)
cfgmode=0
mf_cfgmode(mf_cur)=cfgmode
return
endif
# clear channels
call mfit_set(1)
call DELDATA(1,mf_max)
# uses QHKL from the 1st channel as the starting point
call ADDDATA( 'channel01',1,mf_cur,2)
call BEFORE
# scan through all other channels
do i=2,np
write(chn,2) i
if(chn(1:1).eq. ' ') chn(1:1)= '0'
call ROTA4(res_dat(i_da4)*deg,res_dat(i_qh)) ! add setp in A4
call ADDDATA( 'channel'//chn,1,mf_cur+1,2)
call BEFORE
enddo
call mfit_set(1) ! return to the 1st channel
endif
end
#-----------------------------------------------
SUBROUTINE RUNMC(cond,nev)
# Run MC for all loaded data with NEV events
# All calls of Monte Caro should be made through this subroutine !!
# COND=1 ... run anyway
# COND=0 ... run only when necessary
# if NEV<=0 => prompt for number of events
# Dialog arguments:
# DLGARG(0) ... number of events / 1000
#-----------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
integer*4 cond,nev,i,cur0,n
data n/0/
cur0=mf_cur
do i=1,mf_max
call mfit_set(i)
if (nev.gt.0) n=nev
if (nev.le.0.or.n.eq.0.or.n.gt.maxnev) then
if (cmdmode.eq.1) then
call DLG_DOUBLE( 'Number of events / 1000',dlgarg(0),0,1.d-2,2.d+2)
endif
n=int(dlgarg(0)*1000)
endif
if (n.gt.0) call IFNESS(cond,n)
enddo
if (mf_cur.ne.cur0) call mfit_set(cur0)
end
#-----------------------------------------------------
SUBROUTINE READINIFILE(jobname)
# read initialization file
# CFGNAME = configuration file
# DATAPATH = path to the data files
# OPENFILE = data or RESCAL file to open
# return JOBNAME .. filename of a job file to be executed at the startup
#-----------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'restrax.inc'
character*128 line,s
character*(*) jobname
character*128 inifile
integer*4 ires,is,il,ierr
logical*4 log1
real*8 z
1 format(a)
jobname= ' '
call OPENRESFILE( 'restrax.ini', 'ini',22,0,0,inifile,ires)
if(ires.eq.0) then
do while(ires.eq.0)
read(22,1,end=100,iostat=ires) line
if(line(1:1).ne. '#') then
# start with this configuration file:
call READ_STR( 'CFGNAME',line,cfgname,ierr)
# search data files in this directory
call READ_STR( 'DATAPATH',line,datpath,ierr)
# Set default EXCI module (overrides command-line option)
call READ_STR( 'EXCI',line,excilib,ierr)
# Set startup jobfile
call READ_STR( 'JOB',line,jobname,ierr)
# Set fit tolerance
call READ_R8( 'FITTOL',line,z,ierr)
if (ierr.ge.0) fittol=z
# Set initial weight on steepest descent
call READ_R8( 'FITLAM0',line,z,ierr)
if (ierr.ge.0) fitlam0=z
# Open this file first
call READ_STR( 'OPENFILE',line,s,ierr)
#/ for OPENFILE, test whether it is data file or a *.res file
if (ierr.ge.0) then
is=1
call FINDPAR(s,1,is,il)
log1=(il.gt.4)
if (log1) log1=(log1.and.s(is+il-4:is+il-1).eq. '.res') ! stupid construct, required by Absoft debugger ...
if (log1) then ! RESCAL file
datname= ' '
rescal_name=s(is:is+il-1)
else ! data file
datname=s(is:is+il-1)
rescal_name= ' '
endif
endif
endif
enddo
100 close(22)
endif
end
#-----------------------------------------------------------------------------
SUBROUTINE RESINIT
#// initialize RESTRAX
#// include all actions necessary to allocate memory,
#// initialize variables, print LOGO etc..
#-----------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
INCLUDE 'res_grf.inc'
INCLUDE 'config.inc'
INCLUDE 'restrax.inc'
INCLUDE 'restrax_cmd.inc'
integer*4 i,m,j,il,is,ires
character*60 s1
real*4 GASDEV1
character*10 prompt
character*128 s,jobname
integer*4 iseed,irnd
common /rndgen/ iseed,irnd
integer*4 iargc
call MKUPCASE(sysname)
if (sysname(1:7).eq. 'WINDOWS') then
pathdel= '\'
endif
# REAL*8 DUM4(4),DUM6(6),dum61(6)
# normalizing constant for intensities
znorm=2.d+6*hsqov2m*sqrt(24/pi) ! monitor efficiency included
# initialize LINP
prompt= 'ResTrax'
call LINPSET(res_nvar+res_ncmd,prompt,res_nam,res_hlp)
call LINPSETIO(sinp,sout,smes)
# Handle command-line options
m=iargc()
# Derive path to default setup files from executable pathname
call getarg(0,s)
i=index(s,pathdel// 'bin'//pathdel)
if (i.gt.1) then
cfgpath=s(1:i)// 'setup'//pathdel
else
cfgpath= 'setup'//pathdel
endif
# I=IGETARG(1,S,30)
do i=1,m
call getarg(i,s)
if (s(1:2).eq. '-s') then
read(s(3:30),*) j
if(j.ne.0) then
iseed=abs(j)
write(sout,*) 'SEED=',iseed
endif
endif
if (s(1:3).eq. '-gs') then
read(s(4:30),*,err=9) grfsave
write(sout,*) 'GRFSAVE=',grfsave
goto 10
9 grfsave=0
endif
10 if (s(1:2).eq. '-t') then
read(s(3:30),*) j
call RAN1SEED(iseed)
write(sout,*) 'Test of the random number generator:'
call RAN1TEST(j,1000000*j)
goend=1
endif
if (s(1:4).eq. '-sil') then
read(s(5:30),*,err=12) dsilent
write(sout,*) 'SILENT=',dsilent
endif
12 if (s(1:5).eq. '-ran1') then
irnd=1
write(sout,*) 'Numerical Recipes RAN1 generator'
endif
if (s(1:5).eq. '-rand') then
irnd=2
write(sout,*) 'System random number generator'
endif
if (s(1:5).eq. '-help'.or.s(1:1).eq. '?') then
do j=1,9
write(sout,*) hlpopt(j)
enddo
goend=1
endif
if (s(1:5).eq. '-dir=') then
call BOUNDS(s,is,il)
is=is+5
il=il-5
if(il.gt.0) then
respath=s(is:is+il-1)
write(sout,*) 'dir='//s(is:is+il-1)
endif
endif
if (s(1:6).eq. '-exci=') then
call BOUNDS(s,is,il)
is=is+6
il=il-6
if(il.gt.0) then
excilib=s(is:is+il-1)
write(sout,*) 'exci='//s(is:is+il-1)
endif
endif
enddo
if (goend.ne.0) goto 30
silent=dsilent ! set default silence leveel
call RAN1SEED(iseed) ! Initialize random number generator
call LOGO ! print LOGO
call SETRESPATH(respath) ! set default path for configuration
call READINIFILE(jobname) ! read restrax.ini file
# CALL EXCI(-999)
# CALL EXCI(-999,DUM4,DUM6,DUM61) ! only test call to EXCI, set title etc...
# CALL INITEXCI(1)
call UNITS(cunit) ! set units for energy
call SETPATH( ' ') ! Ask for the data path
call SETCFG(cfgname,0) ! Read the configuration file, but without running TRAX etc...
#// Read environment variables:
#// PGPLOT default device
call getenv( 'PGPLOT_DEV',s1)
if(s1(1:1).ne. ' ') then
call SELGRFDEV(s1,1)
endif
#// Optional plotting command, used e.g. to invoke GSView for plotting restrax.ps
call getenv( 'RESTRAX_PLOT',s1)
if(s1(1:1).ne. ' ') then
extplot=s1
else
extplot= ' '
endif
mf_cur=1
mf_max=1
#// initialize array with gaussian random numbers
do i=1,maxd
erhist(i)=GASDEV1(0.,3.)
enddo
#// Do initial partitioning of the histogram
call HISTINIT
#// load default exci library (without initialization)
call SETEXCI(excilib,0)
#// load parameters or a data file
call ILLNameParse(datname,1) ! if name is an integer, convert to nnnnnn format (ILL convention)
call OPENFILE( ' ',ires)
if (ires.le.0.and.rescal_name.ne. ' ') then ! try the filename defined in restrax.ini
call OPENFILE(rescal_name,ires)
endif
if (ires.gt.0) then
if (jobname(1:1).ne. ' ') then
call REINP(jobname)
call LINPSETIO(sinp,sout,smes)
endif
return
endif
20 format(/, 'Look first for a valid data filename.',/,
* 'It should be either data in the ILL format or a RESCAL file.',/,
* 'Look in the RESTRAX folder for a template (./demo/*.res)')
write(smes,20)
write(smes,*)
30 goend=1
end
#-----------------------------------------------------------------------------
SUBROUTINE RESEND
#// end of RESTRAX
#// include all actions necessary to deallocate memory etc...
#-----------------------------------------------------------------------------
implicit none
INCLUDE 'const.inc'
INCLUDE 'inout.inc'
call REINP( ' ') ! reset input to STDIN
call REOUT( ' ') ! reset output to STDOUT
call NESSEND ! NESSEND must be called to deallocate
write(smes,*) ' -> End of ResTrax'
call releaseexci ! release EXCI module
stop
end