src/resgraph3.f

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

Source module last modified on Thu, 6 Apr 2006, 11: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 3 subroutines (plotting primitives)
#////  PLOT ELEMENTS, SETTING VIEWPORTS, DEFINE PLOT SCALES, FORMAT LABELS,
#////  SAVE RESULTS 
#////
#//////////////////////////////////////////////////////////////////////

#*************************  PLOT ELEMENTS ********************************

#-------------------------------------------------------------
      SUBROUTINE PLOT2D(port,a,nx,ny,ndx,ndy,logscale)
# plots a gray-scale map of the array A to the viewport PORT
# assumes positive (>=0) values in A matrix, negative ones are ignored
#-------------------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
      INCLUDE 'res_grf.inc'

      integer*4 i,j,nx,ny,ndx,ndy 
      real*4 a(ndx,ndy),cx,dx,cy,dy,tr(6),amin,amax
      record/VIEWSET/ port
      real*4 logscale  ! base for logarithmic scale
      real*4 aux(mimax,mimax),l0,range

      dx=(port.wx2-port.wx1)/nx
      cx=(port.wx2+port.wx1)/nx

      dy=(port.wy2-port.wy1)/ny
      cy=(port.wy2+port.wy1)/ny

      tr(1)=-0.5*dx+port.wx1
      tr(2)=dx
      tr(3)=0.
      tr(4)=-0.5*dy+port.wy1
      tr(5)=0.
      tr(6)=dy
      amin=1.d20
      amax=-1.d20
      do 20 i=1,nx
      do 20 j=1,ny
          if (a(i,j).gt.amax) amax=a(i,j)
          if (a(i,j).lt.amin) amin=a(i,j)
20    continue
      if (logscale.le.0.0) then
        call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)
        call pgwindow(port.wx1,port.wx2,port.wy1,port.wy2)
        call pggray(a,ndx,ndy,1,nx,1,ny,amax*1.1,0.,tr)
      else if (amax.gt.0.and.amin.lt.amax) then
        l0=log(10.)
        if (amin.le.0) amin=1.e-20
        range=log(amax/amin)/l0
        if (range.gt.16.0) then
          amin=amax*exp(-16.0*l0) 
          range=16.0
        endif  
        if (range.lt.1.0) then
          amin=amax*exp(-1.0*l0) 
          range=1.0
        endif  
        do i=1,nx
        do j=1,ny
          if (a(i,j).gt.0) then
            aux(i,j)=log(a(i,j)/amin)/l0
          else
            aux(i,j)=0.d0
          endif              
        enddo
        enddo        
        call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)
        call pgwindow(port.wx1,port.wx2,port.wy1,port.wy2)
        call pggray(aux,mimax,mimax,1,nx,1,ny,range*1.1,0.,tr)
      endif  

      return
      end
#-------------------------------------------------------------
      SUBROUTINE PLOTCMAP(port,a,ndx,ndy,c,nc)
#   plots a contour map of the array A to the viewport PORT
#-------------------------------------------------------------
      implicit none

      INCLUDE 'res_grf.inc'

      integer*4 nx,ny,ndx,ndy,nc 
      real*4 a(ndx,ndy),cx,dx,cy,dy,tr(6),c(nc)
      record/VIEWSET/ port

      nx=ndx
      ny=ndy

      dx=(port.wx2-port.wx1)/nx
      cx=(port.wx2+port.wx1)/nx

      dy=(port.wy2-port.wy1)/ny
      cy=(port.wy2+port.wy1)/ny

      tr(1)=-0.5*dx+port.wx1
      tr(2)=dx
      tr(3)=0.
      tr(4)=-0.5*dy+port.wy1
      tr(5)=0.
      tr(6)=dy
      call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)
      call pgwindow(port.wx1,port.wx2,port.wy1,port.wy2)
      call pgcont (a,ndx,ndy,1,nx,1,ny, c, nc, tr)

      return
      end

#---------------------------------------------------------------
      SUBROUTINE PLOTFRAME(port,icl,ils,zch,iax)
