src/resgraph1.f

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

Source module last modified on Tue, 9 May 2006, 20:05;
HTML image of Fortran source automatically generated by for2html on Mon, 29 May 2006, 15:06.


#//////////////////////////////////////////////////////////////////////
#////
#////  R E S T R A X   4.8
#////
#////  Graphics output subroutines (PGPlot library required)
#////  LEVEL 1 subroutines:
#////  RESOLUTION FUNCTIONS, SCAN DATA 
#//////////////////////////////////////////////////////////////////////
        
#************************* RESOLUTION FUNCTIONS ********************************

#----------------------------------------------------------------------
      SUBROUTINE SAVEGRFDATA(port,x,y,n,comment)
#     save data for the graph
#----------------------------------------------------------------------
      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'res_grf.inc'
      
      integer*4 n, i_io,i,ires    
#      REAL*4 X(*),Y(*)
      real*4 x(n),y(n)
      character*72 fname
      character*(*) comment
      record /VIEWSET/ port
      
1     format(a)
4     format(1x,4(2x,g13.5))
12    format(a50)
13    format( ' output filename: ',$)
      
      if (toprint.eq.1) return   ! no file output when printing
      i_io=22
      fname= ' '
      call DLG_FILESAVE(fname, ' ', 'dat',1,0,ires,fname)
      if (ires.gt.0) then
        open(unit=i_io,file=fname,err=999,status= 'Unknown')
        write(i_io,1,err=998) port.head
        write(i_io,*) trim(comment)
        write(i_io,1,err=998) port.xtit// '  '//port.ytit
        do i=1,n
          write(i_io,4,err=998) x(i),y(i)
        enddo  
998     close(i_io)
      endif
      return
      
999   write(smes,*)  'Cannot open output  file as unit ',i_io      
      return      
      end

#----------------------------------------------------------------------      
      SUBROUTINE VIEWSCAN
# Plot R(Q,w) in given direction (integrate over others)
# ITASK is taken from GRFARG(1) array item:
#   ITASK=1 ... QH
#   ITASK=2 ... QK
#   ITASK=3 ... QL
#   ITASK=4 ... EN
#   ITASK=5 ... KF
#   ITASK=6 ... Delta(Phi)=spin-echo phase shift
#----------------------------------------------------------------------      
      implicit none

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

      real*8 ymax
      integer*4  i,itask
      integer*4 nn      
      real*4 xx(65),yy(65)
      real*8 suma,center,fwhm,wspread
      record /VIEWSET/ port
      character*20 ind(6)
      character*55 comment
      data ind/ 'Q\dX\u [\A\u-1\d] ', 'Q\dY\u [\A\u-1\d] ',
     * 'Q\dZ\u [\A\u-1\d] ', '\gDE  meV ', 'k\df\u [\A\u-1\d] ', '\gD\gf'/

14    format( 'fwhm: ',g10.3, ' spread: ',g10.3, ' center: ',g10.3)
15    format( 'SE resolution: ',g10.4, ' [ueV]')      

      itask=nint(grfarg(1))
      if (itask.le.0.or.itask.gt.6) itask=4  ! do scan along kf as default

      nn=25
      call RES_SCAN(itask,xx,yy,nn,suma,center,fwhm,wspread)
      
      port.xtit=ind(itask)
      port.ytit= ' counts '
      port.head= ' R(Q,\gw) profile '
      if (itask.eq.6) port.head= ' Precession phase distribution '
            
      ymax=0
      do i=1,nn
          ymax=max(ymax,1.d0*yy(i))
      enddo
      if (ymax.lt.1e-10) ymax=1e-10
      comment= ' '
      write(sout,14) fwhm,wspread,center
      write(comment,14) fwhm,wspread,center
      port.dx1=0.1
      port.dx2=0.9
      port.dy1=0.3
      port.dy2=0.8
      port.wx1=xx(1)
      port.wx2=xx(nn)
      port.wy1=0.
      port.wy2=ymax*1.1
      call PlotCurve(port,xx,yy,nn,2,4,1)
      call pgslw(2)
      call pgsch(1.0)      
      call pgsci(1) 
      call pgmtext( 'B',10.,-0.05,0.0,comment)
      call pgiden
#      CALL SAVEGRFDATA(PORT,XX,YY,NN,comment)
      if (itask.eq.6) then
        write(sout,15) fwhm*hbar/stp.tauf*1000
      endif
      end

