Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QParticle Class Reference

#include <G4QParticle.hh>

Public Member Functions

 G4QParticle ()
 
 G4QParticle (G4bool f, G4int theQCode)
 
 G4QParticle (G4int thePDG)
 
 G4QParticle (const G4QParticle &right)
 
 G4QParticle (G4QParticle *right)
 
 ~G4QParticle ()
 
const G4QParticleoperator= (const G4QParticle &right)
 
G4bool operator== (const G4QParticle &rhs) const
 
G4bool operator!= (const G4QParticle &rhs) const
 
G4QPDGCode GetQPDG () const
 
G4int GetPDGCode () const
 
G4int GetQCode () const
 
G4int GetSpin () const
 
G4int GetCharge () const
 
G4int GetStrange () const
 
G4int GetBaryNum () const
 
G4QContent GetQContent ()
 
G4QDecayChanVector GetDecayVector ()
 
G4double GetMass ()
 
G4double GetWidth ()
 
G4QDecayChanVector InitDecayVector (G4int Q)
 
void InitPDGParticle (G4int thePDGCode)
 
void InitQParticle (G4int theQCode)
 
G4double MinMassOfFragm ()
 

Detailed Description

Definition at line 45 of file G4QParticle.hh.

Constructor & Destructor Documentation

◆ G4QParticle() [1/5]

G4QParticle::G4QParticle ( )

Definition at line 45 of file G4QParticle.cc.

46{
47#ifdef debug
48 G4cout<<"G4QParticle::Constr: Default constructor is called"<<G4endl;
49#endif
50}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ G4QParticle() [2/5]

G4QParticle::G4QParticle ( G4bool  f,
G4int  theQCode 
)

Definition at line 59 of file G4QParticle.cc.

60{
61 aQPDG = G4QPDGCode(f,theQCode);
62#ifdef debug
63 G4cout<<"G4QParticle::Constr: PDG="<<aQPDG.GetPDGCode()<<G4endl;
64#endif
65 aQuarkCont = aQPDG.GetQuarkContent();
66 aDecay = InitDecayVector(theQCode);
67}
G4QContent GetQuarkContent() const
Definition: G4QPDGCode.cc:2057
G4int GetPDGCode() const
Definition: G4QPDGCode.hh:326
G4QDecayChanVector InitDecayVector(G4int Q)
Definition: G4QParticle.cc:153

◆ G4QParticle() [3/5]

G4QParticle::G4QParticle ( G4int  thePDG)

Definition at line 52 of file G4QParticle.cc.

53{
54 aQPDG = G4QPDGCode(thePDG);
55 aQuarkCont = aQPDG.GetQuarkContent();
56 aDecay = InitDecayVector(aQPDG.GetQCode());
57}
G4int GetQCode() const
Definition: G4QPDGCode.hh:327

◆ G4QParticle() [4/5]

G4QParticle::G4QParticle ( const G4QParticle right)

Definition at line 69 of file G4QParticle.cc.

70{
71 aQPDG = right.aQPDG;
72 //aDecay (Vector)
73 G4int nD = right.aDecay.size();
74 if(nD) for(G4int id=0; id<nD; id++)
75 {
76 G4QDecayChan* curD = new G4QDecayChan(right.aDecay[id]);
77 aDecay.push_back(curD);
78 }
79
80 aQuarkCont = right.aQuarkCont;
81}
int G4int
Definition: G4Types.hh:66

◆ G4QParticle() [5/5]

G4QParticle::G4QParticle ( G4QParticle right)

Definition at line 83 of file G4QParticle.cc.

84{
85 aQPDG = right->aQPDG;
86 //aDecay (Vector)
87 G4int nD = right->aDecay.size();
88 if(nD) for(G4int id=0; id<nD; id++)
89 {
90 G4QDecayChan* curD = new G4QDecayChan(right->aDecay[id]);
91 aDecay.push_back(curD);
92 }
93
94 aQuarkCont = right->aQuarkCont;
95}

◆ ~G4QParticle()

G4QParticle::~G4QParticle ( )

Definition at line 97 of file G4QParticle.cc.

98{
99 G4int nDC=aDecay.size();
100 //G4cout<<"G4QParticle::Destructor: Before nDC="<<nDC<<G4endl; // TMP
101 if(nDC) std::for_each(aDecay.begin(), aDecay.end(), DeleteQDecayChan());
102 //G4cout<<"G4QParticle::Destructor: After"<<G4endl; // TMP
103 aDecay.clear();
104}

Member Function Documentation

◆ GetBaryNum()

G4int G4QParticle::GetBaryNum ( ) const
inline

Definition at line 107 of file G4QParticle.hh.

107{return aQuarkCont.GetBaryonNumber();}
G4int GetBaryonNumber() const
Definition: G4QContent.cc:1182

◆ GetCharge()

G4int G4QParticle::GetCharge ( ) const
inline

Definition at line 105 of file G4QParticle.hh.

105{return aQuarkCont.GetCharge();}
G4int GetCharge() const
Definition: G4QContent.cc:1159

◆ GetDecayVector()

G4QDecayChanVector G4QParticle::GetDecayVector ( )
inline

Definition at line 109 of file G4QParticle.hh.

