source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/body/src/G4QParticle.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 38.6 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26//
27// $Id: G4QParticle.cc,v 1.36 2010/06/25 14:03:44 mkossov Exp $
28// GEANT4 tag $Name: hadr-chips-V09-03-08 $
29//
30//      ---------------- G4QParticle ----------------
31//             by Mikhail Kossov, Sept 1999.
32//      class for Particles in the CHIPS Model
33// -------------------------------------------------------------------
34//       1         2         3         4         5         6         7         8         9
35//34567890123456789012345678901234567890123456789012345678901234567890123456789012345678901
36// -----------------------------------------------------------------------
37// Short description: The G4QParticle is a part of the CHIPS World. It is
38// characterized by the quark content, spin, mass, width and a vector of
39// the decay channels (G4QDecayCannelVector).
40// -----------------------------------------------------------------------
41//#define debug
42//#define pdebug
43
44#include "G4QParticleVector.hh"
45
46G4QParticle::G4QParticle()
47{
48#ifdef debug
49  G4cout<<"G4QParticle::Constr: Default constructor is called"<<G4endl;
50#endif
51}
52
53G4QParticle::G4QParticle(G4int thePDG)
54{
55  aQPDG      = G4QPDGCode(thePDG);
56  aQuarkCont = aQPDG.GetQuarkContent();
57  aDecay     = InitDecayVector(aQPDG.GetQCode());
58}
59
60G4QParticle::G4QParticle(G4bool f, G4int theQCode)
61{
62  aQPDG      = G4QPDGCode(f,theQCode);
63#ifdef debug
64  G4cout<<"G4QParticle::Constr: PDG="<<aQPDG.GetPDGCode()<<G4endl;
65#endif
66  aQuarkCont = aQPDG.GetQuarkContent();
67  aDecay     = InitDecayVector(theQCode);
68}
69
70G4QParticle::G4QParticle(const G4QParticle& right)
71{
72  aQPDG                = right.aQPDG;
73  //aDecay (Vector)
74  G4int nD             = right.aDecay.size(); 
75  if(nD) for(G4int id=0; id<nD; id++)
76  {
77    G4QDecayChan* curD = new G4QDecayChan(right.aDecay[id]);
78    aDecay.push_back(curD);
79  }
80
81  aQuarkCont           = right.aQuarkCont;
82}
83
84G4QParticle::G4QParticle(G4QParticle* right)
85{
86  aQPDG                = right->aQPDG;
87  //aDecay (Vector)
88  G4int nD             = right->aDecay.size();
89  if(nD) for(G4int id=0; id<nD; id++)
90  {
91    G4QDecayChan* curD = new G4QDecayChan(right->aDecay[id]);
92    aDecay.push_back(curD);
93  }
94
95  aQuarkCont           = right->aQuarkCont;
96}
97
98G4QParticle::~G4QParticle() 
99{
100  G4int nDC=aDecay.size();
101  //G4cout<<"G4QParticle::Destructor: Before nDC="<<nDC<<G4endl; // TMP
102  if(nDC) std::for_each(aDecay.begin(), aDecay.end(), DeleteQDecayChan());
103  //G4cout<<"G4QParticle::Destructor: After"<<G4endl; // TMP
104  aDecay.clear();
105}
106
107// Assignment operator
108const G4QParticle& G4QParticle::operator=(const G4QParticle &right)
109{
110  if(this != &right)                          // Beware of self assignment
111  {
112    aQPDG                = right.aQPDG;
113    //aDecay (Vector)
114    G4int iD             = aDecay.size();
115    if(iD) for(G4int jd=0; jd<iD; jd++) delete aDecay[jd];
116    aDecay.clear();
117    G4int nD             = right.aDecay.size();
118    if(nD) for(G4int id=0; id<nD; id++)
119    {
120      G4QDecayChan* curD = new G4QDecayChan(right.aDecay[id]);
121      aDecay.push_back(curD);
122    }
123
124    aQuarkCont           = right.aQuarkCont;
125  }
126  return *this;
127}
128
129// Standard output for QParticle
130std::ostream& operator<<(std::ostream& lhs, G4QParticle& rhs)
131//       =========================================
132{
133  G4QPDGCode rhsQPDG = rhs.GetQPDG();
134  lhs << G4endl << "Particle with PDG=" << rhsQPDG << ", Spin=" << rhs.GetSpin()
135      << ", mass=" << rhs.GetMass() << ", width=" << rhs.GetWidth() << G4endl;
136  lhs<<" Quark Content of the Particle="<<rhs.GetQContent()<<", Decay Channels:"<<G4endl;
137  G4QDecayChanVector DCV = rhs.GetDecayVector();
138  G4int n = DCV.size();
139  for (int i=0; i<n; i++)
140  {
141    lhs << DCV[i]->GetDecayChanLimit() << "PDG codes";
142    G4QPDGCodeVector PCV=DCV[i]->GetVecOfSecHadrons();
143    G4int m = PCV.size();
144    for (int j=0; j<m; j++)
145    {
146      if(!j) lhs << ":";
147      else   lhs << ",";
148      lhs << PCV[j]->GetPDGCode() ;
149    }
150  }
151  return lhs;
152}
153
154// Initialize the PDG-Particle by QCode @@ Can be improved, using PDG.DATA file
155G4QDecayChanVector G4QParticle::InitDecayVector(G4int nQ)
156//    ===================================================
157{
158  //static G4int nP = 486;                  // Up to A=80
159  //static const G4int nP = 494;              // Up to A=80(?) "Isonuclear revision"
160  static const G4int nP = 512; // A<57 "Leptons/Hypernuclei" G4QCHIPSWorld::GetParticles(!)
161  //static G4QDecayChanVector* DecayDB = new G4QDecayChanVector[nP];
162  static G4QDecayChanVector DecayDB[nP];
163  static int limit= 0;
164  if(nQ>=limit && nQ<nP)
165  {
166    //*** Secondary PDG-particles should be ordered in a Channel by increasing width***!***
167    //** Channels should be ordered by increasing minimum mass of the secondary particles**
168    //if(limit<=  0 && nQ>=  0)DecayDB[  0] = 0;    // e-     (11)
169    //if(limit<=  1 && nQ>=  1)DecayDB[  1] = 0;    // nu_e   (12)
170    //if(limit<=  2 && nQ>=  2)DecayDB[  2] = 0;    // mu-    (13)
171    //if(limit<=  3 && nQ>=  3)DecayDB[  3] = 0;    // mu_e   (14)
172    //if(limit<=  4 && nQ>=  4)DecayDB[  4] = 0;    // tau-   (15)
173    //if(limit<=  5 && nQ>=  5)DecayDB[  5] = 0;    // nu_tau (16)
174    //if(limit<=  6 && nQ>=  6)DecayDB[  6] = 0;    // gamma  (22)
175    if(limit<=  7 && nQ>=  7)                       // Z0     (23)
176    {
177      DecayDB[  7].push_back(new G4QDecayChan(.036,  11, -11));
178      DecayDB[  7].push_back(new G4QDecayChan(.073,  13, -13));
179      DecayDB[  7].push_back(new G4QDecayChan(.107,  15, -15));
180      DecayDB[  7].push_back(new G4QDecayChan(.174,  12, -12)); // @@ Fake invisible decays
181      DecayDB[  7].push_back(new G4QDecayChan(.240,  14, -14));
182      DecayDB[  7].push_back(new G4QDecayChan(.307,  16, -16));
183      DecayDB[  7].push_back(new G4QDecayChan(.400,2112,-2112)); // @@ Fake Hadronic decays
184      DecayDB[  7].push_back(new G4QDecayChan(.500,2212,-2212)); // @@ Need heavy quarks
185      DecayDB[  7].push_back(new G4QDecayChan(.600,2212,-2212, 111));
186      DecayDB[  7].push_back(new G4QDecayChan(.700,2112,-2112, 111));
187      DecayDB[  7].push_back(new G4QDecayChan(.800,2212,-2112,-211));
188      DecayDB[  7].push_back(new G4QDecayChan(.990,2112,-2212, 211));
189      DecayDB[  7].push_back(new G4QDecayChan(1.00,2112,-3122, 111));
190    }
191    if(limit<=  8 && nQ>=  8)                       // W-     (24) @@ Update HadronicDecays
192    {
193      DecayDB[  8].push_back(new G4QDecayChan(.107,  11, -12));
194      DecayDB[  8].push_back(new G4QDecayChan(.214,  13, -14));
195      DecayDB[  8].push_back(new G4QDecayChan(.321,  15, -16));
196      DecayDB[  8].push_back(new G4QDecayChan(.421,2112,-2212)); // @@ Fake Hadronic decays
197      DecayDB[  8].push_back(new G4QDecayChan(.521,2112,-2112,-211));
198      DecayDB[  8].push_back(new G4QDecayChan(.621,2212,-2212,-211));
199      DecayDB[  8].push_back(new G4QDecayChan(.721,3122,-3122,-211));
200      DecayDB[  8].push_back(new G4QDecayChan(.821,2112,-2212, 111));
201      DecayDB[  8].push_back(new G4QDecayChan(.921,3122,-2212, 111));
202      DecayDB[  8].push_back(new G4QDecayChan(1.00,2112,-3122,-211));
203    }
204    //if(limit<=  9 && nQ>=  9)DecayDB[  9] = 0;    // H0     (25)
205    //if(limit<= 10 && nQ>= 10)DecayDB[ 10] = 0;    // H-     (37)
206    if(limit<= 11 && nQ>= 11)    // Low sigma=pi,pi S-wave : f_0 (800)
207    {
208      DecayDB[ 11].push_back(new G4QDecayChan(.333,211,-211));
209      DecayDB[ 11].push_back(new G4QDecayChan(1.00,111, 111));
210    }
211    if(limit<= 12 && nQ>= 12)    // Midle Regeon-Pomeron   : f_0 (980)
212    {
213      DecayDB[ 12].push_back(new G4QDecayChan(.333,211,-211));
214      DecayDB[ 12].push_back(new G4QDecayChan(1.00,111, 111));
215    }
216    if(limit<= 13 && nQ>= 13)    // High Regeon-Pomeron    : f_0 (1500)
217    {
218      DecayDB[ 13].push_back(new G4QDecayChan(.019,221, 331));
219      DecayDB[ 13].push_back(new G4QDecayChan(.070,221, 221));
220      DecayDB[ 13].push_back(new G4QDecayChan(.113,311,-311));
221      DecayDB[ 13].push_back(new G4QDecayChan(.156,321,-321));
222      DecayDB[ 13].push_back(new G4QDecayChan(.578,211,-211)); //@@ include 4pi decays
223      DecayDB[ 13].push_back(new G4QDecayChan(1.00,111, 111));
224    }
225    //if(limit<= 14 && nQ>= 14)DecayDB[ 14].push_back(new G4QDecayChan(1.00,22,22));//Pi0
226    //if(limit<= 15 && nQ>= 15)DecayDB[ 15] = 0;    // Pi +
227    if(limit<= 16 && nQ>= 16)    // eta
228    {
229      DecayDB[ 16].push_back(new G4QDecayChan(.226,211,-211,111));
230      DecayDB[ 16].push_back(new G4QDecayChan(.551,111, 111,111));
231      DecayDB[ 16].push_back(new G4QDecayChan(.598,211,-211, 22));
232      DecayDB[ 16].push_back(new G4QDecayChan(.606, 11, -11, 22)); //@@ .002 (pi+)(pi-)2gam
233      DecayDB[ 16].push_back(new G4QDecayChan(1.00, 22,  22));
234    }
235    //if(limit<= 17 && nQ>= 17)    // K 0 (K_short - probab 1/2) @@@@@@@@@@@@
236    //{
237    //  DecayDB[ 17].push_back(new G4QDecayChan(.6861,211,-211));
238    //  DecayDB[ 17].push_back(new G4QDecayChan(1.00, 111, 111));
239    //}
240    //if(limit<= 18 && nQ>= 18)DecayDB[  8] = 0;    // K +
241    if(limit<= 19 && nQ>= 19)    // eta'
242    {
243      DecayDB[ 19].push_back(new G4QDecayChan(.443,211,-211,221));
244      DecayDB[ 19].push_back(new G4QDecayChan(.652,111, 111,221));
245      DecayDB[ 19].push_back(new G4QDecayChan(.947, 22, 223));
246      DecayDB[ 19].push_back(new G4QDecayChan(.949,111, 111,111));
247      DecayDB[ 19].push_back(new G4QDecayChan(.979, 22, 113));
248      DecayDB[ 19].push_back(new G4QDecayChan(1.00, 22,  22));
249    }
250    //if(limit<= 20 && nQ>= 20)DecayDB[ 20] = 0;    // n
251    //if(limit<= 21 && nQ>= 21)DecayDB[ 21] = 0;    // p
252    //if(limit<= 22 && nQ>= 22)    // Lambda ===>>> all week decays are closed at this time
253    //{
254    //  DecayDB[ 22].push_back(new G4QDecayChan(.640,2212,-211));
255    //  DecayDB[ 22].push_back(new G4QDecayChan(1.00,2112, 111));
256    //}
257    //if(limit<= 23 &&nQ>=23)DecayDB[23].push_back(new G4QDecayChan(1.,2112,-211));//Sigma-
258    if(limit<= 24 &&nQ>=24)DecayDB[24].push_back(new G4QDecayChan(1.,3122,22));//Sigma0(EM)
259    //if(limit<= 25 && nQ>= 25)    // Sigma +
260    //{
261    //  DecayDB[ 25].push_back(new G4QDecayChan(.484,2112, 211));
262    //  DecayDB[ 25].push_back(new G4QDecayChan(1.00,2212, 111));
263    //}
264    //if(limit<= 26 && nQ>=26)DecayDB[26].push_back(new G4QDecayChan(1.,3122,-211));// Ksi-
265    //if(limit<= 27 && nQ>=27)DecayDB[27].push_back(new G4QDecayChan(1.,3122, 111));// Ksi0
266    if(limit<= 28 && nQ>= 28)DecayDB[ 28].push_back(new G4QDecayChan(1., 211,-211));// rho0
267    if(limit<= 29 && nQ>= 29)DecayDB[ 29].push_back(new G4QDecayChan(1., 211, 111));// rho+
268    if(limit<= 30 && nQ>= 30)    // omega
269    {
270      DecayDB[ 30].push_back(new G4QDecayChan(.891, 211,-211,111));
271      DecayDB[ 30].push_back(new G4QDecayChan(.908, 211,-211));
272      DecayDB[ 30].push_back(new G4QDecayChan(.997,  22, 111));
273      DecayDB[ 30].push_back(new G4QDecayChan(.998,  11, -11, 111)); //@@NeedsMoreAccurate
274      DecayDB[ 30].push_back(new G4QDecayChan(.998,  13, -13, 111));
275      DecayDB[ 30].push_back(new G4QDecayChan(1.00,  22, 221));
276    }
277    if(limit<= 31 && nQ>= 31)    // K* 0
278    {
279      DecayDB[ 31].push_back(new G4QDecayChan(.667,-211, 321));
280      DecayDB[ 31].push_back(new G4QDecayChan(1.00, 111, 311));
281    }
282    if(limit<= 32 && nQ>= 32)    // K* +
283    {
284      DecayDB[ 32].push_back(new G4QDecayChan(.667, 211, 311));
285      DecayDB[ 32].push_back(new G4QDecayChan(1.00, 111, 321));
286    }
287    if(limit<= 33 && nQ>= 33)    // phi
288    {
289      DecayDB[ 33].push_back(new G4QDecayChan(.491, 311,-311));
290      DecayDB[ 33].push_back(new G4QDecayChan(.831, 321,-321));
291      DecayDB[ 33].push_back(new G4QDecayChan(.844,  22, 221));
292      DecayDB[ 33].push_back(new G4QDecayChan(.846,  22, 111));
293      DecayDB[ 33].push_back(new G4QDecayChan(.897, 211,-213));
294      DecayDB[ 33].push_back(new G4QDecayChan(.948,-211, 213));
295      DecayDB[ 33].push_back(new G4QDecayChan(1.00, 111, 113));
296    }
297    if(limit<= 34 && nQ>= 34)DecayDB[34].push_back(new G4QDecayChan(1.,2112,-211));//Delta-
298    if(limit<= 35 && nQ>= 35)    // Delta 0
299    {
300      DecayDB[ 35].push_back(new G4QDecayChan(.333,2212,-211));
301      DecayDB[ 35].push_back(new G4QDecayChan(1.00,2112, 111));
302    }
303    if(limit<= 36 && nQ>= 36)    // Delta +
304    {
305      DecayDB[ 36].push_back(new G4QDecayChan(.333,2112, 211));
306      DecayDB[ 36].push_back(new G4QDecayChan(1.00,2212, 111));
307    }
308    if(limit<= 37 && nQ>= 37)DecayDB[37].push_back(new G4QDecayChan(1.,2212,211));//Delta++
309    if(limit<= 38 && nQ>= 38)    // Lambda* (1520)
310    {
311      DecayDB[ 38].push_back(new G4QDecayChan(.225,3112,-311));
312      DecayDB[ 38].push_back(new G4QDecayChan(.450,3222,-321));
313      DecayDB[ 38].push_back(new G4QDecayChan(.453,3112,211,111));
314      DecayDB[ 38].push_back(new G4QDecayChan(.456,3212,211,-211));
315      DecayDB[ 38].push_back(new G4QDecayChan(.459,3212,111,111));
316      DecayDB[ 38].push_back(new G4QDecayChan(.462,3222,-211,111));
317      DecayDB[ 38].push_back(new G4QDecayChan(.512,3122,211,-211));
318      DecayDB[ 38].push_back(new G4QDecayChan(.562,3122,111,111));
319      DecayDB[ 38].push_back(new G4QDecayChan(.702,3222,-211));
320      DecayDB[ 38].push_back(new G4QDecayChan(.842,3212, 111));
321      DecayDB[ 38].push_back(new G4QDecayChan(.982,3112, 211));
322      DecayDB[ 38].push_back(new G4QDecayChan(1.00,3122,  22));
323    }
324    if(limit<= 39 && nQ>= 39)    // Sigma* -
325    {
326      DecayDB[ 39].push_back(new G4QDecayChan(.060,3112, 111));
327      DecayDB[ 39].push_back(new G4QDecayChan(.120,3212,-211));
328      DecayDB[ 39].push_back(new G4QDecayChan(1.00,3122,-211));
329    }
330    if(limit<= 40 && nQ>= 40)    // Sigma* 0
331    {
332      DecayDB[ 40].push_back(new G4QDecayChan(.040,3112, 211));
333      DecayDB[ 40].push_back(new G4QDecayChan(.080,3222,-211));
334      DecayDB[ 40].push_back(new G4QDecayChan(.120,3212, 111));
335      DecayDB[ 40].push_back(new G4QDecayChan(1.00,3122, 111));
336    }
337    if(limit<= 41 && nQ>= 41)    // Sigma* +
338    {
339      DecayDB[ 41].push_back(new G4QDecayChan(.060,3212, 211));
340      DecayDB[ 41].push_back(new G4QDecayChan(.120,3222, 111));
341      DecayDB[ 41].push_back(new G4QDecayChan(1.00,3122, 211));
342    }
343    if(limit<= 42 && nQ>= 42)    // Ksi* -
344    {
345      DecayDB[ 42].push_back(new G4QDecayChan(.667,3322,-211));
346      DecayDB[ 42].push_back(new G4QDecayChan(1.00,3312, 111));
347    }
348    if(limit<= 43 && nQ>= 43)    // Ksi* 0
349    {
350      DecayDB[ 43].push_back(new G4QDecayChan(.667,3312, 211));
351      DecayDB[ 43].push_back(new G4QDecayChan(1.00,3322, 111));
352    }
353    //if(limit<= 44 && nQ>= 44)    // OMEGA - (Weak)
354    //{
355    //  DecayDB[ 44].push_back(new G4QDecayChan(.678,3122, 321));
356    //  DecayDB[ 44].push_back(new G4QDecayChan(.914,3322,-211));
357    //  DecayDB[ 44].push_back(new G4QDecayChan(1.00,3312, 111));
358    //}
359    if(limit<= 45 && nQ>= 45)    // a_2 0
360    {
361      DecayDB[ 45].push_back(new G4QDecayChan(.070, 211,-211,223));
362      DecayDB[ 45].push_back(new G4QDecayChan(.106, 111, 111,223));
363      DecayDB[ 45].push_back(new G4QDecayChan(.131, 321,-321));
364      DecayDB[ 45].push_back(new G4QDecayChan(.156, 311,-311));
365      DecayDB[ 45].push_back(new G4QDecayChan(.301, 111, 221));
366      DecayDB[ 45].push_back(new G4QDecayChan(.534,-211, 213));
367      DecayDB[ 45].push_back(new G4QDecayChan(.767, 211,-213));
368      DecayDB[ 45].push_back(new G4QDecayChan(1.00, 111, 113));
369    }
370    if(limit<= 46 && nQ>= 46)    // a_2 +
371    {
372      DecayDB[ 46].push_back(new G4QDecayChan(.106,111,211,223));
373      DecayDB[ 46].push_back(new G4QDecayChan(.156, 321,-311));
374      DecayDB[ 46].push_back(new G4QDecayChan(.301, 211, 221));
375      DecayDB[ 46].push_back(new G4QDecayChan(.651, 211, 113));
376      DecayDB[ 46].push_back(new G4QDecayChan(1.00, 111, 213));
377    }
378    if(limit<= 47 && nQ>= 47)    // f_2 0
379    {
380      DecayDB[ 47].push_back(new G4QDecayChan(.005, 221, 221));
381      DecayDB[ 47].push_back(new G4QDecayChan(.028, 311,-311));
382      DecayDB[ 47].push_back(new G4QDecayChan(.051, 321,-321));
383      DecayDB[ 47].push_back(new G4QDecayChan(.123, 111, 113));
384      DecayDB[ 47].push_back(new G4QDecayChan(.126, 111, 221));
385      DecayDB[ 47].push_back(new G4QDecayChan(.152, 211,-211,113));
386      DecayDB[ 47].push_back(new G4QDecayChan(.717, 211,-211));
387      DecayDB[ 47].push_back(new G4QDecayChan(1.00, 111, 111));
388    }
389    if(limit<= 48 && nQ>= 48)    // K_2 0
390    {
391      DecayDB[ 48].push_back(new G4QDecayChan(.028, 311, 223));
392      DecayDB[ 48].push_back(new G4QDecayChan(.074, 211,-211,313));
393      DecayDB[ 48].push_back(new G4QDecayChan(.143,111,-211,323));
394      DecayDB[ 48].push_back(new G4QDecayChan(.166,111, 111,313));
395      DecayDB[ 48].push_back(new G4QDecayChan(.190,-211, 323));
396      DecayDB[ 48].push_back(new G4QDecayChan(.314, 111, 313));
397      DecayDB[ 48].push_back(new G4QDecayChan(.357, 311, 113));
398      DecayDB[ 48].push_back(new G4QDecayChan(.500, 321,-213));
399      DecayDB[ 48].push_back(new G4QDecayChan(.750,-211, 321));
400      DecayDB[ 48].push_back(new G4QDecayChan(1.00, 111, 311));
401    }
402    if(limit<= 49 && nQ>= 49)    // K_2 +
403    {
404      DecayDB[ 49].push_back(new G4QDecayChan(.028, 321, 223));
405      DecayDB[ 49].push_back(new G4QDecayChan(.074,211,-211,323));
406      DecayDB[ 49].push_back(new G4QDecayChan(.143,111, 211,313));
407      DecayDB[ 49].push_back(new G4QDecayChan(.166,111, 111,323));
408      DecayDB[ 49].push_back(new G4QDecayChan(.190, 211, 313));
409      DecayDB[ 49].push_back(new G4QDecayChan(.314, 111, 323));
410      DecayDB[ 49].push_back(new G4QDecayChan(.357, 311, 213));
411      DecayDB[ 49].push_back(new G4QDecayChan(.500, 321, 113));
412      DecayDB[ 49].push_back(new G4QDecayChan(.750, 211, 311));
413      DecayDB[ 49].push_back(new G4QDecayChan(1.00, 111, 321));
414    }
415    if(limit<= 50 && nQ>= 50)    // f_2' 0
416    {
417      DecayDB[ 50].push_back(new G4QDecayChan(.103, 221, 221));
418      DecayDB[ 50].push_back(new G4QDecayChan(.547, 311,-311));
419      DecayDB[ 50].push_back(new G4QDecayChan(.991, 321,-321));
420      DecayDB[ 50].push_back(new G4QDecayChan(.997, 211,-211));
421      DecayDB[ 50].push_back(new G4QDecayChan(1.00, 111, 111));
422    }
423    if(limit<= 51 && nQ>= 51)    // N_5/2 0
424    {
425      DecayDB[ 51].push_back(new G4QDecayChan(.040, 211, 1114));
426      DecayDB[ 51].push_back(new G4QDecayChan(.080, 111, 2114));
427      DecayDB[ 51].push_back(new G4QDecayChan(.120,-211, 2214));
428      DecayDB[ 51].push_back(new G4QDecayChan(.180, 2112, 113));
429      DecayDB[ 51].push_back(new G4QDecayChan(.210, 2212,-213));
430      DecayDB[ 51].push_back(new G4QDecayChan(.340, 2112, 110));
431      DecayDB[ 51].push_back(new G4QDecayChan(.780, 2212,-211));
432      DecayDB[ 51].push_back(new G4QDecayChan(1.00, 2112, 111));
433    }
434    if(limit<= 52 && nQ>= 52)    // N_5/2 +
435    {
436      DecayDB[ 52].push_back(new G4QDecayChan(.040,-211, 2224));
437      DecayDB[ 52].push_back(new G4QDecayChan(.080, 211, 2114));
438      DecayDB[ 52].push_back(new G4QDecayChan(.120, 111, 2214));
439      DecayDB[ 52].push_back(new G4QDecayChan(.180, 2112, 213));
440      DecayDB[ 52].push_back(new G4QDecayChan(.210, 2212, 113));
441      DecayDB[ 52].push_back(new G4QDecayChan(.340, 2212, 229));
442      DecayDB[ 52].push_back(new G4QDecayChan(.780, 2112, 211));
443      DecayDB[ 52].push_back(new G4QDecayChan(1.00, 2212, 111));
444    }
445    if(limit<= 53 && nQ>= 53)    // LAMBDA_5/2
446    {
447      DecayDB[ 53].push_back(new G4QDecayChan(.350, 2112,-311));
448      DecayDB[ 53].push_back(new G4QDecayChan(.700, 2212,-321));
449      DecayDB[ 53].push_back(new G4QDecayChan(.740, 211, 3114));
450      DecayDB[ 53].push_back(new G4QDecayChan(.780,-211, 3224));
451      DecayDB[ 53].push_back(new G4QDecayChan(.820, 111, 3214));
452      DecayDB[ 53].push_back(new G4QDecayChan(.880, 3112, 211));
453      DecayDB[ 53].push_back(new G4QDecayChan(.940, 3222,-211));
454      DecayDB[ 53].push_back(new G4QDecayChan(1.00, 3212, 111));
455    }
456    if(limit<= 54 && nQ>= 54)    // SIGMA_5/2 -
457    {
458      DecayDB[ 54].push_back(new G4QDecayChan(.600, 2112,-321));
459      DecayDB[ 54].push_back(new G4QDecayChan(.660,-211, 3214));
460      DecayDB[ 54].push_back(new G4QDecayChan(.720, 111, 3114));
461      DecayDB[ 54].push_back(new G4QDecayChan(.810, 3212,-211));
462      DecayDB[ 54].push_back(new G4QDecayChan(.900, 3112, 111));
463      DecayDB[ 54].push_back(new G4QDecayChan(1.00, 3122,-211));
464    }
465    if(limit<= 55 && nQ>= 55)    // SIGMA_5/2 0
466    {
467      DecayDB[ 55].push_back(new G4QDecayChan(.300, 2112,-311));
468      DecayDB[ 55].push_back(new G4QDecayChan(.600, 2212,-321));
469      DecayDB[ 55].push_back(new G4QDecayChan(.640, 211, 3114));
470      DecayDB[ 55].push_back(new G4QDecayChan(.680,-211, 3224));
471      DecayDB[ 55].push_back(new G4QDecayChan(.720, 111, 3214));
472      DecayDB[ 55].push_back(new G4QDecayChan(.780, 3112, 211));
473      DecayDB[ 55].push_back(new G4QDecayChan(.840, 3222,-211));
474      DecayDB[ 55].push_back(new G4QDecayChan(.900, 3212, 111));
475      DecayDB[ 55].push_back(new G4QDecayChan(1.00, 3122, 111));
476    }
477    if(limit<= 56 && nQ>= 56)    // SIGMA_5/2 +
478    {
479      DecayDB[ 56].push_back(new G4QDecayChan(.600, 2212,-311));
480      DecayDB[ 56].push_back(new G4QDecayChan(.660, 211, 3214));
481      DecayDB[ 56].push_back(new G4QDecayChan(.720, 111, 3224));
482      DecayDB[ 56].push_back(new G4QDecayChan(.810, 3212, 211));
483      DecayDB[ 56].push_back(new G4QDecayChan(.900, 3222, 111));
484      DecayDB[ 56].push_back(new G4QDecayChan(1.00, 3122, 211));
485    }
486    if(limit<= 57 && nQ>= 57)    // KSI_5/2 -
487    {
488      DecayDB[ 57].push_back(new G4QDecayChan(.400, 3112,-311));
489      DecayDB[ 57].push_back(new G4QDecayChan(.800, 3212,-321));
490      DecayDB[ 57].push_back(new G4QDecayChan(1.00, 3122,-321));
491    }
492    if(limit<= 58 && nQ>= 58)    // KSI_5/2 0
493    {
494      DecayDB[ 58].push_back(new G4QDecayChan(.400, 3212,-311));
495      DecayDB[ 58].push_back(new G4QDecayChan(.800, 3222,-321));
496      DecayDB[ 58].push_back(new G4QDecayChan(1.00, 3122,-311));
497    }
498    if(limit<= 59 && nQ>= 59)    // rho_3 0
499    {
500      DecayDB[ 59].push_back(new G4QDecayChan(.019,311,-313));
501      DecayDB[ 59].push_back(new G4QDecayChan(.038,321,-323));
502      DecayDB[ 59].push_back(new G4QDecayChan(.046,311,-311));
503      DecayDB[ 59].push_back(new G4QDecayChan(.054,321,-321));
504      DecayDB[ 59].push_back(new G4QDecayChan(.224,111, 223));
505      DecayDB[ 59].push_back(new G4QDecayChan(.404,111,-211,213));
506      DecayDB[ 59].push_back(new G4QDecayChan(.584,111, 211,-213));
507      DecayDB[ 59].push_back(new G4QDecayChan(.764,111, 111,113));
508      DecayDB[ 59].push_back(new G4QDecayChan(1.00,211,-211));
509    }
510    if(limit<= 60 && nQ>= 60)    // rho_3 +
511    {
512      DecayDB[ 60].push_back(new G4QDecayChan(.019, 321,-313));
513      DecayDB[ 60].push_back(new G4QDecayChan(.038,-311, 323));
514      DecayDB[ 60].push_back(new G4QDecayChan(.054, 321,-311));
515      DecayDB[ 60].push_back(new G4QDecayChan(.224, 211, 223));
516      DecayDB[ 60].push_back(new G4QDecayChan(.404,211,-211,213));
517      DecayDB[ 60].push_back(new G4QDecayChan(.584,211,211,-213));
518      DecayDB[ 60].push_back(new G4QDecayChan(.764,211,111,113));
519      DecayDB[ 60].push_back(new G4QDecayChan(1.00, 211, 111));
520    }
521    if(limit<= 61 && nQ>= 61)    // omega_3
522    {
523      DecayDB[ 61].push_back(new G4QDecayChan(.020,211,-211,223));
524      DecayDB[ 61].push_back(new G4QDecayChan(.040,111, 111,223));
525      DecayDB[ 61].push_back(new G4QDecayChan(.060, 211,-213));
526      DecayDB[ 61].push_back(new G4QDecayChan(.080,-211, 213));
527      DecayDB[ 61].push_back(new G4QDecayChan(1.00, 111, 113));
528    }
529    if(limit<= 62 && nQ>= 62)    // K_3 0
530    {
531      DecayDB[ 62].push_back(new G4QDecayChan(.030, 111, 315));
532      DecayDB[ 62].push_back(new G4QDecayChan(.060,-211, 325));
533      DecayDB[ 62].push_back(new G4QDecayChan(.340, 311, 331));
534      DecayDB[ 62].push_back(new G4QDecayChan(.400, 111, 313));
535      DecayDB[ 62].push_back(new G4QDecayChan(.520,-211, 323));
536      DecayDB[ 62].push_back(new G4QDecayChan(.620, 311, 113));
537      DecayDB[ 62].push_back(new G4QDecayChan(.820, 321,-213));
538      DecayDB[ 62].push_back(new G4QDecayChan(.940,-211, 321));
539      DecayDB[ 62].push_back(new G4QDecayChan(1.00, 111, 311));
540    }
541    if(limit<= 63 && nQ>= 63)    // K_3 +
542    {
543      DecayDB[ 63].push_back(new G4QDecayChan(.030, 211, 315));
544      DecayDB[ 63].push_back(new G4QDecayChan(.060, 111, 325));
545      DecayDB[ 63].push_back(new G4QDecayChan(.340, 321, 331));
546      DecayDB[ 63].push_back(new G4QDecayChan(.400, 211, 313));
547      DecayDB[ 63].push_back(new G4QDecayChan(.520, 111, 323));
548      DecayDB[ 63].push_back(new G4QDecayChan(.620, 311, 213));
549      DecayDB[ 63].push_back(new G4QDecayChan(.820, 321, 113));
550      DecayDB[ 63].push_back(new G4QDecayChan(.940, 211, 311));
551      DecayDB[ 63].push_back(new G4QDecayChan(1.00, 111, 321));
552    }
553    if(limit<= 64 && nQ>= 64)    // phi_3
554    {
555      DecayDB[ 64].push_back(new G4QDecayChan(.250, 321,-321));
556      DecayDB[ 64].push_back(new G4QDecayChan(.500, 311,-311));
557      DecayDB[ 64].push_back(new G4QDecayChan(.625, 321,-323));
558      DecayDB[ 64].push_back(new G4QDecayChan(.750,-321, 323));
559      DecayDB[ 64].push_back(new G4QDecayChan(.875, 311,-313));
560      DecayDB[ 64].push_back(new G4QDecayChan(1.00,-311, 313));
561    }
562    if(limit<= 65 && nQ>= 65)    // DELTA_7/2 -
563    {
564      DecayDB[ 65].push_back(new G4QDecayChan(.200, 2112,-213 ));
565      DecayDB[ 65].push_back(new G4QDecayChan(.320,-211 , 2114));
566      DecayDB[ 65].push_back(new G4QDecayChan(.500, 111 , 1114));
567      DecayDB[ 65].push_back(new G4QDecayChan(1.00, 2112,-211 ));
568    }
569    if(limit<= 66 && nQ>= 66)    // DELTA_7/2 0
570    {
571      DecayDB[ 66].push_back(new G4QDecayChan(.133, 2112, 113 ));
572      DecayDB[ 66].push_back(new G4QDecayChan(.200, 2212,-213 ));
573      DecayDB[ 66].push_back(new G4QDecayChan(.360,-211 , 2214));
574      DecayDB[ 66].push_back(new G4QDecayChan(.480, 211 , 1114));
575      DecayDB[ 66].push_back(new G4QDecayChan(.500, 111 , 2114));
576      DecayDB[ 66].push_back(new G4QDecayChan(.666, 2212,-211 ));
577      DecayDB[ 66].push_back(new G4QDecayChan(1.00, 2112, 111 ));
578    }
579    if(limit<= 67 && nQ>= 67)    // DELTA_7/2 +
580    {
581      DecayDB[ 67].push_back(new G4QDecayChan(.133, 2112, 213 ));
582      DecayDB[ 67].push_back(new G4QDecayChan(.200, 2212, 113 ));
583      DecayDB[ 67].push_back(new G4QDecayChan(.360,-211 , 2224));
584      DecayDB[ 67].push_back(new G4QDecayChan(.480, 211 , 2114));
585      DecayDB[ 67].push_back(new G4QDecayChan(.500, 111 , 2214));
586      DecayDB[ 67].push_back(new G4QDecayChan(.666, 2112, 211 ));
587      DecayDB[ 67].push_back(new G4QDecayChan(1.00, 2212, 111 ));
588    }
589    if(limit<= 68 && nQ>= 68)    // DELTA_7/2 ++
590    {
591      DecayDB[ 68].push_back(new G4QDecayChan(.200, 2212, 213 ));
592      DecayDB[ 68].push_back(new G4QDecayChan(.320, 211 , 2214));
593      DecayDB[ 68].push_back(new G4QDecayChan(.500, 111 , 2224));
594      DecayDB[ 68].push_back(new G4QDecayChan(1.00, 2212, 211 ));
595    }
596    if(limit<= 69 && nQ>= 69)     // LAMBDA_7/2
597    {
598      DecayDB[ 69].push_back(new G4QDecayChan(.160, 3122, 223 ));
599      DecayDB[ 69].push_back(new G4QDecayChan(.260, 2112,-313 ));
600      DecayDB[ 69].push_back(new G4QDecayChan(.360, 2212,-323 ));
601      DecayDB[ 69].push_back(new G4QDecayChan(.400, 3312, 321 ));
602      DecayDB[ 69].push_back(new G4QDecayChan(.440, 3322, 311 ));
603      DecayDB[ 69].push_back(new G4QDecayChan(.480, 3122, 221 ));
604      DecayDB[ 69].push_back(new G4QDecayChan(.520, 2112,-311 ));
605      DecayDB[ 69].push_back(new G4QDecayChan(.560, 2212,-321 ));
606      DecayDB[ 69].push_back(new G4QDecayChan(.600, 3112, 211 ));
607      DecayDB[ 69].push_back(new G4QDecayChan(.800, 3222,-211 ));
608      DecayDB[ 69].push_back(new G4QDecayChan(1.00, 3212, 111 ));
609    }
610    if(limit<= 70 && nQ>= 70)    // SIGMA_7/2 -
611    {
612      DecayDB[ 70].push_back(new G4QDecayChan(.030, 2112,-323 ));
613      DecayDB[ 70].push_back(new G4QDecayChan(.165,-311 , 1114));
614      DecayDB[ 70].push_back(new G4QDecayChan(.210,-321 , 2114));
615      DecayDB[ 70].push_back(new G4QDecayChan(.390,-211 , 3124));
616      DecayDB[ 70].push_back(new G4QDecayChan(.450,-211 , 3214));
617      DecayDB[ 70].push_back(new G4QDecayChan(.510, 111 , 3114));
618      DecayDB[ 70].push_back(new G4QDecayChan(.540, 311 , 3314));
619      DecayDB[ 70].push_back(new G4QDecayChan(.570,-211 , 3212));
620      DecayDB[ 70].push_back(new G4QDecayChan(.600, 111 , 3112));
621      DecayDB[ 70].push_back(new G4QDecayChan(.780,-321 , 2112));
622      DecayDB[ 70].push_back(new G4QDecayChan(1.00,-211 , 3122));
623    }
624    if(limit<= 71 && nQ>= 71)    // SIGMA_7/2 0
625    {
626      DecayDB[ 71].push_back(new G4QDecayChan(.015, 2112,-313 ));
627      DecayDB[ 71].push_back(new G4QDecayChan(.030, 2212,-321 ));
628      DecayDB[ 71].push_back(new G4QDecayChan(.120,-311 , 2114));
629      DecayDB[ 71].push_back(new G4QDecayChan(.210,-321 , 2214));
630      DecayDB[ 71].push_back(new G4QDecayChan(.390, 111 , 3124));
631      DecayDB[ 71].push_back(new G4QDecayChan(.450,-211 , 3224));
632      DecayDB[ 71].push_back(new G4QDecayChan(.510, 211 , 3114));
633      DecayDB[ 71].push_back(new G4QDecayChan(.525, 311 , 3324));
634      DecayDB[ 71].push_back(new G4QDecayChan(.540, 321 , 3314));
635      DecayDB[ 71].push_back(new G4QDecayChan(.570,-211 , 3222));
636      DecayDB[ 71].push_back(new G4QDecayChan(.600, 211 , 3112));
637      DecayDB[ 71].push_back(new G4QDecayChan(.690,-311 , 2112));
638      DecayDB[ 71].push_back(new G4QDecayChan(.780,-321 , 2212));
639      DecayDB[ 71].push_back(new G4QDecayChan(1.00, 111 , 3122));
640    }
641    if(limit<= 72 && nQ>= 72)    // SIGMA_7/2 +
642    {
643      DecayDB[ 72].push_back(new G4QDecayChan(.030, 2212,-313 ));
644      DecayDB[ 72].push_back(new G4QDecayChan(.165,-321 , 2224));
645      DecayDB[ 72].push_back(new G4QDecayChan(.210,-311 , 2214));
646      DecayDB[ 72].push_back(new G4QDecayChan(.390, 211 , 3124));
647      DecayDB[ 72].push_back(new G4QDecayChan(.450, 211 , 3214));
648      DecayDB[ 72].push_back(new G4QDecayChan(.510, 111 , 3224));
649      DecayDB[ 72].push_back(new G4QDecayChan(.540, 321 , 3324));
650      DecayDB[ 72].push_back(new G4QDecayChan(.570, 211 , 3212));
651      DecayDB[ 72].push_back(new G4QDecayChan(.600, 111 , 3222));
652      DecayDB[ 72].push_back(new G4QDecayChan(.780,-311 , 2212));
653      DecayDB[ 72].push_back(new G4QDecayChan(1.00, 211 , 3122));
654    }
655    if(limit<= 73 && nQ>= 73)    // KSI_7/2 -
656    {
657      DecayDB[ 73].push_back(new G4QDecayChan(.400, 3112,-311));
658      DecayDB[ 73].push_back(new G4QDecayChan(.800, 3212,-321));
659      DecayDB[ 73].push_back(new G4QDecayChan(1.00, 3122,-321));
660    }
661    if(limit<= 74 && nQ>= 74)    // KSI_7/2 0
662    {
663      DecayDB[ 74].push_back(new G4QDecayChan(.400, 3212,-311));
664      DecayDB[ 74].push_back(new G4QDecayChan(.800, 3222,-321));
665      DecayDB[ 74].push_back(new G4QDecayChan(1.00, 3122,-311));
666    }
667    if(limit<= 75 && nQ>= 75)    // OMEGA_7/2 -
668    {
669      DecayDB[ 75].push_back(new G4QDecayChan(.250,-311 , 3314));
670      DecayDB[ 75].push_back(new G4QDecayChan(.500,-321 , 3324));
671      DecayDB[ 75].push_back(new G4QDecayChan(.750, 3312,-313 ));
672      DecayDB[ 75].push_back(new G4QDecayChan(1.00, 3322,-323 ));
673    }
674    if(limit<= 76 && nQ>= 76)    // a_4 0
675    {
676      DecayDB[ 76].push_back(new G4QDecayChan(.200, 311,-311));
677      DecayDB[ 76].push_back(new G4QDecayChan(.400, 321,-321));
678      DecayDB[ 76].push_back(new G4QDecayChan(.600, 111, 221));
679      DecayDB[ 76].push_back(new G4QDecayChan(1.00, 111, 211,-211));
680    }
681    if(limit<= 77 && nQ>= 77)    // a_4 +
682    {
683      DecayDB[ 77].push_back(new G4QDecayChan(.400, 321,-311));
684      DecayDB[ 77].push_back(new G4QDecayChan(.600, 211, 221));
685      DecayDB[ 77].push_back(new G4QDecayChan(.800, 211, 211,-211));
686      DecayDB[ 77].push_back(new G4QDecayChan(1.00, 211, 111, 111));
687    }
688    if(limit<= 78 && nQ>= 78)    // f_4 0
689    {
690      DecayDB[ 78].push_back(new G4QDecayChan(.020, 333, 333));
691      DecayDB[ 78].push_back(new G4QDecayChan(.340, 223, 223));
692      DecayDB[ 78].push_back(new G4QDecayChan(.350, 221, 221));
693      DecayDB[ 78].push_back(new G4QDecayChan(.360, 311,-311));
694      DecayDB[ 78].push_back(new G4QDecayChan(.370, 321,-321));
695      DecayDB[ 78].push_back(new G4QDecayChan(.670, 213,-213));
696      DecayDB[ 78].push_back(new G4QDecayChan(.820, 113, 113));
697      DecayDB[ 78].push_back(new G4QDecayChan(.940, 211,-211));
698      DecayDB[ 78].push_back(new G4QDecayChan(1.00, 111, 111));
699    }
700    if(limit<= 79 && nQ>= 79)    // K_4 0
701    {
702      DecayDB[ 79].push_back(new G4QDecayChan(.060, 333, 313));
703      DecayDB[ 79].push_back(new G4QDecayChan(.260, 223, 313));
704      DecayDB[ 79].push_back(new G4QDecayChan(.380, 313, 113));
705      DecayDB[ 79].push_back(new G4QDecayChan(.400, 323,-213));
706      DecayDB[ 79].push_back(new G4QDecayChan(.500, 223, 311));
707      DecayDB[ 79].push_back(new G4QDecayChan(.550, 311, 113));
708      DecayDB[ 79].push_back(new G4QDecayChan(.600, 321,-213));
709      DecayDB[ 79].push_back(new G4QDecayChan(.700, 311, 113));
710      DecayDB[ 79].push_back(new G4QDecayChan(.800, 321,-213));
711      DecayDB[ 79].push_back(new G4QDecayChan(.900, 311, 111));
712      DecayDB[ 79].push_back(new G4QDecayChan(1.00, 321,-211));
713    }
714    if(limit<= 80 && nQ>= 80)    // K_4 +
715    {
716      DecayDB[ 80].push_back(new G4QDecayChan(.060, 333, 323));
717      DecayDB[ 80].push_back(new G4QDecayChan(.260, 223, 323));
718      DecayDB[ 80].push_back(new G4QDecayChan(.380, 313, 213));
719      DecayDB[ 80].push_back(new G4QDecayChan(.400, 323, 113));
720      DecayDB[ 80].push_back(new G4QDecayChan(.500, 223, 321));
721      DecayDB[ 80].push_back(new G4QDecayChan(.550, 321, 113));
722      DecayDB[ 80].push_back(new G4QDecayChan(.600, 311, 213));
723      DecayDB[ 80].push_back(new G4QDecayChan(.700, 311, 211));
724      DecayDB[ 80].push_back(new G4QDecayChan(.800, 321, 111));
725      DecayDB[ 80].push_back(new G4QDecayChan(.900, 311, 211));
726      DecayDB[ 80].push_back(new G4QDecayChan(1.00, 321, 111));
727    }
728    if(limit<=81&&nQ>=81)DecayDB[81].push_back(new G4QDecayChan(1., 333,333));//phi_4(2300)
729    if(limit<=82&&nQ>=82)DecayDB[82].push_back(new G4QDecayChan(1.,2212, 2224));//pDelta++
730    if(limit<=83&&nQ>=83)DecayDB[83].push_back(new G4QDecayChan(1.,2112, 1114));//nDelta-
731    if(limit<=84&&nQ>=84)DecayDB[84].push_back(new G4QDecayChan(1.,2224, 2224));//D++D++
732    if(limit<=85&&nQ>=85)DecayDB[85].push_back(new G4QDecayChan(1.,1114, 1114));//Del-Del-
733    if(limit<=86&&nQ>=86)DecayDB[86].push_back(new G4QDecayChan(1.,2212,2212,2224));//ppD++
734    if(limit<=87&&nQ>=87)DecayDB[87].push_back(new G4QDecayChan(1.,2112,2112,1114));//nnD-
735    if(limit<=88&&nQ>=88)DecayDB[88].push_back(new G4QDecayChan(1.,2212,2224,2224));//p2D++
736    if(limit<=89&&nQ>=89)DecayDB[89].push_back(new G4QDecayChan(1.,2112,1114,1114));//nD-D-
737    //if(limit<=90&&nQ>=90)DecayDB[90] = 0; // neutron (n) as a quark exchange fragment
738    //if(limit<=91&&nQ>=91)DecayDB[91] = 0; // proton (p)  as a quark exchange fragment
739    //if(limit<=92&&nQ>=92)DecayDB[92] = 0; // neutron (L/Sigma0) as aQuarkExchangeFragment
740    //if(limit<=93&&nQ>=93)DecayDB[93] = 0; // neutron (Sigma-) as a quarkExchangeFragment
741    //if(limit<=94&&nQ>=94)DecayDB[94] = 0; // neutron (Sigma+) as a quarkExchangeFragment
742    //if(limit<=95&&nQ>=95)DecayDB[95] = 0; // neutron (Xi-) as a quark exchange fragment
743    //if(limit<=96&&nQ>=96)DecayDB[96] = 0; // neutron (Xi0) as a quark exchange fragment
744    //if(limit<=97&&nQ>=97)DecayDB[97] = 0; // neutron (Omega-) as a quarkExchangeFragment
745    if(limit<=98&&nQ>=98)DecayDB[98].push_back(new G4QDecayChan(1.,2112, 2112)); //nn (src)
746    if(limit<=99&&nQ>=99)DecayDB[99].push_back(new G4QDecayChan(1.,2212, 2112)); //d/pn(?)
747    if(limit<=100&&nQ>=100)DecayDB[100].push_back(new G4QDecayChan(1.,2212,2212));//pp(src)
748    if(limit<=101&&nQ>=101)DecayDB[101].push_back(new G4QDecayChan(1.,3122,2112));//Ln
749    if(limit<=102&&nQ>=102)DecayDB[102].push_back(new G4QDecayChan(1.,3122,2212));//Lp
750    if(limit<=104&&nQ>=104)DecayDB[104].push_back(new G4QDecayChan(1.,3112,2112));//nSig-
751    if(limit<=103&&nQ>=103)DecayDB[103].push_back(new G4QDecayChan(1.,3122,3122));//LL
752    if(limit<=105&&nQ>=105)DecayDB[105].push_back(new G4QDecayChan(1.,3222,2212));//pSig+
753    //if(limit<=106&&nQ>=106)DecayDB[106] = 0; // t
754    //if(limit<=107&&nQ>=107)DecayDB[107] = 0; // He3
755    //Lnn=>Lambda+n+n decay to avoid the final state Hypernucleus
756    if(limit<=108&&nQ>=108)DecayDB[108].push_back(new G4QDecayChan(1.,3122,2112,2112));
757    if(limit<=109&&nQ>=109)DecayDB[109].push_back(new G4QDecayChan(1.,3122,90001001));// Ld
758    //Lpp=>Lambda+p+p decay to avoid the final state Hypernucleus
759    if(limit<=110&&nQ>=110)DecayDB[110].push_back(new G4QDecayChan(1.,3122,2212,2212));
760    //LLn=>Lambda+Lambda+n decay to avoid the final state Hypernucleus
761    if(limit<=111&&nQ>=111)DecayDB[111].push_back(new G4QDecayChan(1.,3122,3122,2112));
762    //LLp=>Lambda+Lambda+p decay to avoid the final state Hypernucleus
763    if(limit<=112&&nQ>=112)DecayDB[112].push_back(new G4QDecayChan(1.,3122,3122,2212));
764    // nnSigma-=>n+n+Sigma- decay to avoid the final state Hypernucleus
765    if(limit<=113&&nQ>=113)DecayDB[113].push_back(new G4QDecayChan(1.,2112,2112,3112));
766    // ------- Nuclear fragments
767    //if(limit<=114 && nQ>=114)
768    //{
769    //  if(limit<114) limit=101;
770    //  for (int i=limit; i<nQ; i++) DecayDB[i] = 0;
771    //}
772    //Update the limit
773    limit=nQ+1;
774#ifdef debug
775    G4cout<<"G4QParticle::InitDecayVector: limit is set to "<<limit<<G4endl;
776#endif
777  }
778  //if(!nQ)G4cout<<"G4QParticle::InitDecayVector:Q=0,nD="<<DecayDB[abs(nQ)].size()<<G4endl;
779  return DecayDB[std::abs(nQ)];
780}
781
782// Initialize the Particle by a Q Code
783void G4QParticle::InitQParticle(G4int theQCode)
784//    =============================================
785{
786  aQPDG.InitByQCode(theQCode);
787  aQuarkCont = aQPDG.GetQuarkContent();
788  aDecay     = InitDecayVector(theQCode);
789  //if(!theQCode)G4cout<<"G4QPar::InitQP:PDG="<<GetPDGCode()<<",n="<<aDecay.size()<<G4endl;
790}
791
792// Initialize the Particle by a PDG Code
793void G4QParticle::InitPDGParticle(G4int thePDGCode)
794//    =============================================
795{
796  aQPDG      = G4QPDGCode(thePDGCode);
797  aQuarkCont = aQPDG.GetQuarkContent();
798  aDecay     = InitDecayVector(aQPDG.GetQCode());
799}
Note: See TracBrowser for help on using the repository browser.