src/res_grf.f

Fortran project RESTRAX, source module src/res_grf.f.

Source module last modified on Sat, 6 May 2006, 15:54;
HTML image of Fortran source automatically generated by for2html on Mon, 29 May 2006, 15:06.


#//////////////////////////////////////////////////////////////////////
#////  $Id: res_grf.f,v 1.6 2006/05/06 13:54:58 saroun Exp $
#////  R E S T R A X   4.6
#////
#////  Graphics output subroutines (PGPlot library required)
#////  Top level subroutines
#////
#////
#//////////////////////////////////////////////////////////////////////

#-------------------------------------------------
      SUBROUTINE SELGRFDEV(sarg,iquiet)
#  Select an existing graphics device
#  if SARG given, use it as requested device string,
#  otherwise start dialog
#  no info if IQUIET>0
#-------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'res_grf.inc'
      
      character*(*) sarg
      integer*4 nd,nf,iquiet
      parameter(nd=32)
      integer*4 isel,i,ic,iq,iss,ils,lls,isd,ild,lld,ipos,ila,lla
      character*128 s
      character*8 dtype(nd),sd
      
1     format(a)
2     format( 'current graphics device: ',a,2x, '(',i3, ')')      
3     format( 'selected graphics device: ',a,2x, '(',i3, ')')      
      isel=0  ! selection index
      
      if (iquiet.le.0) write(sout,2) trim(devstr),devid
      if (sarg.ne. ' ') then  ! non-interactive mode
        ic=19  ! only one attempt
        iq=1  ! quiet, non-interactive
        call BOUNDS(sarg,iss,ils)
        lls=iss+ils-1
        s=sarg(iss:lls)
        iss=1
        ils=lls
      else   ! interactive mode
        s=devstr
        call BOUNDS(s,iss,ils)
        lls=iss+ils-1
        iq=0  ! interactive
        ic=0  ! attempts counter
      endif  
              
#// make up to 20 attempts to select a valid device      
      do while (isel.le.0.and.ic.lt.20)
        s=s(iss:lls)
        call LISTGRFDEV(dtype,nd,nf,iquiet)  ! list available drivers
        if (iq.le.0) call DLG_STRING( 'select device: ',s,1)        
        call BOUNDS(s,iss,ils)
        lls=iss+ils-1  ! last non-space character  
            
        call LASTSUBSTR(s, '/',ipos)  ! find device type substring
        if (ipos.gt.0) then
          ila=min(lls-ipos+1,8)  ! max. 8 characters for device type
          lla=ipos+ila-1
          call MKUPCASE(s(ipos:lla))
# search for matching device
          isel=0
          i=1
          do while (isel.le.0.and.i.le.nf)
            sd=dtype(i)  ! get upper case version of the i-th device type
            call MKUPCASE(sd)
            call BOUNDS(sd,isd,ild)  ! calculate SD string limits
            lld=isd+ild-1
#       write(*,*) '<',S(IPOS:LLA),'><',SD(ISD:LLD),'>'
            if (index(sd(isd:lld),s(ipos:lla)).eq.1) isel=i  ! check match
#           IF (ILA.EQ.ILD.AND.S(IPOS:LLA).EQ.SD(ISD:LLD)) ISEL=I ! check match
            i=i+1
          enddo
        else
          isel=0  ! not a valid device name
        endif     
        if (isel.le.0) ic=ic+1
      enddo
      if (isel.gt.0) then
         if (ipos.gt.iss) then
           devstr=s(iss:ipos-1)//sd(isd:lld)
         else
           devstr=sd(isd:lld)
         endif
         devid=isel
      else
         if (iquiet.le.0) write(smes,*) "device not available: ",trim(s)
      endif       
      write(sout,3) trim(devstr),devid
      end


#-------------------------------------------------
      SUBROUTINE LISTGRFDEV(dtype,nd,nf,iq)
#  Print an indexed list of available devices on sout
#  return device types in DTYPE array
#  return number of found devoces in NF
#  IQ .. quiet if IQ>0
#-------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      integer*4 nd,nf,iq
      integer*4 i,tlen,dlen,inter,is,il
      character*8 dtype(nd)
      character*64 descr
      character*3 ci(0:1)
      data ci / 'no ', 'yes'/