#----------------------------------------------------------------------      
      SUBROUTINE RES_SCAN(itask,xx,yy,nn,suma,center,fwhm,wspread)
# Makes a scan through R(Q,w) in given direction (integrate over others)
#   ITASK=1 ... QH
#   ITASK=2 ... QK
#   ITASK=3 ... QL
#   ITASK=4 ... EN
#   ITASK=5 ... KF
#   ITASK=6 ... Delta(Phi)=spin-echo phase shift
# OUTPUT:
#   XX(NN),YY(NN) ... arrays with scan variable and intensity
#----------------------------------------------------------------------      
      implicit none

      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'restrax.inc'
      
      integer*4 nn      
      real*4 xx(nn),yy(nn)

      real*4 xmin(4),xmax(4)
      real*8 x(4),p,z,dx,kf
      integer*4  ix,i,j,nh0,nph,itask
      real*8 suma,center,fwhm,wspread,phi

14    format( 'fwhm: ',g10.3, ' spread: ',g10.3, ' center: ',g10.3)
15    format( 'SE resolution: ',g10.4, ' [ueV]')


      call KSTACK_N(nph,mf_cur)
      if(nph.le.0) then
           write(smes,*)  ' M.C. events not generated !'           
           return
      endif

      nh0=(nn+1)/2
      do i=1,nn
          xx(i)=(i-nh0)
          yy(i)=0.
      end do      
      
      if (itask.ge.1.and.itask.le.4) then
         ix=itask
         call GETMSCALE(xmin,xmax)
         dx=(xmax(ix)-xmin(ix))/nn
         if (dx.le.1e-10) return
         do i=1,nph
             call GETQE(i,mf_cur,x,p)
             z=x(ix)-xmin(ix)
             j=int(z/dx)+1
             if((j.ge.1).and.(j.le.nn)) then
                yy(j)=yy(j)+p
             endif
         end do
      else if (itask.eq.5) then
         ix=1
         xmin(ix)=+1e+10
         xmax(ix)=-1e+10
         do i=1,nph
             call KSTACK_KF(i,mf_cur,kf,p)
             xmin(ix)=min(1.d0*xmin(ix),kf)
             xmax(ix)=max(1.d0*xmax(ix),kf)
         end do
         dx=(xmax(ix)-xmin(ix))/nn
         xmin(ix)=xmin(ix)-dx
         xmax(ix)=xmax(ix)+dx
         dx=(xmax(ix)-xmin(ix))/nn
         if (dx.le.1e-10) return
         do i=1,nph
             call KSTACK_KF(i,mf_cur,kf,p)
             z=kf-xmin(ix)
             j=int(z/dx)+1
             if((j.ge.1).and.(j.le.nn)) then
                yy(j)=yy(j)+p
             endif
         end do
      else if (itask.eq.6) then
         ix=1
         xmin(ix)=+1e+10
         xmax(ix)=-1e+10
         do i=1,nph
             call KSTACK_PHI(i,mf_cur,phi,p)
             xmin(ix)=min(1.d0*xmin(ix),phi)
             xmax(ix)=max(1.d0*xmax(ix),phi)
         end do
         dx=(xmax(ix)-xmin(ix))/nn
         xmin(ix)=xmin(ix)-dx
         xmax(ix)=xmax(ix)+dx
         dx=(xmax(ix)-xmin(ix))/nn
         if (dx.le.1e-10) return
         do i=1,nph
             call KSTACK_PHI(i,mf_cur,phi,p)
             z=phi-xmin(ix)
             j=int(z/dx)+1
             if((j.ge.1).and.(j.le.nn)) then
                yy(j)=yy(j)+p
             endif
         end do
      endif

      do i=1,nn
          xx(i)=(i-1+0.5)*dx+xmin(ix)
      enddo
      call GETPEAKPARAM(xx,yy,nn,suma,center,fwhm,wspread)

      end


#----------------------------------------------------------------------      
      SUBROUTINE RES_IMAGE(id)
