source: PSPA/madxPSPA/src/plot.f90 @ 430

Last change on this file since 430 was 430, checked in by touze, 11 years ago

import madx-5.01.00

File size: 92.6 KB
Line 
1subroutine pecat1(rb, ra, rd)
2
3  implicit none
4
5  !----------------------------------------------------------------------*
6  ! purpose:
7  !   concatenate two transport maps
8  ! input:
9  !   rb(6,6) second map in beam line order.
10  !   ra(6,6) first map in beam line order.
11  ! output:
12  !   rd(6,6) result map.
13  !----------------------------------------------------------------------*
14
15  !--- type definition of the routine arguments
16
17  double precision rb(6,6), ra(6,6), rd(6,6)
18
19  !--- type definitions of local variables
20
21  integer k, j
22
23  !--- Routine body
24
25  do k = 1, 6
26     do j = 1, 6
27        rd(j,k) = rb(j,1) * ra(1,k) + rb(j,2) * ra(2,k)               &
28             + rb(j,3) * ra(3,k) + rb(j,4) * ra(4,k)                           &
29             + rb(j,5) * ra(5,k) + rb(j,6) * ra(6,k)
30     enddo
31  enddo
32end subroutine pecat1
33!***********************************************************************
34
35subroutine pecurv (ncc, spname, annh, usex, sych, ippar,          &
36     np, xval, yval, window, actwin, ierr)
37
38  use plotfi
39  implicit none
40
41  !----------------------------------------------------------------------*
42  ! Purpose:
43  !   Plot one curve
44  ! Input:
45  !   NCC      (integer)  current curve count (1,2, etc.)
46  !   SPNAME    (char)     curve annotation string
47  !   ANNH     (real)     character height
48  !   USEX     (real)     user character height expansion
49  !   SYCH     (real)     symbol character height
50  !   IPPAR    (integer)  array containing the plot parameters
51  !   NP       (integer)  no. of points to plot
52  !   XVAL     (real)     x values
53  !   YVAL     (real)     y values
54  !   WINDOW   (real)     array containing the window to use
55  !   ACTWIN   (real)     active (inside frame) window
56  ! Output:
57  !   IERR     (integer)  0 if OK, else GXPLOT error
58  !
59  ! calls the routines peiact and pegacn in this file.
60  ! calls the routines gxsave, gxswnd, gxpnbl, jsln, jsplci,  gxpl,
61  !                    gxplt1, gvpl, jsmk, gxpmsw, jspmci, jschh, jstxal,
62  !                    gxqrvp, jschxp, gxtx, jstxci, gxrest
63  !                    defined in file gxx11.F.
64  !
65  ! it is called by the routine peplot in this file.
66  !----------------------------------------------------------------------*
67
68  !--- type definition of the routine arguments
69
70  integer  ncc, np, ierr, ippar(*)
71  character spname*(*)
72  real annh, sych, xval(*), yval(*), window(*)
73  real actwin(*), usex
74
75  !--- type definitions of local variables
76
77  real xpos, ypos, xwpos, ywpos, xf, wsclx, wscly
78  real xaux, yaux
79  real rsave(20), act(4), xpl(2), ypl(2), winnor(4)
80  real xreal(maxpnt), yreal(maxpnt)
81  real xmd
82  integer isymb, icolr, ilb, ispli, kft, klt, kact,                 &
83       kf, kl, npt, j, iecub
84  integer   istyl, ipbar, isave(20)
85  character symloc*1
86
87  double precision get_value !hbu
88  logical zero_suppr !hbu
89  logical marker_plot !hg
90
91  !--- Initialisation of local variables
92
93  data winnor /0., 1., 0., 1./
94  save act
95
96  iecub = 0
97
98  zero_suppr=get_value('plot ','zero_suppr ').ne.0 !hbu
99  marker_plot = get_value('plot ','marker_plot ').ne.0 !hg
100  !--- Output initialisation
101
102  ierr = 0
103
104  !--- Routine body
105
106  !--- save GKS settings
107
108  call gxsave (isave, rsave, ierr)
109  call gxswnd (window)
110  wsclx = 1. / (window(2) - window(1))
111  wscly = 1. / (window(4) - window(3))
112  xmd = 1.e-8 * (window(2) - window(1))**2
113  if (ncc .eq. 1)  then
114
115     !--- first curve in frame - reset label position array, get act.
116     !    window in NDC
117
118     act(1) = (actwin(1) - window(1)) * wsclx
119     act(2) = (actwin(2) - window(1)) * wsclx
120     act(3) = (actwin(3) - window(3)) * wscly
121     act(4) = (actwin(4) - window(3)) * wscly
122  endif
123  istyl = ippar(1)
124  ipbar = ippar(3)
125  isymb = ippar(4)
126  icolr = ippar(5)
127  if(icolr .eq. 100)  icolr = mod(ncc-1,4) + 1
128  icolr = max(1, min(icolr, 7))
129  ilb = -1
130  if (istyl .ne. 0)  then
131
132     !--- polyline requested
133
134     if (np .lt. 2) goto 999
135     if(istyl .eq. 100)  istyl = mod(ncc-1,4) + 1
136     ispli = ippar(2)
137
138     !--- get first and last blank in annotation
139
140     call gxpnbl (spname, kft, klt)
141     if (kft .ne. 0) then
142
143        !--- annotation exists
144
145        ilb    = 0
146     endif
147
148     !--- set line style
149
150     call jsln (max (1, min (4, istyl)))
151
152     !--- set line colour
153
154     call jsplci(icolr)
155     kact = 1
15610   continue
157
158     !--- get first and last point inside
159
160     call peiact(kact, np, xval, yval, actwin, kf, kl)
161
162     !--- quit if no points inside
163
164     if (kf .eq. 0) goto 40
165     kf = max(1, kf - 1)
166     kl = min(np, kl + 1)
167     npt      = 1
168     xreal(1) = xval(kf)
169     yreal(1) = yval(kf)
170     do j = kf + 1, kl
171
172        !--- avoid identical points
173
174        if ( (marker_plot.and..not.zero_suppr)   &
175             .or. (xreal(npt) - xval(j))**2 +      &
176             (yreal(npt) - yval(j))**2 .gt. xmd   .and.                &
177             ( yval(j).ne. 0 .or. .not.zero_suppr) ) then ! hbu optionally skip 0 points
178           npt        = npt + 1
179           xreal(npt) = xval(j)
180           yreal(npt) = yval(j)
181        endif
182        if ((j .eq. kl .and. npt .ge. 2) .or. npt .eq. maxpnt) then
183
184           !--- plot - get first curve annotation position
185
186           if (ilb .eq. 0)                                             &
187                call pegacn(ncc, window, act, xreal, yreal, npt, usex,            &
188                xwpos, xpos, ypos, ilb)
189           if (interf .eq. 0 .or. npt .eq. 2 .or. ispli .ne. 0) then
190
191              !--- no spline
192
193              call gxpl (npt, xreal, yreal, actwin)
194              if (ilb .gt. 0) then
195
196                 !--- get y pos. on curve for label
197
198                 ywpos = yreal(ilb - 1) + (yreal(ilb) - yreal(ilb-1))    &
199                      * (xwpos  - xreal(ilb - 1))                                       &
200                      / (xreal(ilb) - xreal(ilb - 1))
201                 ilb = -2
202              endif
203           else
204
205              !--- spline
206
207              call gxplt1 (npt, xreal, yreal, actwin)
208              if (ilb .gt. 0) ilb   = -2
209
210              !--- get y pos. on curve for label
211
212           endif
213           xreal(1) = xreal(npt)
214           yreal(1) = yreal(npt)
215           npt = 1
216        endif
217     enddo
218     if (kl .lt. np)  then
219        kact = kl + 1
220        goto 10
221     endif
222  else
223
224     !--- no polyline
225
226     if (np .eq. 0) goto 999
227  endif
228
229  !--- plot symbols or bars if requested
230
231  if (ipbar .ne. 0)  then
232     call jsln (1)
233
234     !--- set line colour
235
236     call jsplci(icolr)
237     do j = 1, np
238        xpl(1) = xval(j)
239        xpl(2) = xval(j)
240        ypl(1) = yval(j)
241        ypl(2) = actwin(3)
242        call gvpl (2, xpl, ypl)
243     enddo
244  endif
24540 continue
246  if (isymb .ne. 0)  then
247     if (isymb .le. 5)  then
248        call jsmk (isymb)
249        call gxpmsw (np, xval, yval, actwin)
250     elseif (isymb .eq. 100)  then
251        if (istyl .ne. 0)  then
252
253           !--- use current curve count
254
255           write (symloc, '(I1)')  mod (ncc, 10)
256        endif
257
258        !--- set marker colour
259
260        call jspmci(icolr)
261
262        !--- plot one character symbol
263        !    switch to normalized window
264
265        call gxswnd (winnor)
266
267        !--- set character height
268
269        call jschh (sych)
270
271        !--- text alignment
272
273        call jstxal (2, 3)
274
275        !--- text expansion factor - mind distorted viewports
276
277        call gxqrvp (xf)
278        call jschxp (xf)
279        do j = 1, np
280           if (isymb .eq. 100 .and. istyl .eq. 0)  then
281
282              !--- use current point number
283
284              write (symloc, '(I1)')  mod (j, 10)
285           endif
286           xaux = wsclx * (xval(j) - window(1))
287           yaux = wscly * (yval(j) - window(3))
288           if (xaux .gt. act(1) .and. xaux .lt. act(2)                 &
289                .and. yaux .gt. act(3) .and. yaux .lt. act(4))                    &
290                call gxtx (xaux, yaux, symloc)
291        enddo
292     endif
293  endif
294  if (ilb .eq. -2)  then
295
296     !--- plot annotation
297     !    switch to normalized window
298
299     call gxswnd (winnor)
300
301     !--- set character height
302
303     call jschh (annh)
304     !--- text alignment
305
306     call jstxal (2, 5)
307
308     !--- text expansion factor - mind distorted viewports
309
310     call gxqrvp (xf)
311     call jschxp (xf)
312
313     !--- set marker colour
314
315     call jstxci(icolr)
316
317     !--- plot annotation string
318
319     call gxtx (xpos, ypos, spname(kft:klt))
320
321     !--- connect to curve
322
323     xpl(1) = xpos
324     xpl(2) = xpos
325     ypl(1) = ypos
326     ypl(2) = (ywpos - window(3)) * wscly
327     if (ypl(2) .gt. ypl(1))  ypl(1) = ypl(1) + .02
328
329     !--- set dotted line
330
331     call jsln (3)
332
333     !--- set line colour
334
335     call jsplci(icolr)
336
337     !--- plot line
338
339     call gxpl (2, xpl, ypl, act)
340  endif
341
342  !--- restore
343
344  call gxrest (isave, rsave)
345999 continue
346end subroutine pecurv
347
348!***********************************************************************
349
350subroutine pefill(ierr)
351
352  use plotfi
353  use plot_bfi
354  use plot_cfi
355  use plot_mathfi
356  implicit none
357
358  !----------------------------------------------------------------------*
359  ! Purpose:
360  !   fill plot arrays with coordinate values, set up machine plot
361  !
362  ! Output:  ierr  (int)     =0: OK, >0: error
363  !
364  ! calls the functions double_from_table_row, advance_to_pos
365  ! table_length and restart_sequ defined in file madxn.c.
366  ! calls the function lastnb defined in file util.F.
367  ! calls the routines peintp in this file.
368  !
369  ! it is called by the routine exec_plot in file madxn.c after pesopt.
370  !----------------------------------------------------------------------*
371
372  integer ierr
373  integer i,j,k,l,new1,new2,currtyp,crow,pltyp,pos_flag
374  integer ilist(mtype), p(mxplot), proc_n(2, mxcurv)
375  integer nint
376  double precision fact, d_val, d_val1, get_value
377  double precision currpos, currleng, currtilt, currk1l, currk1sl,  &
378       currk2l, currk2sl, currk3l, currk3sl
379  real tval, step, mystep
380  logical machp, rselect, marker_plot, range_plot
381  character*120 msg
382
383  !--- definitions of function primitives
384
385  integer double_from_table_row, restart_sequ,advance_to_pos
386  integer lastnb, table_length
387  integer advance_node, get_option
388  double precision node_value
389
390  !--- codes see in peschm
391
392  data ilist /                                                      &
393       0, 21, 1, 0, 2, 10, 12, 8, 0, 9,                                  &
394       6, 0, 0, 0, 0,  0,  4, 4, 4, 0,                                   &
395       0, 0, 0, 0, 0,  0, 14, 0, 0, 0,                                   &
396       20 * 0 /
397
398  ! Initialize marker_plot logical
399
400  marker_plot=get_value('plot ','marker_plot ').ne.zero
401  range_plot=get_value('plot ','range_plot ').ne.zero
402
403  !--- Output initialisation
404
405  ierr = 0
406
407  !--- Initialisation of variables in common peaddi
408
409  nelmach = 0
410  do j = 1, mxcurv
411     nqval(j) = 0
412  enddo
413  do j = 1 , maxseql
414     ieltyp(j) = 0
415  enddo
416
417  !--- Initialisation of variables in common peaddr
418
419  do i = 1 , maxseql
420     estart(i) = 0.0
421     eend(i) = 0.0
422     do j = 1 , mxcurv
423        qhval(i,j) = 0.0
424        qvval(i,j) = 0.0
425     enddo
426  enddo
427
428  !--- Initialisation of variables in common peotcl
429
430  dpp_flag = .false.
431
432  !--- Initialisation of local variables
433  currtilt =0
434  currk1l = 0
435  currk1sl =0
436  currk2l = 0
437  currk2sl =0
438  currk3l = 0
439  currk3sl =0
440
441  new1 = 0
442  new2 = 0
443  currtyp = 0
444  nint = 0
445  step = 0
446  crow = 0
447  pltyp = 0
448  pos_flag = 3
449  machp = itbv .ne. 0
450  do i=1,mxplot
451     p(i)=0
452  enddo
453
454  !--- save process flags, proc_flag may be modified in peintp
455
456  do i = 1, mxcurv
457     do j = 1, 2
458        proc_n(j,i) = proc_flag(j,i)
459     enddo
460  enddo
461
462  !--- Routine body
463
464  !--- No interpolation if centre option is set
465
466  if (get_option('centre ') .ne. 0) then
467     if (interf .eq. 1) then
468        write (msg, 910)
469        call aawarn('PLOT: ',msg)
470     endif
471     interf = 0
472     pos_flag = 2
473  endif
474  k = double_from_table_row(tabname, horname, 1, d_val)
475  if (k .lt. 0)  then
476     if (k .eq. -1)  then
477        print *, 'Warning: table ', tabname, ' not found'
478     elseif (k .eq. -2)  then
479        print *, 'Warning: hor. variable ', horname,                  &
480             ' not in table ', tabname
481     else
482        print *, 'Warning: table ', tabname, ' is empty'
483     endif
484     ierr = 1
485     goto 999
486  endif
487  if (horname .eq. 'dpp') dpp_flag = .true.
488  !rdemaria 1/6/2005: if the horizzontal name is deltap don't plot
489  !momentum offset.
490  if(horname .eq. 'deltap') dpp_flag=.true.
491  rselect = machp .and. hrange(2) .gt. hrange(1)
492  do l = 1, nivvar
493     k = double_from_table_row(tabname, sname(l), 1, d_val)
494     if (k .lt. 0)  then
495        print *, 'Warning: vertical variable: ',                      &
496             sname(l)(:lastnb(sname(l))), ' not in table ',                    &
497             tabname
498        ierr = 1
499        goto 999
500     endif
501     if (sname(l) .eq. 'dpp') dpp_flag = .true.
502  enddo
503  if (rselect .or. range_plot)  then
504
505     !-------adjust element range to horizontal range
506
507     new1 = nrrang(1)
508     new2 = nrrang(2)
509     crow = nrrang(1)
510     do j = nrrang(1), nrrang(2)
511        k = double_from_table_row(tabname, horname, j, d_val)
512        tval = d_val
513        if (tval .lt. hrange(1)) new1 = j
514        if (tval .lt. hrange(2)) new2 = j
515     enddo
516     nrrang(1) = new1
517     if (nrrang(2) .gt. new2+2) nrrang(2) = new2 + 2
518  endif
519  if (itbv .eq. 0 .and. .not. range_plot)  then
520     nrrang(1) = 1
521     nrrang(2) = table_length(tabname)
522  endif
523  if (nrrang(1) .eq. 0) nrrang(1) = 1
524
525  !--- get interpolation interval size
526  if (machp)  then
527     k = double_from_table_row(tabname, horname, nrrang(1), d_val)
528     k = double_from_table_row(tabname, horname, nrrang(2), d_val1)
529     step = (d_val1 - d_val) / (maxpnt / 2)
530  endif
531  fact = pos_flag - 1
532  j = restart_sequ()
533
534  do j = nrrang(1), nrrang(2)
535     crow = j
536     if (itbv .eq. 1 .and. advance_to_pos(tabname, j) .eq. 0)        &
537          goto 10
538     k = double_from_table_row(tabname, horname, j, currpos)
539
540     if (itbv .eq. 1)  then
541        currtyp = node_value('mad8_type ')
542        if(currtyp.eq.39) currtyp=15
543        if(currtyp.eq.38) currtyp=24
544        if (currtyp .le. mtype) pltyp = ilist(currtyp)
545
546        !--- get element parameters & build up plytp (to be used by the routine peschm)
547
548        currleng = node_value('l ')
549        if (currleng .gt. 0.d0 .and. currtyp .gt. 1                   &
550             .and. currtyp .lt. 8) then
551           currtilt = node_value('tilt ')
552           k = double_from_table_row(tabname, 'k1l ' , j, currk1l)
553           currk1l = currk1l/currleng
554           k = double_from_table_row(tabname, 'k1sl ', j, currk1sl)
555           currk1sl = currk1sl/currleng
556           k = double_from_table_row(tabname, 'k2l ' , j, currk2l)
557           currk2l = currk2l/currleng
558           k = double_from_table_row(tabname, 'k2sl ', j, currk2sl)
559           currk2sl = currk2sl/currleng
560           k = double_from_table_row(tabname, 'k3l ' , j, currk3l)
561           currk3l = currk3l/currleng
562           k = double_from_table_row(tabname, 'k3sl ', j, currk3sl)
563           currk3sl = currk3sl/currleng
564        endif
565
566        !--- sbend, rbend
567
568        if (mod(pltyp,20) .eq. 1 .and. currtilt .ne. zero)            &
569             pltyp = pltyp + 6
570
571        !--- quad
572
573        if (pltyp .eq. 2 .and. min(currk1l, currk1sl) .lt. zero)      &
574             pltyp = 3
575
576        !--- sext
577
578        if (pltyp .eq. 10 .and. min(currk2l, currk2sl) .lt. zero)     &
579             pltyp = 11
580
581        !--- oct
582
583        if (pltyp .eq. 12 .and. min(currk3l, currk3sl) .lt. zero)     &
584             pltyp = 13
585
586        !--- Compute the start & end position
587
588        if (machp)  then
589           nelmach = nelmach + 1
590           estart(nelmach) = currpos - fact * half * currleng
591           eend(nelmach)   = estart(nelmach) + currleng
592           ieltyp(nelmach) = pltyp
593        endif
594
595        !--- Interpolation if required
596
597        if (machp .and. j .gt. nrrang(1) .and. currleng .gt. zero     &
598             .and. interf .gt. 0 .and. step .gt. zero .and. .not. ptc_flag)    &
599             then
600
601           nint = currleng / step
602           if (nint .lt. 2) nint = 2
603           call peintp(crow, nint, proc_n, currleng, ierr)
604
605           if (ierr .eq. 1)  then
606              ierr = 0
607              print *, 'Warning: plot buffer full, plot truncated'
608              goto 100
609           elseif (ierr .ne. 0)  then
610              goto 999
611           endif
612        endif
613     endif
614
615     mystep=0.1d0 * step
616     do l = 1, nivvar
617        if (nqval(l) .eq. maxseql)  then
618           print *, 'Warning: plot buffer full, plot truncated'
619           goto 100
620        elseif (nqval(l) .eq. 0)  then
621           nqval(l) = nqval(l) + 1
622           qhval(nqval(l),l) = currpos
623           k = double_from_table_row(tabname, sname(l), j, d_val)
624           k = p(l)
625           qvval(nqval(l),l) = d_val
626           if (proc_flag(1,l) .eq. 1) then
627              qvval(nqval(l),l) = sqrt(abs(qvval(nqval(l),l)))
628           endif
629        elseif (itbv .eq. 0 .or. currpos - qhval(nqval(l),l)          &
630             .gt. mystep .or. (marker_plot .and. currtyp .eq. 25)) then
631           nqval(l) = nqval(l) + 1
632           qhval(nqval(l),l) = currpos
633           k = double_from_table_row(tabname, sname(l), j, d_val)
634           k = p(l)
635           qvval(nqval(l),l) = d_val
636           if (proc_flag(1,l) .eq. 1) then
637              qvval(nqval(l),l) = sqrt(abs(qvval(nqval(l),l)))
638           endif
639        endif
640     enddo
641
642     k = advance_node()
643
644     if (itbv .eq. 1 .and. k .eq. 0) goto 10
645  enddo
64610 continue
647100 continue
648  fpmach = machp .and. nelmach .gt. 0 .and. noline .eq. 0
649999 continue
650910 format('Interpolation is not compatible with the',                &
651       ' Twiss centre option')
652end subroutine pefill
653!***********************************************************************
654
655subroutine pegacn(ncc, window, act, xreal, yreal, np, usex,       &
656     xwpos, xpos, ypos, ilb)
657
658  use plotfi
659  implicit none
660
661  !----------------------------------------------------------------------*
662  ! Purpose:
663  !   Find suitable position for the curve annotation
664  ! Input:
665  !   NCC      (integer)  current curve count (1,2, etc.)
666  !   WINDOW   (real)     array containing the window to use
667  !   ACT      (real)     window in NDC
668  !   XREAL    (real)     x values of curve
669  !   YREAL    (real)     y values of curve
670  !   NP       (integer)  no. of points to plot
671  !   USEX     (real)     user character height expansion
672  ! Output:
673  !   XWPOS    (real)     x position of label in world coords.
674  !   XPOS     (real)     x pos. of label in NDC
675  !   YPOS     (real)     y pos. of label in NDC
676  !   ILB      (integer)  number of point behind label, or 0 if no
677  !                       label possible
678  !
679  ! calls the utility routine iucomp in this file.
680  !
681  ! it is called by the routine pecurv in this file.
682  !----------------------------------------------------------------------*
683
684  !--- type definition of the routine arguments
685
686  integer ncc, np, ilb
687  real window(*), act(*), xreal(*), yreal(*)
688  real usex, xwpos, xpos, ypos
689
690  !--- type definitions of local variables
691
692  real ywpos, xmax, xmin, xdiff, ydiff, d, t, eps
693  real xdiag(2,2), ydiag(2,2)
694  real xadd, yadd
695  integer i, iapos, iposx, iposy, j, iy
696  integer iucomp, kapos(mposx, mposy)
697
698  !--- Initialisation of local variables
699
700  save kapos
701
702  !--- Output initialisation
703
704  ilb = 0
705
706  !--- Routine body
707
708  !--- reset position array if first curve in frame
709
710  if (ncc .eq. 1)  then
711     do i = 1, mposx
712        do j = 1, mposy
713           kapos(i,j) = 0
714        enddo
715     enddo
716  endif
717  xdiff = window(2) - window(1)
718  ydiff = window(4) - window(3)
719  eps = 1.e-6 * max(xdiff, ydiff)
720  xmax  = xreal(1)
721  xmin  = xmax
722  do i = 2, np
723     xmin = min(xmin, xreal(i))
724     xmax = max(xmax, xreal(i))
725  enddo
72620 continue
727
728  !--- find first unoccupied position
729
730  iapos  = iucomp(0, kapos, mpost)
731  if (iapos .eq. 0)  then
732     ilb = 0
733  else
734     iposx  = mod (iapos-1, mposx) + 1
735     iposy  = (iapos-1) / mposx + 1
736     kapos(iposx,iposy) = -1
737
738     !--- annot. pos. in NDC
739
740     xpos = act(1) +                                                 &
741          0.125 * usex * (iposx - .5) * (act(2) - act(1))
742     ypos = act(4) -                                                 &
743          usex * (0.05 * (act(4) - act(3)) + 0.03 * (iposy - 1))
744
745     !---- annot. position in world coord.
746
747     xwpos = window(1) + xpos * xdiff
748
749     !--- get next if outside x values of curve
750
751     if (xwpos .le. xmin .or. xwpos .gt. xmax) goto 20
752     ywpos = window(3) + ypos * ydiff
753
754     !--- get endpoint of both diagonals of box
755
756     xadd = 0.0625 * xdiff
757     yadd = 0.03  * ydiff
758     xdiag(1,1) = xwpos - xadd
759     xdiag(2,1) = xwpos + xadd
760     xdiag(1,2) = xwpos - xadd
761     xdiag(2,2) = xwpos + xadd
762     ydiag(1,1) = ywpos
763     ydiag(2,1) = ywpos + yadd
764     ydiag(1,2) = ywpos + yadd
765     ydiag(2,2) = ywpos
766
767     !--- make sure no part of curve cuts these lines (curve approx. by
768     !    straight line segments)
769
770     do i = 2, np
771        if (xwpos .gt. xreal(i-1) .and. xwpos .le. xreal(i)) ilb = i
772        do j = 1, 2
773           d = (xdiag(2,j) - xdiag(1,j)) * (yreal(i-1) - yreal(i)) -   &
774                (ydiag(2,j) - ydiag(1,j)) * (xreal(i-1) - xreal(i))
775           if (abs(d) .lt. eps) goto 30
776           t = (xreal(i-1) - xdiag(1,j)) * (yreal(i-1) - yreal(i)) -   &
777                (yreal(i-1) - ydiag(1,j)) * (xreal(i-1) - xreal(i))
778           t = t / d
779           if (t .lt. 0. .or. t .gt. 1.) goto 30
780           t = (xdiag(2,j) - xdiag(1,j)) * (yreal(i-1) - ydiag(1,j)) - &
781                (ydiag(2,j) - ydiag(1,j)) * (xreal(i-1) - xdiag(1,j))
782           t = t / d
783           if (t .ge. 0. .and. t .le. 1.) goto 20
78430         continue
785        enddo
786     enddo
787  endif
788  if (ilb .gt. 0)  then
789     do iy = 1, mposy
790        kapos(iposx, iy) = 1
791     enddo
792  endif
793  do i = 1, mposx
794     do j = 1, mposy
795        kapos(i,j) = max(0, kapos(i,j))
796     enddo
797  enddo
798end subroutine pegacn
799!***********************************************************************
800
801subroutine pegaxn (nax, vax, sax, ns)
802
803  implicit none
804
805  !----------------------------------------------------------------------*
806  ! Purpose:
807  !   Returns compound vertical axis annotation
808  !
809  !--- Input
810  !   NAX       (integer) no. of vert. var. names in VAX
811  !   VAX          (char) vert. var. names
812  !---Output
813  !   SAX          (char) remaining (possibly truncated) names
814  !   NS        (integer) no. of names in SAX
815  ! calls the routines gxpnbl in file gxx11.F.
816  !
817  ! it is called by the routine pemima in this file.
818  !----------------------------------------------------------------------*
819
820  !--- type definition of the routine arguments
821
822  character * 16 vax(*), sax(*)
823  integer nax, ns
824
825  !--- type definitions of local variables
826
827  character * 16 scut, saloc
828  integer i, k, k1, k2, j, k1f, k2f
829
830  !--- Initialisation of local variables
831
832  ns = 0
833  sax(1) = ' '
834
835  !--- Routine body
836
837  if (nax .gt. 0)  then
838     do  i = 1, nax
839        saloc = vax(i)
840        call gxpnbl(saloc, k1, k2)
841        if (k2 .gt. 1 .and. index('XY', saloc(k2:k2)) .ne. 0)  then
842           scut = saloc(:k2-1)
843           do j = 1, ns
844              if (scut .eq. sax(j))  goto 10
845           enddo
846           do j = i + 1, nax
847              call gxpnbl(vax(j), k1f, k2f)
848              if (k2 .eq. k2f)  then
849                 if (index('XY', vax(j)(k2:k2)) .ne. 0)  then
850                    if (saloc(:k2-1) .eq. vax(j)(:k2-1))  then
851                       saloc = scut
852                       do k = 1, ns
853                          if (saloc .eq. sax(k))  goto 10
854                       enddo
855                    endif
856                 endif
857              endif
858           enddo
859        endif
860        ns      = ns + 1
861        sax(ns) = saloc
86210      continue
863     enddo
864  endif
865end subroutine pegaxn
866!***********************************************************************
867
868subroutine pegetn (iflag, svar, it, ipflg, sovar, reqann)
869
870  use plotfi
871  use plot_bfi
872  implicit none
873
874  !----------------------------------------------------------------------*
875  ! Purpose:                                                             *
876  !   Finds variable, dependent variables, axis and curve annotations    *
877  ! Input:                                                               *
878  !   IFLAG    (integer)  0 for dependent variables and process flag,    *
879  !                       1 for axis, 2 for curve, 3 for trunc. name,    *
880  !                       4 to print the axis names on IQLOG             *
881  !   SVAR        (char)  variable to be looked up.                      *
882  !   IT          (int)   table number (see PLGTBS).                     *
883  ! Output:                                                              *
884  !   IPFLG(1) (integer)  process flag: 0 as is, 1 take root, else call  *
885  !                       function PLPVAL                                *
886  !   IPFLG(2) (integer)  interpol. flag: 0 spline, else call            *
887  !                       function PEINTP                                *
888  !   SOVAR       (char)  array of (up to MXDEP) dependent variables     *
889  !   reqann      (char)  requested annotation                           *
890  !                                                                      *
891  ! calls the utility routine pupnbl in this file.                       *
892  !                                                                      *
893  ! it is called by the routines pemima, peplot and pesopt in this file. *
894  !----------------------------------------------------------------------*
895
896  !--- type definition of the routine arguments
897
898  integer iflag, it, ipflg(2)
899  character         svar*(mcnam), sovar(*)*(mcnam), reqann * (*)
900
901  !--- type definitions of local variables
902
903  character         svlabl(mnvar)*(mxlabl)
904  character         svanno(mnvar)*(mxlabl)
905  character         svname(mnvar)*(mcnam)
906  !--- strings:
907  !   SVLABL   plot prescriptions for variables on axis labels
908  !   SVANNO   plot prescriptions for variables in annotations
909  !   SVNAME   names of variables known to the program
910  integer i, iref, k1, k2, k1f, k2f, j
911  integer              iproc(mnvar,3), intpo(mnvar)
912  integer              ivdep(mnvar,mxdep,3)
913
914  !--- Initialisation of local variables
915
916  data (svname(j), j = 1, 32) /                                     &
917       's', 'size', 'deltap',                                            &
918       'qs', 'x', 'y', 'xsize', 'ysize',                                 &
919       'dt', 'xn', 'yn', 'pxn', 'pyn',                                   &
920       'gammatr', 'xrms', 'yrms',                                        &
921       'xmax', 'ymax', 'bxmax',                                          &
922       'bymax', 'dxmax', 'dymax',                                        &
923       'tn', 't', 'turns', 'particle', 'alfa',                           &
924       'ptn', 'wt', 'phit',                                              &
925       'rbxmax',                                                         &
926       'rbymax' /
927  data (svname(j), j = 33, mnvar) /                                 &
928       'betx', 'rbetx',                                                  &
929       'alfx', 'mux', 'dx',                                              &
930       'dpx', 'qx', 'px', 'wx',                                          &
931       'phix', 'dmux',                                                   &
932       'ddx', 'ddpx', 'iwx',                                             &
933       'xix',                                                            &
934       'bety', 'rbety',                                                  &
935       'alfy', 'muy', 'dy',                                              &
936       'dpy', 'qy', 'py', 'wy',                                          &
937       'phiy', 'dmuy',                                                   &
938       'ddy', 'ddpy', 'iwy',                                             &
939       'xiy', 'xns', 'pxns', 'wxs',                                      &
940       'yns', 'pyns', 'wys',                                             &
941       'energy', 'spintune',                                             &
942       'poltotal', 'poldiffx', 'poldiffy', 'poldiffs' /
943
944  data (svlabl(j), j = 1, 32) /                                     &
945       's (m)', 'n<G>s<G> (mm)', '<G>d<G><?>E<?>/p<?>0<?>c',             &
946       'Q<?>s<?>', 'x (m)', 'y (m)', 'n<G>s<G> (mm)', 'n<G>s<G> (mm)',   &
947       'ct (m)', 'x<?>n<?>', 'y<?>n<?>', 'p<?>xn<?>', 'p<?>yn<?>',       &
948       '<G>g<G><?>tr<?>', 'X<?>rms<?> (m)', 'Y<?>rms<?> (m)',            &
949       'X<?>max<?> (m)', 'Y<?>max<?> (m)', '<G>b<G><?>x_max<?> (m)',     &
950       '<G>b<G><?>y_max<?> (m)', 'D<?>x_max<?> (m)', 'D<?>y_max<?> (m)', &
951       't<?>n<?>', 'ct (m)', 'turns', 'particle', '<G>a<G>',             &
952       'p<?>t_n<?>', 'W<?>t<?>', '<G>F<G><?>t<?> (rad/2<G>p<G>)',        &
953       '<G>b<G><?>x_max<?><!>1/2<!> (m<!>1/2<!>)',                       &
954       '<G>b<G><?>y_max<?><!>1/2<!> (m<!>1/2<!>)' /
955  data (svlabl(j), j = 33, mnvar) /                                 &
956       '<G>b<G><?>x<?> (m)', '<G>b<G><?>x<?><!>1/2<!> (m<!>1/2<!>)',     &
957       '<G>a<G><?>x<?>', '<G>m<G><?>x<?> (rad/2<G>p<G>)', 'D<?>x<?> (m)',&
958       'D<?>px<?>', 'Q<?>x<?>', 'p<?>x<?>/p<?>0<?>', 'W<?>x<?>',         &
959       '<G>F<G><?>x<?> (rad/2<G>p<G>)', 'd<G>m<G><?>x<?>/d<G>d<G>',      &
960       'dD<?>x<?>/d<G>d<D> (m)', 'dD<?>px<?>/d<G>d<G>', 'W<?>x<?> (m)',  &
961       'XI<?>x<?>',                                                      &
962       '<G>b<G><?>y<?> (m)', '<G>b<G><?>y<?><!>1/2<!> (m<!>1/2<!>)',     &
963       '<G>a<G><?>y<?>', '<G>m<G><?>y<?> (rad/2<G>p<G>)', 'D<?>y<?> (m)',&
964       'D<?>py<?>', 'Q<?>y<?>', 'p<?>y<?>/p<?>0<?>', 'W<?>y<?>',         &
965       '<G>F<G><?>y<?> (rad/2<G>p<G>)', 'd<G>m<G><?>y<?>/d<G>d<G>',      &
966       'dD<?>y<?>/d<G>d<D> (m)', 'dD<?>py<?>/d<G>d<G>', 'W<?>y<?> (m)',  &
967       'XI<?>y<?>', 'x<?>ns<?>', 'p<?>x_ns<?>', 'W<?>xs<?>',             &
968       'y<?>ns<?>', 'p<?>y_ns<?>', 'W<?>ys<?>',                          &
969       'E[GeV]', 'spintune',                                             &
970       'polarization','polarization','polarization','polarization' /
971
972  data (svanno(j), j = 1, 32) /                                     &
973       's', 'n<G>s<G>', '<G>d<G>',                                       &
974       'Q<?>s<?>', 'x', 'y', 'n<G>s<G><?>x<?>', 'n<G>s<G><?>y<?>',       &
975       'ct', 'x<?>n<?>', 'y<?>n<?>', 'p<?>xn<?>', 'p<?>yn<?>',           &
976       '<G>g<G><?>tr<?>', 'X<?>rms<?>', 'Y<?>rms<?>',                    &
977       'X<?>max<?>', 'Y<?>max<?>', '<G>b<G><?>x_max<?>',                 &
978       '<G>b<G><?>y_max<?>', 'D<?>x_max<?>', 'D<?>y_max<?>',             &
979       't<?>n<?>', 't', 'turns', 'particle', '<G>a<G>',                  &
980       'p<?>t_n<?>', 'W<?>t<?>', '<G>F<G><?>t<?>',                       &
981       '<G>b<G><?>x_max<?><!>1/2<!>',                                    &
982       '<G>b<G><?>y_max<?><!>1/2<!>' /
983  data (svanno(j), j = 33, mnvar) /                                 &
984       '<G>b<G><?>x<?>', '<G>b<G><?>x<?><!>1/2<!>',                      &
985       '<G>a<G><?>x<?>', '<G>m<G><?>x<?>', 'D<?>x<?>',                   &
986       'D<?>px<?>', 'Q<?>x<?>', 'p<?>x<?>', 'W<?>x<?>',                  &
987       '<G>F<G><?>x<?>', '<G>m<G><?>x<?>''',                             &
988       'D<?>x<?>''', 'D<?>px<?>''', 'W<?>x<?>',                          &
989       'XI<?>x<?>',                                                      &
990       '<G>b<G><?>y<?>', '<G>b<G><?>y<?><!>1/2<!>',                      &
991       '<G>a<G><?>y<?>', '<G>m<G><?>y<?>', 'D<?>y<?>',                   &
992       'D<?>py<?>', 'Q<?>y<?>', 'p<?>y<?>', 'W<?>y<?>',                  &
993       '<G>F<G><?>y<?>', '<G>m<G><?>y<?>''',                             &
994       'D<?>y<?>''', 'D<?>py<?>''', 'W<?>y<?>',                          &
995       'XI<?>y<?>', 'x<?>ns<?>', 'p<?>x_ns<?>', 'W<?>xs<?>',             &
996       'y<?>ns<?>', 'p<?>y_ns<?>', 'W<?>ys<?>',                          &
997       'E', ' ',                                                         &
998       'p<?>tot<?>', 'p<?>diff_x<?>', 'p<?>diff_y<?>', 'p<?>diff_s<?>'/
999
1000  data (iproc(j,1), j = 1, 32) /                                    &
1001       0, 0, 0,                                                          &
1002       0, 0, 0, 0, 0,                                                    &
1003       0, 2, 3, 4, 5,                                                    &
1004       0, 0, 0,                                                          &
1005       0, 0, 0,                                                          &
1006       0, 0, 0,                                                          &
1007       6, 0, 0, 0, 0,                                                    &
1008       7, 8, 9,                                                          &
1009       1,                                                                &
1010       1 /
1011
1012  data (iproc(j,1), j = 33, mnvar) /                                &
1013       0, 1,                                                             &
1014       0, 0, 0,                                                          &
1015       0, 0, 0, 0,                                                       &
1016       0, 0,                                                             &
1017       0, 0, 0,                                                          &
1018       0,                                                                &
1019       0, 1,                                                             &
1020       0, 0, 0,                                                          &
1021       0, 0, 0, 0,                                                       &
1022       0, 0,                                                             &
1023       0, 0, 0,                                                          &
1024       0, 14, 15, 16,                                                    &
1025       17, 18, 19,                                                       &
1026       0, 0,                                                             &
1027       0, 0, 0, 0 /
1028
1029  data (iproc(j,2), j = 1, 32) /                                    &
1030       0, 0, 0,                                                          &
1031       0, 0, 0, 0, 0,                                                    &
1032       0, 2, 3, 4, 5,                                                    &
1033       0, 0, 0,                                                          &
1034       0, 0, 0,                                                          &
1035       0, 0, 0,                                                          &
1036       6, 0, 0, 0, 0,                                                    &
1037       7, 8, 9,                                                          &
1038       1,                                                                &
1039       1 /
1040
1041  data (iproc(j,2), j = 33, mnvar) /                                &
1042       0, 1,                                                             &
1043       0, 0, 0,                                                          &
1044       0, 0, 0, 10,                                                      &
1045       11, 0,                                                            &
1046       0, 0, 0,                                                          &
1047       0,                                                                &
1048       0, 1,                                                             &
1049       0, 0, 0,                                                          &
1050       0, 0, 0, 12,                                                      &
1051       13, 0,                                                            &
1052       0, 0, 0,                                                          &
1053       0, 14, 15, 16,                                                    &
1054       17, 18, 19,                                                       &
1055       0, 0,                                                             &
1056       0, 0, 0, 0 /
1057
1058  data (iproc(j,3), j = 1, 32) /                                    &
1059       0, 0, 0,                                                          &
1060       0, 0, 0, 0, 0,                                                    &
1061       0, 2, 3, 4, 5,                                                    &
1062       0, 0, 0,                                                          &
1063       0, 0, 0,                                                          &
1064       0, 0, 0,                                                          &
1065       6, 0, 0, 0, 0,                                                    &
1066       7, 8, 9,                                                          &
1067       1,                                                                &
1068       1 /
1069
1070  data (iproc(j,3), j = 33, mnvar) /                                &
1071       0, 1,                                                             &
1072       0, 0, 0,                                                          &
1073       0, 0, 0, 10,                                                      &
1074       11, 0,                                                            &
1075       0, 0, 0,                                                          &
1076       0,                                                                &
1077       0, 1,                                                             &
1078       0, 0, 0,                                                          &
1079       0, 0, 0, 12,                                                      &
1080       13, 0,                                                            &
1081       0, 0, 0,                                                          &
1082       0, 14, 15, 16,                                                    &
1083       17, 18, 19,                                                       &
1084       0, 0,                                                             &
1085       0, 0, 0, 0 /
1086
1087  data (intpo(j), j = 1, 32) / 32 * 0 /
1088  !--- in INTPO, n+100 means: take SQRT of var. n
1089  data (intpo(j), j = 33, mnvar) /                                  &
1090       1, 101,                                                           &
1091       2, 3, 4,                                                          &
1092       5, 0, 0, 0,                                                       &
1093       0, 0,                                                             &
1094       0, 0, 0,                                                          &
1095       0,                                                                &
1096       6, 106,                                                           &
1097       7, 8, 9,                                                          &
1098       10, 0, 0, 0,                                                      &
1099       0, 0,                                                             &
1100       0, 0, 0,                                                          &
1101       0, 0, 0, 0,                                                       &
1102       0, 0, 0,                                                          &
1103       0, 0,                                                             &
1104       0, 0, 0, 0 /
1105
1106  data (ivdep(j,1,1), j = 1, 32) /                                  &
1107       1, 2, 3,                                                          &
1108       4, 5, 6, 7, 8,                                                    &
1109       9, 5, 6, 5, 6,                                                    &
1110       14, 15, 16,                                                       &
1111       17, 18, 19,                                                       &
1112       20, 21, 22,                                                       &
1113       24, 24, 25, 26, 27,                                               &
1114       3, 3, 3,                                                          &
1115       19,                                                               &
1116       20 /
1117  data (ivdep(j,2,1), j = 1, 32) /                                  &
1118       0, 0, 0,                                                          &
1119       0, 0, 0, 0, 0,                                                    &
1120       0, 0, 0, 40, 55,                                                  &
1121       0, 0, 0,                                                          &
1122       0, 0, 0,                                                          &
1123       0, 0, 0,                                                          &
1124       0, 0, 0, 0, 0,                                                    &
1125       0, 24, 24,                                                        &
1126       0,                                                                &
1127       0 /
1128  data (ivdep(j,1,1), j = 33, mnvar) /                              &
1129       33, 33,                                                           &
1130       35, 36, 37,                                                       &
1131       38, 39, 40, 41,                                                   &
1132       42, 43,                                                           &
1133       44, 45, 46,                                                       &
1134       47,                                                               &
1135       48, 48,                                                           &
1136       50, 51, 52,                                                       &
1137       53, 54, 55, 56,                                                   &
1138       57, 58,                                                           &
1139       59, 60, 61,                                                       &
1140       62, 5, 5, 5,                                                      &
1141       6, 6, 6,                                                          &
1142       69, 70,                                                           &
1143       71, 72, 73, 74 /
1144  data (ivdep(j,2,1), j = 33, mnvar) /                              &
1145       0, 0,                                                             &
1146       0, 0, 0,                                                          &
1147       0, 0, 0, 0,                                                       &
1148       0, 0,                                                             &
1149       0, 0, 0,                                                          &
1150       0,                                                                &
1151       0, 0,                                                             &
1152       0, 0, 0,                                                          &
1153       0, 0, 0, 0,                                                       &
1154       0, 0,                                                             &
1155       0, 0, 0,                                                          &
1156       0, 0, 40, 40,                                                     &
1157       0, 55, 55,                                                        &
1158       0, 0,                                                             &
1159       0, 0, 0, 0 /
1160
1161  data (ivdep(j,1,2), j = 1, 32) /                                  &
1162       1, 2, 3,                                                          &
1163       4, 5, 6, 7, 8,                                                    &
1164       9, 5, 6, 5, 6,                                                    &
1165       14, 15, 16,                                                       &
1166       17, 18, 19,                                                       &
1167       20, 21, 22,                                                       &
1168       24, 24, 25, 26, 27,                                               &
1169       3, 3, 3,                                                          &
1170       19,                                                               &
1171       20 /
1172  data (ivdep(j,2,2), j = 1, 32) /                                  &
1173       0, 0, 0,                                                          &
1174       0, 0, 0, 0, 0,                                                    &
1175       0, 0, 0, 40, 55,                                                  &
1176       0, 0, 0,                                                          &
1177       0, 0, 0,                                                          &
1178       0, 0, 0,                                                          &
1179       0, 0, 0, 0, 0,                                                    &
1180       0, 24, 24,                                                        &
1181       0,                                                                &
1182       0 /
1183  data (ivdep(j,1,2), j = 33, mnvar) /                              &
1184       33, 33,                                                           &
1185       35, 36, 37,                                                       &
1186       38, 39, 40, 5,                                                    &
1187       5, 43,                                                            &
1188       44, 45, 46,                                                       &
1189       47,                                                               &
1190       48, 48,                                                           &
1191       50, 51, 52,                                                       &
1192       53, 54, 55, 6,                                                    &
1193       6, 58,                                                            &
1194       59, 60, 61,                                                       &
1195       62, 5, 5, 5,                                                      &
1196       6, 6, 6,                                                          &
1197       69, 70,                                                           &
1198       71, 72, 73, 74 /
1199  data (ivdep(j,2,2), j = 33, mnvar) /                              &
1200       0, 0,                                                             &
1201       0, 0, 0,                                                          &
1202       0, 0, 0, 40,                                                      &
1203       40, 0,                                                            &
1204       0, 0, 0,                                                          &
1205       0,                                                                &
1206       0, 0,                                                             &
1207       0, 0, 0,                                                          &
1208       0, 0, 0, 55,                                                      &
1209       55, 0,                                                            &
1210       0, 0, 0,                                                          &
1211       0, 0, 40, 40,                                                     &
1212       0, 55, 55,                                                        &
1213       0, 0,                                                             &
1214       0, 0, 0, 0 /
1215
1216  data (ivdep(j,1,3), j = 1, 32) /                                  &
1217       1, 2, 3,                                                          &
1218       4, 5, 6, 7, 8,                                                    &
1219       9, 5, 6, 5, 6,                                                    &
1220       14, 15, 16,                                                       &
1221       17, 18, 19,                                                       &
1222       20, 21, 22,                                                       &
1223       24, 24, 25, 26, 27,                                               &
1224       3, 3, 3,                                                          &
1225       19,                                                               &
1226       20 /
1227  data (ivdep(j,2,3), j = 1, 32) /                                  &
1228       0, 0, 0,                                                          &
1229       0, 0, 0, 0, 0,                                                    &
1230       0, 0, 0, 40, 55,                                                  &
1231       0, 0, 0,                                                          &
1232       0, 0, 0,                                                          &
1233       0, 0, 0,                                                          &
1234       0, 0, 0, 0, 0,                                                    &
1235       0, 24, 24,                                                        &
1236       0,                                                                &
1237       0 /
1238  data (ivdep(j,1,3), j = 33, mnvar) /                              &
1239       33, 33,                                                           &
1240       35, 36, 37,                                                       &
1241       38, 39, 40, 5,                                                    &
1242       5, 43,                                                            &
1243       44, 45, 46,                                                       &
1244       47,                                                               &
1245       48, 48,                                                           &
1246       50, 51, 52,                                                       &
1247       53, 54, 55, 6,                                                    &
1248       6, 58,                                                            &
1249       59, 60, 61,                                                       &
1250       62, 5, 5, 5,                                                      &
1251       6, 6, 6,                                                          &
1252       69, 70,                                                           &
1253       71, 72, 73, 74 /
1254  data (ivdep(j,2,3), j = 33, mnvar) /                              &
1255       0, 0,                                                             &
1256       0, 0, 0,                                                          &
1257       0, 0, 0, 40,                                                      &
1258       40, 0,                                                            &
1259       0, 0, 0,                                                          &
1260       0,                                                                &
1261       0, 0,                                                             &
1262       0, 0, 0,                                                          &
1263       0, 0, 0, 55,                                                      &
1264       55, 0,                                                            &
1265       0, 0, 0,                                                          &
1266       0, 0, 40, 40,                                                     &
1267       0, 55, 55,                                                        &
1268       0, 0,                                                             &
1269       0, 0, 0, 0 /
1270
1271  !--- Routine body
1272
1273  if (it .le. 0 .or. it .gt. 3)  then
1274     sovar(1) = svar
1275     sovar(2) = ' '
1276     reqann   = svar
1277     ipflg(1) = 0
1278     ipflg(2) = 0
1279     goto 999
1280  endif
1281  sovar(1) = ' '
1282  reqann = svar
1283
1284  !--- search in list of known variables
1285
1286  do  iref = 1, mnvar
1287     if (svar .eq. svname(iref))  goto 9
1288  enddo
1289  call pupnbl(svar, k1, k2)
1290  do  iref = 1, mnvar
1291     call pupnbl(svname(iref), k1f, k2f)
1292     if (k2 + 1 .eq. k2f)  then
1293        if (index('xy', svname(iref)(k2f:k2f)) .ne. 0)  then
1294           if (svar(:k2) .eq. svname(iref)(:k2))  goto 9
1295        endif
1296     endif
1297  enddo
1298  goto 999
12999 continue
1300  if (iflag .eq. 0)  then
1301     reqann = svname(iref)
1302     ipflg(1) = iproc(iref,it)
1303     ipflg(2) = intpo(iref)
1304     do  j = 1, mxdep
1305        if (ivdep(iref,j,it) .eq. 0)  then
1306           sovar(j) = ' '
1307        else
1308           sovar(j) = svname(ivdep(iref,j,it))
1309        endif
1310     enddo
1311  elseif (iflag .eq. 1) then
1312     reqann = svlabl(iref)
1313     if (svar .ne. svname(iref))  then
1314
1315        !--- incomplete match
1316        !    replace x or y in name by blank
1317
1318        call pupnbl(reqann, k1, k2)
1319        do  i = 2, k2
1320           if (index('XYxy', reqann(i:i)) .ne. 0)  then
1321              reqann(i:i) = ' '
1322           endif
1323        enddo
1324     endif
1325  elseif (iflag .eq. 2) then
1326     reqann = svanno(iref)
1327  elseif (iflag .eq. 3) then
1328     if (svar .eq. svname(iref))  then
1329        reqann = svname(iref)
1330     else
1331        reqann = svname(iref)(:k2)
1332     endif
1333  else
1334     reqann = svar
1335  endif
1336999 end subroutine pegetn
1337
1338  !***********************************************************************
1339
1340subroutine peiact(kact, np, x, y, ac, kf, kl)
1341
1342  implicit none
1343
1344  !----------------------------------------------------------------------*
1345  ! Purpose:
1346  !   Return first and last point of curve inside active window
1347  ! Input:
1348  !   KACT        (int)   starting point for check
1349  !   NP          (int)   number of points in XVAL, YVAL
1350  !   X           (real)  x values
1351  !   Y           (real)  y values
1352  !   AC          (real)  active window in WC
1353  ! Output:
1354  !   KF          (int)   first point inside, or 0
1355  !   KL          (int)   last  point inside, or 0
1356  !
1357  !----------------------------------------------------------------------*
1358
1359  !--- type definition of the routine arguments
1360
1361  integer kact, np, kf, kl
1362  real x(*), y(*), ac(4)
1363
1364  !--- type definitions of local variables
1365
1366  integer i
1367  real toleps, xtol, ytol
1368  parameter (toleps = 1.e-5)
1369
1370  !--- Routine body
1371
1372  xtol = toleps * (ac(2) - ac(1))
1373  ytol = toleps * (ac(4) - ac(3))
1374  kf = 0
1375  kl = 0
1376  do i = kact, np
1377     if(x(i) + xtol .lt. ac(1)) goto 10
1378     if(x(i) - xtol .gt. ac(2)) goto 10
1379     if(y(i) + ytol .lt. ac(3)) goto 10
1380     if(y(i) - ytol .gt. ac(4)) goto 10
1381     kf = i
1382     goto 20
138310   continue
1384  enddo
1385
1386  !--- no point inside
1387
1388  goto 999
138920 continue
1390  do i = kf, np
1391     if(x(i) + xtol .lt. ac(1)) goto 40
1392     if(x(i) - xtol .gt. ac(2)) goto 40
1393     if(y(i) + ytol .lt. ac(3)) goto 40
1394     if(y(i) - ytol .gt. ac(4)) goto 40
1395  enddo
139640 kl = i - 1
1397999 continue
1398end subroutine peiact
1399!***********************************************************************
1400
1401subroutine peintp(crow, nint, proc, length, ierr)
1402
1403  use plotfi
1404  use plot_bfi
1405  use plot_mathfi
1406  implicit none
1407
1408  !----------------------------------------------------------------------*
1409  !     purpose:
1410  !     interpolate variables plotted against s
1411  !     input:
1412  !     crow        (int)   table row number at start of element
1413  !     nint        (int)   number of interpolation intervals
1414  !     type        (int)   (local) element type
1415  !     proc        (int)   original process flags
1416  !     step        (d.p.)  max. dist. between two successive hor. values
1417  !     length      (d.p.)  element length
1418  !     output:
1419  !     ierr        (int)   0 if ok, else > 0
1420  !     the results are stored in qhval and qvval
1421  !
1422  !     calls the functions double_from_table_row and get_value
1423  !     defined in file madxn.c.
1424  !     calls the function peelma in this file.
1425  !
1426  !     it is called by the routine pefill in this file.
1427  !----------------------------------------------------------------------*
1428
1429  !---  type definition of the routine arguments
1430
1431  integer crow, nint, proc(2,*), ierr
1432  double precision length
1433
1434  !---  type definitions of local variables
1435
1436  double precision tw1(mintpl)
1437  double precision ex, ey
1438  double precision xn1, pxn1, yn1, pyn1
1439  double precision s_elem, s_incr, s, gamx, gamy
1440  integer i, j, k, ipc
1441
1442  !---  definitions of function primitives
1443
1444  integer double_from_table_row
1445  integer interpolate_node, reset_interpolation
1446  integer embedded_twiss
1447  double precision get_value
1448
1449  !---  Output initialisation
1450
1451  ierr = 0
1452
1453  !---  Routine body
1454
1455  if (crow .eq. 1) then
1456     k = double_from_table_row(tabname, 'x ', 1, tw1(11))
1457     k = double_from_table_row(tabname, 'px ', 1, tw1(12))
1458     k = double_from_table_row(tabname, 'betx ', 1, tw1(1))
1459     k = double_from_table_row(tabname, 'alfx ', 1, tw1(2))
1460     k = double_from_table_row(tabname, 'mux ', 1, tw1(3))
1461     k = double_from_table_row(tabname, 'dx ', 1, tw1(4))
1462     k = double_from_table_row(tabname, 'dpx ', 1, tw1(5))
1463     k = double_from_table_row(tabname, 'y ', 1, tw1(13))
1464     k = double_from_table_row(tabname, 'py ', 1, tw1(14))
1465     k = double_from_table_row(tabname, 'bety ', 1, tw1(6))
1466     k = double_from_table_row(tabname, 'alfy ', 1, tw1(7))
1467     k = double_from_table_row(tabname, 'muy ', 1, tw1(8))
1468     k = double_from_table_row(tabname, 'dy ', 1, tw1(9))
1469     k = double_from_table_row(tabname, 'dpy ', 1, tw1(10))
1470     k = double_from_table_row(tabname, 's ', 1, s_incr)
1471     s = 0.0
1472     ex = get_value('beam ','ex ')
1473     ey = get_value('beam ','ey ')
1474
1475     !---  xn, pxn, yn, pyn
1476
1477     if (ex * tw1(1).eq. zero)  then
1478        xn1 = zero
1479     else
1480        xn1 = tw1(11) / sqrt(ex*abs(tw1(1)))
1481     endif
1482     if (ey * tw1(6) .eq. zero)  then
1483        yn1 = zero
1484     else
1485        yn1 = tw1(13) / sqrt(ey*abs(tw1(6)))
1486     endif
1487     pxn1 = tw1(12) * gamx
1488     pyn1 = tw1(14) * gamy
1489     tw1(15) = xn1
1490     tw1(16) = pxn1
1491     tw1(17) = yn1
1492     tw1(18) = pyn1
1493
1494     !---  loop over variables, interpolate those with codes
1495
1496     do j = 1, nivvar
1497        ipc = mod(proc(2,j), 100)
1498        if (ipc .gt. 0)  then
1499           if (nqval(j) .eq. maxseql)  then
1500              ierr = 1
1501              goto 999
1502           endif
1503           ipparm(2,j) = 1
1504           nqval(j) = nqval(j) + 1
1505           qhval(nqval(j),j) = s
1506           if (proc(1,j) .gt. 0)  then
1507              qvval(nqval(j), j) = sqrt(abs(tw1(ipc)))
1508           else
1509              qvval(nqval(j), j) = tw1(ipc)
1510           endif
1511        else
1512           ipparm(2,j) = 0
1513        endif
1514     enddo
1515     goto 999
1516  endif
1517  if (length .eq. zero)  goto 999
1518  k = double_from_table_row(tabname, horname, crow - 1, s_elem)
1519
1520  !---  set flag for correct interpolation
1521
1522
1523  !---  get intermediate s values, and interpolate twiss parameters
1524
1525  k = interpolate_node(nint)
1526  k = embedded_twiss()
1527  do i = 1, nint
1528     k = double_from_table_row('embedded_twiss_table ',                  &
1529          'x ', i, tw1(11))
1530     k = double_from_table_row('embedded_twiss_table ',                  &
1531          'px ', i, tw1(12))
1532     k = double_from_table_row('embedded_twiss_table ',                  &
1533          'betx ', i, tw1(1))
1534     k = double_from_table_row('embedded_twiss_table ',                  &
1535          'alfx ', i, tw1(2))
1536     k = double_from_table_row('embedded_twiss_table ',                  &
1537          'mux ', i, tw1(3))
1538     k = double_from_table_row('embedded_twiss_table ',                  &
1539          'dx ', i, tw1(4))
1540     k = double_from_table_row('embedded_twiss_table ',                  &
1541          'dpx ', i, tw1(5))
1542     k = double_from_table_row('embedded_twiss_table ',                  &
1543          'y ', i, tw1(13))
1544     k = double_from_table_row('embedded_twiss_table ',                  &
1545          'py ', i, tw1(14))
1546     k = double_from_table_row('embedded_twiss_table ',                  &
1547          'bety ', i, tw1(6))
1548     k = double_from_table_row('embedded_twiss_table ',                  &
1549          'alfy ', i, tw1(7))
1550     k = double_from_table_row('embedded_twiss_table ',                  &
1551          'muy ', i, tw1(8))
1552     k = double_from_table_row('embedded_twiss_table ',                  &
1553          'dy ', i, tw1(9))
1554     k = double_from_table_row('embedded_twiss_table ',                  &
1555          'dpy ', i, tw1(10))
1556     k = double_from_table_row('embedded_twiss_table ',                  &
1557          's ', i, s_incr)
1558     s = s_elem + s_incr
1559     ex = get_value('beam ','ex ')
1560     ey = get_value('beam ','ey ')
1561
1562     !---  xn, pxn, yn, pyn
1563
1564     if (ex * tw1(1).eq. zero)  then
1565        xn1 = zero
1566     else
1567        xn1 = tw1(11) / sqrt(ex*abs(tw1(1)))
1568     endif
1569     if (ey * tw1(6) .eq. zero)  then
1570        yn1 = zero
1571     else
1572        yn1 = tw1(13) / sqrt(ey*abs(tw1(6)))
1573     endif
1574     if (tw1(1) .ne. zero)  then
1575        gamx = (one + tw1(2)**2) / tw1(1)
1576     else
1577        gamx = zero
1578     endif
1579     if (tw1(6) .ne. zero)  then
1580        gamy = (one + tw1(7)**2) / tw1(6)
1581     else
1582        gamy = zero
1583     endif
1584     pxn1 = tw1(12) * gamx
1585     pyn1 = tw1(14) * gamy
1586     tw1(15) = xn1
1587     tw1(16) = pxn1
1588     tw1(17) = yn1
1589     tw1(18) = pyn1
1590
1591     !---  loop over variables, interpolate those with codes
1592
1593     do j = 1, nivvar
1594        ipc = mod(proc(2,j), 100)
1595        if (ipc .gt. 0)  then
1596           if (nqval(j) .eq. maxseql)  then
1597              ierr = 1
1598              goto 999
1599           endif
1600           ipparm(2,j) = 1
1601           nqval(j) = nqval(j) + 1
1602           qhval(nqval(j),j) = s
1603           if (proc(1,j) .gt. 0)  then
1604              qvval(nqval(j), j) = sqrt(abs(tw1(ipc)))
1605           else
1606              qvval(nqval(j), j) = tw1(ipc)
1607           endif
1608        else
1609           ipparm(2,j) = 0
1610        endif
1611     enddo
1612  enddo
1613  k = reset_interpolation(nint)
1614
1615999 end subroutine peintp
1616  !***********************************************************************
1617
1618subroutine pemima
1619
1620  use plotfi
1621  use plot_bfi
1622  implicit none
1623
1624  !----------------------------------------------------------------------*
1625  ! Purpose:                                                             *
1626  !   Constrain axis reference, find minima and maxima of coordinates,   *
1627  !   construct axis labels                                              *
1628  ! calls the routine gxpnbl in file gxx11.F.                            *
1629  ! calls the routines pegetn and pegaxn in this file.                   *
1630  !                                                                      *
1631  ! it is called by the routine exec_plot in file madxn.c after pefill.  *
1632  !----------------------------------------------------------------------*
1633
1634  !--- type definitions of local variables
1635
1636  integer i, j, k, iv, idum(2), ns, k1, k2, i1, i2, it(4)
1637  character * (mtitl)  s
1638  character * (mxlabl) slab
1639  character * (mcnam) sdum(mxcurv), saxis(mxcurv), vaxis(mxcurv,4)
1640
1641  !--- Initialisation of variables in common peaddi
1642
1643  numax = 0
1644
1645  !--- Initialisation of local variables
1646
1647  do i = 1, 4
1648     it(i) = 0
1649  enddo
1650
1651  !--- Routine body
1652
1653  do  j = 1, nivvar
1654     do i = 1, numax
1655        if (it(i) .eq. naxref(j))  goto 10
1656     enddo
1657     if (numax .eq. 4)  then
1658        naxref(j) = it(4)
1659     else
1660        numax = numax + 1
1661        it(numax) = naxref(j)
1662     endif
166310   continue
1664  enddo
1665  do i = 1, 4
1666     do j = 1, numax - 1
1667        if (it(j) .gt. it(j+1))  then
1668           k = it(j)
1669           it(j) = it(j+1)
1670           it(j+1) = k
1671        endif
1672     enddo
1673  enddo
1674  do j = 1, nivvar
1675     do i = 1, numax
1676        if (naxref(j) .eq. it(i))  then
1677           naxref(j) = i
1678           goto 50
1679        endif
1680     enddo
168150   continue
1682  enddo
1683  do j = 1, nivvar
1684     k = naxref(j)
1685     do i = 1, nqval(j)
1686        hmima(1) = min(hmima(1), qhval(i,j))
1687        hmima(2) = max(hmima(2), qhval(i,j))
1688        vmima(1,k) = min(vmima(1,k), qvval(i,j))
1689        vmima(2,k) = max(vmima(2,k), qvval(i,j))
1690     enddo
1691  enddo
1692
1693  !--- get axis annotation
1694
1695  do j = 1, nivvar
1696     k = naxref(j)
1697     nvvar(k) = nvvar(k) + 1
1698     vaxis(nvvar(k),k) = slabl(j)
1699  enddo
1700  do iv = 1, 4
1701     if (nvvar(iv) .gt. 0)  then
1702        if (nvvar(iv) .eq. 1)  then
1703           call pegetn (1, vaxis(1,iv), itbv, idum, sdum, slab)
1704           ns = 1
1705        else
1706           call pegaxn (nvvar(iv), vaxis(1,iv), saxis, ns)
1707           call pegetn (1, saxis(1), itbv, idum, sdum, slab)
1708        endif
1709        call gxpnbl (slab, k1, k2)
1710        s  = '<#>' // slab
1711        k2 = k2 + 3
1712        do i = 2, ns
1713           call pegetn (1, saxis(i), itbv, idum, sdum, slab)
1714           call gxpnbl (slab, i1, i2)
1715           if (index(s(:k2),slab(:i2)) .eq. 0)  then
1716              s(k2 + 1:) = ', ' // slab(:i2)
1717              k2 = k2 + i2 + 2
1718           endif
1719        enddo
1720        axlabel(iv) = s
1721     endif
1722  enddo
1723end subroutine pemima
1724!***********************************************************************
1725
1726subroutine peplot
1727
1728  use plotfi
1729  use plot_bfi
1730  use plot_mathfi
1731  implicit none
1732
1733  !----------------------------------------------------------------------*
1734  ! Purpose:
1735  !   Plot all types of graphs from MAD.
1736  !   Uses GXPLOT with underlying X-Windows (PostScript)
1737  !
1738  ! calls the function plot_option and the routine comm_para
1739  ! defined in file madxn.c.
1740  ! calls the routines pegetn, pecurv and peschm in this file.
1741  ! calls the routines gxsdef, gxsvar, jslwsc, gxsvpt, gxpnbl, gxqaxs,
1742  !                    gxsaxs, iaxseq, gxqcrv, gxscrv, axlabel, gxfrm1,
1743  !                    gxswnd  defined in file gxx11.F.
1744  !
1745  ! it is called by the routine plotit in this file.
1746  !----------------------------------------------------------------------*
1747
1748  !--- type definitions of local variables
1749  !--- strings:
1750  !   SVAR     buffer for variable names etc.
1751  !   SLOCN    local name buffer (without leading "_")
1752  !   STEMP    temporary buffer for titles
1753  !   STEXT    buffer for labels etc.
1754  !   SFORM    format buffer
1755  !--- reals:
1756  !   PRMACH fraction of viewport taken by machine plot
1757  !   SYMCH  preset symbol character height
1758
1759  character svar*(mcnam)
1760  character *(mxlabl) slocn, slname
1761  character stemp*(mtitl), stext*300, sform*20
1762  character sdum(mxdep)*(mcnam)
1763  character * 1 ssymb
1764  character*80 ch
1765  double precision deltap
1766  real prmach, symch, tmpval, yvtop, fdum, chh
1767  real vpt(4), window(4,4), actwin(4,4), range(2), xax(2), yax(8)
1768  integer ipar(50), nptval(4), ipxval(4), ipyval(4), icvref(4)
1769  integer lastnb, iaxseq(4)
1770  integer idum, k1dum, k2dum, k3dum, i, npar, ivvar, nvax, ivax,    &
1771       ierr, vdum(2), j, k
1772
1773  !--- definitions of function primitives
1774
1775  double precision plot_option
1776  integer double_from_table_row, double_from_table_header
1777  integer table_column_exists, table_header_exists
1778
1779
1780  !--- Initialisation of local variables
1781
1782  data prmach /0.1/, symch /0.01/
1783  data iaxseq / 1, 4, 2, 3 /
1784
1785  ssymb = ' '
1786  fdum = 0.0D0
1787  chh = 0.0D0
1788  do i = 1 , 4
1789     vpt(i) = 0.0D0
1790     nptval(i) = 0
1791     ipxval(i) = 0
1792     ipyval(i) = 0
1793     icvref(i) = 0
1794     do j = 1 , 4
1795        window(i,j) = 0.0D0
1796        actwin(i,j) = 0.0D0
1797     enddo
1798  enddo
1799  do i = 1 , 2
1800     range(i) = 0.0D0
1801     xax(i) = 0.0D0
1802  enddo
1803  do i = 1 , 8
1804     yax(i) = 0.0D0
1805  enddo
1806  do i = 1 , 50
1807     ipar(i) = 0
1808  enddo
1809
1810  !--- Routine body
1811
1812  !--- Acquire deltap
1813  deltap=zero
1814  if(.not.ptc_flag) then
1815    if(tabname.eq."summ") then
1816      if (table_column_exists('summ ', 'deltap ').ne.0) then
1817        k = double_from_table_row('summ ', 'deltap ', 1, deltap)
1818      endif
1819    else if (table_header_exists(tabname, 'deltap ').ne.0) then
1820        k = double_from_table_header(tabname, 'deltap ', deltap)
1821    endif
1822  else if (table_header_exists(tabname, 'deltap ').ne.0) then
1823        k = double_from_table_header(tabname, 'deltap ', deltap)
1824  endif
1825
1826  deltap = deltap*100.0
1827  write(ch,'(f9.4)') deltap
1828  if(dpp_flag) then
1829     slocn = ' <#>'
1830  else
1831     slocn = 'Momentum offset = '//ch(:7)//' %'//'&'//'<#>'
1832  endif
1833
1834  !--- reset axis and curve defaults
1835
1836  call gxsdef ('AXIS', 0)
1837  call gxsdef ('CURVE', 0)
1838
1839  !--- set "new line" character (change default = '/')
1840
1841  call gxsvar ('SDEFNL', idum, fdum, '&')
1842
1843  !--- set top of viewport - leave space to plot machine if required
1844
1845  if (fpmach)  then
1846     yvtop = 1. - prmach
1847  else
1848     yvtop = 1.
1849  endif
1850
1851  !--- set line width scale factor
1852
1853  tmpval = plot_option('lwidth ')
1854  if (tmpval .eq. 0.) tmpval = 1.
1855  call jslwsc (tmpval)
1856
1857  !--- loop over frames
1858  !--- set viewport
1859
1860  vpt(1) = 0.
1861  vpt(2) = 1.
1862  vpt(3) = 0.
1863  vpt(4) = yvtop
1864  call gxsvpt (vpt)
1865
1866  !--- find variable name in list
1867
1868  svar = horname
1869  slname = "_"
1870  call pegetn (1, svar, itbv, vdum, sdum, slname)
1871  call gxpnbl(slname, k1dum, k2dum)
1872  k3dum = 31
1873  if (dpp_flag) k3dum = 4
1874  slocn = ' '
1875  if(k1dum.gt.0.and.k2dum.gt.0) then
1876     do idum = k1dum, k2dum
1877        if (slname(idum:idum) .ne. '_') then
1878           k3dum = k3dum + 1
1879           slocn(k3dum:k3dum) = slname(idum:idum)
1880        endif
1881     enddo
1882  endif
1883
1884  !--- prepare horizontal axis
1885
1886  do i = 1, 4, 3
1887     call gxqaxs ('X', i, npar, ipar, range, stext, sform)
1888
1889     !--- set character sizes for labels and text including user requests
1890
1891     ipar(7) = max (mlsize * qlscl + .01, 1.1)
1892     ipar(13) = max( mtsize * qtscl + .01, 1.1)
1893
1894     !--- text left adjusted
1895
1896     ipar(10) = 1
1897
1898     !--- font
1899
1900     ipar(11) = plot_option('font ')
1901     if (ipar(11) .eq. 0) ipar(11) = 1
1902
1903     !--- axis ref. number
1904
1905     ipar(21) = 1
1906
1907     !--- range centre etc.
1908
1909     if (hrange(1) .lt. hrange(2)) then
1910
1911        !--- use range as is
1912
1913        ipar(23) = 1
1914        range(1) = hrange(1)
1915        range(2) = hrange(2)
1916
1917        !--- set min. and max. for horizontal axis
1918
1919        xax(1) = hrange(1)
1920        xax(2) = hrange(2)
1921     else
1922        xax(1) = hmima(1)
1923        xax(2) = hmima(2)
1924     endif
1925     if (i .eq. 1) then
1926
1927        !--- bottom title
1928
1929        stext = slocn(:lastnb(slocn))
1930     else
1931
1932        !--- suppress labels on upper axis
1933
1934        ipar(3) = 0
1935
1936        !--- ticks below axis
1937
1938        ipar(4) = 1
1939
1940        !--- top title
1941
1942        stext = toptitle
1943     endif
1944
1945     !--- set axis parameters
1946
1947     call gxsaxs ('X', i, npar, ipar, range, stext, sform)
1948  enddo
1949  do nvax = 1, numax
1950
1951     !--- set curve parameters for frame call
1952
1953     ivax = iaxseq(nvax)
1954     call gxqcrv (nvax, npar, ipar, ssymb)
1955     ipar(2) = ivax
1956     call gxscrv (nvax, npar, ipar, ' ')
1957     call gxqaxs ('Y', ivax, npar, ipar, range, stext, sform)
1958
1959     !--- set character sizes for labels and text including user requests
1960
1961     ipar(7) = max (mlsize * qlscl + .01, 1.1)
1962     ipar(13) = max (mtsize * qtscl + .01, 1.1)
1963
1964     !--- right adjusted label
1965
1966     ipar(10) = 3
1967
1968     !--- font
1969
1970     ipar(11) = plot_option('font ')
1971     if (ipar(11) .eq. 0) ipar(11) = 1
1972
1973     !--- range centre etc
1974
1975     if (vrange(1,nvax) .lt. vrange(2,nvax)) then
1976
1977        !--- use range as is
1978
1979        ipar(23) = 1
1980        range(1) = vrange(1,nvax)
1981        range(2) = vrange(2,nvax)
1982
1983        !--- store y values for frame scaling
1984
1985        yax(2 * nvax - 1) = vrange(1,nvax)
1986        yax(2 * nvax) = vrange(2,nvax)
1987     else
1988
1989        !--- store y values for frame scaling
1990
1991        yax(2 * nvax - 1) = vmima(1,nvax)
1992        yax(2 * nvax) = vmima(2,nvax)
1993     endif
1994
1995     !--- get axis annotation
1996
1997     slocn = axlabel(nvax)
1998     stemp = ' '
1999     call gxpnbl(slocn, k1dum, k2dum)
2000     k3dum = 0
2001     do idum = k1dum, k2dum
2002        if (slocn(idum:idum) .ne. '_') then
2003           k3dum = k3dum + 1
2004           stemp(k3dum:k3dum) = slocn(idum:idum)
2005        endif
2006     enddo
2007     if (nvax .eq. 1) then
2008        stext = '&' // stemp
2009     else
2010        stext = stemp
2011     endif
2012     call gxsaxs ('Y', ivax, npar, ipar, range, stext, sform)
2013     nptval(nvax) = 2
2014     ipxval(nvax) = 1
2015     ipyval(nvax) = 2 * nvax - 1
2016     icvref(nvax) = nvax
2017  enddo
2018
2019  !--- if only one y axis, plot right axis with ticks only
2020
2021  if (numax .eq. 1) then
2022     ivax = 4
2023     call gxqaxs ('Y', ivax, npar, ipar, range, stext, sform)
2024     ipar(3) = 0
2025     ipar(4) = 1
2026     ipar(21) = 1
2027     call gxsaxs ('Y', ivax, npar, ipar, range, stext, sform)
2028  endif
2029
2030  !--- plot frame, keep windows for curves + clipping
2031
2032  call gxfrm1 (numax, nptval, ipxval, ipyval, icvref, xax, yax,     &
2033       window, actwin, ierr)
2034  if (ierr .ne. 0) goto 120
2035
2036  !--- now loop over vertical variables for real curve plotting
2037
2038  do ivvar = 1, nivvar
2039     nvax = naxref(ivvar)
2040     ivax = iaxseq(nvax)
2041
2042     !--- find variable name in list for annotation
2043
2044     svar = slabl(ivvar)
2045     slname = "_"
2046     call pegetn (2, svar, itbv, vdum, sdum, slname)
2047     call gxpnbl(slname, k1dum, k2dum)
2048     k3dum = 0
2049     slocn = ' '
2050     if(k1dum.gt.0.and.k2dum.gt.0) then
2051        do idum = k1dum, k2dum
2052           if (slname(idum:idum) .ne. '_') then
2053              k3dum = k3dum + 1
2054              slocn(k3dum:k3dum) = slname(idum:idum)
2055           endif
2056        enddo
2057     endif
2058
2059     !--- character height including user request
2060
2061     chh = 0.001 * masize * qascl
2062
2063     !--- call curve plot routine with simple arrays and flags
2064
2065     call pecurv (ivvar, slocn, chh, qascl,                          &
2066          symch * qsscl, ipparm(1,ivvar), nqval(ivvar), qhval(1,ivvar),     &
2067          qvval(1,ivvar), window(1,nvax), actwin(1,nvax), ierr)
2068     if (ierr .ne. 0) goto 150
2069  enddo
2070  if (fpmach)  then
2071     vpt(1) = 0.
2072     vpt(2) = 1.
2073     vpt(3) = yvtop
2074     vpt(4) = 1.
2075     call gxsvpt (vpt)
2076     window(3,1) = -1.
2077     window(4,1) = 1.
2078     call gxswnd (window)
2079     call peschm (nelmach, ieltyp, xax, estart, eend, actwin)
2080  endif
2081  goto 999
2082
2083120 continue
2084
2085  !--- curve for vert. var. missing
2086
2087150 continue
2088999 continue
2089end subroutine peplot
2090!***********************************************************************
2091
2092subroutine peschm (nel, ityp, hr, es, ee, actwin)
2093
2094  use plotfi
2095  implicit none
2096
2097  !----------------------------------------------------------------------*
2098  ! Purpose:
2099  !   Plot schema
2100  ! Input:
2101  !   nel      (integer)  no. of elements
2102  !   ityp     (integer)  array with element types:
2103  !                       0: drift                                       *
2104  !                       1: sbend, zero tilt                            *
2105  !                       2: focussing quad                              *
2106  !                       3: defocussing quad                            *
2107  !                       4: monitor                                     *
2108  !                       5: collimator                                  *
2109  !                       6: electrostatic separator                     *
2110  !                       7: sbend, non-zero tilt                        *
2111  !                       8: multipole                                   *
2112  !                       9: RF cavity                                   *
2113  !                       10: positive sext                              *
2114  !                       11: negative sext                              *
2115  !                       12: positive oct                               *
2116  !                       13: negative oct                               *
2117  !                       14: lcavity                                    *
2118  !                       21: rbend, zero tilt                           *
2119  !                       27: rbend, non-zero tilt                       *
2120  !   hr          (real)  horizontal range (lower and upper)
2121  !   es          (real)  array with element start position
2122  !   ee          (real)  array with element end position
2123  !   actwin      (real)  active window for curve plot (array of 4)
2124  !
2125  ! calls the routines jsln, gvpl defined in file gx11.F
2126  ! it is called by the routine peplot in this file.
2127  !----------------------------------------------------------------------*
2128
2129  !--- type definition of the routine arguments
2130
2131  integer nel, ityp(*)
2132  real hr(2), es(*), ee(*), actwin(4)
2133
2134  !--- type definitions of local variables
2135
2136  integer i, it, j, j_nodrift, im1
2137  integer npst(mobj), npnd(mobj), npsl(msize), i_nodrift(maxseql)
2138  real ell, shapex(msize), shapey(msize)
2139  real txp(2), typ(2), typz(2)
2140
2141  !--- Initialisation of local variables
2142
2143  data npst   / 1,  6, 11, 16, 21,                                  &
2144       33, 43, 48,                                                       &
2145       50,                                                               &
2146       64, 69, 74, 79, 84 /
2147  data npnd   / 5, 10, 15, 20, 32,                                  &
2148       42, 47, 49,                                                       &
2149       63,                                                               &
2150       68, 73, 78, 83, 88 /
2151  data npsl   /5 * 1, 5 * 1, 5 * 1, 5 * 3, 5 * 1, 0, 4 * 1, 0, 1,   &
2152       1, 1, 0, 1, 0, 1, 0, 1, 0, 1, 5 * 1, 2 * 1,                       &
2153       6 * 1, 0, 5 * 1, 0, 1,                                            &
2154       5 * 1, 5 * 1, 5 * 1, 5 * 1, 5 * 1 /
2155  data typz   / 2 * 0. /
2156  data shapex /0., 1., 1., 0., 0.,                                  &
2157       0., 1., 1., 0., 0.,                                               &
2158       0., 1., 1., 0., 0.,                                               &
2159       0., 1., 1., 0., 0.,                                               &
2160       0., 1., 1., 0., 0., 0., 1., 1., 0., 0., 0., 1.,                   &
2161       0., 1., 0.5, 0.5, 0., 1., 0.5, 0.5, 0., 1.,                       &
2162       0., 1., 1., 0., 0.,                                               &
2163       0., 0.,                                                           &
2164       0., 0.25, 0.25, 0.75, 0.75, 1.,                                   &
2165       0., 0.25, 0.25, 0.75, 0.75, 1., 0., 1.,                           &
2166       0., 1., 1., 0., 0.,                                               &
2167       0., 1., 1., 0., 0.,                                               &
2168       0., 1., 1., 0., 0.,                                               &
2169       0., 1., 1., 0., 0.,                                               &
2170       0., 1., 1., 0., 0. /
2171  data shapey /0.6, 0.6, -0.6, -0.6, 0.6,                           &
2172       0., 0., 0.8, 0.8, 0.,                                             &
2173       0., 0., -0.8, -0.8, 0.,                                           &
2174       0.6, 0.6, -0.6, -0.6, 0.6,                                        &
2175       0.8, 0.8, 0.4, 0.4, 0.8, -0.8, -0.8, -0.4, -0.4, -0.8, 0., 0.,    &
2176       0.4, 0.4, 0.8, 0.4, -0.4, -0.4, -0.8, -0.4, 0., 0.,               &
2177       0.5, 0.5, -0.5, -0.5, 0.5,                                        &
2178       0.5, -0.5,                                                        &
2179       0.2, 0.2, 0.8, 0.8, 0.2, 0.2,                                     &
2180       -0.2, -0.2, -0.8, -0.8, -0.2, -0.2, 0., 0.,                       &
2181       0., 0., 0.5, 0.5, 0.,                                             &
2182       0., 0., -0.5, -0.5, 0.,                                           &
2183       0., 0., 0.25, 0.25, 0.,                                           &
2184       0., 0., -0.25, -0.25, 0.,                                         &
2185       0.2, 0.2, -0.2, -0.2, 0.2 /
2186
2187  !--- Routine body
2188
2189  j_nodrift = 0
2190  im1 = 0
2191
2192  !--- set line style to solid
2193
2194  do i = 1, nel
2195     call jsln(1)
2196     it = mod(ityp(i), 20)
2197     if (it .eq. 0) goto 10
2198     j_nodrift = j_nodrift + 1
2199     i_nodrift(j_nodrift) = i
2200     if (j_nodrift .gt. 1) im1 = i_nodrift(j_nodrift-1)
2201     ell = ee(i) - es(i)
2202     if (j_nodrift .eq. 1) then
2203        if(es(i) .gt. hr(1))  then
2204           txp(1) = hr(1)
2205           txp(2) = es(i)
2206           call gvpl (2, txp, typz)
2207        endif
2208     else
2209        if (ee(im1) .lt. es(i))  then
2210           txp(1) = ee(im1)
2211           txp(2) = es(i)
2212           call gvpl (2, txp, typz)
2213        endif
2214     endif
2215     if (es(i) .gt. actwin(2)) goto 50
2216     if (ee(i) .ge. actwin(1)) then
2217        txp(1) = es(i) + shapex(npst(it)) * ell
2218        typ(1) = shapey(npst(it))
2219        do  j = npst(it)+1, npnd(it)
2220           txp(2) = es(i) + shapex(j) * ell
2221           typ(2) = shapey(j)
2222           if (npsl(j) .gt. 0)  then
2223              call jsln(npsl(j))
2224              call gvpl(2, txp, typ)
2225           endif
2226           txp(1) = txp(2)
2227           typ(1) = typ(2)
2228        enddo
2229     endif
223010   continue
2231  enddo
223250 continue
2233  call jsln(1)
2234  j = i_nodrift(j_nodrift)
2235  if (ee(j) .lt. hr(2))  then
2236     txp(1) = ee(j)
2237     txp(2) = hr(2)
2238     call gvpl (2, txp, typz)
2239  endif
2240end subroutine peschm
2241!***********************************************************************
2242
2243subroutine pesopt(ierr)
2244
2245  use plotfi
2246  use plot_bfi
2247  implicit none
2248
2249  !----------------------------------------------------------------------*
2250  ! Purpose:
2251  !   Stores plot options and values, checks
2252  !
2253  ! Output:  ierr  (int)     =0: OK, >0: error
2254  !
2255  ! calls the function plot_option and the routines comm_para,
2256  ! get_title, get_version and table_range defined in file madxn.c.
2257  ! calls the utility routine pesplit and the routine pegetn in this file.
2258  !
2259  ! it is called by the routine exec_plot in file madxn.c
2260  !----------------------------------------------------------------------*
2261
2262  integer ierr
2263
2264  integer i, j, k, notitle, noversi, nivaxs,inter_setplot
2265  character * (mcnam) sdum(1)
2266  integer nint, ndble, int_arr(szcompar), char_l(szcompar)
2267  integer plot_style(szcompar),plot_symbol(szcompar)
2268  double precision d_arr(szcompar)
2269  double precision plot_option
2270  character * (szchara) char_a, version
2271
2272  ierr = 0
2273
2274  !--- Initialisation of variables in common peaddi
2275
2276  itbv = 0
2277  nivvar = 0
2278  interf = 0
2279  noline = 0
2280  do i = 1, 4
2281     nvvar(i) = 0
2282  enddo
2283  do i = 1 , 2
2284     do j = 1 , mxcurv
2285        proc_flag(i,j)= 0
2286     enddo
2287  enddo
2288  do i = 1 , mpparm
2289     do j = 1 , mxcurv
2290        ipparm(i,j) = 0
2291     enddo
2292  enddo
2293  do i = 1, mxcurv
2294     naxref(i) = 0
2295  enddo
2296  nrrang(1) = 0
2297  nrrang(2) = 0
2298
2299  !--- Initialisation of variables in common peaddr
2300
2301  qascl = 0.0
2302  qlscl = 0.0
2303  qsscl = 0.0
2304  qtscl = 0.0
2305  hrange(1) = 0.0
2306  hrange(2) = 0.0
2307  hmima(1) = 1.e20
2308  hmima(2) = -1.e20
2309  do j = 1 , 4
2310     vrange(1,j) = 0.0
2311     vrange(2,j) = 0.0
2312     vmima(1,j) = 1.e20
2313     vmima(2,j) = -1.e20
2314  enddo
2315
2316  !--- Initialisation of variables in common peaddc
2317
2318  horname = ' '
2319  tabname = ' '
2320  toptitle = ' '
2321  do i = 1 , mxcurv
2322     sname(i) = ' '
2323  enddo
2324
2325  !--- Initialisation of variables in common peotcl
2326
2327  fpmach = .false.
2328  ptc_flag = .false.
2329
2330  !--- Initialisation of local variables
2331
2332  nivaxs = 0
2333  notitle = 0
2334  noversi = 0
2335  do i = 1, szcompar
2336     int_arr(i) = 0
2337     char_l(i) = 0
2338     d_arr(i)=0.0d0
2339  enddo
2340  char_a = ' '
2341  sdum(1) = ' '
2342
2343  !--- Routine body
2344
2345  !--- ptc flag setting
2346
2347  call comm_para('ptc ', nint, ndble, k, int_arr, d_arr,            &
2348       char_a, char_l)
2349  if (nint .gt. 0 .and. int_arr(1) .eq. 1) ptc_flag = .true.
2350
2351  !--- Get notitle
2352
2353  call comm_para('notitle ', nint, ndble, k, int_arr, d_arr,        &
2354       char_a, char_l)
2355  if (nint .gt. 0) notitle = int_arr(1)
2356
2357  !--- get noversi
2358
2359  call comm_para('noversion ', nint, ndble, k, int_arr, d_arr,      &
2360       char_a, char_l)
2361  if (nint .gt. 0) noversi = int_arr(1)
2362
2363  !--- if ptc flag is on look for the ptc_table
2364
2365  if(ptc_flag) then
2366     call comm_para('ptc_table ', nint, ndble, k, int_arr, d_arr,    &
2367          char_a, char_l)
2368     if (k .gt. 0) tabname = char_a
2369  else
2370
2371     !--- else normal twiss treatment : any table - for hor = s plot machine
2372
2373     call comm_para('table ', nint, ndble, k, int_arr, d_arr,        &
2374          char_a, char_l)
2375     if (k .gt. 0) tabname = char_a
2376  endif
2377
2378  !--- Horizontal variable
2379
2380  char_a = ' '
2381  call comm_para( 'haxis ', nint, ndble, k, int_arr, d_arr,         &
2382       char_a, char_l)
2383  if (k .eq. 0)  then
2384     print *, 'no horizontal variable'
2385     ierr = 1
2386     goto 999
2387  else
2388     horname = char_a
2389  endif
2390  itbv = 0
2391  if (horname .eq. 's')  itbv = 1
2392
2393  !--- Prepare title
2394
2395  if (notitle .eq. 0)  then
2396     char_a = ' '
2397     call comm_para('title ', nint, ndble, k, int_arr, d_arr,        &
2398          char_a, char_l)
2399     if (k .eq. 0) then
2400        call get_title(char_a, k)
2401     else
2402        k = char_l(1)
2403     endif
2404     if (noversi .eq. 0) then
2405        call get_version(version, j)
2406        if (k .gt. 0)  then
2407           toptitle = char_a(:k) // '<#>' // version(:j)
2408        else
2409           toptitle = '<#>' // version(:j)
2410        endif
2411     else
2412        if (k .gt. 0) toptitle = char_a(:k)
2413     endif
2414  endif
2415
2416  qascl = plot_option('ascale ')
2417  qlscl = plot_option('lscale ')
2418  qsscl = plot_option('sscale ')
2419  qtscl = plot_option('rscale ')
2420
2421  char_a = ' '
2422  call comm_para('range ', nint, ndble, k, int_arr, d_arr, char_a, char_l)
2423  call table_range(tabname, char_a, nrrang)
2424  if (nrrang(1) .eq. 0 .and. nrrang(2) .eq. 0)  then
2425     print *, 'unknown table or illegal range, skipped'
2426     ierr = 1
2427     goto 999
2428  endif
2429
2430  char_a = ' '
2431
2432  call comm_para('noline ', nint, ndble, k, int_arr, d_arr,         &
2433       char_a, char_l)
2434  if (nint .gt. 0) noline = int_arr(1)
2435
2436  call comm_para('hmin ', nint, ndble, k, int_arr, d_arr,           &
2437       char_a, char_l)
2438  if (ndble .gt. 0) hrange(1) = d_arr(1)
2439
2440  call comm_para('hmax ', nint, ndble, k, int_arr, d_arr,           &
2441       char_a, char_l)
2442  if (ndble .gt. 0) hrange(2) = d_arr(1)
2443
2444  call comm_para('vmin ', nint, ndble, k, int_arr, d_arr,           &
2445       char_a, char_l)
2446
2447  do i = 1, ndble
2448     vrange(1,i) = d_arr(i)
2449  enddo
2450
2451  call comm_para('vmax ', nint, ndble, k, int_arr, d_arr,           &
2452       char_a, char_l)
2453
2454  do i = 1, ndble
2455     vrange(2,i) = d_arr(i)
2456  enddo
2457
2458  !--- Check that STYLE & SYMBOL are both non zero
2459
2460  call comm_para('style ', nint, ndble, k, plot_style, d_arr,       &
2461       char_a, char_l)
2462  call comm_para('symbol ', nint, ndble, k, plot_symbol, d_arr,     &
2463       char_a, char_l)
2464  if (plot_style(1) + plot_symbol(1) .eq. 0) then
2465     print *,'Warning: style & symbol attributes will make plot invisible. Thus style is set to 1.'
2466     plot_style(1) = 1
2467  endif
2468  ipparm(1,1) = plot_style(1)
2469  ipparm(4,1) = plot_symbol(1)
2470
2471  char_a = ' '
2472  call comm_para('bars ', nint, ndble, k, ipparm(3,1), d_arr,       &
2473       char_a, char_l)
2474
2475  call comm_para('colour ', nint, ndble, k, ipparm(5,1), d_arr,     &
2476       char_a, char_l)
2477
2478  !--- if ptc_flag is on, no interpolation and check only ptc-related attributes
2479
2480  if (ptc_flag .and. itbv .eq. 0) goto 999
2481
2482  !--- Interpolation is not possible for ptc twiss variables
2483
2484  if (.not. ptc_flag) then
2485     call comm_para('spline ', nint,ndble,k,ipparm(2,1),d_arr,       &
2486          char_a,char_l)
2487     if (i .eq. 1) print *,'SPLINE attribute is obsolete, no action taken, use interpolate attribute instead.'
2488
2489     ipparm(2,1) = 0
2490     inter_setplot = plot_option('interpolate ')
2491     if (inter_setplot .eq. 0) then
2492        call comm_para('interpolate ', nint, ndble, k, ipparm(2,1),   &
2493             d_arr,char_a, char_l)
2494     else
2495        ipparm(2,1) = inter_setplot
2496     endif
2497  endif
2498
2499  !--- Continue fetching variables to be plotted
2500
2501  interf = ipparm(2,1)
2502
2503  char_a = ' '
2504  call comm_para('vaxis ', nint, ndble, k, int_arr, d_arr,          &
2505       char_a, char_l)
2506  if (k .gt. 0)  then
2507     nivaxs = 1
2508     nivvar = min(k, mxcurv)
2509     call pesplit(k, char_a, char_l, slabl)
2510     do j = 1, nivvar
2511        naxref(j) = 1
2512     enddo
2513  else
2514     char_a = ' '
2515     call comm_para('vaxis1 ', nint, ndble, k, int_arr, d_arr,       &
2516          char_a, char_l)
2517     if (k .gt. 0)  then
2518        if (nivvar+k .gt. mxcurv) then
2519           print *, 'Warning: # vertical variables cut at = ', nivvar
2520
2521           goto 110
2522        endif
2523        nivaxs = nivaxs + 1
2524        call pesplit(k, char_a, char_l, slabl(nivvar+1))
2525        do j = 1, k
2526           nivvar = nivvar + 1
2527           naxref(nivvar) = 1
2528        enddo
2529     endif
2530     char_a = ' '
2531     call comm_para('vaxis2 ', nint, ndble, k, int_arr, d_arr,       &
2532          char_a, char_l)
2533     if (k .gt. 0)  then
2534        if (nivvar+k .gt. mxcurv) then
2535           print *, 'Warning: # vertical variables cut at = ', nivvar
2536           goto 110
2537        endif
2538        nivaxs = nivaxs + 1
2539        call pesplit(k, char_a, char_l, slabl(nivvar+1))
2540        do j = 1, k
2541           nivvar = nivvar + 1
2542           naxref(nivvar) = 2
2543        enddo
2544     endif
2545     char_a = ' '
2546     call comm_para('vaxis3 ', nint, ndble, k, int_arr, d_arr,       &
2547          char_a, char_l)
2548     if (k .gt. 0)  then
2549        if (nivvar+k .gt. mxcurv) then
2550           print *, 'Warning: # vertical variables cut at = ', nivvar
2551           goto 110
2552        endif
2553        nivaxs = nivaxs + 1
2554        call pesplit(k, char_a, char_l, slabl(nivvar+1))
2555        do j = 1, k
2556           nivvar = nivvar + 1
2557           naxref(nivvar) = 3
2558        enddo
2559     endif
2560     char_a = ' '
2561     call comm_para('vaxis4 ', nint, ndble, k, int_arr, d_arr,       &
2562          char_a, char_l)
2563     if (k .gt. 0)  then
2564        if (nivvar+k .gt. mxcurv) then
2565           print *, 'Warning: # vertical variables cut at = ', nivvar
2566           goto 110
2567        endif
2568        nivaxs = nivaxs + 1
2569        call pesplit(k, char_a, char_l, slabl(nivvar+1))
2570        do j = 1, k
2571           nivvar = nivvar + 1
2572           naxref(nivvar) = 4
2573        enddo
2574     endif
2575  endif
2576  if (nivvar .eq. 0)  then
2577     print *, 'Warning: no vertical plot variables, plot skipped'
2578     ierr = 1
2579     goto 999
2580  endif
2581  do i = 2, mxcurv
2582     ipparm(1,i) = ipparm(1,1)
2583     ipparm(2,i) = ipparm(2,1)
2584     ipparm(3,i) = ipparm(3,1)
2585     ipparm(4,i) = ipparm(4,1)
2586     ipparm(5,i) = ipparm(5,1)
2587  enddo
2588110 continue
2589
2590  do j = 1, nivvar
2591     call pegetn (0, slabl(j), itbv, proc_flag(1,j), sname(j),       &
2592          sdum(1))
2593     if (slabl(j) .eq. 'rbetx')  then
2594        sname(j) = 'betx'
2595        proc_flag(1,j) = 1
2596     else if (slabl(j) .eq. 'rbety')  then
2597        sname(j) = 'bety'
2598        proc_flag(1,j) = 1
2599     else
2600        sname(j) = slabl(j)
2601        proc_flag(1,j) = 0
2602     endif
2603  enddo
2604999 end subroutine pesopt
2605  !***********************************************************************
2606
2607subroutine pesplit(n_str, char_a, char_l, char_buff)
2608
2609  implicit none
2610
2611  !----------------------------------------------------------------------*
2612  !
2613  !   Utility routine
2614  !   Purpose: splits a string in several sub-strings
2615  !
2616  !--- Input
2617  !   n_str      number of sub-strings
2618  !   char_l     number of characters of each sub-string
2619  !   char_a     character string
2620  !--- Output
2621  !   cahr_buf   sub_strings
2622  !
2623  !----------------------------------------------------------------------*
2624
2625  !--- type definition of the routine arguments
2626
2627  integer n_str, char_l(*)
2628  character*(*) char_a, char_buff(*)
2629
2630  !--- type definitions of local variables
2631
2632  integer i, k, l
2633
2634  !--- Initialisation of local variables
2635
2636  k = 0
2637
2638  !--- Routine body
2639
2640  do i = 1, n_str
2641     l = char_l(i)
2642     char_buff(i) = char_a(k+1:k+l)
2643     k = k+l
2644  enddo
2645end subroutine pesplit
2646!***********************************************************************
2647
2648subroutine plginit
2649
2650  use plotfi
2651  use plot_bfi
2652  implicit none
2653
2654  !----------------------------------------------------------------------*
2655  ! Purpose:
2656  !   Overall initialization
2657  !
2658  ! calls the function plot_option and the routine comm_para
2659  ! defined in file madxn.c.
2660  ! calls the function intrac defined in file madxu.c.
2661  ! calls the routines gxtint, gxsvar, gxasku, gxinit, gxclos
2662  ! defined in file gxx11.F
2663  !
2664  ! it is called by the routine plotit in this file.
2665  !----------------------------------------------------------------------*
2666
2667  !--- type definitions of local variables
2668
2669  integer ipseps, iset, nint, ndble, k, int_arr(100), char_l(100)
2670  double precision d_arr(100)
2671  real tmpval
2672  character * 40 char_a
2673  character * 8 tmp_a
2674
2675  !--- definitions of function primitives
2676
2677  double precision plot_option
2678  logical intrac
2679
2680  !--- Initialisation of local variables
2681
2682  data iset / 0 /
2683
2684  !--- Routine body
2685
2686  call gxtint
2687  call gxsvar ('INUNIT', 5, 0., ' ')
2688  call gxsvar ('IOUNIT', 6, 0., ' ')
2689  char_a = ' '
2690  tmp_a = 'file '
2691  call comm_para(tmp_a, nint, ndble, k, int_arr, d_arr,             &
2692       char_a, char_l)
2693  if (k .gt. 0) then
2694     plfnam = char_a(:char_l(1))
2695  else
2696     plfnam = 'madx'
2697  endif
2698  ipseps = plot_option('post ')
2699  if (ipseps .eq. 0 .and. .not. intrac())  then
2700     ipseps = 2
2701  endif
2702  if (iset .eq. 0 .and. ipseps .ne. 0) then
2703     iset = 1
2704     call gxsvar ('SMETNM', 0, 0., plfnam)
2705     call gxsvar('IPSEPS', ipseps, 0., ' ')
2706  endif
2707  if (intrac())  then
2708
2709     !--- set wait time to 1 sec.
2710
2711     call gxsvar ('WTTIME', 0, 1., ' ')
2712     call gxasku
2713  endif
2714
2715  !--- reduce window size (only X11)
2716
2717  call gxsvar('NYPIX', 670, 0., ' ')
2718
2719  !--- set bounding box (only X11)
2720
2721  tmpval=plot_option('xsize ')
2722  call gxsvar('XMETAF', 0, tmpval, ' ')
2723  tmpval=plot_option('ysize ')
2724  call gxsvar('YMETAF', 0, tmpval, ' ')
2725
2726  !--- inhibit initial X-Window (only X11)
2727
2728  call gxsvar('ITSEOP', 1, 0., ' ')
2729  call gxinit
2730  call gxclos
2731end subroutine plginit
2732!***********************************************************************
2733
2734subroutine plotit(initfl)
2735
2736  use plotfi
2737  use plot_bfi
2738  implicit none
2739
2740  !----------------------------------------------------------------------*
2741  ! Purpose:
2742  !   Plots on screen and/or file
2743  !
2744  ! calls the routines gxsvar, gxterm, gxinit, gxopen, gxwait, gxclrw and
2745  ! gxclos in file gxx11.F.
2746  ! calls  the routines plginit and peplot in this file.
2747  !
2748  ! it is called by the routine exec_plot in file madxn.c after pemima.
2749  !----------------------------------------------------------------------*
2750
2751  !--- type definition of the routine arguments
2752
2753  integer initfl
2754
2755  !--- type definitions of local variables
2756
2757  character * (mfile) plpnam
2758  integer plot_No
2759
2760  !--- Initialisation of local variables
2761
2762  save plpnam
2763  save plot_No
2764
2765  !--- Routine body
2766  !
2767  if (initfl .eq. 0)  then
2768
2769     !--- overall initialization
2770     plot_No = 0
2771     call plginit
2772     plpnam = plfnam
2773  endif
2774  plot_No = plot_No + 1
2775  print *,"plot number = ",plot_No
2776  if (plpnam .ne. plfnam)  then
2777     call gxsvar ('SMETNM', 0, 0., plfnam)
2778
2779     !--- close current .ps file if any
2780
2781     call gxterm
2782     plpnam = plfnam
2783     call gxinit
2784  endif
2785
2786  call gxopen
2787  call peplot
2788  call gxwait
2789  call gxclrw
2790  call gxclos
2791end subroutine plotit
2792!***********************************************************************
2793
2794subroutine pupnbl(string,ifirst,ilast)
2795
2796  implicit none
2797
2798  !----------------------------------------------------------------------*
2799  !
2800  !   Utility routine
2801  !   Purpose: returns position of first and last non-blank in STRING
2802  !
2803  !--- Input
2804  !   string     character string
2805  !--- Output
2806  !   ifirst     first non-blank in string, or 0 if only blanks
2807  !   ilast      last non-blank
2808  !
2809  !   Author: H. Grote / CERN                        date: June 16, 1987
2810  !                                             last mod: Sept. 13, 2001
2811  !
2812  !----------------------------------------------------------------------*
2813
2814  !--- type definition of the routine arguments
2815
2816  character *(*)  string
2817  integer ifirst, ilast
2818
2819  !--- type definitions of local variables
2820
2821  integer i
2822
2823  !--- Output initialisation
2824
2825  ifirst=0
2826  ilast=0
2827
2828  !--- Routine body
2829
2830  do i=1,len(string)
2831     if(string(i:i).ne.' ') then
2832        ifirst=i
2833        goto 20
2834     endif
2835  enddo
2836  goto 999
283720 continue
2838  do i=len(string),1,-1
2839     if(string(i:i).ne.' ') then
2840        ilast=i
2841        goto 999
2842     endif
2843  enddo
2844999 end subroutine pupnbl
2845  !***********************************************************************
2846
2847integer function iucomp(comp, arr, n)
2848
2849  implicit none
2850
2851  !----------------------------------------------------------------------*
2852  ! Utility function
2853  ! Purpose:
2854  !   Find first occurrence of integer in integer array
2855  !
2856  !---Input:
2857  !   comp        integer being looked up
2858  !   arr         integer array being searched
2859  !   n           length of arr
2860  !
2861  !  returns 0 if not found, else position in arr
2862  !----------------------------------------------------------------------*
2863
2864  !--- type definition of the routine arguments
2865
2866  integer comp, arr(*), n
2867
2868  !--- type definitions of local variables
2869
2870  integer j
2871
2872  !--- Output initialisation
2873
2874  iucomp = 0
2875
2876  !--- Routine body
2877
2878  do j = 1, n
2879     if (comp .eq. arr(j))  then
2880        iucomp = j
2881        return
2882     endif
2883  enddo
2884end function iucomp
Note: See TracBrowser for help on using the repository browser.