1     format( '    name    [int]  description ',/,79( '-'))
2     format(a10,2x, '[',a3, ']',2x,a)
      
      if (iq.le.0) write(sout,1)
      call pgqndt(nf)
      do i=1,nf
        call pgqdt(i, dtype(i), tlen, descr, dlen, inter)
        call BOUNDS(dtype(i)(1:tlen),is,il)
        if (iq.le.0) then
           write(sout,2) dtype(i)(is:is+il-1),trim(ci(inter)),descr(1:dlen)
        endif
      enddo
      end
      
#-------------------------------------------------
      SUBROUTINE INITGRF(iq)
#  Initialization of the PGPlot graphics device.
#  IQ=0 ... output to the current device
#  IQ=1 ... output to the PostScript file "restrax.ps"
#  IQ=2 ... prompts for another output device
#-------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
#      INCLUDE 'inout.inc'
      INCLUDE 'res_grf.inc'

      character*60 dst,aux
      integer*4 iq,idevstr,pgopen,ifirst
      save idevstr,ifirst

      if(iq.eq.2) call SELGRFDEV( ' ',0)
      dst=devstr
      if(iq.eq.1) dst= '"restrax.ps"/vps'
     
201   idevstr=pgopen(trim(dst)) 
#      write(*,*) IDEVSTR
      if (idevstr.le.0) then
          write(*,*)  'PGOPEN not succesful, set to /NULL'
          call SELGRFDEV( '/NULL',1)
          call pgqinf( 'DEV/TYPE',devstr,devid)
          return
      end if
      
#// set portrait window at the beginning
      if (ifirst.eq.0) then
        aux=dst
        call MKUPCASE(aux)
        if (index(trim(aux), '/XSERV').gt.0) then
          call pgpap(4.0,1.4)
          ifirst=1
        endif
      endif
#// open the page
      call pgpage
      if(iq.eq.1) then
         call pgslw(4)  ! set wider lines for printing
      else         
         call pgslw(2)
      endif 
#// keep device info up to date
      if(iq.ne.1) call pgqinf( 'DEV/TYPE',devstr,devid)
      end

#----------------------------------------------------------------------
      SUBROUTINE PLOTOUT
# Top-level subroutine for handling graphics output in RESTRAX.
# Uses command line arguments stored in common RET,NOS array (inout.inc)
# TOPRINT>0  ... print the last output. Repeat the last command with the output 
#                device redirected to a PostScript file and prints. 
#                Plotting subroutines should use RET(11..40) to store interactively 
#                entered arguments for subsequent printing.
# GRFARG(0)=ICOM   ... command ID
# Calls appropriate subroutine accorrding to the ICOM value:
#   -3     ... RES_IMAGE           
#   -4     ... RES_IMAGEALL        
#   -5     ... AB_IMAGE(ig_FCRES)  
#   -6     ... AB_IMAGE(ig_SQMAP)
#   -1..3  ... PAGE1(ICOM)
#    4..8  ... PAGE2(ICOM-4) (single data) or  PLOT_MDATA (multiple data)     
#   14..18 ... PLOT_MRES(ICOM-14) 
#    9     ... VIEWSCAN(NINT(RET(3)))
# CALLED BY: PLOT command by RESTRAX
#----------------------------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'

      integer*4 icom

      
      if (swplot.eq.0) then
        write(sout,*)  'graphics output switched off'
        return
      endif
      
      icom=nint(grfarg(0))

#// open PGPLOT device
      call INITGRF(toprint)

#// Dispatch plot jobs according to ICOM
      if(icom.eq.-3) then
          call RES_IMAGE(mf_cur)    ! R(Q,E) projection
      else if(icom.eq.-4) then
          call RES_IMAGE(0)         ! R(Q,E) projection merged for all data sets
      else if(icom.eq.-5) then
          call AB_IMAGE(ig_fcres)   ! R(Q,E) for flat-cone
      else if(icom.eq.-6) then
          call AB_IMAGE(ig_sqmap)   ! S(Q) map in E=const.
      else if(icom.eq.-7) then
          call AB_IMAGE(ig_fcdata)   ! data for flat-cone
      else if(icom.le.3) then
          call PAGE1(icom)         ! 4 projections in C&N space + parameters
      else if(icom.eq.4) then
          call PAGE2(0)  ! plot both res.function and scan profile       
      else if(icom.eq.5) then
          call PLOT_MDATA  ! plot scan profiles for multiple data 
      endif      
      if(icom.eq.9) call VIEWSCAN   ! scan through R(Q,w)
      if(icom.eq.14.or.icom.eq.15) call PLOT_MRES(1)  ! plot res. function ellipsoid
      if(icom.eq.16) call PLOT_MRES(0)  ! plot res. function - image  