#  plots frame, axes and titles of the viewport PORT with given
#  collor index (ICL), line style (ILS) and character size (ZCH)
#  If IAX=0, axes at (0,0) are not plotted.
#---------------------------------------------------------------
      implicit none

      INCLUDE 'res_grf.inc'
      
      integer*4 icl,ils,iax
      real*4 zch

      record/VIEWSET/ port

      call pgsch(zch)
      call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)
      call pgwindow(port.wx1,port.wx2,port.wy1,port.wy2)
      call pgsci(icl)
      call pgsls(ils)
      if(iax.eq.0) then
         call pgbox( 'BCNST',0.0,0, 'BCNST',0.0,0)
      else
         call pgbox( 'ABCNST',0.0,0, 'ABCNST',0.0,0)
      endif
      call pglab(port.xtit,port.ytit,port.head)
      call pgsch(1.0)

      return
      end
      
#------------------------------------------------------------------
      SUBROUTINE PLOTELL(port,icl,ils,a,u,is)
#  plots projection (IS=0) or section (IS=1) of the 4-dim. ellipsoid
#  defined by the matrix A(4,4) and the centre position U(4).
#  The projection plane is specified by the IX and IY fields of
#  the PORT record.
#  color index (ICL), line style (ILS), character size (ZCH)
#-------------------------------------------------------------------
      implicit none
      
      integer*4 icl,ils,is,np,i
      parameter(np=100)
      real*8 a(4,4),b(2,2),u(4)
      real*8 d,d1,det,step,cf,z
      real*4 x(np),f(np)


      INCLUDE 'res_grf.inc'

      record/VIEWSET/ port

      call REDUCE42(a,b,port.ix,port.iy,is)

      call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)
      call pgwindow(port.wx1,port.wx2,port.wy1,port.wy2)
      call pgsci(icl)
      call pgsls(ils)


#      if(PORT.IX.EQ.1.AND.PORT.IY.EQ.2) then
#10     format('ELL: ',3(G12.5,1x))
#      write(*,10) B(1,1),B(1,2),B(2,2)
#      write(*,10) U(PORT.IX),U(PORT.IY)
#      endif

      cf=2.0d0*log(2.0d0)
      d1= sqrt(cf*b(2,2)/(b(2,2)*b(1,1)-b(1,2)**2))
      step=2.*d1/(np-1)
      do 1 i=1,np
        d=i-1
        z=-d1+d*step
        x(i)=z+u(port.ix)
        det=(b(1,2)*z)**2-b(2,2)*(b(1,1)*z**2-cf)
        f(i)=(-b(1,2)*z+sqrt(abs(det)))/b(2,2)+u(port.iy)
1     continue
      call pgline(np,x,f)
      do 2 i=1,np
        x(i)=-x(i)+2.*u(port.ix)
        f(i)=-f(i)+2.*u(port.iy)
2     continue
      call pgline(np,x,f)
      call pgsci(1)
      call pgsls(1)
      return
      end


#-----------------------------------------------------
      SUBROUTINE PLOTLINE(port,icl,ils,x,y,np)
#  plots line in physical coordinates given by vectors
#  X(NP) and Y(NP) to the PORT.
#  ICL .. collor,  ILS .. line style
#-----------------------------------------------------
      implicit none

      INCLUDE 'res_grf.inc'
      
      integer*4 icl,ils,np
      real*4 x(np),y(np)
      record/VIEWSET/ port

      call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)
      call pgwindow(port.wx1,port.wx2,port.wy1,port.wy2)
      call pgsci(icl)
      call pgsls(ils)

      call pgline(np,x,y)

      return
      end


#--------------------------------------------------------------------
      SUBROUTINE PlotCurve(port,xfx,fy,nf,ic,ip,iline)
#  Plot a curve on PORT.
#   PORT       ... plotting viewport
#   XFX        ... x-values
#   FY         ... y-values
#   NF         ... number of points
#   IC         ... color
#   IP         ... point style
#   ILINE      ... line style
# CALLS:
# CALLED BY: VIEWSCAN
#--------------------------------------------------------------------
      implicit none

      INCLUDE 'res_grf.inc'

      integer*4 nf
      record /VIEWSET/ port

      real*4 xfx(nf),fy(nf)
      real*8 centre,range
      integer*4 i,ic,ip,iline,j1,j2,jr

#//   define viewport:

      if (port.wx1.eq.port.wx2) then
        centre=(xfx(1)+xfx(nf))/2.
        i=1
        do 41 while ((fy(i).eq.0).and.(i.lt.nf))
41          i=i+1
        if(i.gt.1) i=i-1
        j1=i
        i=nf
        do 42 while ((fy(i).eq.0).and.(i.gt.1))
