File: RHadrons.cc

package info (click to toggle)
pythia8 8.1.65-1
  • links: PTS, VCS
  • area: main
  • in suites: wheezy
  • size: 22,660 kB
  • sloc: cpp: 59,593; xml: 30,509; php: 6,649; sh: 796; makefile: 353; ansic: 33
file content (1437 lines) | stat: -rw-r--r-- 50,399 bytes parent folder | download
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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
// RHadrons.cc is a part of the PYTHIA event generator.
// Copyright (C) 2012 Torbjorn Sjostrand.
// PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
// Please respect the MCnet Guidelines, see GUIDELINES for details.

// Function definitions (not found in the header) for the RHadrons class.

#include "RHadrons.h"

namespace Pythia8 {
 
//==========================================================================

// The RHadrons class.

//--------------------------------------------------------------------------

// Constants: could be changed here if desired, but normally should not.
// These are of technical nature, as described for each.

const int RHadrons::IDRHADSB[14] = {  1000512, 1000522, 1000532,
  1000542, 1000552, 1005113, 1005211, 1005213, 1005223, 1005311,
  1005313, 1005321, 1005323, 1005333 };

const int RHadrons::IDRHADST[14] = {  1000612, 1000622, 1000632,
  1000642, 1000652, 1006113, 1006211, 1006213, 1006223, 1006311,
  1006313, 1006321, 1006323, 1006333 };

// Gluino code and list of gluino R-hadron codes.
const int RHadrons::IDRHADGO[38] = {  1000993, 1009113, 1009213, 
  1009223, 1009313, 1009323, 1009333, 1009413, 1009423, 1009433, 
  1009443, 1009513, 1009523, 1009533, 1009543, 1009553, 1091114, 
  1092114, 1092214, 1092224, 1093114, 1093214, 1093224, 1093314, 
  1093324, 1093334, 1094114, 1094214, 1094224, 1094314, 1094324, 
  1094334, 1095114, 1095214, 1095224, 1095314, 1095324, 1095334 };

// Allow a few tries for flavour selection since failure is allowed.
const int RHadrons::NTRYMAX = 10;

// Safety margin (in GeV) when constructing kinematics of system.
const double RHadrons::MSAFETY = 0.1;

// Maximal energy to borrow for gluon to insert on leg in to junction.
const double RHadrons::EGBORROWMAX = 4.;

//--------------------------------------------------------------------------

// Main routine to initialize R-hadron handling.

bool RHadrons::init( Info* infoPtrIn, Settings& settings, 
  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn) {

  // Store input pointers for future use. 
  infoPtr          = infoPtrIn;
  particleDataPtr  = particleDataPtrIn;  
  rndmPtr          = rndmPtrIn;

  // Flags and parameters related to R-hadron formation and decay.
  allowRH          = settings.flag("RHadrons:allow");
  maxWidthRH       = settings.parm("RHadrons:maxWidth");
  idRSb            = settings.mode("RHadrons:idSbottom");
  idRSt            = settings.mode("RHadrons:idStop");
  idRGo            = settings.mode("RHadrons:idGluino");
  setMassesRH      = settings.flag("RHadrons:setMasses");
  probGluinoballRH = settings.parm("RHadrons:probGluinoball");
  mOffsetCloudRH   = settings.parm("RHadrons:mOffsetCloud");
  mCollapseRH      = settings.parm("RHadrons:mCollapse");
  diquarkSpin1RH   = settings.parm("RHadrons:diquarkSpin1"); 

  // Check whether sbottom, stop or gluino should form R-hadrons. 
  allowRSb         = allowRH && idRSb > 0 
    && (particleDataPtr->mWidth(idRSb) < maxWidthRH);
  allowRSt         = allowRH && idRSt > 0 
    && (particleDataPtr->mWidth(idRSt) < maxWidthRH);
  allowRGo         = allowRH && idRGo > 0 
    && (particleDataPtr->mWidth(idRGo) < maxWidthRH);
  allowSomeR       = allowRSb || allowRSt || allowRGo;

  // Set nominal masses of sbottom R-mesons and R-baryons.
  if (allowRSb) {
    m0Sb = particleDataPtr->m0(idRSb);
    if (setMassesRH) {
      for (int i = 0; i < 14; ++i) {
        int idR = IDRHADSB[i]; 
        double m0RHad = m0Sb + mOffsetCloudRH;
        m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
        if (i > 4) 
        m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
        particleDataPtr->m0( idR, m0RHad);    
      }
    }

    // Set widths and lifetimes of sbottom R-hadrons.
    double mWidthRHad = particleDataPtr->mWidth(idRSb);
    double tau0RHad   = particleDataPtr->tau0(  idRSb);
    for (int i = 0; i < 14; ++i) {
      particleDataPtr->mWidth( IDRHADSB[i], mWidthRHad);
      particleDataPtr->tau0(   IDRHADSB[i],   tau0RHad);
    }
  }

  // Set nominal masses of stop R-mesons and R-baryons.
  if (allowRSt) {
    m0St = particleDataPtr->m0(idRSt);
    if (setMassesRH) {
      for (int i = 0; i < 14; ++i) {
        int idR = IDRHADST[i]; 
        double m0RHad = m0St + mOffsetCloudRH;
        m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
        if (i > 4) 
        m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
        particleDataPtr->m0( idR, m0RHad);    
      }
    }

    // Set widths and lifetimes of stop R-hadrons.
    double mWidthRHad = particleDataPtr->mWidth(idRSt);
    double tau0RHad   = particleDataPtr->tau0(  idRSt);
    for (int i = 0; i < 14; ++i) {
      particleDataPtr->mWidth( IDRHADST[i], mWidthRHad);
      particleDataPtr->tau0(   IDRHADST[i],   tau0RHad);
    }
  }

  // Set nominal masses of gluino R-glueballs, R-mesons and R-baryons.
  if (allowRGo) {
    m0Go = particleDataPtr->m0(idRGo);
    if (setMassesRH) {
      particleDataPtr->m0( IDRHADGO[0], m0Go + 2. * mOffsetCloudRH 
        + particleDataPtr->constituentMass(21) );
      for (int i = 1; i < 38; ++i) {
        int idR = IDRHADGO[i]; 
        double m0RHad = m0Go + 2. * mOffsetCloudRH;
        m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
        m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
        if (i > 15) 
        m0RHad += particleDataPtr->constituentMass( (idR%10000)/1000 );
        particleDataPtr->m0( idR, m0RHad);    
      }
    }

    // Set widths and lifetimes of gluino R-hadrons.
    double mWidthRHad = particleDataPtr->mWidth(idRGo);
    double tau0RHad   = particleDataPtr->tau0(  idRGo);
    for (int i = 0; i < 38; ++i) {
      particleDataPtr->mWidth( IDRHADGO[i], mWidthRHad);
      particleDataPtr->tau0(   IDRHADGO[i],   tau0RHad);
    }
  }   

  // Done. 
  return true;

}
//--------------------------------------------------------------------------

// Check if a given particle can produce R-hadrons.

bool RHadrons::givesRHadron( int id) {
  if (allowRSb && abs(id) == idRSb) return true;
  if (allowRSt && abs(id) == idRSt) return true;
  if (allowRGo && id == idRGo) return true;
  return false;

}

//--------------------------------------------------------------------------

// Produce R-hadrons by fragmenting them off from existing strings.

bool RHadrons::produce( ColConfig& colConfig, Event& event) {

  // Check whether some sparticles are allowed to hadronize.
  if (!allowSomeR) return true;

  // Reset arrays and values.
  iBefRHad.resize(0);
  iCreRHad.resize(0);
  iRHadron.resize(0);
  iAftRHad.resize(0);
  isTriplet.resize(0);
  nRHad = 0;

  // Loop over event and identify hadronizing sparticles.
  for (int i = 0; i < event.size(); ++i) 
   if (event[i].isFinal() && givesRHadron(event[i].id())) { 
    iBefRHad.push_back(i);
    iCreRHad.push_back(i);
    iRHadron.push_back(0);
    iAftRHad.push_back(0);
    isTriplet.push_back(true);
  } 
  nRHad = iRHadron.size();
  
  // Done if no hadronizing sparticles.
  if (nRHad == 0) return true;

  // Max two R-hadrons. Randomize order of processing.
  if (nRHad > 2) {
     infoPtr->errorMsg("Error in RHadrons::produce: "
       "cannot handle more than two R-hadrons");
     return false;
  }
  if (nRHad > 1 && rndmPtr->flat() > 0.5) swap( iBefRHad[0], iBefRHad[1]);

  // Split a system with both a sparticle and a junction.
  iBef = iBefRHad[0];  
  iSys = colConfig.findSinglet( iBef);
  systemPtr = &colConfig[iSys];
  if (systemPtr->hasJunction && !splitOffJunction( colConfig, event)) {
    infoPtr->errorMsg("Error in RHadrons::produce: "
      "cannot handle system with junction");
    return false;
  }
  if (nRHad == 2) {
    iBef = iBefRHad[1];  
    iSys = colConfig.findSinglet( iBefRHad[1]);
    systemPtr = &colConfig[iSys];
    if (systemPtr->hasJunction && !splitOffJunction( colConfig, event)) {
      infoPtr->errorMsg("Error in RHadrons::produce: "
        "cannot handle system with junction");
      return false;
    }
  }

  // Open up a closed gluon/gluino loop.
  iBef = iBefRHad[0];  
  iSys = colConfig.findSinglet( iBef);
  systemPtr = &colConfig[iSys];
  if (systemPtr->isClosed && !openClosedLoop( colConfig, event)) {
    infoPtr->errorMsg("Error in RHadrons::produce: "
      "cannot open up closed gluon/gluino loop");
    return false;
  }
  if (nRHad == 2) {
    iBef = iBefRHad[1];  
    iSys = colConfig.findSinglet( iBefRHad[1]);
    systemPtr = &colConfig[iSys];
    if (systemPtr->isClosed && !openClosedLoop( colConfig, event)) {
      infoPtr->errorMsg("Error in RHadrons::produce: "
        "cannot open up closed gluon/gluino loop");
      return false;
    }
  }

  // Split up a colour singlet system that should give two R-hadrons.
  if (nRHad == 2) {
    int iSys1 = colConfig.findSinglet( iBefRHad[0]);
    int iSys2 = colConfig.findSinglet( iBefRHad[1]);
    if (iSys2 == iSys1) { 
      iSys = iSys1;
      systemPtr = &colConfig[iSys];
      if ( !splitSystem( colConfig, event) ) { 
        infoPtr->errorMsg("Error in RHadrons::produce: "
          "failed to handle two sparticles in same system");
        return false;
      }
    } 
  }
    
  // Loop over R-hadrons to be formed. Find its colour singlet system.
  for (iRHad = 0; iRHad < nRHad; ++iRHad) {
    iBef = iBefRHad[iRHad];  
    iSys = colConfig.findSinglet( iBef);
    if (iSys < 0) {
      infoPtr->errorMsg("Error in RHadrons::produce: "
        "sparticle not in any colour singlet");
      return false;
    }
    systemPtr = &colConfig[iSys];

    // For now don't handle systems involving junctions or loops.
    if (systemPtr->hasJunction) {
      infoPtr->errorMsg("Error in RHadrons::produce: "
        "cannot handle system with junction");
      return false;
    }
    if (systemPtr->isClosed) {
      infoPtr->errorMsg("Error in RHadrons::produce: "
        "cannot handle closed colour loop");
      return false;
    }

    // Handle formation of R-hadron separately for gluino and squark.
    if (event[iBef].id() == idRGo) isTriplet[iRHad] = false;
    bool formed = (isTriplet[iRHad]) ? produceSquark( colConfig, event)
                                     : produceGluino( colConfig, event);
    if (!formed) return false;

  // End of loop over R-hadrons. Done.
  } 
  return true;

}

//--------------------------------------------------------------------------

// Decay R-hadrons by resolving them into string systems and letting the
// heavy unstable particle decay as normal.

bool RHadrons::decay( Event& event) {

  // Loop over R-hadrons to decay. 
  for (iRHad = 0; iRHad < nRHad; ++iRHad) {
    int    iRNow  = iRHadron[iRHad]; 
    int    iRBef  = iBefRHad[iRHad];
    int    idRHad = event[iRNow].id();
    double mRHad  = event[iRNow].m();
    double mRBef  = event[iRBef].m();
    int    iR0    = 0;
    int    iR2    = 0; 

    // Find flavour content of squark or gluino R-hadron.
    pair<int,int> idPair = (isTriplet[iRHad]) 
      ? fromIdWithSquark( idRHad) : fromIdWithGluino( idRHad);
    int id1 = idPair.first;
    int id2 = idPair.second;

    // Sharing of momentum: the squark/gluino should be restored
    // to original mass, but error if negative-mass spectators.
    double fracR = mRBef / mRHad;
    if (fracR >= 1.) {
      infoPtr->errorMsg("Error in RHadrons::decay: "
          "too low R-hadron mass for decay");
      return false;
    }

    // Squark: new colour needed in the breakup.
    if (isTriplet[iRHad]) {
      int colNew = event.nextColTag();
      int col    = (event[iRBef].col() != 0) ? colNew : 0;
      int acol   = (col == 0) ? colNew : 0; 
      
      // Store the constituents of a squark R-hadron.
      iR0 = event.append( id1, 106, iRNow, 0, 0, 0, col, acol,
        fracR * event[iRNow].p(), fracR * mRHad, 0.);
      iR2 = event.append( id2, 106, iRNow, 0, 0, 0, acol, col, 
        (1. - fracR) * event[iRNow].p(), (1. - fracR) * mRHad, 0.);

    // Gluino: set mass sharing between two spectators.
    } else {
      double m1Eff  = particleDataPtr->constituentMass(id1) + mOffsetCloudRH;
      double m2Eff  = particleDataPtr->constituentMass(id2) + mOffsetCloudRH;
      double frac1 = (1. - fracR) * m1Eff / ( m1Eff + m2Eff); 
      double frac2 = (1. - fracR) * m2Eff / ( m1Eff + m2Eff); 
   
      // Two new colours needed in the breakups.
      int col1 = event.nextColTag();
      int col2 = event.nextColTag();

      // Store the constituents of a gluino R-hadron.
      iR0 = event.append( idRGo, 106, iRNow, 0, 0, 0, col2, col1,
        fracR * event[iRNow].p(), fracR * mRHad, 0.);
      event.append( id1, 106, iRNow, 0, 0, 0, col1, 0, 
        frac1 * event[iRNow].p(), frac1 * mRHad, 0.);
      iR2 = event.append( id2, 106, iRNow, 0, 0, 0, 0, col2, 
        frac2 * event[iRNow].p(), frac2 * mRHad, 0.);
    }

    // Mark R-hadron as decayed and update history.
    event[iRNow].statusNeg();
    event[iRNow].daughters( iR0, iR2);
    iAftRHad[iRHad] = iR0;

    // Set secondary vertex for decay products, but no lifetime.
    Vec4 vDec = event[iRNow].vProd() + event[iRNow].tau()
              * event[iR0].p() / event[iR0].m();
    for (int iRd = iR0; iRd <= iR2; ++iRd) event[iRd].vProd( vDec);

  // End loop over R-hadron decays, based on velocity of squark.
  
  }

  // Done.
  return true;

}

//--------------------------------------------------------------------------

// Split a system that contains both a sparticle and a junction. 

bool RHadrons::splitOffJunction( ColConfig& colConfig, Event& event) {

  // Classify system into three legs, and find sparticle location.
  vector<int> leg1, leg2, leg3;
  int iLegSP = 0;
  int iPosSP = 0;
  int iLeg = 0;
  int iPos = 0;
  for (int i = 0; i < systemPtr->size(); ++i) {
    ++iPos;
    int iP = systemPtr->iParton[i];
    if (iP < 0) {
      ++iLeg;
      iPos = -1;
    } else if (iLeg == 1) leg1.push_back( iP);
    else if   (iLeg == 2) leg2.push_back( iP);
    else if   (iLeg == 3) leg3.push_back( iP);
    if (iP == iBef) {
      iLegSP = iLeg;
      iPosSP = iPos;
    }
  }
  if (iLegSP == 0) return false;
  
  // Swap so leg 1 contains sparticle. If not innermost sparticle then
  // skip for now, and wait for this other (gluino!) to be split off. 
  if      (iLegSP == 2) swap(leg2, leg1);
  else if (iLegSP == 3) swap(leg3, leg1); 
  for (int i = 0; i < iPosSP; ++i)
    if (event[leg1[i]].id() != 21) return true;
 
  // No gluon between sparticle and junction: find kinetic energy of system.
  if (iPosSP == 0) {
    Vec4 pSP  = event[iBef].p();
    Vec4 pRec = 0.;
    for (int i = 0; i < int(leg2.size()); ++i) pRec += event[leg2[i]].p();
    for (int i = 0; i < int(leg3.size()); ++i) pRec += event[leg3[i]].p();
    double mSys  = (pSP + pRec).mCalc();
    double mSP   = pSP.mCalc();
    double mRec  = pRec.mCalc();
    double eKin  = mSys - mSP - mRec;

    // Insert new gluon that borrows part of kinetic energy.
    double mNewG  = EGBORROWMAX * eKin / (EGBORROWMAX + eKin) ;
    Vec4   pNewG  = (mNewG / mSys) * (pSP + pRec);
    int    newCol = event.nextColTag();
    bool   isCol  = (event[leg1.back()].col() > 0);
    int    colNG  = (isCol)? newCol :  event[iBef].acol();
    int    acolNG = (isCol) ? event[iBef].col() : newCol;
    int    iNewG  = event.append( 21, 101, iBef, 0, 0, 0, colNG, acolNG, 
      pNewG, mNewG, 0.);
    leg1.insert( leg1.begin(), iNewG);
    ++iPosSP;

    // Boosts for remainder systems that gave up energy.
    double mNewSys = mSys - mNewG;
    double pAbsOld = 0.5 * sqrtpos( pow2(mSys*mSys - mSP*mSP - mRec*mRec)
                   - pow2(2. * mSP * mRec) ) / mSys;
    double pAbsNew = 0.5 * sqrtpos( pow2(mNewSys*mNewSys - mSP*mSP - mRec*mRec)
                   - pow2(2. * mSP * mRec) ) / mNewSys;
    RotBstMatrix MtoCM, MfromCM, MSP, MRec;
    MtoCM.toCMframe( pSP, pRec);
    MfromCM = MtoCM;
    MfromCM.invert(); 
    MSP = MtoCM;
    MSP.bst( 0., 0., -pAbsOld / sqrt(mSP * mSP + pAbsOld * pAbsOld));
    MSP.bst( 0., 0.,  pAbsNew / sqrt(mSP * mSP + pAbsNew * pAbsNew)); 
    MSP.rotbst( MfromCM);
    MRec = MtoCM;
    MRec.bst( 0., 0.,  pAbsOld / sqrt(mRec * mRec + pAbsOld * pAbsOld));
    MRec.bst( 0., 0., -pAbsNew / sqrt(mRec * mRec + pAbsNew * pAbsNew)); 
    MRec.rotbst( MfromCM);

    // Copy down recoling partons and boost their momenta.
    int iNewSP  = event.copy( iBef, 101);
    event[iNewSP].rotbst( MSP);
    leg1[iPosSP]   = iNewSP;
    if (iBefRHad[0] == iBef) iBefRHad[0] = iNewSP;
    else if (nRHad > 1 && iBefRHad[1] == iBef) iBefRHad[1] = iNewSP;
    iBef = iNewSP;
    for (int i = 0; i < int(leg2.size()); ++i) {
      int iCopy = event.copy( leg2[i], 101);  
      event[iCopy].rotbst( MRec);
      if (iBefRHad[0] == leg2[i]) iBefRHad[0] = iCopy;
      else if (nRHad > 1 && iBefRHad[1] == leg2[i]) iBefRHad[1] = iCopy;
      leg2[i] = iCopy;
    }   
    for (int i = 0; i < int(leg3.size()); ++i) {
      int iCopy = event.copy( leg3[i], 101);  
      event[iCopy].rotbst( MRec);
      if (iBefRHad[0] == leg3[i]) iBefRHad[0] = iCopy;
      else if (nRHad > 1 && iBefRHad[1] == leg3[i]) iBefRHad[1] = iCopy;
      leg3[i]   = iCopy;
    }
 
  // Now always at least one gluon between sparticle and junction.  
  }

  // Find gluon with largest energy in sparticle rest frame.
  int    iGspl = 0;
  double eGspl = event[leg1[0]].p() * event[iBef].p();
  for (int i = 1; i < iPosSP; ++i) {
    double eGnow = event[leg1[i]].p() * event[iBef].p();
    if (eGnow > eGspl) {
      iGspl = i;
      eGspl = eGnow;
    }
  }
  int iG = leg1[iGspl];
   
  // Split this gluon into a collinear quark.antiquark pair.
  int idNewQ = flavSelPtr->pickLightQ(); 
  int iNewQ  = event.append(  idNewQ, 101, iG, 0, 0, 0, event[iG].col(), 0, 
    0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
  int iNewQb = event.append( -idNewQ, 101, iG, 0, 0, 0, 0, event[iG].acol(), 
    0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
  event[iG].statusNeg();
  event[iG].daughters( iNewQ, iNewQb); 
  if (event[leg1.back()].col() == 0) swap( iNewQ, iNewQb);

  // Collect two new systems after split.
  vector<int> iNewSys1, iNewSys2;
  iNewSys1.push_back( iNewQb);
  for (int i = iGspl + 1; i < int(leg1.size()); ++i)
    iNewSys1.push_back( leg1[i]);
  iNewSys2.push_back( -10);
  for (int i = 0; i < iGspl; ++i) iNewSys2.push_back( leg1[i]);
  iNewSys2.push_back( iNewQ);
  iNewSys2.push_back( -11);
  for (int i = 0; i < int(leg2.size()); ++i) iNewSys2.push_back( leg2[i]);
  iNewSys2.push_back( -12);
  for (int i = 0; i < int(leg3.size()); ++i) iNewSys2.push_back( leg3[i]);

  // Remove old system and insert two new ones.
  colConfig.erase(iSys);
  colConfig.insert( iNewSys1, event);
  colConfig.insert( iNewSys2, event);   

  // Done. 
  return true;

}

//--------------------------------------------------------------------------

// Open up a closed gluon/gluino loop.
  
bool RHadrons::openClosedLoop( ColConfig& colConfig, Event& event) {

  // Find gluon with largest energy in gluino rest frame.
  int    iGspl = -1;
  double eGspl = 0.;
  for (int i = 0; i < systemPtr->size(); ++i) {
    int  iGNow = systemPtr->iParton[i];
    if (event[iGNow].id() == 21) {
      double eGnow = event[iGNow].p() * event[iBef].p();
      if (eGnow > eGspl) {
        iGspl = i;
        eGspl = eGnow;
      }
    }
  }
  if (iGspl == -1) return false;
   
  // Split this gluon into a collinear quark.antiquark pair.
  int iG     = systemPtr->iParton[iGspl];
  int idNewQ = flavSelPtr->pickLightQ(); 
  int iNewQ  = event.append(  idNewQ, 101, iG, 0, 0, 0, event[iG].col(), 0, 
    0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
  int iNewQb = event.append( -idNewQ, 101, iG, 0, 0, 0, 0, event[iG].acol(), 
    0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
  event[iG].statusNeg();
  event[iG].daughters( iNewQ, iNewQb); 
   
  // Order partons in new system. Note order of colour flow.
  int iNext = iGspl + 1;
  if (iNext == systemPtr->size()) iNext = 0; 
  if (event[ systemPtr->iParton[iNext]].acol() != event[iNewQ].col())
    swap( iNewQ, iNewQb);        
  vector<int> iNewSys;
  iNewSys.push_back( iNewQ);
  for (int i = iGspl + 1; i < systemPtr->size(); ++i)
    iNewSys.push_back( systemPtr->iParton[i]);
  for (int i = 0; i < iGspl; ++i)
    iNewSys.push_back( systemPtr->iParton[i]);
  iNewSys.push_back( iNewQb);  

  // Erase the old system and insert the new one instead.
  colConfig.erase(iSys);
  colConfig.insert( iNewSys, event);

  // Done. 
  return true;

}

//--------------------------------------------------------------------------

// Split a single colour singlet that contains two sparticles.
// To fix : if nLeg >= 3 && mMin large handle as in nLeg == 1??

bool RHadrons::splitSystem( ColConfig& colConfig, Event& event) {

  // First and second R-hadron mother. Number of legs in between.
  int iFirst  = -1;
  int iSecond = -1;
  for (int i = 0; i < int(systemPtr->size()); ++i) {
    int  iTmp   = systemPtr->iParton[i];
    if ( givesRHadron( event[iTmp].id()) ) { 
      if (iFirst == -1) iFirst  = i;
      else              iSecond = i;
    }
  }
  int nLeg = iSecond - iFirst;

  // New flavour pair for breaking the string, and its mass.
  int    idNewQ = flavSelPtr->pickLightQ();
  double mNewQ  = particleDataPtr->constituentMass( idNewQ);
  vector<int> iNewSys1, iNewSys2;

  // If sparticles are neighbours then need new q-qbar pair in between.
  if (nLeg == 1) {

    // The location of the two sparticles.
    int i1Old = systemPtr->iParton[iFirst];
    int i2Old = systemPtr->iParton[iSecond];

    // Take a fraction of momentum to give breakup quark pair.
    double mHat = (event[i1Old].p() + event[i2Old].p()).mCalc();
    double mMax = mHat - event[i1Old].m() - event[i2Old].m(); 
    if (mMax < 2. * (mNewQ + MSAFETY)) return false;
    double mEff = min( 2. * (mNewQ + mOffsetCloudRH), mMax - 2. * MSAFETY);
    double frac = mEff / mHat;
    Vec4   pEff = frac * (event[i1Old].p() + event[i2Old].p());
  
    // New kinematics by (1) go to same mHat with bigger masses, and 
    // (2) rescale back down to original masses and reduced mHat.
    Vec4 p1New, p2New;
    if ( !newKin( event[i1Old].p(), event[i2Old].p(), 
      event[i1Old].m() / (1. - frac), event[i2Old].m() / (1. - frac), 
      p1New, p2New) ) return false; 
    p1New *= 1. - frac;
    p2New *= 1. - frac;

    // Fill in new partons after branching, with correct colour/flavour sign.
    int i1New, i2New, i3New, i4New;
    int newCol = event.nextColTag();
    i1New = event.copy( i1Old, 101);
    if (event[i2Old].acol() == event[i1Old].col()) {
      i3New = event.append( -idNewQ, 101, i1Old, 0, 0, 0, 
        0, event[i2Old].acol(), 0.5 * pEff, 0.5 * mEff, 0.);
      i2New = event.copy( i2Old, 101);
      event[i2New].acol( newCol);
      i4New = event.append(  idNewQ, 101, i2Old, 0, 0, 0, 
        newCol, 0, 0.5 * pEff, 0.5 * mEff, 0.); 
    } else {
      i3New = event.append(  idNewQ, 101, i1Old, 0, 0, 0, 
        event[i2Old].col(), 0, 0.5 * pEff, 0.5 * mEff, 0.);
      i2New = event.copy( i2Old, 101);
      event[i2New].col( newCol);
      i4New = event.append( -idNewQ, 101, i2Old, 0, 0, 0, 
        0, newCol, 0.5 * pEff, 0.5 * mEff, 0.); 
    }

    // Modify momenta and history. For iBotCopyId tracing asymmetric 
    // bookkeeping where one sparticle is radiator and the other recoiler.
    event[i1New].p( p1New);
    event[i2New].p( p2New);
    event[i1Old].daughters( i1New, i3New);
    event[i1New].mother2( 0);
    event[i2Old].daughters( i2New, i4New);
    event[i2New].mother2( 0);
    iBefRHad[0] = i1New;
    iBefRHad[1] = i2New;
 
    // Partons in the two new systems.
    for (int i = 0; i < iFirst; ++i) 
      iNewSys1.push_back( systemPtr->iParton[i]);
    iNewSys1.push_back( i1New);
    iNewSys1.push_back( i3New);
    iNewSys2.push_back( i4New);
    iNewSys2.push_back( i2New);
    for (int i = iSecond + 1; i < int(systemPtr->size()); ++i) 
      iNewSys2.push_back( systemPtr->iParton[i]);

  // If one gluon between sparticles then split it into a
  // collinear q-qbar pair.
  } else if (nLeg == 2) {

    // Gluon to be split and its two daughters.
    int iOld  = systemPtr->iParton[iFirst + 1];
    int i1New = event.append(  idNewQ, 101, iOld, 0, 0, 0, 
      event[iOld].col(), 0, 0.5 * event[iOld].p(), 
      0.5 * event[iOld].m(), 0.);
    int i2New = event.append( -idNewQ, 101, iOld, 0, 0, 0, 
      0, event[iOld].acol(), 0.5 * event[iOld].p(), 
      0.5 * event[iOld].m(), 0.);
    event[iOld].statusNeg();
    event[iOld].daughters( i1New, i2New);
     
    // Partons in the two new systems.
    if (event[systemPtr->iParton[iFirst]].col() == event[i2New].acol()) 
      swap( i1New, i2New);
    for (int i = 0; i <= iFirst; ++i) 
      iNewSys1.push_back( systemPtr->iParton[i]);
    iNewSys1.push_back( i1New);
    iNewSys2.push_back( i2New);
    for (int i = iSecond; i < int(systemPtr->size()); ++i) 
      iNewSys2.push_back( systemPtr->iParton[i]);

  // If several gluons between sparticles then find lowest-mass gluon pair
  // and replace it by a q-qbar pair.
  } else {

    // Find lowest-mass gluon pair and adjust effective quark mass.
    int    iMin  = 0;
    int    i1Old = 0;
    int    i2Old = 0;
    double mMin  = 1e20;
    for (int i = iFirst + 1; i < iSecond - 1; ++i) { 
      int    i1Tmp = systemPtr->iParton[i];
      int    i2Tmp = systemPtr->iParton[i + 1];
      double mTmp  = (event[i1Tmp].p() + event[i2Tmp].p()).mCalc();
      if (mTmp < mMin) {
        iMin  = i;
        i1Old = i1Tmp;
        i2Old = i2Tmp;
        mMin  = mTmp;
      }
    }
    double mEff = min( mNewQ + mOffsetCloudRH, 0.4 * mMin);

    // New kinematics  by sharing mHat appropriately.
    Vec4 p1New, p2New;
    if ( !newKin( event[i1Old].p(), event[i2Old].p(), 
      mEff, mEff, p1New, p2New, false) ) return false; 

    // Insert new partons and update history.
    int i1New, i2New;
    if (event[systemPtr->iParton[0]].acol() == 0) {
      i1New = event.append( -idNewQ, 101, i1Old, 0, 0, 0, 
        0, event[i1Old].acol(), p1New, mEff, 0.);
      i2New = event.append(  idNewQ, 101, i2Old, 0, 0, 0, 
        event[i2Old].col(), 0, p2New, mEff, 0.);
    } else {
      i1New = event.append(  idNewQ, 101, i1Old, 0, 0, 0, 
        event[i1Old].col(), 0, p1New, mEff, 0.);
      i2New = event.append( -idNewQ, 101, i2Old, 0, 0, 0, 
        0, event[i2Old].acol(), p2New, mEff, 0.);
    } 
    event[i1Old].statusNeg();
    event[i2Old].statusNeg();
    event[i1Old].daughters( i1New, 0);
    event[i2Old].daughters( i2New, 0);
     
    // Partons in the two new systems.
    for (int i = 0; i < iMin; ++i) 
      iNewSys1.push_back( systemPtr->iParton[i]);
    iNewSys1.push_back( i1New);
    iNewSys2.push_back( i2New);
    for (int i = iMin + 2; i < int(systemPtr->size()); ++i) 
      iNewSys2.push_back( systemPtr->iParton[i]);
  }

  // Erase the old system and insert the two new ones instead.
  colConfig.erase(iSys);
  colConfig.insert( iNewSys1, event);
  colConfig.insert( iNewSys2, event);   

  // Done. 
  return true;

}

//--------------------------------------------------------------------------

// Produce a R-hadron from a squark and another string end.

bool RHadrons::produceSquark( ColConfig& colConfig, Event& event) {

  // Initial values.
  int    nBody  = 0;
  int    iRNow  = 0;
  int    iNewQ  = 0;
  int    iNewL  = 0;

  // Check at which end of the string the squark is located.
  int    idAbsTop = event[ systemPtr->iParton[0] ].idAbs(); 
  bool   sqAtTop  = (allowRSb && idAbsTop == idRSb) 
                 || (allowRSt && idAbsTop == idRSt);

  // Copy down system consecutively from squark end.
  int    iBeg = event.size();
  iCreRHad[iRHad] = iBeg;
  if (sqAtTop) for (int i = 0; i < systemPtr->size(); ++i) 
    event.copy( systemPtr->iParton[i], 102);
  else         for (int i = systemPtr->size() - 1; i >= 0; --i) 
    event.copy( systemPtr->iParton[i], 102);
  int    iEnd = event.size() - 1; 

  // Input flavours. From now on H = heavy and L = light string ends. 
  int    idOldH = event[iBeg].id(); 
  int    idOldL = event[iEnd].id();

  // Pick new flavour to form R-hadron. 
  FlavContainer flavOld( idOldH%10);
  int    idNewQ = flavSelPtr->pick(flavOld).id;
  int    idRHad = toIdWithSquark( idOldH, idNewQ);
  if (idRHad == 0) {
     infoPtr->errorMsg("Error in RHadrons::produceSquark: "
       "cannot form R-hadron code");
     return false;
  }

  // Target mass of R-hadron and z value of fragmentation function.
  double mRHad  = particleDataPtr->m0(idRHad) + event[iBeg].m() 
                - ( (abs(idOldH) == idRSb) ? m0Sb : m0St );
  double z      = zSelPtr->zFrag( idOldH, idNewQ, mRHad*mRHad);

  // Basic kinematics of string piece where break is to occur.
  Vec4   pOldH  = event[iBeg].p();
  int    iOldL  = iBeg + 1;
  Vec4   pOldL  = event[iOldL].p();
  double mOldL  = event[iOldL].m();
  double mNewH  = mRHad / z;
  double sSys   = (pOldH + pOldL).m2Calc();
  double sRem   = (1. - z) * (sSys - mNewH*mNewH);
  double sMin   = pow2(mOldL + mCollapseRH); 

  // If too low remaining mass in system then add one more parton to it. 
  while ( ( sRem < sMin || sSys < pow2(mNewH + mOldL + MSAFETY) )
    && iOldL < iEnd ) {    
    ++iOldL;
    pOldL      += event[iOldL].p();
    mOldL       = event[iOldL].m();
    sSys        = (pOldH + pOldL).m2Calc();
    sRem        = (1. - z) * (sSys - mNewH*mNewH);
    sMin        = pow2(mOldL + mCollapseRH); 
  }

  // If enough mass then split off R-hadron and reduced system.
  if ( sRem > sMin && sSys > pow2(mNewH + mOldL + MSAFETY) ) {
    Vec4 pNewH, pNewL;
    if ( !newKin( pOldH, pOldL, mNewH, mOldL, pNewH, pNewL) ) {
      infoPtr->errorMsg("Error in RHadrons::produceSquark: "
       "failed to construct kinematics with reduced system");
      return false;
    }

    // Insert R-hadron with new momentum. 
    iRNow       = event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
      z * pNewH, mRHad, 0.);
 
    // Reduced system with new string endpoint and modified recoiler. 
    idNewQ      = -idNewQ;
    bool hasCol = (idNewQ > 0 && idNewQ < 10) || idNewQ < -10;
    int  col    = (hasCol) ? event[iOldL].acol() : 0;
    int  acol   = (hasCol) ? 0 : event[iOldL].col(); 
    iNewQ       = event.append( idNewQ, 105, iBeg, iOldL, 0, 0, col, acol,
      (1. - z) * pNewH, (1. - z) * mNewH, 0.);
    iNewL       = event.copy( iOldL, 105);
    event[iNewL].mothers( iBeg, iOldL);
    event[iNewL].p( pNewL);

    // Done with processing of split to R-hadron and reduced system.
    nBody = 3;
  }

  // Two-body final state: form light hadron from endpoint and new flavour.
  if (nBody == 0) {
    FlavContainer flav1( idOldL);
    FlavContainer flav2( -idNewQ);
    int iTry   = 0;
    int idNewL = flavSelPtr->combine( flav1, flav2);
    while (++iTry < NTRYMAX && idNewL == 0) 
      idNewL = flavSelPtr->combine( flav1, flav2);    
    if (idNewL == 0) {
       infoPtr->errorMsg("Error in RHadrons::produceSquark: "
         "cannot form light hadron code");
       return false;
    }
    double mNewL = particleDataPtr->mass( idNewL); 
    
    // Check that energy enough for light hadron and R-hadron.
    if ( sSys > pow2(mRHad + mNewL + MSAFETY) ) { 
      Vec4 pRHad, pNewL;
      if ( !newKin( pOldH, pOldL, mRHad, mNewL, pRHad, pNewL) ) {
        infoPtr->errorMsg("Error in RHadrons::produceSquark: "
         "failed to construct kinematics for two-hadron decay");
        return false;
      }

      // Insert R-hadron and light hadron. 
      iRNow = event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
        pRHad, mRHad, 0.);
      event.append( idNewL, 105, iBeg, iOldL, 0, 0, 0, 0, 
        pNewL, mNewL, 0.);
   
      // Done for two-body case.
      nBody = 2;
    }
  }

  // Final case: for very low mass collapse to one single R-hadron.  
  if (nBody == 0) { 
    idRHad = toIdWithSquark( idOldH, idOldL);
    if (idRHad == 0) {
       infoPtr->errorMsg("Error in RHadrons::produceSquark: "
         "cannot form R-hadron code");
       return false;
    }

    // Insert R-hadron with new momentum. 
    iRNow = event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
      systemPtr->pSum, systemPtr->mass, 0.);

    // Done with one-body case.
    nBody = 1;
  }
      
  // History bookkeeping: the R-hadron and processed partons. 
  iRHadron[iRHad] = iRNow;
  int iLast = event.size() - 1;
  for (int i = iBeg; i <= iOldL; ++i) {
    event[i].statusNeg(); 
    event[i].daughters( iRNow, iLast);
  }  

  // Remove old string system and, if needed, insert a new one.
  colConfig.erase(iSys);
  if (nBody == 3) {
    vector<int> iNewSys;
    iNewSys.push_back( iNewQ);
    iNewSys.push_back( iNewL);
    for ( int i = iOldL + 1; i <= iEnd; ++i) iNewSys.push_back( i);
    colConfig.insert( iNewSys, event);
  }     

  // Copy lifetime and vertex from sparticle to R-hadron.
  event[iRNow].tau( event[iBef].tau() );
  if (event[iBef].hasVertex()) event[iRNow].vProd( event[iBef].vProd() );
 
  // Done with production of a R-hadron from a squark.  
  return true;

}

//--------------------------------------------------------------------------

// Produce a R-hadron from a gluino and two string ends (or a gluon).

bool RHadrons::produceGluino( ColConfig& colConfig, Event& event) {

  // Initial values.
  int    iGlui   = 0; 
  int    idSave  = 0; 
  int    idQLeap = 0;
  bool   isDiq1  = false;
  vector<int> iSide1, iSide2, iSideTmp, iNewSys1, iNewSys2;
  Vec4   pGlui, pSide1, pSide2;

  // Decide whether to produce a gluinoball or not.
  bool isGBall = (rndmPtr->flat() < probGluinoballRH);
         
  // Extract one string system on either side of the gluino.
  for (int i = 0; i < systemPtr->size(); ++i) {
    int  iTmp   = systemPtr->iParton[i];
    if (iGlui == 0 && event[ iTmp ].id() == idRGo) {
      iGlui     = iTmp;
      pGlui     = event[ iTmp ].p();
    } else if (iGlui == 0) {
      iSideTmp.push_back( iTmp);
      pSide1   += event[ iTmp ].p();
    } else {
      iSide2.push_back( iTmp);
      pSide2   += event[ iTmp ].p();
    }
  }
      
  // Order sides from gluino outwards and with lowest relative mass first.
  for (int i = int(iSideTmp.size()) - 1; i >= 0; --i) 
    iSide1.push_back( iSideTmp[i]);
  double m1H    = (pGlui + pSide1).mCalc() - event[ iSide1.back() ].m();
  double m2H    = (pGlui + pSide2).mCalc() - event[ iSide2.back() ].m();
  if (m2H < m1H) {
    swap( iSide1, iSide2);
    swap( pSide1, pSide2);
  }

  // Begin loop over two sides of gluino, with lowest relative mass first.
  for (int iSide = 1; iSide < 3; ++iSide) {

    // Begin processing of lower-mass string side.
    int    idRHad = 0;
    double mRHad  = 0.;
    int    nBody  = 0;
    int    iRNow  = 0;
    int    iNewQ  = 0;
    int    iNewL  = 0;
    int    statusRHad = 0;

    // Copy down system consecutively from gluino end.
    int    iBeg   = event.size();
    if (iSide == 1) {
      iCreRHad[iRHad] = iBeg;
      event.copy( iGlui, 102);
      for (int i = 0; i < int(iSide1.size()); ++i) 
        event.copy( iSide1[i], 102);
    } else {
      event.copy( iGlui, 102);
      for (int i = 0; i < int(iSide2.size()); ++i) 
        event.copy( iSide2[i], 102);
    }
    int    iEnd   = event.size() - 1; 

    // Pick new flavour to help form R-hadron. Initial colour values.
    int    idOldL = event[iEnd].id();
    int    idNewQ = 21;
    if (!isGBall) {
      do {
        FlavContainer flavOld( idOldL);
        idNewQ = -flavSelPtr->pick(flavOld).id;
      } while (iSide == 2 && isDiq1 && abs(idNewQ) > 10);
      if (iSide == 1) isDiq1 = (abs(idNewQ) > 10);
    }
    bool   hasCol = (event[iEnd].col() > 0);
    int    colR   = 0;
    int    acolR  = 0;

    // Target intermediary mass or R-hadron mass.
    if (iSide == 1) {
      idSave      = idNewQ;
      idRHad      = (hasCol) ? 1009002 : -1009002;
      if (hasCol) colR  = event[iBeg].col();
      else        acolR = event[iBeg].acol();
      statusRHad  = 103;
      double mNewQ = particleDataPtr->constituentMass( idNewQ);
      if (isGBall) mNewQ *= 0.5;
      mRHad       = event[iBeg].m() + mOffsetCloudRH + mNewQ;
    } else {
      idRHad      = toIdWithGluino( idSave, idNewQ);
      if (idRHad == 0) {
         infoPtr->errorMsg("Error in RHadrons::produceGluino: "
           "cannot form R-hadron code");
         return false;
      }
      statusRHad  = 104;
      mRHad       = particleDataPtr->m0(idRHad) + event[iBeg].m() - m0Go;
    }
      
    // z value of fragmentation function.
    double z      = zSelPtr->zFrag( idRGo, idNewQ, mRHad*mRHad);

    // Basic kinematics of string piece where break is to occur.
    Vec4   pOldH  = event[iBeg].p();
    int    iOldL  = iBeg + 1;
    Vec4   pOldL  = event[iOldL].p();
    double mOldL  = event[iOldL].m();
    double mNewH  = mRHad / z;
    double sSys   = (pOldH + pOldL).m2Calc();
    double sRem   = (1. - z) * (sSys - mNewH*mNewH);
    double sMin   = pow2(mOldL + mCollapseRH); 

    // If too low remaining mass in system then add one more parton to it. 
    while ( ( sRem < sMin || sSys < pow2(mNewH + mOldL + MSAFETY) )
      && iOldL < iEnd ) {    
      ++iOldL;
      pOldL      += event[iOldL].p();
      mOldL       = event[iOldL].m();
      sSys        = (pOldH + pOldL).m2Calc();
      sRem        = (1. - z) * (sSys - mNewH*mNewH);
      sMin        = pow2(mOldL + mCollapseRH); 
    }

    // If enough mass then split off R-hadron and reduced system.
    if ( sRem > sMin && sSys > pow2(mNewH + mOldL + MSAFETY) ) {
      Vec4 pNewH, pNewL;
      if ( !newKin( pOldH, pOldL, mNewH, mOldL, pNewH, pNewL) ) {
        infoPtr->errorMsg("Error in RHadrons::produceGluino: "
         "failed to construct kinematics with reduced system");
        return false;
      }

      // Insert intermediary or R-hadron with new momentum, less colour.
      iRNow       = event.append( idRHad, statusRHad, iBeg, iOldL, 
        0, 0, colR, acolR, z * pNewH, mRHad, 0.);
 
      // Reduced system with new string endpoint and modified recoiler. 
      if (!isGBall) idNewQ = -idNewQ;
      int  colN   = (hasCol) ? 0 : event[iOldL].acol();
      int  acolN  = (hasCol) ? event[iOldL].col() : 0; 
      iNewQ       = event.append( idNewQ, 105, iBeg, iOldL, 0, 0, 
        colN, acolN, (1. - z) * pNewH, (1. - z) * mNewH, 0.);
      iNewL       = event.copy( iOldL, 105);
      event[iNewL].mothers( iBeg, iOldL);
      event[iNewL].p( pNewL);

      // Collect partons of new string system.
      if (iSide == 1) {
        iNewSys1.push_back( iNewQ);
        iNewSys1.push_back( iNewL);
        for ( int i = iOldL + 1; i <= iEnd; ++i) iNewSys1.push_back( i);
      } else {
        iNewSys2.push_back( iNewQ);
        iNewSys2.push_back( iNewL);
        for ( int i = iOldL + 1; i <= iEnd; ++i) iNewSys2.push_back( i);
      }     

      // Done with processing of split to R-hadron and reduced system.
      nBody = 3;
    }

    // For side-1 low-mass glueball system reabsorb full momentum. 
    if (nBody == 0 && isGBall && iSide == 1) { 
      idQLeap    = event[iOldL].id();
      Vec4 pNewH = event[iBeg].p() + pOldL;

      // Insert intermediary R-hadron with new momentum, less colour.
      iRNow      = event.append( idRHad, statusRHad, iBeg, iEnd, 
        0, 0, colR, acolR, pNewH, pNewH.mCalc(), 0.);
      nBody = 1;
    }
      
    // Two-body final state: form light hadron from endpoint and new flavour.
    // Also for side 2 if gluinoball formation gives quark from side 1.
    if (nBody == 0 && (!isGBall || (iSide == 2 && idQLeap != 0) )) {
      if (isGBall) idNewQ = -idQLeap;
      FlavContainer flav1( idOldL);
      FlavContainer flav2( -idNewQ);
      int iTry   = 0;
      int idNewL = flavSelPtr->combine( flav1, flav2);
      while (++iTry < NTRYMAX && idNewL == 0) 
        idNewL = flavSelPtr->combine( flav1, flav2);    
      if (idNewL == 0) {
         infoPtr->errorMsg("Error in RHadrons::produceGluino: "
           "cannot form light hadron code");
         return false;
      }
      double mNewL = particleDataPtr->mass( idNewL); 
    
      // Check that energy enough for light hadron and R-hadron.
      if ( sSys > pow2(mRHad + mNewL + MSAFETY) ) { 
        Vec4 pRHad, pNewL;
        if ( !newKin( pOldH, pOldL, mRHad, mNewL, pRHad, pNewL) ) {
          infoPtr->errorMsg("Error in RHadrons::produceGluino: "
           "failed to construct kinematics for two-hadron decay");
          return false;
        }

        // Insert intermediary or R-hadron and light hadron. 
        iRNow = event.append( idRHad, statusRHad, iBeg, iOldL, 0, 0,
          colR, acolR, pRHad, mRHad, 0.);
        event.append( idNewL, 105, iBeg, iOldL, 0, 0, 0, 0, 
          pNewL, mNewL, 0.);
   
        // Done for two-body case.
        nBody   = 2;
        // For this case gluinoball should be handled as normal flavour.
        isGBall = false;
      }
    }

    // Final case: for very low mass collapse to one single R-hadron.  
    if (nBody == 0 && (!isGBall || (iSide == 2 && idQLeap != 0) )) { 
      if (isGBall) idSave = idQLeap;
      if (iSide == 1) idSave = idOldL;
      else            idRHad = toIdWithGluino( idSave, idOldL);
      if (idRHad == 0) {
         infoPtr->errorMsg("Error in RHadrons::produceGluino: "
           "cannot form R-hadron code");
         return false;
      }

      // Insert R-hadron with new momentum. 
      iRNow = event.append( idRHad, statusRHad, iBeg, iOldL, 0, 0, 
        colR, acolR, pOldH + pOldL, (pOldH + pOldL).mCalc(), 0.);

      // Done with one-body case.
      // Even if hoped-for, it was not possible to create a gluinoball.
      isGBall = false;
    }
      
    // History bookkeeping: the processed partons. 
    int iLast = event.size() - 1;
    for (int i = iBeg; i <= iOldL; ++i) {
      event[i].statusNeg(); 
      event[i].daughters( iRNow, iLast);
    }  

    // End of loop over two sides of the gluino.
    iGlui   = iRNow;
  }

  // History bookkeeping: insert R-hadron; delete old string system. 
  if (iGlui == 0) {
     infoPtr->errorMsg("Error in RHadrons::produceGluino: "
           "could not handle gluinoball kinematics");
     return false;
  }
  iRHadron[iRHad] = iGlui;
  colConfig.erase(iSys);

  // Non-gluinoball: insert (up to) two new string systems.
  if (!isGBall) {
    if (iNewSys1.size() > 0) colConfig.insert( iNewSys1, event);
    if (iNewSys2.size() > 0) colConfig.insert( iNewSys2, event);

  // Gluinoball with enough energy in first string: 
  // join two temporary gluons into one. 
  } else if (idQLeap == 0) {
    int iG1   = iNewSys1[0];
    int iG2   = iNewSys2[0];
    int colG  = event[iG1].col()  + event[iG2].col();  
    int acolG = event[iG1].acol() + event[iG2].acol();  
    Vec4 pG   = event[iG1].p()    + event[iG2].p(); 
    int iG12  = event.append( 21, 105, iG1, iG2, 0, 0, colG, acolG, 
      pG, pG.mCalc(), 0.);

    // Temporary gluons no longer needed, but new colour to avoid warnings.
    event[iG1].id( 21);
    event[iG2].id( 21);
    event[iG1].statusNeg();
    event[iG2].statusNeg();
    int colBridge = event.nextColTag();
    if (event[iG1].col() == 0) {
      event[iG1].col(  colBridge);
      event[iG2].acol( colBridge);
    } else {
      event[iG1].acol( colBridge);
      event[iG2].col(  colBridge);
    } 

    // Form the remnant system from which the R-hadron has been split off. 
    vector<int> iNewSys12;
    for (int i = int(iNewSys1.size()) - 1; i > 0; --i) 
      iNewSys12.push_back( iNewSys1[i]);
    iNewSys12.push_back( iG12);
    for (int i = 1; i < int(iNewSys2.size()); ++i) 
      iNewSys12.push_back( iNewSys2[i]);
    colConfig.insert( iNewSys12, event); 

  // Gluinoball where side 1 was fully eaten, and its flavour content
  // should leap over to the other side, to a gluon there.
  } else {
    int iG2   = iNewSys2[0];
    event[iG2].id( idQLeap);
    colConfig.insert( iNewSys2, event);
  }

  // Copy lifetime and vertex from sparticle to R-hadron.
  event[iGlui].tau( event[iBef].tau() );
  if (event[iBef].hasVertex()) event[iGlui].vProd( event[iBef].vProd() );
 
  // Done with production of a R-hadron from a gluino.  
  return true;

}

//--------------------------------------------------------------------------

// Form a R-hadron code from a squark and a (di)quark code.
// First argument is the (anti)squark, second the (anti)(di)quark.

int RHadrons::toIdWithSquark( int id1, int id2) {

  // Check that physical combination; return 0 if not.
  int id1Abs = abs(id1);
  int id2Abs = abs(id2);
  if (id2Abs < 10 && id1 > 0 && id2 > 0) return 0;
  if (id2Abs < 10 && id1 < 0 && id2 < 0) return 0;
  if (id2Abs > 10 && id1 > 0 && id2 < 0) return 0;
  if (id2Abs > 10 && id1 < 0 && id2 > 0) return 0;

  // Form R-hadron code. Flip sign for antisquark. 
  bool isSt = (id1Abs == idRSt);
  int idRHad = 1000000;
  if (id2Abs < 10) idRHad += ((isSt) ? 600 : 500) + 10 * id2Abs + 2;
  else idRHad += ((isSt) ? 6000 : 5000) + 10 * (id2Abs/100) + id2Abs%10;
  if (id1 < 0) idRHad = -idRHad;

  // Done.
  return idRHad;

}

//--------------------------------------------------------------------------

// Resolve a R-hadron code into a squark and a (di)quark.
// Return a pair, where first is the squark and the second the (di)quark.

pair<int,int> RHadrons::fromIdWithSquark( int idRHad) {

  // Find squark flavour content. 
  int idLight = (abs(idRHad) - 1000000) / 10;
  int idSq    = (idLight < 100) ? idLight/10 : idLight/100;
  int id1     = (idSq == 6) ? idRSt : idRSb;
  if (idRHad < 0) id1 = -id1;

  // Find light (di)quark flavour content. 
  int id2     =  (idLight < 100) ? idLight%10 : idLight%100;
  if (id2 > 10) id2 = 100 * id2 + abs(idRHad)%10;
  if ((id2 < 10 && idRHad > 0) || (id2 > 10 && idRHad < 0)) id2 = -id2;

  // Done.
  return make_pair( id1, id2);

}   
 
//--------------------------------------------------------------------------

// Form a R-hadron code from two quark/diquark endpoints and a gluino.

int RHadrons::toIdWithGluino( int id1, int id2) {

  // Check that physical combination; return 0 if not. Handle gluinoball.
  int id1Abs = abs(id1);
  int id2Abs = abs(id2);
  if (id1Abs == 21 && id2Abs == 21) return 1000993;
  int idMax  = max( id1Abs, id2Abs);
  int idMin  = min( id1Abs, id2Abs);
  if (idMin > 10) return 0;
  if (idMax > 10 && id1 > 0 && id2 < 0) return 0;
  if (idMax > 10 && id1 < 0 && id2 > 0) return 0;
  if (idMax < 10 && id1 > 0 && id2 > 0) return 0;
  if (idMax < 10 && id1 < 0 && id2 < 0) return 0;

  // Form R-meson code. Flip sign for antiparticle.
  int idRHad = 0;
  if (idMax < 10) {
    idRHad = 1009003 + 100 * idMax + 10 * idMin;
    if (idMin != idMax && idMax%2 == 1) {
      if (id1Abs == idMax && id1 > 0) idRHad = -idRHad;
      if (id2Abs == idMax && id2 > 0) idRHad = -idRHad;
    }
    if (idMin != idMax && idMax%2 == 0) {
      if (id1Abs == idMax && id1 < 0) idRHad = -idRHad;
      if (id2Abs == idMax && id2 < 0) idRHad = -idRHad;
    }

  // Form R-baryon code. Flip sign for antiparticle.
  } else {
    int idA = idMax/1000;
    int idB = (idMax/100)%10;
    int idC = idMin;
    if (idC > idB) swap( idB, idC);
    if (idB > idA) swap( idA, idB);    
    if (idC > idB) swap( idB, idC);
    idRHad  = 1090004 + 1000 * idA + 100 * idB + 10 * idC;
    if (id1 < 0) idRHad = -idRHad;
  }

  // Done.
  return idRHad;

}

//--------------------------------------------------------------------------

// Resolve a R-hadron code into two quark/diquark endpoints (and a gluino).
// Return a pair, where first carries colour and second anticolour.

pair<int,int> RHadrons::fromIdWithGluino( int idRHad) {

  // Find light flavour content of R-hadron.
  int idLight = (abs(idRHad) - 1000000) / 10; 
  int id1, id2, idTmp, idA, idB, idC; 

  // Gluinoballs: split g into d dbar or u ubar.
  if (idLight < 100) {
    id1 = (rndmPtr->flat() < 0.5) ? 1 : 2;
    id2 = -id1;

  // Gluino-meson: split into q + qbar.
  } else if (idLight < 1000) {
    id1 = (idLight / 10) % 10;  
    id2 = -(idLight % 10);
    // Flip signs when first quark of down-type.
    if (id1%2 == 1) {
      idTmp = id1;
      id1   = -id2;
      id2   = -idTmp;
    }

  // Gluino-baryon: split to q + qq (diquark). 
  // Pick diquark at random, except if c or b involved.
  } else {
    idA = (idLight / 100) % 10;
    idB = (idLight / 10) % 10;
    idC = idLight % 10;
    double rndmQ = 3. * rndmPtr->flat();
    if (idA > 3) rndmQ = 0.5;
    if (rndmQ < 1.) {
      id1 = idA;
      id2 = 1000 * idB + 100 * idC + 3;
      if (idB != idC && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2; 
    } else if (rndmQ < 2.) {
      id1 = idB;
      id2 = 1000 * idA + 100 * idC + 3;
      if (idA != idC && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2; 
    } else {
      id1 = idC;
      id2 = 1000 * idA + 100 * idB +3;
      if (idA != idB && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2; 
    }
  }

  // Flip signs for anti-R-hadron.
  if (idRHad < 0) {
    idTmp = id1;
    id1   = -id2;
    id2   = -idTmp;
  }

  // Done.
  return make_pair( id1, id2);

}   
 
//--------------------------------------------------------------------------

// Construct modified four-vectors to match modified masses:
// minimal reshuffling of momentum along common axis. 
// Note that last two arguments are used to provide output!

bool RHadrons::newKin( Vec4 pOld1, Vec4 pOld2, double mNew1, double mNew2,
  Vec4& pNew1, Vec4& pNew2, bool checkMargin) {

  // Squared masses in initial and final kinematics.
  double sSum   = (pOld1 + pOld2).m2Calc();
  double sOld1  = pOld1.m2Calc();
  double sOld2  = pOld2.m2Calc();
  double sNew1  = mNew1 * mNew1;
  double sNew2  = mNew2 * mNew2;

  // Check that kinematically possible.
  if (checkMargin && pow2(mNew1 + mNew2 + MSAFETY) > sSum) return false;

  // Transfer coefficients to give four-vectors with the new masses.
  double lamOld = sqrt( pow2(sSum - sOld1 - sOld2) - 4. * sOld1 * sOld2 );
  double lamNew = sqrt( pow2(sSum - sNew1 - sNew2) - 4. * sNew1 * sNew2 );   
  double move1  = (lamNew * (sSum - sOld1 + sOld2) 
                -  lamOld * (sSum - sNew1 + sNew2)) / (2. * sSum * lamOld);
  double move2  = (lamNew * (sSum + sOld1 - sOld2) 
                -  lamOld * (sSum + sNew1 - sNew2)) / (2. * sSum * lamOld);
  
  // Construct final vectors. Done.
  pNew1 = (1. + move1) * pOld1 - move2 * pOld2;
  pNew2 = (1. + move2) * pOld2 - move1 * pOld1;
  return true;

}   
 
//==========================================================================

} // end namespace Pythia8