109{return aDecay;}

Referenced by G4Quasmon::DecayQHadron(), and operator<<().

◆ GetMass()

G4double G4QParticle::GetMass ( )
inline

Definition at line 110 of file G4QParticle.hh.

110{return aQPDG.GetMass();}
G4double GetMass()
Definition: G4QPDGCode.cc:693

Referenced by MinMassOfFragm(), and operator<<().

◆ GetPDGCode()

G4int G4QParticle::GetPDGCode ( ) const
inline

Definition at line 103 of file G4QParticle.hh.

103{return aQPDG.GetPDGCode();}

◆ GetQCode()

G4int G4QParticle::GetQCode ( ) const
inline

Definition at line 102 of file G4QParticle.hh.

102{return aQPDG.GetQCode();}

◆ GetQContent()

G4QContent G4QParticle::GetQContent ( )
inline

Definition at line 108 of file G4QParticle.hh.

108{return aQuarkCont;}

Referenced by operator<<().

◆ GetQPDG()

G4QPDGCode G4QParticle::GetQPDG ( ) const
inline

Definition at line 101 of file G4QParticle.hh.

101{return aQPDG;}

Referenced by operator<<().

◆ GetSpin()

G4int G4QParticle::GetSpin ( ) const
inline

Definition at line 104 of file G4QParticle.hh.

104{return aQPDG.GetSpin();}
G4int GetSpin() const
Definition: G4QPDGCode.hh:330

Referenced by operator<<().

◆ GetStrange()

G4int G4QParticle::GetStrange ( ) const
inline

Definition at line 106 of file G4QParticle.hh.

106{return aQuarkCont.GetStrangeness();}
G4int GetStrangeness() const
Definition: G4QContent.hh:184

◆ GetWidth()

G4double G4QParticle::GetWidth ( )
inline

Definition at line 111 of file G4QParticle.hh.

111{return aQPDG.GetWidth();}
G4double GetWidth()
Definition: G4QPDGCode.cc:740

Referenced by MinMassOfFragm(), and operator<<().

◆ InitDecayVector()

G4QDecayChanVector G4QParticle::InitDecayVector ( G4int  Q)

Definition at line 153 of file G4QParticle.cc.

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

Referenced by G4QParticle(), InitPDGParticle(), and InitQParticle().

◆ InitPDGParticle()

void G4QParticle::InitPDGParticle ( G4int  thePDGCode)

Definition at line 789 of file G4QParticle.cc.

790{
791 aQPDG = G4QPDGCode(thePDGCode);
792 aQuarkCont = aQPDG.GetQuarkContent();
793 aDecay = InitDecayVector(aQPDG.GetQCode());
794}

◆ InitQParticle()

void G4QParticle::InitQParticle ( G4int  theQCode)

Definition at line 780 of file G4QParticle.cc.

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}
void InitByQCode(G4int QCode)
Definition: G4QPDGCode.hh:356

◆ MinMassOfFragm()

G4double G4QParticle::MinMassOfFragm ( )
inline

Definition at line 113 of file G4QParticle.hh.

114{
115 G4int nCh=aDecay.size();
116 G4double mass=GetMass();
117 G4double min=mass;
118 if(nCh)
119 {
120 min=aDecay[0]->GetMinMass();
121 if(nCh>1) for(G4int j=1; j<nCh; j++)
122 {
123 G4double next=aDecay[j]->GetMinMass();
124 if(next<min) min=next;
125 }
126 }
127 G4double w=GetWidth();
128 G4double lim=mass+.001;
129 if(w) lim-=1.5*w;
130 if(min<lim) min=lim;
131 return min;
132}
double G4double
Definition: G4Types.hh:64
G4double GetMass()
Definition: G4QParticle.hh:110
G4double GetWidth()
Definition: G4QParticle.hh:111

Referenced by G4Quasmon::DecayQHadron(), G4QFragmentation::Fragment(), G4QIonIonCollision::Fragment(), and G4QHadron::RandomizeMass().

◆ operator!=()

G4bool G4QParticle::operator!= ( const G4QParticle rhs) const
inline

Definition at line 99 of file G4QParticle.hh.

99{return this!=&rhs;}

◆ operator=()

const G4QParticle & G4QParticle::operator= ( const G4QParticle right)

Definition at line 107 of file G4QParticle.cc.

108{
109 if(this != &right) // Beware of self assignment
110 {
111 aQPDG = right.aQPDG;
112 //aDecay (Vector)
113 G4int iD = aDecay.size();
114 if(iD) for(G4int jd=0; jd<iD; jd++) delete aDecay[jd];
115 aDecay.clear();
116 G4int nD = right.aDecay.size();
117 if(nD) for(G4int id=0; id<nD; id++)
118 {
119 G4QDecayChan* curD = new G4QDecayChan(right.aDecay[id]);
120 aDecay.push_back(curD);
121 }
122
123 aQuarkCont = right.aQuarkCont;
124 }
125 return *this;
126}

◆ operator==()

G4bool G4QParticle::operator== ( const G4QParticle rhs) const
inline

Definition at line 98 of file G4QParticle.hh.

98{return this==&rhs;}

The documentation for this class was generated from the following files: