IAP GITLAB

Skip to content
Snippets Groups Projects
make22.f 81.1 KiB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000
c$Id: make22.f,v 1.36 2003/05/02 11:06:54 weber Exp $
C####C##1#########2#########3#########4#########5#########6#########7##
      subroutine make22(iio,e,ii1,iiz1,mm1,xfac1,ii2,iiz2,mm2,xfac2)
c
cinput iio   : label for exit-channel
cinput e     : $\sqrt{s}$ of process
cinput ii1   : ID of incoming particle 1
cinput iiz1  : $2\cdot I_3$ of incoming particle 1
cinput mm1   : mass of incoming particle 1
cinput xfac1 : scaling factor for preformed hadron
cinput ii2   : ID of incoming particle 2
cinput iiz2  : $2\cdot I_3$ of incoming particle 2
cinput mm2   : mass of incoming particle 2
cinput xfac2 : scaling factor for preformed hadron
c
c  output:   exit channel via common-blocks in {\tt newpart.f}
c
c {\tt make22} generates the final state for all scatterings and
c decays. Due to the diverse nature of the interactions handled
c many special cases have to be taken care of. The label {\tt iio}
c matches in most cases the respective label in subroutine {\tt crossx},
c which returns the respective partial cross sections.
c
C####C##1#########2#########3#########4#########5#########6#########7##
      implicit none

      include 'coms.f'
      include 'options.f'
      include 'newpart.f'
      include 'comres.f'



      integer io,iio,i1,i2,i3,i4,i,k,i1p,i1m,i2p,i2m
      integer ifqrk1,ifqrk2,ifdiq1,ifdiq2,ifqrk3,ifqrk4,ifdiq3,ifdiq4
      integer iq1(2),iq2(3),kq1,kq2,kqq1,kqq2,iii,iff(3)
      real*8 sig,e,m1,m2,m3,m4,gam,mm1,mm2,m1m,m2m,xfac1,xfac2

      integer iz1,iz2,iz3,iz4,errflg,icnt,ntry,ib1,ib2,i1old,i2old

      logical bit


c...functions
      integer isoit,whichres
      real*8 getmass,fmsr,ranf,pcms,massit,widit,mminit