42         i=i-1
        if(i.lt.nf) i=i+1
        j2=i
        jr=j2-j1+5
        jr=max(jr,20)
        j1=max(j1-jr/2,1)
        j2=min(j2+jr/2,nf)
        range=abs(xfx(j2)-xfx(j1))
        port.wx1=centre-range/2.
        port.wx2=centre+range/2.
      endif
      if (port.wy1.eq.port.wy2) then
         port.wy1=+1e+20
         port.wy2=-1e+20
         do i=1,nf
           port.wy2=max(port.wy2,fy(i))
           port.wy1=min(port.wy1,fy(i))
         end do
         port.wy1=port.wy1*1.1
         port.wy2=port.wy2*1.1
      endif

      if (port.wy1.eq.port.wy2) return
      if (port.wx1.eq.port.wx2) return
      call CLRPORT(port)
      call PLOTFRAME(port,1,1,1.3,0)
      if(iline.eq.0) then
          call pgsls(1)
      else
          call pgsls(iline)
      endif
      call pgsci(ic)
      if (ip.ne.0) call pgpoint(nf,xfx,fy,ip)
      if (iline.ne.0) call pgline(nf,xfx,fy)
      call pgsci(1)
      call pgsls(1)

      end

#*************************  SETTING VIEWPORTS ********************************

#-------------------------------------------------------------
      SUBROUTINE CLRPORT(port)
#     clear viewport
#-------------------------------------------------------------
      implicit none 

      INCLUDE 'res_grf.inc'

      record/VIEWSET/ port

      call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)
      call pgwindow(0.,1.,0.,1.)
      call pgsci(0)
      call pgrect(0.,1.,0.,1.)
      call pgsci(1)
      end


#---------------------------------------------------------
      SUBROUTINE MK_PORTS(icom)
#   Defines viewport parameters for the four projections
#   of the ellipsoid defined by the matrix in C&N system:
#   ICOM=2,3   ...  ANESS(4,4) from the common /NESSIF/
#   ICOM=1     ...  A(4,4) from the common /MATRIX/
#---------------------------------------------------------
      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'
      
      integer*4 plane(2,4),i,icom
      character*17 ind(4)
      real*4 xmin(4),xmax(4),qm

      data ind/ 'Q\dX\u [\A\u-1\d] ', 'Q\dY\u [\A\u-1\d] ',
     * 'Q\dZ\u [\A\u-1\d] ', '   '/

      ind(4)=  '   \gDE '//cunit


#///  Specify the four planes for drawing projections

      plane(1,1)=1
      plane(2,1)=4

      plane(1,2)=1
      plane(2,2)=2

      plane(1,3)=2
      plane(2,3)=3

      plane(1,4)=2
      plane(2,4)=4


#///  Specify 4 viewports  :

      vport(4).dx1=0.1
      vport(4).dx2=0.47
      vport(4).dy1=0.33
      vport(4).dy2=0.58

      vport(3).dx1=vport(4).dx1+0.5
      vport(3).dx2=vport(4).dx2+0.5
      vport(3).dy1=vport(4).dy1
      vport(3).dy2=vport(4).dy2

      vport(2).dx1=vport(4).dx1+0.5
      vport(2).dx2=vport(4).dx2+0.5
      vport(2).dy1=vport(4).dy1+0.35
      vport(2).dy2=vport(4).dy2+0.35

      vport(1).dx1=vport(4).dx1
      vport(1).dx2=vport(4).dx2
      vport(1).dy1=vport(4).dy1+0.35
      vport(1).dy2=vport(4).dy2+0.35

#//// Set scales:


      if(icom.eq.3) then
         call GETMSCALE(xmin,xmax)
         do i=1,4
            vport(i).wx1=xmin(plane(1,i))
            vport(i).wx2=xmax(plane(1,i))
            vport(i).wy1=xmin(plane(2,i))
            vport(i).wy2=xmax(plane(2,i))
         end do
      endif


      do 10 i=1,4
       if(icom.eq.2) then
         call GETSCALE(aness,plane(1,i),plane(2,i),0,
     1   vport(i).wx1,vport(i).wx2,vport(i).wy1,vport(i).wy2)
       else if(icom.le.1) then
         call GETSCALE(atrax,plane(1,i),plane(2,i),0,
     1   vport(i).wx1,vport(i).wx2,vport(i).wy1,vport(i).wy2)
       endif
       vport(i).xtit=ind(plane(1,i))
       vport(i).ytit=ind(plane(2,i))
       vport(i).head= ' '
       vport(i).ix=plane(1,i)
       vport(i).iy=plane(2,i)
10    continue

      qm=vport(2).wx2
      qm=max(qm,vport(2).wy2)
#      QM=MAX(QM,VPORT(3).WY2)
      do 12 i=1,4
       if((plane(1,i).ne.4).and.(plane(1,i).ne.3)) then
         vport(i).wx2=qm
         vport(i).wx1=-qm
       endif
       if((plane(2,i).ne.4).and.(plane(2,i).ne.3)) then
         vport(i).wy2=qm
         vport(i).wy1=-qm
       endif
12    continue

      return
      end


#----------------------------------------------------------------------
      SUBROUTINE SETPORTCELL(port,nrow,ncol,n)
#   Set viewport dimensions and position for given layout and cell index
# NCOL ... number of coulmns
# NROW ... number of rows
# N .. cell index (from top left to bottom right)
#----------------------------------------------------------------------
      implicit none
      INCLUDE 'res_grf.inc'
      
      record /VIEWSET/ port
      integer*4 nrow,ncol,n
      real*4 aspect,height,width,m1,z
      integer*4 icol,irow
      real*4 mar(6)    ! margins: left,right,top,bottom, H-spacing, V-spacing in page fraction
      data mar /0.08,0.02,0.02,0.10,0.08,0.07/


      m1=mar(1)
      if(ncol.eq.1) m1=2.*mar(1)
      aspect=(1.*nrow)/ncol
      irow=int(0.999*n/ncol)+1
      icol=n-(irow-1)*ncol
      width=(1.-m1-mar(2)-mar(5)*(ncol-1))/ncol
      port.dx1=m1+(icol-1)*(width+mar(5))
      port.dx2=port.dx1+width
      z=1.
      if (nrow.le.2) z=1.5
      height=(1.-mar(3)-mar(4)-z*mar(6)*(nrow-1))/nrow
      height=min(height,aspect*width)
      port.dy2=1-mar(3)-(irow-1)*(height+mar(6)*z)
      port.dy1=port.dy2-height
      end



#*************************  DEFINE PLOT SCALES ********************************


#----------------------------------------------------------
      SUBROUTINE GETSCALE(a,ix,iy,is,x1,x2,y1,y2)
# Calculates physical scales of a graph for ellipsoid projections 
#   A       ... resolution matrix
#   IX,IY   ... defines projection plane
#   IS<>0   ... use ellipsoid section instead of projection
#   X1..Y2  ... output: plot scales
# CALLED BY: PlotResol,PlotResQE
#----------------------------------------------------------
      implicit none
      
      real*8 a(4,4),b(2,2),dd(2),cf,fact1,fact2
      integer*4 code(2),expon(2),ix,iy,is,i
      real*4 x1,x2,y1,y2

      call REDUCE42(a,b,ix,iy,is)

#///  Scales of the axes are calculated from the ellipsoid limits:

      cf=2.0*log(2.0d0)
      dd(1)= sqrt(abs(cf*b(2,2)/(b(2,2)*b(1,1)-b(1,2)**2)))
      dd(2)= sqrt(abs(cf*b(1,1)/(b(1,1)*b(2,2)-b(1,2)**2)))
#      do i=1,2
#        write(*,*) (B(i,j),j=1,2)
#      enddo
#      write(*,*) DD(1),DD(2),IX,IY

      do 10,i=1,2
        expon(i)=abs(log10(dd(i)))
        if (dd(i).lt.1) expon(i)=-expon(i)-1
        if(dd(i)/10.**expon(i).lt.3) then
          code(i)=1
        else if(dd(i)/10.**expon(i).lt.6) then
          code(i)=2
        else
          code(i)=5
        endif
10    continue

      fact1=code(1)*10.0**expon(1)
      fact2=code(2)*10.0**expon(2)
      x1=-3.0*fact1
      y1=-3.0*fact2
      x2=-x1
      y2=-y1
      if (y1.eq.y2) y2=y2+1e-15
      if (x1.eq.x2) x2=x2+1e-15

      return
      end

#--------------------------------------------------------------------
      SUBROUTINE GETMSCALE(xmin,xmax)
# Calculates physical limits of a graph for the Monte Carlo events.
# Works with current dataset.
#   XMIN(4), XMAX(4)  ... plot limits in C&N coordinates
# CALLED BY: VIEWSCAN,MK_PORTS
#--------------------------------------------------------------------
      implicit none

      INCLUDE 'const.inc'      
      INCLUDE 'restrax.inc'
      
      real*8 dd(4)      !,dd1(4)
      integer*4 code(4),expon(4),i,j,ncnt
      real*4 xmin(4),xmax(4)
      real*8 e1(4),fact(4),p

      call KSTACK_N(ncnt,mf_cur)
      if (ncnt.le.0) return

      do i=1,4
        dd(i)=0
#        dd1(j)=0
      end do
      do i=1,ncnt
         call GETQE(i,mf_cur,e1,p)
         do j=1,4
           dd(j)=max(abs(e1(j)),dd(j))

#           dd1(j)=min(e1(j),dd1(j))
         end do
      enddo

      do 10,i=1,4
        expon(i)=abs(log10(dd(i)))
        if (dd(i).lt.1) expon(i)=-expon(i)-1
        if(dd(i)/10.**expon(i).lt.2) then
          code(i)=2
        else if(dd(i)/10.**expon(i).lt.4) then
          code(i)=4
        else if(dd(i)/10.**expon(i).lt.6) then
          code(i)=6
        else if(dd(i)/10.**expon(i).lt.8) then
          code(i)=8
        else
          code(i)=10
        endif
        fact(i)=code(i)*10.0**expon(i)
        xmin(i)=-fact(i)
        xmax(i)=+fact(i)
        if (xmax(i).eq.xmin(i)) xmax(i)=xmax(i)+1e-15
10    continue

      return
      end

#----------------------------------------------------------------------
      SUBROUTINE GETDATASCALE(n,centre,range,xstep,x0,indx)
# Get optimum scale for data X-axis of N-th dataset
#   N        ... Index of dataset
#   CENTRE  ... X-axis centre
#   RANGE   ... X-axis range
#   XSTEP   ... mean step (used to construct phys. values from XHIST etc..)
#   X0      ... scan centre
#   INDX    ... INDX = 1..4 for QH,QK,QL,E respectively
#               denotes the most rapidly varying coordinate
# CALLED BY: PLOTCELL
#----------------------------------------------------------------------
      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'rescal.inc'
      INCLUDE 'restrax.inc'

      integer*4 nf,np,ib,n,ibh,indx
      real*4 fx(nhi*4),fy(nhi*4)
      real*8 centre,range,xstep,x0
      integer*4 i,j1,j2,jr,k
      
      np=npt(n)-npt(n-1)    ! number of points incurrent data set
      if(np.gt.nhi*4) np=nhi*4
      
#// get QHKL,E component which varies most rapidly            
#// if spectrum is loaded, get x-axis points from it
      if (np.gt.0) then
        ib=npt(n-1)+1         ! base index for the incurrent data set
        xstep=0.
        do i=1,5
          if(abs(dqe0(i,n)).gt.xstep) then
            k=i 
          endif
          xstep=abs(dqe0(i,n))
        enddo
#        XSTEP=DQE0(K,N)
#        X0=QE0(K,N)
        xstep=mf_par(i_dqh+k-1,n)
        x0=mf_par(i_qh+k-1,n)
        if (k.eq.5) x0=0
        
        centre=(spx(np+ib-1)+spx(ib))*xstep/2. +x0
        range=abs(spx(np+ib-1)-spx(ib))*2*xstep
#// otherwise use step size and histogram content
      else 
        nf=nhist(n)-nhist(n-1)    ! number of points in current histogram
        if (nf.le.0) goto 99      ! nothing to plot -> exit
        ibh=nhist(n-1)+1          ! base index for current histogram
        do i=1,nf
          fx(i)=xhist(i+ibh-1)
          fy(i)=rhist(i+ibh-1)
        end do      
        xstep=0.
        do i=1,4
          if(abs(mf_par(i_dqh+i-1,n)).gt.xstep) then
            k=i 
          endif
          xstep=abs(mf_par(i_dqh+i-1,n))
        enddo
        xstep=mf_par(i_dqh+k-1,n)     
        x0=mf_par(i_qh+k-1,n)      
        if (k.eq.5) x0=0
        do i=1,nf
          fx(i)=fx(i)*xstep + x0
        enddo        
        centre=(fx(1)+fx(nf))/2.
        i=1
        do 41 while ((fy(i).eq.0).and.(i.lt.nf))
41          i=i+1
        if(i.gt.1) i=i-1
        j1=i
        i=nf
        do 42 while ((fy(i).eq.0).and.(i.gt.1))
42          i=i-1
        if(i.lt.nf) i=i+1
        j2=i
        jr=j2-j1+5      ! number of non-zero points + 5
        jr=max(jr,20)   ! show at least 20 points
        j1=max(j1-jr/2,1)
        j2=min(j2+jr/2,nf)
        range=abs(fx(j2)-fx(j1))
      endif
      
      indx=k
      return

99    centre=0
      range=1
      xstep=0.1
      x0=-0.5
      indx=0
      end


#---------------------------------------------------------      
      SUBROUTINE FCONE_RANGE(xmin,xmax,ymin,ymax,eqscale)
# define automatically viewport scales for flat-cone scan
# if EQSCALE, then make aspect ratio=1
#---------------------------------------------------------            
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'rescal.inc'     
      INCLUDE 'restrax.inc
      INCLUDE 'res_grf.inc'
      
      real*8 xmin,xmax,ymin,ymax
      logical*4 eqscale
      real*4 qi(3) 
      real*8 qf1(3),qf2(3)
      real*8 x0,x1,x2,y0,y1,y2,xmi,xma,ymi,yma,da3,axn,bxn
      integer*4 id,j,np
      real*8 sc(17),dum,dx,dy
      real*8 ROUNDSCALE
      real*8 ax(3),bx(3),QxQ 
      equivalence (ax(1),mf_par(i_ax,1))
      equivalence (bx(1),mf_par(i_bx,1))
      
10    format(a,6(1x,g10.4))

      do j=1,17
        sc(j)=1+j*0.5  ! possible rounded scale limits
      enddo
      
      np=nhist(1)
      da3=(np-1)/2*mf_par(i_da3,1)
      call QNORM(ax,axn,dum)
      call QNORM(bx,bxn,dum)
      
      xmi=1d10
      xma=-1d10       
      ymi=1d10       
      yma=-1d10       
      do id=1,mf_max
        do j=1,3
          qi(j)=mf_par(i_qh+j-1,id)
        enddo  
        call ROTA3(qi,da3,qf1)
        call ROTA3(qi,-da3,qf2)
        x0=QxQ(mf_par(i_qh,id),ax)/axn
        x1=QxQ(qf1,ax)/axn
        x2=QxQ(qf2,ax)/axn
        y0=QxQ(mf_par(i_qh,id),bx)/bxn
        y1=QxQ(qf1,bx)/bxn
        y2=QxQ(qf2,bx)/bxn
        xmi=min(xmi,x0,x1,x2)
        ymi=min(ymi,y0,y1,y2)
        xma=max(xma,x0,x1,x2)
        yma=max(yma,y0,y1,y2)
      enddo
      dx=xma-xmi
      dy=yma-ymi
      
      xmin=ROUNDSCALE(xmi-dx/20.,-1,sc,17)
      xmax=ROUNDSCALE(xma+dx/20.,1,sc,17)
      ymin=ROUNDSCALE(ymi-dy/20.,-1,sc,17)
      ymax=ROUNDSCALE(yma+dy/20.,1,sc,17)
      if (xmin.gt.0.and.xmin.lt.0.1*dx) xmin=0
      if (ymin.gt.0.and.ymin.lt.0.1*dy) ymin=0
      if (xmax.lt.0.and.xmax.gt.-0.1*dx) xmax=0
      if (ymax.lt.0.and.ymax.gt.-0.1*dy) ymax=0
      if (abs(xmin).lt.0.05*dx) xmin=0
      if (abs(ymin).lt.0.05*dx) ymin=0
      if (abs(xmax).lt.0.05*dx) xmax=0
      if (abs(ymax).lt.0.05*dx) ymax=0
      if (eqscale) then
        xma=(xmax-xmin)
        yma=(ymax-ymin)
        if (xma.gt.0) then        
          da3=yma/xma
          if (da3.gt.1) then
             xmax=xmin+yma          
          else
             ymax=ymin+xma
          endif   
        endif       
      endif
#      write(*,10) 'rounded: ',XMIN,XMAX,YMIN,YMAX
      
      end
      
#*************************  FORMAT LABELS ********************************

#--------------------------------------------------------------------
      SUBROUTINE GetAxTitle(v,s)
#// format the a label to write dir. V in symbolic format, if possible 
#--------------------------------------------------------------------
      implicit none
      character*3 chind(3)
      character*30 s
      real*8 v(3),z
      integer*4 j,i