#// close PGPLOT device
      call pgclos
      if (extplot.ne. ' ') then
        call DOSHELL(extplot)
      endif

#// Try to send restrax.ps on a printer defined as:
#// 1) 'PGPLOT_ILL_PRINT_CMD'  .... PGPLOT enviroment variable
#// 2) 'PRINTER'               .... system enviroment variable
#// 3)  lpr                    .... system standard printer 

      if (toprint.eq.1) then
         call PRINTFILE( 'restrax.ps')
         toprint=0
      endif
      end


#----------------------------------------------------------------------
      SUBROUTINE PAGE1(icom)
# Creates the page 1: 4 projections of the resolution ellipsoid + param. list
# for current dataset
#   ICOM=0..1 ... plots only the resol. ellipsoid from TRAX
#   ICOM=2    ...  plots also external ellipsoid from M.C.
#   ICOM=3    ...  plots also 2D-image from M.C.
#   ICOM=-1   ...  plots also a direction of the scan made by a position-sensitive 
#                  detector if it is used instead of a conventional one.
# CALLS: FILARRAY, PSDSCAN
# CALLED BY: PLOTOUT(-1..3) 
#----------------------------------------------------------------------
      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'rescal.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'

      integer*4 nima,icom,k,ierr,iuser,inow
      parameter(nima=64)

      real*8 zero(4),psd(4,2)
      real*4 aima(nima,nima),rx(2),ry(2)
      character*20 user,now
      character*12 fname
      common /error/ierr
      data zero/0.,0.,0.,0./

10    format( 'ness_2d',i1, '.mat')

      if ((icom.eq.2).or.(icom.eq.3)) then
        if(aness(1,1).eq.0) goto 998
      else
        if(atrax(1,1).eq.0) return
      endif

      if((icom.eq.2.or.icom.eq.3).and.(aness(1,1).eq.0)) goto 998

      call MK_PORTS(icom)              ! prepares 4 viewports

#/// writes the headline:
      call pgvport(0.1,0.97,0.9,1.0)
      call pgwindow(0.0,1.0,0.0,1.0)
      call pgsch(1.3)
      call pgscf(2)
      call pgptext(0.5,0.75,0.0,0.5,
     * 'Projections of the Resolution Function')
      call pgsch(1.0)
      call pgscf(1)
      call pgqinf( 'USER',user,iuser)
      call pgqinf( 'NOW',now,inow)
      call pgiden


#/// **** Begin of the cycle in which the 4 projections are plotted ***

      do 200 k=1,4

      call PLOTFRAME(vport(k),1,1,0.75,1)
      if (icom.eq.3) then
        call FILARRAY(vport(k),0,0.d0,0.d0,aima,nima,nima,0)          
        if (grfsave.ge.2) then
          write(fname,10) k
          call WriteMap(fname,aima,nima,vport(k).xtit,vport(k).ytit,
     &        vport(k).wx1,vport(k).wx2,vport(k).wy1,vport(k).wy2,
     &        mf_par(i_qh,mf_cur),mf_par(i_en,mf_cur))
        endif

        call PLOT2D(vport(k),aima,nima,nima,nima,nima,0.)
        call PLOTFRAME(vport(k),1,1,0.75,0)
      endif
      if (atrax(1,1).ne.0) then
        call PLOTELL(vport(k),1,1,atrax,zero,0)
        call PLOTELL(vport(k),1,2,atrax,zero,1)
      endif
      if ((icom.eq.2).and.(aness(1,1).ne.0)) then
        call PLOTELL(vport(k),2,1,aness,amean,0)
        call PLOTELL(vport(k),2,2,aness,amean,1)
      endif
200   continue

      if (icom.eq.-1) then
        call PSDSCAN(psd(1,1),psd(1,2),psd(2,1),psd(2,2),psd(4,1),
     1               psd(4,2))
        psd(3,1)=0
        psd(3,2)=0
        do k=1,4
           rx(1)=psd(vport(k).ix,1)
           rx(2)=psd(vport(k).ix,2)
           ry(1)=psd(vport(k).iy,1)
           ry(2)=psd(vport(k).iy,2)
           call PLOTLINE(vport(k),2,1,rx,ry,2)
        end do
      endif