# Plot R(Q,E) for dataset(s) projected on a specified plane in 
# (h,k,l,E) space. Asks for the projection plane interactively.
# INPUT:
#   ID    ... dataset index, if=0 then merge all data sets
#----------------------------------------------------------------------            
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'res_grf.inc'
     
      integer*4 id,nima
      parameter(nima=128)
      real*4 aima(nima,nima)     
      character*16 labels
      character*50 comment
      integer*4 ix,iy,i,j,i_io,is,il,is1,il1
      real*8 xmax,ymax,xmin,ymin
  
      data labels/ 'h:k:l:dE [meV]'/
      
9     format(128(1x,g10.4))
13    format( 'scale (',g12.5, ',',g12.5, ',',g12.5, ',',g12.5, ')')

#// graph attributes are taken from GRFARG array above index=40
      ix=nint(grfarg(1))
      iy=nint(grfarg(2))
      xmin=grfarg(3)
      xmax=grfarg(4)
      ymin=grfarg(5)
      ymax=grfarg(6)
      comment=grfstr
        
      call MRES_ALL(id,ix,iy,xmin,xmax,ymin,ymax,comment,aima,nima)
           
      if (toprint.eq.1) return  ! dont save ASCII data when printing
      if (grfsave.ge.2) then 
        i_io=25
        close(i_io)
        open(unit=i_io,file= 'res_2d.dat',err=999,status= 'Unknown')
        call FINDSTRPAR(labels, ':',ix,is,il)
        call FINDSTRPAR(labels, ':',iy,is1,il1)
        write(i_io,*) 'projection ('//labels(is:is+il-1)// ', '
     &  //labels(is1:is1+il1-1)// ')'
        write(i_io,13) xmin,xmax,ymin,ymax                  
        do i=1,nima
          write(i_io,9) (aima(i,j),j=1,nima)
        enddo
        close(i_io) 
        return      
999     write(*,*)  'Cannot open file as unit ',i_io      
        return
      endif           

      end

#--------------------------------------------------------------------     
      SUBROUTINE AB_IMAGE(icom)
# Plot projections on the scattering plane (for flatcone arrangement) 
# Asks for scale and comment string, calculates range and calls AB_MAP
#   ICOM  ... as in AB_MAP
#--------------------------------------------------------------------            
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'lattice.inc'
      INCLUDE 'rescal.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'
     
      
      integer*4 nima
      parameter(nima=128)
      real*4 aima(nima,nima)     
      character*50 s,fname
      character*128 comment
      integer*4 i,j,i_io,icom
      real*8 xmax,ymax,xmin,ymin,scale,aux(4),dx,dy
      data scale /1.d0/

1     format(128(1x,g10.4))
2     format( 'scale (',g10.4, ',',g10.4, '),(',g10.4, ',',g10.4, ')')

      comment=grfstr          
      scale=grfarg(1)
      if (scale.lt.0.01.or.scale.gt.100) scale=1.
      if (icom.eq.ig_sqmap) then 
        call M4xV4_3(mabr,res_dat(i_qh),aux)
        xmin=aux(1)-scale/2.
        xmax=aux(1)+scale/2.
        ymin=aux(2)-scale/2.
        ymax=aux(2)+scale/2.
      else  
#      write(*,*) 'FCONE_RANGE'  
        call FCONE_RANGE(xmin,xmax,ymin,ymax,.true.) 
        dx=(scale-1.0)*(xmax-xmin)/2.
        dy=(scale-1.0)*(ymax-ymin)/2.       
        xmin=xmin-dx
        xmax=xmax+dx
        ymin=ymin-dy
        ymax=ymax+dy
#      write(*,*) 'FCONE_RANGE done'  
      endif
      write(*,2)  xmin,xmax,ymin,ymax         
22    call AB_MAP(xmin,xmax,ymin,ymax,comment,aima,nima,icom)
      
      if (grfsave.ge.2) then 
#// write map to a default ASCII file
        if (icom.eq.ig_fcres) then
           fname= 'fcres_2d.dat'
        else  if (icom.eq.ig_fcdata) then 
           fname= 'fcdata_2d.dat'
        else  if (icom.eq.ig_sqmap) then 
           fname= 'sqmap_2d.dat'
        endif                  
        i_io=25
        close(i_io)
        open(unit=i_io,file=fname,err=999,status= 'Unknown')
        call FORMAT_HKL(mf_par(i_ax,1),s,50)
        write(i_io,*)  'x-axis :',s
        call FORMAT_HKL(mf_par(i_bx,1),s,50)
        write(i_io,*)  'y-axis :',s
        write(i_io,2) xmin,xmax,ymin,ymax      
        write(i_io,*) comment      
        do j=1,nima
          write(i_io,1) (aima(i,j),i=1,nima)
        enddo
        close(i_io)           
        return      