101   format( '[',f6.2, ' ',f6.2, ' ',f6.2, '] / rel.')
104   format( '[ ',a3, '  ',a3, '  ',a3, ' ] '

      j=0
      z=0
      do i=1,3
        if (v(i).ne.0.and.v(i).ne.z) then
          j=j+1
          if(z.eq.0) z=v(i)
        endif  
      end do
      if (j.eq.1) then
        do i=1,3
          if(v(i).ne.0) then 
            chind(i)= '\gc'
          else
            chind(i)= ' 0 '
          endif
        end do
        write(s,104) (chind(i),i=1,3)        
      else
        write(s,101) v
      endif      
    
      end


#----------------------------------------      
      SUBROUTINE FORMAT_HKL(v,s,ns)
#----------------------------------------            
      implicit none
      character*(*) s
      character*128 s1,s2
      real*8 v(3)
      integer*4 i,is,il,ns
      logical*4 log
10    format( '[',3(1x,i3), ' ]')      
20    format( '[',3(1x,f6.1), ' ]')      
      
      
      log=.true.
      do i=1,3
        log=(log.and.(abs(v(i)-nint(v(i))).lt.1e-2))
      enddo
      if (log) then
#        write(*,10) (NINT(V(I)),I=1,3)
        write(s1,10)  (nint(v(i)),i=1,3)
      else
#        WRITE(*,20) V
        write(s1,20) v
      endif    
      is=1
      call BOUNDS(s1,is,il)
      i=1
      do while (i.gt.0)  ! remove extra spaces
        i=index(s1(is:is+il-1), '  ')
        if (i.gt.0) then
          s2=s1(is:i)//s1(i+2:is+il-1)          
          s1=s2  
          call BOUNDS(s1,is,il)        
        endif
      enddo  
      if (il.gt.ns) il=ns
      s=s1(is:is+il-1)
      end
      
#----------------------------------------------------------------------
      SUBROUTINE LEGFIT(port,size,offset)
#  Show fit parameters in the area defined by PORT
#  PORT     ... plotting viewport
#  SIZE     ... character size
#  OFFSET   ... left margin (in page width units)
#----------------------------------------------------------------------
      implicit none

      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'
      INCLUDE 'restrax.inc'
      INCLUDE 'res_grf.inc'
      INCLUDE 'exciimp.inc'
            
      record /MODEL/ rm
      
      real*8 eps
      real*4 size,offset
      parameter (eps=1e-20) 

      record /VIEWSET/ port

      character*60 leg1
      integer*4 i,j,k
      real*4 rx

101   format( '\gx\u2\d   =',g11.4)
102   format(a, '=',g11.4)

      call getmodel(rm)

      call pgsch(size)

      call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)

      rx=0.04
      if (jfit.gt.0) then
         write(leg1,101) chisqr
         call pgmtext( 'T',-offset,rx,0.0,leg1)
      endif
      j=0
      do i=1,nfpar
         k=index(rm.parname(i), ' ')
         if (k.le.0) k=10
         if (rm.fixparam(i).gt.0) then
           j=j+1
           write(leg1,102) rm.parname(i)(1:k), fpar(i)
           call pgmtext( 'T',-offset-1.5*j,rx,0.0,leg1)
         endif  
      end do
      call pgsch(1.0)
      call pgsci(1)

      end
      
#----------------------------------------------------------------------
      SUBROUTINE LEGFILE(port)
# Show filenames below the lower-left corner of the PORT
# CALLED BY: PAGE2
#----------------------------------------------------------------------
      implicit none

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

      record /VIEWSET/ port

      character*60 leg1,leg2,leg3

106   format( 'DAT: ',a20)
110   format( 'CFG: ',a20)
111   format( 'RES: ',a20)

      call pgsch(0.8)

      call pgvport(port.dx1,port.dx2,0,port.dy2)
      write(leg1,110) cfgname
      if (datname.ne. ' ') then 
         write(leg2,106) datname
      else
         write(leg2,106) rescal_name
      endif 
      write(leg3,111) resname
      call pgmtext( 'B',-0.5,0.,0.,leg1(1:26))
      call pgmtext( 'B',-2.0,0.,0.,leg2(1:26))
      call pgmtext( 'B',-3.5,0.,0.,leg3(1:26))
      call pgvport(port.dx1,port.dx2,port.dy1,port.dy2)
      call pgsch(1.0)
      call pgsci(1)

      end
      
#*************************  SAVE RESULTS ********************************

#------------------------------------------------------------------
      SUBROUTINE SAVEELL(fname,port,a,u,is)
