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

Last change on this file since 967 was 962, checked in by garnier, 15 years ago

update processes

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