c...vacuum quantumnumber(s) for special string decay (don't touch)
      real*8 valint(1)
      common /values/ valint

c...string
      real*8 b1,b2,ba1(3),ba2(3)
      integer j,l,ii1,ii2,iiz1,iiz2,iexopt,iddum,ibar,jbar
      logical fboost,switips



      integer ident(2,mprt)
      real*8 part(9,mprt),ms1,ms2,msmin1,msmin2,tau,esum


      bit=.true.
      switips=.false.
      io=mod(iio,200)
      i1=ii1
      i2=ii2
      iz1=iiz1
      iz2=iiz2
      m1=mm1
      m2=mm2
      icnt=0
      ntry=0
      ibar=0

      if(iabs(i1).lt.minmes)ibar=ibar+isign(1,i1)
      if(iabs(i2).lt.minmes)ibar=ibar+isign(1,i2)

c in case of a MB-reaction, sort particles (but keep track of
c any id-switch with the 'switips'-flag)
      if(iabs(i1).ge.minmes.and.iabs(i1).le.maxmes.and.
     &   iabs(i2).ge.minbar.and.iabs(i2).le.maxbar)then
        call swpizm(i1,iz1,m1,i2,iz2,m2)
        switips=.true.
      endif

      if(i1+i2.eq.0.and.iz1+iz2.eq.0.and.CTOption(20).ne.0.and.
     .    io.gt.20)goto 27 !e+e-

      if(io.lt.0)goto(100,100,100,100,100,100,100,29)-io

c      if(i1+i2.gt.2)write(6,*)'make22:',i1,i2
ctp060202 1007 continue
      goto(1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,9,17,9,17,20,
     ,9,13,23,14,12,26,27,15,14,100,29,100,14,15,14,36,36,13)io

      write(6,*)'make22: unknown channel requested io:',io
      write(6,*)'  ',e,i1,iz1,m1,i2,iz2,m2
      stop

 1    continue
c...pp->ND
      i3=minnuc
      i4=mindel
      m3=massit(i3)

      call getmas(massit(i4),widit(i4),i4,isoit(i4),mminit(i4),
     .     e-m3,m3,m4)

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 1008

 2    continue
c...pp->pp*
      i3=minnuc
      m3=massit(i3)
      if(bit)call getres(io,e,minnuc+1,maxnuc,i4)
      call getmas(massit(i4),widit(i4),i4,isoit(i4),mminit(i4),
     .     e-m3,m3,m4)

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 1008

 3    continue
c...pp->ND*
      i3=minnuc
      m3=massit(i3)
      if(bit)call getres(io,e,mindel+1,maxdel,i4)
      call getmas(massit(i4),widit(i4),i4,isoit(i4),mminit(i4),
     .     e-m3,m3,m4)

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 1008

 4    continue
c...pp->DD
      i3=mindel
      i4=mindel
      call getmas(massit(i3),widit(i3),i3,isoit(i3),
     .     mminit(i4),e-mminit(i4),mminit(i4),m3)
      call getmas(massit(i4),widit(i4),i4,isoit(i4),mminit(i4),
     .     e-m3,m3,m4)

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 1008

 5    continue
c...pp->DN*, DN* with factor 4/3
      i3=mindel
      if(bit)call getres(io,e,minnuc+1,maxnuc,i4)
      call getmas(massit(i3),widit(i3),i3,isoit(i3),
     .     mminit(i3),e-mminit(i4),mminit(i4),m3)
      call getmas(massit(i4),widit(i4),i4,isoit(i4),mminit(i4),
     .     e-m3,m3,m4)

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 1008

 6    continue
c...pp->DD*, DN* with factor 4/3
      i3=mindel
      if(bit)call getres(io,e,mindel+1,maxdel,i4)
      call getmas(massit(i3),widit(i3),i3,isoit(i3),
     .     mminit(i3),e-mminit(i4),mminit(i4),m3)
      call getmas(massit(i4),widit(i4),i4,isoit(i4),mminit(i4),
     .     e-m3,m3,m4)

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 1008

 7    continue
c...pp -> generate N*N*,N*D*,D*D*
      m3=fmsr(mminit(mindel),e-mresmin)
      m4=fmsr(mminit(mindel),e-m3)
      if(bit)i3=whichres(m3,3)
      if(bit)i4=whichres(m4,3)
      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 1008

 8    continue
c...ND->DD
      i3=mindel
      i4=mindel
      call getmas(massit(i3),widit(i3),i3,isoit(i3),
     .     mminit(i4),e-mminit(i4),mminit(i4),m3)
      call getmas(massit(i4),widit(i4),i4,isoit(i4),mminit(i4),
     .     e-m3,m3,m4)

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 1008

 9    continue
      write(6,*)'make22: channel no.',io,'not implemented.'
      stop

 10   continue
c...MB->B',MM->M* annihilations
      call anndec(0,m1,i1,iz1,m2,i2,iz2,e,sig,gam)

      goto 2002
c      return

 11   continue
c...MM->M'
      call anndec(0,m1,i1,iz1,m2,i2,iz2,e,sig,gam)

      goto 2002
c      return


 12   continue
      write(6,*)'make22: channel no.',io,'should correspond to',
     ,'total cross section, i.e. make22 should not be called.'
      stop

 13   continue
c...elastic scattering
c
      if(switips)then
         call swpizm(i1,iz1,m1,i2,iz2,m2)
         switips=.false.
      endif

      call setizm(i1,iz1,m1,i2,iz2,m2,i3,iz3,m3,i4,iz4,m4)

      if(mminit(i4)+mminit(i3).gt.e)then
ctp060926         write(6,*)'make22(el):threshold violated'
ctp060926         write(6,*)'m3:',m3,mminit(i3)
ctp060926         write(6,*)'m4:',m4,mminit(i4)
      else
c
c
        if(m3.lt.mminit(i3))m3=mminit(i3)
        if(m4.lt.mminit(i4))m4=mminit(i4)
      end if

      goto 2001

 14   continue
c...inelastic scattering (aqm for nonstrange resonances)
c   for arbitrary resonances getinw should be modified
c   arbitrary particles are excited (whithout charge exchange)

c...get minimal masses
      call getirg(i1,i1m,i1p)
      m1m=mminit(i1m)
      call getirg(i2,i2m,i2p)
      m2m=mminit(i2m)

      if(e-m1m-m2m.le.0d0)then
c...elastic scattering
        call setizm(i1,iz1,m1,i2,iz2,m2,i3,iz3,m3,i4,iz4,m4)
        goto 2001
c        return
      end if

c...masses
      if(ranf(0).lt.5d-1)then
         m3=fmsr(m1m,e-m2m)
         m4=fmsr(m2m,e-m3)
      else
         m4=fmsr(m2m,e-m1m)
         m3=fmsr(m1m,e-m4)
      end if

      if(e-m3-m4.lt.0d0)goto 13

c...itype3
      if(i1m.lt.i1p.and.widit(i1m).lt.1d-3.and.
     .    m3.le.massit(i1m+1)-0.5*widit(i1m+1))then
c...lowest itype of this kind of particles is stable
c & mass .lt. approximate minimal mass of lowest itype + 1
        m3=massit(i1m)
        i3=i1m
      else if(i1m.eq.i1p.and.widit(i1m).lt.1d-3) then
c...class with only one narrow particle
        i3=i1
        m3=massit(i3)
      else
c...itypes of this kind are all unstable
        call whichi(i3,i1m,i1p,m3)
      end if

c...itype4
      if(i2m.lt.i2p.and.widit(i2m).lt.1d-3.and.
     .    m4.le.massit(i2m+1)-0.5*widit(i2m+1))then
c...lowest itype of this kind of particles is stable
        m4=massit(i2m)
        i4=i2m
      else if(i2m.eq.i2p.and.widit(i2m).lt.1d-3) then
c...class with only one narrow particle
        i4=i2
        m4=massit(i4)
      else
c...itypes of this kind are all unstable
        call whichi(i4,i2m,i2p,m4)
      end if

c...no charge transfer
      iz3=iz1
      iz4=iz2
      i3=isign(i3,ii1)
      i4=isign(i4,ii2)
      goto 2001
c      return

 15    continue
c...XX -> 2 strings
      if(CTOption(12).ne.0)then
        write(6,*)' *** error(make22): string section is called ',
     .            'while strings are switched off:CTOption(12).ne.0'
        stop
      end if

c
      if(switips)then
         call swpizm(i1,iz1,m1,i2,iz2,m2)
         switips=.false.
      endif


c store old itypes
      i1old=i1
      i2old=i2

 155  continue

c allow for deexcitation: 'excitation' starts from groundstate

      if(iabs(i1).ge.minmes)then
c.. meson resonances may also be deexcitated: assign the particle
c id's of the lowest multiplet (with the same quark content)
        call ityp2id(i1,iz1,iq1(1),iq1(2))
        if(iabs(iq1(1)).gt.iabs(iq1(2)))then
          iddum=iq1(1)
          iq1(1)=iq1(2)
          iq1(2)=iddum
        endif
        iddum=isign(100*iabs(iq1(1))+10*iabs(iq1(2)),iq1(1))
        call id2ityp(iddum,0d0,i1,iz1)
      else
         call getirg(i1,i1m,i1p)
      endif

c same for second particle
      if(iabs(i2).ge.minmes)then
        call ityp2id(i2,iz2,iq1(1),iq1(2))
        if(iabs(iq1(1)).gt.iabs(iq1(2)))then
          iddum=iq1(1)
          iq1(1)=iq1(2)
          iq1(2)=iddum
        endif
        iddum=isign(100*iabs(iq1(1))+10*iabs(iq1(2)),iq1(1))
        call id2ityp(iddum,0d0,i2,iz2)
      else
         call getirg(i2,i2m,i2p)
      endif

c
c BEWARE: now i1 and i2 do not contain anymore the ityps of the ingoing
c particles, but the ityps of the lowest possible states (groundstates).
c This is needed in order for the string-excitation to be able to
c excite ALL states and not only those above the state of the ingoing
c particle.
c If you need the old ityps (i.e. for elastic scattering) then reset them
c via i1old and i2old


c if the incoming masses are changed due to excitation, then the
c new masses should not be lower than:
      m1m=mminit(i1)
      m2m=mminit(i2)

c..discriminate between BB and MB collisions:
      ib1=1
      ib2=1
      if(iabs(i1).ge.minmes)ib1=0
      if(iabs(i2).ge.minmes)ib2=0


c set minimum energy for the two strings
      msmin1=m1m+CTParam(2)
      msmin2=m2m+CTParam(2)

c no (meson) string below 1 gev:
      if(msmin1.lt.1d0)msmin1=1d0
      if(msmin2.lt.1d0)msmin2=1d0

c if energy too low: do elastic collision
      if(m1m+m2m+CTParam(34).ge.e)then
ctp060926        write(*,*)'make22: not enough energy for string exc. ->elastic/
ctp060926     &deexcitation'
ctp060926        write(*,*)' i1, i2, m1, m2, e: ',i1,i2,m1,m2,e
         i1=i1old
         i2=i2old
        goto 13
      endif


c convert to quark-IDs
         call ityp2id(i1,iz1,ifdiq1,ifqrk1)
         call ityp2id(i2,iz2,ifdiq2,ifqrk2)

         iexopt=CTOption(22)
 81      continue
c 100 tries for excitation, otherwise elastic scattering
         ntry=ntry+1
         if(ntry.gt.100)then
ctp060926          write(*,*)'make22: too many tries for string exc. ->elastic/
ctp060926     &           deexcitation'
ctp060926          write(*,*)' i1, i2, m1, m2, e: ',i1,i2,m1,m2,e
         i1=i1old
         i2=i2old
          goto 13
         endif

c string-excitation:
c get string masses ms1,ms2 and the leading quarks
         call STREXCT(IFdiq1,IFqrk1,ib1,M1m,
     &           ifdiq2,ifqrk2,ib2,M2m,E,
     &           iexopt,
     &           ba1,ms1,ba2,ms2,
     &           ifdiq3,ifqrk3,ifdiq4,ifqrk4)
c the boost parameters are now fixed for the masses ms1, ms2. If the
c masses will be changed, set the parameter fboost to "false":
      fboost=.true.

c accept deexcitation of one of the hadrons:
      if(ms1.le.msmin1.and.ms2.ge.msmin2)then
        ms1=massit(i1)
        fboost=.false.
      else if(ms2.le.msmin2.and.ms1.ge.msmin1)then
        ms2=massit(i2)
        fboost=.false.
c don't accept elastic-like (both masses too low):
      else if(ms1.lt.msmin1.and.ms2.lt.msmin2)then
        goto 81
      endif

c in case of deexcitation new masses are necessary
c single diffractive, mass excitation according 1/m
      if(ms1.le.msmin1)then
        ms2=fmsr(msmin2,e-ms1)
        fboost=.false.
      elseif(ms2.le.msmin2)then
        ms1=fmsr(msmin1,e-ms2)
        fboost=.false.
      endif

c quark-quark scattering -> only elastic !
      if(xfac1.lt..999d0.and.xfac2.lt..999d0
     &     .and.ranf(0).lt..25) then
         i1=i1old
         i2=i2old
         goto 13
      endif

c single diffractive if one particle is a quark state,
c mass excitation according 1/m
      if(xfac1.lt..999d0.and.ranf(0).lt..5)then
        ms1=massit(i1)
        ms2=fmsr(msmin2,e-ms1)
        fboost=.false.
      elseif(xfac2.lt..999d0.and.ranf(0).lt..5)then
        ms2=massit(i2)
        ms1=fmsr(msmin1,e-ms2)
        fboost=.false.
      endif

c take care that the particles will be able to decay lateron:
      if(ms1.lt.m1m)then
        ms1=m1m
        fboost=.false.
      endif
      if(ms2.lt.m2m)then
        ms2=m2m
        fboost=.false.
      endif

c avoid energy conservation violation
      if(ms1+ms2.gt.e)goto 155

       if(CTOption(22).ne.1.or..not.fboost)then
c the boost parameters have to be calculated:
         b1=2.*e*pcms(e,ms1,ms2)/(e**2+ms1**2-ms2**2)
         b2=2.*e*pcms(e,ms1,ms2)/(e**2-ms1**2+ms2**2)
         ba1(3)=+b1
         ba2(3)=-b2
         do 151 j=1,2
            ba1(j)=0d0
 151        ba2(j)=0d0
      end if

c now write information to newpart common-blocks
      mstring(1)=ms1
      mstring(2)=ms2

c string#1
      if(ms1.le.msmin1)then
        nstring1=1
        l=1
        do j=1,3
          part(j,l)=0d0
          part(j+5,l)=0d0
        enddo
        part(4,l)=ms1
        part(4+5,l)=0d0
        part(5,l)=ms1
        ident(1,l)=i1
        ident(2,l)=iz1
      else
       call qstring(ifdiq3,ifqrk3,ms1,part,ident,nstring1)
       if(nstring1.eq.0) goto 155
      end if

      esum=0d0

      jbar=0

      do l=1,nstring1
         pnew(5,l)=part(5,l)
         itypnew(l)=ident(1,l)

         if(iabs(ident(1,l)).lt.minmes)jbar=jbar+isign(1,ident(1,l))

         i3new(l)=ident(2,l)
         do j=1,4
           pnew(j,l)=part(j,l)
           xnew(j,l)=part(j+5,l)
         enddo
         pnew(3,l)=part(3,l)
         xnew(3,l)=part(3+5,l)
         call rotbos(0d0,0d0,ba1(1),ba1(2),ba1(3),
     ,    pnew(1,l),pnew(2,l),pnew(3,l),pnew(4,l))
         call rotbos(0d0,0d0,ba1(1),ba1(2),ba1(3),
     ,    xnew(1,l),xnew(2,l),xnew(3,l),xnew(4,l))
         esum=esum+pnew(4,l)
      enddo
      call leadhad(1,nstring1,1)


c string #2
      if(ms2.le.msmin2) then
        nstring2=1
        l=1
        do j=1,3
          part(j,l)=0d0
          part(j+5,l)=0d0
        enddo
        part(4,l)=ms2
        part(4+5,l)=0d0
        part(5,l)=ms2
        ident(1,l)=i2
        ident(2,l)=iz2
      else
       call qstring(ifdiq4,ifqrk4,ms2,part,ident,nstring2)
       if(nstring2.eq.0) goto 155
      end if

      esum=0d0
       do l=1,nstring2
        pnew(5,nstring1+l)=part(5,l)
        itypnew(nstring1+l)=ident(1,l)

        if(iabs(ident(1,l)).lt.minmes)jbar=jbar+isign(1,ident(1,l))

        i3new(nstring1+l)=ident(2,l)
        do j=1,4
          pnew(j,nstring1+l)=part(j,l)
          xnew(j,nstring1+l)=part(j+5,l)
        enddo
        pnew(3,nstring1+l)=-pnew(3,nstring1+l)
        xnew(3,nstring1+l)=-xnew(3,nstring1+l)
        call rotbos(0d0,0d0,ba2(1),ba2(2),ba2(3),
     ,    pnew(1,nstring1+l),pnew(2,nstring1+l),pnew(3,nstring1+l),
     ,     pnew(4,nstring1+l))
        call rotbos(0d0,0d0,ba2(1),ba2(2),ba2(3),
     ,    xnew(1,nstring1+l),xnew(2,nstring1+l),xnew(3,nstring1+l),
     ,     xnew(4,nstring1+l))
        esum=esum+pnew(4,nstring1+l)
      enddo
      call leadhad(nstring1+1,nstring1+nstring2,1)

c error check
ctp060926      if(ibar.ne.jbar)then
ctp060926         write(6,*)' *** (E) no baryon number conservation', ibar,jbar
ctp060926         write(6,*)'     ',i1,i2,ms1,ms2
ctp060926         write(6,'(5i4)')(itypnew(l),l=1,nstring1+nstring2)
ctp060926         end if


      return


ctp060202 718  format(i2,i4,i3,1x,10(f10.4,1x))



 17   continue

      iz3=iz1
      iz4=iz2
      i3=i1
      i4=i2


      m3=m1
      m4=m2

c the following lines MUST be there in order to set nucleons on-shell
c after their first collision
      if(m3.lt.mminit(i3))m3=mminit(i3)
      if(m4.lt.mminit(i4))m4=mminit(i4)


      goto 2001

 20   continue
c...decays
      if(ityptd(1,pslot(1)).eq.0) then
c normal decay
c note: m4,i4 and iz4 are dummies in this call
         call anndec(1,m1,i1,iz1,m4,i4,iz4,e,sig,gam)
      else
c forward time-delay
         pnew(5,1)=fmasstd(1,pslot(1))
         itypnew(1)=ityptd(1,pslot(1))
         i3new(1)=iso3td(1,pslot(1))
         pnew(5,2)=fmasstd(2,pslot(1))
         itypnew(2)=ityptd(2,pslot(1))
         i3new(2)=iso3td(2,pslot(1))
      endif

      if(nexit.eq.2) then
         i3=itypnew(1)
         iz3=i3new(1)
         m3=pnew(5,1)
         i4=itypnew(2)
         iz4=i3new(2)
         m4=pnew(5,2)
         goto 2001
      else
c three or four body decay
         nstring1=1
         nstring2=nexit-1
         do i=1,4
            do j=1,nexit
               pnew(i,j)=0d0
               xnew(i,j)=0d0
            end do
         end do
         mstring(1)=pnew(5,1)
c
         mstring(2)=pnew(5,2)
         do 91 j=3,nexit
            mstring(2)=mstring(2)+pnew(5,j)
 91      continue
c
c now call routine for momentum phase space...
         call nbodydec(e)
         return
      endif

 23   continue
c...annihilation -> string
      if(CTOption(12).ne.0)then
        write(6,*)' *** error(make22): string section is called ',
     .            'while strings are switched off:CTOption(12).ne.0'
        stop
      end if

      ms1=e/2.
      ms2=e/2.
      mstring(1)=ms1
      mstring(2)=ms2

c determine flavour content of b-bbar-system
      call ityp2id(i1,iz1,ifdiq1,ifqrk1)
      call ityp2id(i2,iz2,ifdiq2,ifqrk2)
c...create string 1 out of quark-antiquark pair
      call qstring(ifqrk1,ifqrk2,ms1,part,ident,nstring1)
      esum=0d0

      do k=1,nstring1
        l=k
        do j=1,4
          pnew(j,l)=part(j,l)
          xnew(j,l)=part(j+5,l)
        enddo
        esum=esum+pnew(4,l)
        pnew(5,l)=part(5,l)
        itypnew(l)=ident(1,l)
        i3new(l)=ident(2,l)
        tau=part(9,l)/ (part(4,l)/part(5,l))
      enddo
      call leadhad(1,nstring1,0)

c...create string 2 out of diquark-antidiquark-pair
c   use one quark and one antiquark for the string-ends:
       ifqrk1=int(ifdiq1/1000)
       ifqrk2=int(ifdiq2/1000)
c...store remaining flavour quantum numbers in ctp(26), they will
c   be passed to the 'clustr'-routine:
       ifdiq1=mod(ifdiq1/100,10)
       ifdiq2=mod(ifdiq2/100,10)
       valint(1)=dble(((abs(ifdiq1*10.d0)+abs(ifdiq2*1.d0))
     &            *isign(1,ifdiq1)))
       valint(1)=sign(valint(1),dble(ifdiq1))
       call qstring(ifqrk1,ifqrk2,ms2,part,ident,nstring2)
       valint(1)=0.d0
      esum=0d0

      do l=1,nstring2
        do j=1,4
          pnew(j,nstring1+l)=part(j,l)
          xnew(j,nstring1+l)=part(j+5,l)
        enddo
        esum=esum+pnew(4,nstring1+l)
        pnew(5,nstring1+l)=part(5,l)
        itypnew(nstring1+l)=ident(1,l)
        i3new(nstring1+l)=ident(2,l)
        tau=part(9,nstring1+l)/ (part(4,l)/part(5,l))
      enddo
      call leadhad(nstring1+1,nstring1+nstring2,0)

      return

 26   continue
c...elastic MB scattering (the outgoing particle id's must not be
c   assigned randomly like at label 2001)
c
      if(switips)then
         call swpizm(i1,iz1,m1,i2,iz2,m2)
         switips=.false.
      endif

      call setizm(i1,iz1,m1,i2,iz2,m2,i3,iz3,m3,i4,iz4,m4)

      if(mminit(i4)+mminit(i3).gt.e)then
ctp060926         write(6,*)'make22(el):threshold violated'
ctp060926         write(6,*)'m3:',m3,mminit(i3)
ctp060926         write(6,*)'m4:',m4,mminit(i4)
      else
        if(m3.lt.mminit(i3))m3=mminit(i3)
        if(m4.lt.mminit(i4))m4=mminit(i4)
      end if

c... get momenta & fill newpart, 2 particle exit-channel

      nstring1=1
      nstring2=1
      nexit=2
      do i=1,4
         do j=1,2
            pnew(i,j)=0d0
            xnew(i,j)=0d0
         end do
      end do

c...boost to 2-particle cms

      pnew(3,1)=pcms(e,m3,m4)
      pnew(3,2)=-pcms(e,m3,m4)

      pnew(4,1)=sqrt(m3**2+pnew(3,1)**2)
      pnew(4,2)=sqrt(m4**2+pnew(3,2)**2)

      pnew(5,1)=m3
      mstring(1)=m3
      itypnew(1)=i3
      i3new(1)=iz3

      pnew(5,2)=m4
      mstring(2)=m4
      itypnew(2)=i4
      i3new(2)=iz4

      return

 27   continue
c XX-> 1 string : e+e- , MB
      if(CTOption(12).ne.0)then
        write(6,*)' *** error(make22): string section is called ',
     .            'while strings are switched off:CTOption(12).ne.0'
        stop
      end if

c quark-quark scattering -> only elastic !
      if(xfac1.lt..999d0.and.xfac2.lt..999d0
     &   .and.ranf(0).lt.0.25) goto 26

      ms1=e
      mstring(1)=ms1
      mstring(2)=0d0

c determine flavour content of string
      if(CTOption(20).eq.1)then
c..e+e- annihilation
      if(ranf(0).lt.CTParam(6))then
        ifqrk1=3  ! ssbar
        ifqrk2=-3
      else
        call ityp2id(104,0,ifqrk1,ifqrk2)  ! qqbar
      end if
      else
c...MB annihilation. the quark content must be known:
      call ityp2id(i2,iz2,iq1(1),iq1(2))
      call ityp2id(i1,iz1,ifdiq2,iq2(3))

      if(abs(i1).ge.minmes) then
         iq2(1)=ifdiq2
         iq2(2)=iq2(3)
         iq2(3)=0
      else
         iq2(1)=mod(ifdiq2/100,10)
         iq2(2)=int(ifdiq2/1000)
      endif

      do 312 kq1=1,2
       do 412 kq2=1,3
c..two of the quarks must be able to annihilate:
        if(iq1(kq1)+iq2(kq2).eq.0) then
         kqq1=kq1
         kqq2=kq2
         goto 414
        endif
 412   continue
 312  continue
      goto 26 ! could not create double charged strange baryon string
 414  continue
c.. the 'iff'-quarks constitute the produced (anti-)baron
      iff(1)=iq1(3-kqq1)
      iii=1
      do 512 kq2=1,3
       if (kq2.ne.kqq2) then
        iii=iii+1
        iff(iii)=iq2(kq2)
       endif
 512  continue

      if(abs(i1).lt.minmes) then
         call mquarks(iff,ifqrk1,ifqrk2)
      else
         ifqrk1=iff(1)
         ifqrk2=iff(2)
      endif

      endif

c...create string 1 out of quark-antiquark pair
      call qstring(ifqrk1,ifqrk2,ms1,part,ident,nstring1)
      esum=0d0

c primitive bug-fix:
      if(nstring1.eq.0)then
ctp060926        write(6,*)'make22: iline 27 not completed. ->elastic'
        goto 26
      endif

      do 101 k=1,nstring1
        l=k
c no leading hadron in e+e-
        if(CTOption(20).eq.1)then
          leadfac(l)=0.0d0
        endif
        do 102 j=1,4
          pnew(j,l)=part(j,l)
          xnew(j,l)=part(j+5,l)
 102    continue
        esum=esum+pnew(4,l)
        pnew(5,l)=part(5,l)
        itypnew(l)=ident(1,l)
        i3new(l)=ident(2,l)
        tau=part(9,l)/ (part(4,l)/part(5,l))
 101    continue
      if(CTOption(20).eq.1)then
        call leadhad(1,nstring1,3)
      else
        call leadhad(1,nstring1,1)
      endif

      nstring2=0
      return


 29   continue
c...DD->ND detailed balance
      i3=minnuc
      i4=mindel
      m3=massit(i3)
      m4=getmass(e-m3,0)
      if(iabs(iz1+iz2).gt.isoit(i3)+isoit(i4))then
       iz3=-9
       iz4=-9
      else
         nexit=2
         itot(1)=isoit(i3)
         itot(2)=isoit(i4)
         call isocgk4(isoit(i1),iz1,isoit(i2),iz2,itot,i3new,errflg)
         i3=isign(i3,ii1)
         i4=isign(i4,ii2)
         iz3=i3new(1)
         iz4=i3new(2)
      end if

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 2001


 100  continue
c...??->NN detailed balance inverse channels
c iso3 are assigned in detbal
      nexit=2
      i3=minnuc
      i4=minnuc
      m3=massit(minnuc)
      m4=massit(minnuc)
      itot(1)=isoit(i3)
      itot(2)=isoit(i4)
      call isocgk4(isoit(i1),iz1,isoit(i2),iz2,itot,i3new,errflg)
      i3=isign(i3,ii1)
      i4=isign(i4,ii2)
      iz3=i3new(1)
      iz4=i3new(2)

      if(ranf(0).gt.0.5d0)call swpizm(i3,iz3,m3,i4,iz4,m4)
      goto 2001


 36   continue
c...MB->B',MM->M* annihilations (forward time delay)
      call anndec(0,m1,i1,iz1,m2,i2,iz2,e,sig,gam)
      goto 2003

 1008 continue
c...get isospin-3 components

      nexit=2
      itot(1)=isoit(i3)
      itot(2)=isoit(i4)
      call isocgk4(isoit(i1),iz1,isoit(i2),iz2,itot,i3new,errflg)
      iz3=i3new(1)
      iz4=i3new(2)

ctp060926      if(errflg.ne.0)then
ctp060926        write(6,*)'make22: iso-spin conservation ',
ctp060926     ,            'not possible in isocgk: error-flag=',errflg
ctp060926        write(6,*)'      ',isoit(i1),iz1,isoit(i2),iz2,'>',
ctp060926     ,            isoit(i3),isoit(i4),iz3,iz4,
ctp060926     ,            ' process:',e,i1,m1,i2,m2,'>',i3,m3,i4,m4,'io=',io
ctp060926      end if
      i3=isign(i3,ii1)
      i4=isign(i4,ii2)

 2001 continue


c... get momenta & fill newpart, 2 particle exit-channel

      nstring1=1
      nstring2=1
      nexit=2
      do i=1,4
         do j=1,2
            pnew(i,j)=0d0
            xnew(i,j)=0d0
         end do
      end do

c...boost to 2-particle cms


      if(.not.(ityptd(1,pslot(1)).ne.0.and.CTOption(34).eq.2.and.
     &     iline.eq.20)) then
c normal decay

         pnew(3,1)=pcms(e,m3,m4)
         pnew(3,2)=-pcms(e,m3,m4)

         pnew(4,1)=sqrt(m3**2+pnew(3,1)**2)
         pnew(4,2)=sqrt(m4**2+pnew(3,2)**2)
      else
c forward time delay
         pnew(1,1)=pxtd(1,pslot(1))
         pnew(1,2)=pxtd(2,pslot(1))
         pnew(2,1)=pytd(1,pslot(1))
         pnew(2,2)=pytd(2,pslot(1))
         pnew(3,1)=pztd(1,pslot(1))
         pnew(3,2)=pztd(2,pslot(1))
         pnew(4,1)=p0td(1,pslot(1))
         pnew(4,2)=p0td(2,pslot(1))
      endif

      pnew(5,1)=m3
      mstring(1)=m3
      itypnew(1)=i3
      i3new(1)=iz3

      pnew(5,2)=m4
      mstring(2)=m4