#  plots projection (IS=0) or section (IS=1) of the 4-dim. ellipsoid
#  defined by the matrix A(4,4) and the centre position U(4).
#  The projection plane is specified by the IX and IY fields of
#  the PORT record.
#  collor index (ICL), line style (ILS) ,character size (ZCH)
#-------------------------------------------------------------------
      implicit none
      
      integer*4 is,np,np2,ll,i_io,j,i
      parameter(np=128,np2=2*np)
      real*8 a(4,4),b(2,2),u(4)
      real*8 d,d1,det,step,cf,z
      real*4 x(np2),f(np2)
      character*(*) fname

      INCLUDE 'res_grf.inc'

      record/VIEWSET/ port

      call REDUCE42(a,b,port.ix,port.iy,is)

      cf=2.0d0*log(2.0d0)
      d1= sqrt(cf*b(2,2)/(b(2,2)*b(1,1)-b(1,2)**2))
      step=2.*d1/(np-1)
      do 1 i=1,np
        d=i-1
        z=-d1+d*step
        x(i)=z+u(port.ix)
        det=(b(1,2)*z)**2-b(2,2)*(b(1,1)*z**2-cf)
        f(i)=(-b(1,2)*z+sqrt(abs(det)))/b(2,2)+u(port.iy)
1     continue
      do 2 i=1,np
        x(i+np)=-x(i)+2.*u(port.ix)
        f(i+np)=-f(i)+2.*u(port.iy)
2     continue
      
      ll=len(fname)
      i_io=25
      close(i_io)
       open(unit=i_io,file=fname(1:ll),err=999,status= 'Unknown')
      write(i_io,*) port.xtit// ' '//port.ytit            
9     format(2(g12.5,2x))      
      do j=1,np2
        write(i_io,9) x(j),f(j)
      enddo
      close(i_io)

999      return
      end



      
#----------------------------------------------------------------------
      SUBROUTINE GRLIST
# Plots the list of RESTRAX parameters on PAGE 1, below the R(Q,E) projections
#----------------------------------------------------------------------
      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'inout.inc'

      integer*4 j
      character*120 liststr(15)

1     format(2(a4, ' = ',f10.5,7x))
2     format(3(a4, ' = ',f10.2,7x))
3     format(3(a4, ' = ',f10.0,7x))
4     format((a4, ' = ',f10.5,7x, a4, ' = ',f10.0,7x))
5     format(4(a4, ' = ',f10.2,7x))
6     format(3(a4, ' = ',f10.4,7x))
7     format(4(a4, ' = ',f10.4,7x))
8     format(2(a4, ' = ',f10.2,7x))



       write(liststr(1),1) (res_nam(j),res_dat(j),j=1,2)

       write(liststr(2),2) (res_nam(j),res_dat(j),j=3,5)

       write(liststr(3),3) (res_nam(j),res_dat(j),j=6,8)

       write(liststr(4),4) (res_nam(j),res_dat(j),j=9,10)

       write(liststr(5),5) (res_nam(j),res_dat(j),j=11,14)
       write(liststr(6),5) (res_nam(j),res_dat(j),j=15,18)

       write(liststr(7),6) (res_nam(j),res_dat(j),j=19,21)
       write(liststr(8),6) (res_nam(j),res_dat(j),j=22,24)
       write(liststr(9),6) (res_nam(j),res_dat(j),j=25,27)
       write(liststr(10),6) (res_nam(j),res_dat(j),j=28,30)

       write(liststr(11),7) (res_nam(j),res_dat(j),j=31,34)
       write(liststr(12),7) (res_nam(j),res_dat(j),j=35,38)
       write(liststr(13),7) (res_nam(j),res_dat(j),j=39,42)
       write(liststr(14),7) (res_nam(j),res_dat(j),j=43,46)
       write(liststr(15),8) (res_nam(j),res_dat(j),j=47,48)

      call pgsci(1)
      call pgsch(0.75)

      call pgvport(0.1,0.97,0.22,0.25)
      call pgwindow(0.0,1.0,0.0,1.0)
      call pgbox( 'C',0.0,0, ' ',0.0,0)
      call pgptext(0.0,0.5,0.0,0.0,
      'Instrument Configuration: '//cfgname)

      call pgsch(0.6)
      call pgvport(0.1,0.97,0.02,0.22)
      call pgwindow(0.0,1.0,0.0,1.0)
      call pgbox( 'B',0.0,0, ' ',0.0,0)

      do 55 j=1,15
       call pgptext(0.0,(1.05-j/15.),0.0,0.0,liststr(j))
55    continue
      return
      end