#      return

      call GRLIST   ! list of parameters

      return

998   write(smes,*)  'No Monte Carlo events to plot!'
      return
      end

#----------------------------------------------------------------------
      SUBROUTINE PAGE2(icom)
# Creates the page 2: R(Q,E) + dispersion branches in the upper part 
# and the scan curve bellow
#   ICOM  ... as in PlotResol and PlotResQE
#----------------------------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'

      record /VIEWSET/ port

      integer*4 icom,ierr,np,ib,i,ic
      real*4 xfx(101),fy(101)
      common /error/ ierr
      logical*4 emod
      common /mode/ emod

#/// upper part

      ierr=0

      port.dx1=0.15
      port.dx2=0.93
      port.dy1=0.63
      port.dy2=0.99

      if (emod) then
        call PlotResQE(port,mf_cur,icom,1)
      else
        call PlotResol(port,icom)
      endif
      

#/// lower part, plot whatever is in current histogram
      port.dy1=0.12
      port.dy2=0.52
      if (iand(whathis,1).eq.1) then  ! check that the histogram is available

        np=nhist(mf_cur)-nhist(mf_cur-1)    ! number of points in current histogram
        ib=nhist(mf_cur-1)+1                ! base index for current histogram
        do i=1,np
           fy(i)=rhist(i+ib-1)
           xfx(i)=xhist(i+ib-1)
        end do
        if (swraytr.eq.0) then   ! TRAX result
          ic=1
        else            ! ray-tracing result
          ic=2
        endif      
        call PlotScan(port,xfx,fy,np,ic,4,2)
      endif
      if (ic.eq.1) call pgsci(2)
      if (ic.eq.2) call pgsci(1)
      call LegFit(port,0.8,6.5)
      call LegFile(port)
      call pgiden

      return
998    write(smes,*)  'No Monte Carlo events to plot !'
      return
      end


#----------------------------------------------------------------------
      SUBROUTINE PLOT_MRES(icom)
# Plots multiple res. functions, each in one cell, starting with
# the current channel (mf_cur). Maximum is 20 cells.
#   ICOM  ... as in PlotResQE
#----------------------------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
#      INCLUDE 'inout.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'

      record /VIEWSET/ port
      integer*4 icom,j,nrow,ncol,ncell
            
      ncell=mf_max-mf_cur+1   ! Number of data sets to plot on page
      if (ncell.gt.16) ncell=20  ! no more than 20 cells on a page
      
      ncol=nint(sqrt(ncell*1.))
      nrow=int((0.999*ncell)/ncol)+1            
      if(ncol.gt.2) then
         call pgslw(2)
      else
         call pgslw(3)
      endif         
      do j=1,ncell
          call SETPORTCELL(port,nrow,ncol,j) 
          call PlotResQE(port,j+mf_cur-1,icom,0)
      enddo

      call pgiden

      return
      end

#----------------------------------------------------------------------
      SUBROUTINE PLOT_MDATA
# Plots multiple data on one page (complementary to PLOT_MRES)
#   ICOM=0..3   ... used to set colors
# CALLS: PLOTCELL,AB_IMAGE(flatcone only)
# CALLED BY: PLOTOUT(4..8) if mf_max>1
#----------------------------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'

      record /VIEWSET/ port
      real*4 chsize
      integer*4 j,ic,il,nrow,ncol,ncell
            
      ncell=mf_max+1    ! number of plot cells  
      if (ncell.gt.20) ncell=20   ! no more than 20 cells on a page   
      
      ncol=nint(sqrt(ncell*1.))
      nrow=int((0.999*ncell)/ncol)+1            
      if(ncol.gt.2) then
         call pgslw(2)
      else
         call pgslw(3)
      endif               
      if (iand(whathis,2).eq.0) then   ! bit2=0 => TRAX resol.
         ic=1
         il=2
      else
         ic=2
         il=1
      endif      
      do j=1,ncell-1
          call SETPORTCELL(port,nrow,ncol,j) 
          call PLOTCELL(port,j,ic,0,1)
      enddo
      call SETPORTCELL(port,nrow,ncol,ncell) 
      call CLRPORT(port)
      chsize=max(0.7,2.5*(port.dx2-port.dx1))
      chsize=min(chsize,1.2)
      call pgsci(il)
      call LegFit(port,chsize,0.01)

      call pgiden

      return

      end