999     write(*,*)  'Cannot open file as unit ',i_io      
        return
      endif
      end
      
#--------------------------------------------------------------------
      SUBROUTINE PlotResol(port,icom)
# Plot Resolution function + dispersion branches on the given viewport
# Resolution function is plotted with x-axis parallel to GH,GK,GL and 
# y-axis along energy transfer.      
# The kind of R(Q,E) representation depends on the  SWRAYTR and WHATHIS switches
#  PORT         ... plotting viewport
#  ICOM=1       ... show an ellipse instead of a cloud for ray-tracing resolution
# CALLS: GetProj,GETSCALE,FILARRAY,
# CALLED BY: PAGE2
#--------------------------------------------------------------------
      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc
      INCLUDE 'rescal.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'
      INCLUDE 'exciimp.inc'
      
      record /MODEL/ rm
      
      integer*4 nima,nline
      parameter(nima=64,nline=31)
      
      record /VIEWSET/ port
      character*60 leg1,leg2,comment
      character*3 chind(3)
      real*8 en,de
      real*8  a(4,4),a1(4,4),mcn(4,4)
      real*8  azero(4),zero(4),rzero(4),v6(6),w6(6),qen(4),aux(4)
      real*8 gnr,gna,shx,shy
      real*4  aima(nima,nima)
      real*4  xfx(nline),fy(nline),xl(2),yl(2),x1,x2,y1,y2,rx
      integer*4 icom,i,j,k,ib,np, ierr,ibh,nh1
      real*8 centre,range,z
      logical*4 escan
      real*8 QxQ
      common /error/ ierr

      data zero/0.,0.,0.,0./
      en=res_dat(i_en)
      de=res_dat(i_den)

3     format( 'GMOD = ',f10.2, ' ',a3, '/r.l.u.')
101   format( '[',f6.2, ' ',f6.2, ' ',f6.2, '] / r.l.u.')
102   format( 'Q = [',f7.3, ' ',f7.3, ' ',f7.3, ']')
103   format( 'E = ',f7.2, ' ',a3)
104   format( '[ ',a3, ' , ',a3, ' , ',a3, ' ] / r.l.u.'

            
      ierr=0
#// Check, whether dQ=0  =>  E=const scan
      escan=((abs(delq(1))+abs(delq(2))+abs(delq(3))).eq.0)


#// R(Q,E) from Monte Carlo
      if (swraytr.ne.0) then
         if(aness(1,1).eq.0) goto 998    ! no events 
         call GetProj(aness,amean,a1,rzero,mcn)   ! Get projection in r.l.u.
         call GETSCALE(a1,1,4,0,x1,x2,y1,y2)   ! Get extent
         rzero(4)=rzero(4)+en
#// R(Q,E) from TRAX
      else
         if(atrax(1,1).eq.0) return
         call GetProj(atrax,zero,a,azero,mcn)
         call GETSCALE(a,1,4,0,x1,x2,y1,y2)
         azero(4)=azero(4)+en
      endif


#// Define limits from the extent of the resolution ellipsoid

      call CLRPORT(port)
      port.wx1=2.5*x1
      port.wx2=2.5*x2
      port.wy1=2.5*y1+en
      port.wy2=2.5*y2+en
      port.ix=1
      port.iy=4
      shx=0.
      shy=0.

      np=npt(mf_cur)-npt(mf_cur-1)    ! number of points in current data set
      ib=npt(mf_cur-1)+1              ! base index for the in current data set
      
      nh1=nhist(mf_cur)-nhist(mf_cur-1)
      ibh=nhist(mf_cur-1)+1            
            
#// If a spectrum is loaded, define limits from the data range

      if ((np.gt.0)) then
         centre=(spx(ib)+spx(np+ib-1))/2.
         range=abs(spx(np+ib-1)-spx(ib))
         if(de.ne.0) then
            centre=centre*de+en
            range=range*de
            port.wy1=centre-range/2
            port.wy2=centre+range/2
            shy=en-centre
         else
            call QNORM(delq,gnr,gna)
            centre=centre*gnr
            range=range*gnr
            port.wx1=centre-range/2
            port.wx2=centre+range/2
         endif
#//  Calc. difference btw. data and R(Q,E) centre
         call QNORM(grd,gnr,gna)
         do i=1,3
            aux(i)=qe0(i,mf_cur)
         end do   
#         SHX=(QxQ(QHKL,GRD)-QxQ(AUX,GRD))/GNR
         shy=en-qe0(4,mf_cur)
      else
         range=abs(xhist(nh1+ibh-1)-xhist(ibh))
         if(de.ne.0) then
            range=range*de
            port.wy1=-range/2+en
            port.wy2=+range/2+en
         else
            call QNORM(delq,gnr,gna)
            range=range*gnr
            port.wx1=-range/2
            port.wx2=+range/2
         endif      
      endif
      
      

#// format the x-label to write scan dir. in symbolic format if possible 
      j=0
      z=0
      do i=1,3
        if (grd(i).ne.0.and.grd(i).ne.z) then
          j=j+1
          if(z.eq.0) z=grd(i)
        endif  
      end do
      if (j.eq.1) then
        do i=1,3
          if(grd(i).ne.0) then 
            chind(i)= '\gc'
          else
            chind(i)= ' 0 '
          endif
        end do
        write(port.xtit,104) (chind(i),i=1,3)        
      else
        write(port.xtit,101) grd
      endif      
      port.ytit= 'E '//cunit
      port.head= ' '

      call PLOTFRAME(port,1,1,1.,0)

      xl(1)=port.wx1
      xl(2)=port.wx2
      yl(1)=res_dat(i_gmod)*xl(1)+en
      yl(2)=res_dat(i_gmod)*xl(2)+en

      if(swraytr.ne.0) then
        if (icom.eq.1) then
          call PLOTELL(port,2,1,a1,rzero,0)           ! projection by M.C.
          call PLOTELL(port,2,2,a1,rzero,1)           ! section    by M.C.
        else
          call FILARRAY(port,0,shx,shy,aima,nima,nima,-1,mcn)
          call PLOT2D(port,aima,nima,nima,nima,nima,0.)   ! image by M.C.
          if (grfsave.ge.2) then
            call WriteMap( 'ness_2d.mat',aima,nima,port.xtit,port.ytit,
     &        port.wx1,port.wx2,port.wy1,port.wy2,qhkl,en)
          endif
          call PLOTFRAME(port,1,1,1.,0)
        endif
      else
        call PLOTELL(port,1,1,a,azero,0)              ! projection by TRAX
        call PLOTELL(port,1,2,a,azero,1)              ! section    by TRAX
      endif

      xl(1)=port.wx1                                 ! dispersion branches
      xl(2)=port.wx2
      yl(1)=res_dat(i_gmod)*xl(1)+en
      yl(2)=res_dat(i_gmod)*xl(2)+en
#      write(*,*) 'PlotResol: WHATHIS=',WHATHIS
      if(iand(whathis,4).eq.4) then  ! EXCI dispersion
         call getmodel(rm)
         call QNORM(grd,gnr,gna)
         do i=1,nline
           xfx(i)=(i-(nline-1)/2-1)*(xl(2)-xl(1))/(nline-1)
         end do
#         write(*,*) 'nbr = ',nbr
         do j=1,rm.nbr
           do i=1,nline
              do k=1,3
                 qen(k)=qhkl(k)+grd(k)/gnr*xfx(i)
              end do   
              qen(4)=en
              call  EXCI(-1,qen,v6,w6)
#              IF(W6(J).EQ.0) GOTO 12
              fy(i)=v6(j)
#1        format(2(3x,G12.4))
#              write(*,1) XFX(I),FY(I)
           end do
           if(rm.wen(j).eq.0) then
              call pgsls(1)
           else
              call pgsls(2)    ! dashed line if S(Q,E) <> S(Q)*delta(E)
           endif      
           call pgline(nline,xfx,fy)
           if (grfsave.ge.1) then
             comment= 'dispersion curve'
             call SAVEGRFDATA(port,xfx,fy,nline,comment)
           endif  
12         continue
         end do
      else
         call pgline(2,xl,yl)
      endif
      call pgsls(1)

#//   legend:
      call pgsch(0.8)
      rx=0.05
      write(leg1,102) qhkl
      write(leg2,103) en,cunit(2:4)
      call pgmtext( 'T',-1.6,rx,0.0,leg1(1:29))
      call pgmtext( 'T',-3.2,rx,0.0,leg2)
      if(iand(whathis,4).eq.0) then  ! planar dispersion
         write(leg2,3) res_dat(i_gmod),cunit(2:4)
         call pgmtext( 'T',-4.8,rx,0.0,leg2)
      endif
      call pgsch(1.0)


      return
998   ierr=1
      end
      

#--------------------------------------------------------------------
      SUBROUTINE PlotResQE(port,n,icom,smap)
# Plot Resolution function for N-th dataset 
# Similar to PlotResol, but shows R(Q,E) projected on the scattering plane.
# The kind of R(Q,E) representation depends on the  SWRAYTR and WHATHIS switches
#  INPUUT:
#    PORT       ... plotting viewport
#    N          ... dataset index
#    ICOM=1     ... show an ellipse instead of a cloud for ray-tracing resolution
#    SMAP=1     ... plots also S(Q) contours
# CALLS: FILLARRAY, FILLSQ, GETSCALE
# CALLED BY: PLOT_MRES
#--------------------------------------------------------------------
      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc
      INCLUDE 'lattice.inc'
      INCLUDE 'rescal.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'
      
      integer*4 n,nima,nline,ic
      parameter(nima=64,nline=31)
      integer*4 smap
      record /VIEWSET/ port
      character*60 leg1
      character*60 s
      real*8  rzero(4)
      real*8  shx,shy
      real*4  aima(nima,nima)
      real*4  rx
      real*4  chsize
      integer*4 icom,i,j,ierr,lw,nf,ibh
      real*8 range,qq(4),mat(4,4),aux(4,4),dq(4),mcn(4,4)
      real*8 z1,z2
      real*8 amc(4,4),atr(4,4),amc0(4)
      real*4 x1,y1,x2,y2
      integer*4 nlevels
      parameter(nlevels=5)
      real*4 levels(nlevels)      
      common /error/ ierr
      real*8 a(3),b(3)
      equivalence (a(1),res_dat(i_ax))
      equivalence (b(1),res_dat(i_bx))
      data levels / .1,.3,.5,.7,.9/

3     format( 'GMOD = ',f10.2, ' ',a3, '/r.l.u.')
101   format( '[',f6.2, ' ',f6.2, ' ',f6.2, '] / r.l.u.')
102   format( 'E = ',f7.3, ' ',a)
104   format( '\gc [',f3.1, ' ',f3.1, ' ',f3.1, ']')
105   format( '\gc [',i2, ' ',i2, ' ',i2, ']')

      ierr=0

#      NP=NPT(N)-NPT(N-1)   ! number of points in current data set

# Get res. matrices etc. for N-th dataset
      do i=1,4
         rzero(i)=0.
         amc0(i)=mf_amean(i,n)
         do j=1,4
           amc(i,j)=mf_amc(i,j,n)
           atr(i,j)=mf_a(i,j,n)
         enddo
      enddo
            
#// Get QHKL,E in A,B coordinates
      call M4xV4_3(mabr,mf_par(i_qh,n),qq)

#// Get step size in A,B coordinates
      call M4xV4_3(mabr,mf_par(i_dqh,n),dq)

#// Get matrix which transforms vectors from C&N to AB
      call M4XM4_3(mabr,mf_mrc(1,1,n),mcn)

#// R(Q,E) from Monte Carlo
      if (swraytr.ne.0) then
         if(amc(1,1).eq.0) goto 998    ! no events 
         call CN2RLU_MF(amc,aux,n)  ! convert from C&N to r.l.u.
         call RLU2AB(aux,mat)    ! convert from r.l.u. to A,B
         call M4xV4_3(mcn,amc0,rzero)
#        write(*,*) 'Plot MC'
      else
#// R(Q,E) from TRAX
         if(atr(1,1).eq.0) goto 998
         call CN2RLU_MF(atr,aux,n)  ! convert from C&N to r.l.u.
         call RLU2AB(aux,mat)       ! convert from r.l.u. to A,B
#        write(*,*) 'Plot TRAX'
      endif
      call GETSCALE(mat,1,2,0,x1,x2,y1,y2)   ! Get extent
      
#10    format(a,5(2x,G12.6))
#      do i=1,4
#         write(*,10) 'MAT: ',(MAT(i,j),j=1,4)
#      enddo  
#      write(*,10) 'GETSCALE: ',N,X1,X2,Y1,Y2 
      
      nf=nhist(n)-nhist(n-1)    ! number of points in current histogram
      ibh=nhist(n-1)+1           ! base index for current histogram
      
      if (nf.gt.0.and.(abs(dq(1))+abs(dq(2)).gt.1e-6)) then
#// Define limits from the step size
        z1=xhist(ibh)*dq(1)
        z2=xhist(ibh+nf-1)*dq(1)
        x1=min(z1,z2)   !MIN(Z1,Z2,X1)
        x2=max(z1,z2)   !MAX(Z1,Z2,X2)      
        z1=xhist(ibh)*dq(2)
        z2=xhist(ibh+nf-1)*dq(2)
        y1=min(z1,z2)   !MIN(Z1,Z2,Y1)
        y2=max(z1,z2)   !MAX(Z1,Z2,Y2)
        
        x2=max(x2,y2)
        x1=min(x1,y1)
        port.wx1=x1+rzero(1)+qq(1)
        port.wx2=x2+rzero(1)+qq(1)
        port.wy1=x1+rzero(2)+qq(2)
        port.wy2=x2+rzero(2)+qq(2)
      else
#// Define limits from the extent of the resolution ellipsoid
#/ set equal scale on x and y
        x2=max(x2,y2)
        x1=min(x1,y1)
        range=2.
        if(icom.ne.1.and.swraytr.ne.0) range=3.     ! wider range for ray-tracing clouds
        port.wx1=range*x1+rzero(1)+qq(1)
        port.wx2=range*x2+rzero(1)+qq(1)
        port.wy1=range*x1+rzero(2)+qq(2)
        port.wy2=range*x2+rzero(2)+qq(2)
      endif  
      port.ix=1
      port.iy=2
      port.head= ' '
      call CLRPORT(port)
      shx=0.
      shy=0.            

#/ axes titles
      if (int(a(1)+a(2)+a(3)).eq.int(a(1))+int(a(2))+int(a(3))) then
        write(s,105) (int(a(i)),i=1,3)
      else
        write(s,104) (a(i),i=1,3)
      endif
      port.xtit=s
      if (int(b(1)+b(2)+b(3)).eq.int(b(1))+int(b(2))+int(b(3))) then
        write(s,105) (int(b(i)),i=1,3)
      else
        write(s,104) (b(i),i=1,3)
      endif
      port.ytit=s
      
      call pgqlw(lw)      

#// Set optimum character size

      chsize=max(0.7,2.3*(port.dx2-port.dx1))
      chsize=min(chsize,1.2)
      call pgsch(chsize)
      
#// plot R(Q) as an ellipse or a gray-scale image      
      call pgslw(lw)

      if (icom.ne.1.and.swraytr.ne.0) then
         call FILARRAY(port,n,shx,shy,aima,nima,nima,1,mcn)
         call PLOT2D(port,aima,nima,nima,nima,nima,0.)   ! image by M.C.
         call PLOTFRAME(port,1,1,chsize,0)
      else    ! ellipsoid   
        ic=1
        if(swraytr.ne.0) ic=2
        rzero(1)=(port.wx1+port.wx2)/2.
        rzero(2)=(port.wy1+port.wy2)/2.
        call PLOTELL(port,ic,1,mat,rzero,0)           ! projection 
        call PLOTELL(port,ic,2,mat,rzero,1)           ! section    
        call PLOTFRAME(port,1,1,chsize,0)
      endif

      if(smap.eq.1) then
#// MRAB transforms vectors from AB to r.l.u. 
        call FILLSQ(port,aima,nima,nima,mrab)
        call pgslw(1)
        call pgsci(4)
        call PLOTCMAP(port,aima,nima,nima,levels,nlevels)
      endif
      
      call pgslw(lw)
      call pgsci(1)      
      call pgsls(1)

#//   legend:

50    call pgsch(chsize)
      rx=0.05
      write(leg1,102) mf_par(i_en,n),cunit
      call pgmtext( 'T',-1.6,rx,0.0,leg1(1:25))
      call pgsch(1.0)


      return
998   ierr=1
      end