Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VDecayChannel.cc
Go to the documentation of this file.
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// G4VDecayChannel
27//
28// Author: H.Kurashige, 27 July 1996
29// --------------------------------------------------------------------
30
32#include "G4SystemOfUnits.hh"
33#include "G4ParticleTable.hh"
34#include "G4DecayTable.hh"
35#include "G4DecayProducts.hh"
36#include "G4VDecayChannel.hh"
37#include "G4AutoLock.hh"
38#include "Randomize.hh"
39
41
42// --------------------------------------------------------------------
44 : parent_polarization()
45{
46 // set pointer to G4ParticleTable (static and singleton object)
48}
49
50// --------------------------------------------------------------------
52 : kinematics_name(aName), parent_polarization(), verboseLevel(verbose)
53{
54 // set pointer to G4ParticleTable (static and singleton object)
56}
57
58// --------------------------------------------------------------------
60 const G4String& theParentName,
61 G4double theBR,
62 G4int theNumberOfDaughters,
63 const G4String& theDaughterName1,
64 const G4String& theDaughterName2,
65 const G4String& theDaughterName3,
66 const G4String& theDaughterName4 )
67 : kinematics_name(aName), rbranch(theBR), parent_polarization(),
68 numberOfDaughters(theNumberOfDaughters)
69{
70 // set pointer to G4ParticleTable (static and singleton object)
72
73 // parent name
74 parent_name = new G4String(theParentName);
75
76 // cleate array
78 for (G4int index=0; index<numberOfDaughters; ++index)
79 {
80 daughters_name[index] = nullptr;
81 }
82
83 // daughters' name
84 if (numberOfDaughters>0) daughters_name[0] = new G4String(theDaughterName1);
85 if (numberOfDaughters>1) daughters_name[1] = new G4String(theDaughterName2);
86 if (numberOfDaughters>2) daughters_name[2] = new G4String(theDaughterName3);
87 if (numberOfDaughters>3) daughters_name[3] = new G4String(theDaughterName4);
88
89 if (rbranch <0. ) rbranch = 0.0;
90 else if (rbranch >1.0 ) rbranch = 1.0;
91}
92
93// --------------------------------------------------------------------
95{
98 rbranch = right.rbranch;
99 rangeMass = right.rangeMass;
100
101 // copy parent name
102 parent_name = new G4String(*right.parent_name);
103
104 // create array
106
107 daughters_name = nullptr;
108 if ( numberOfDaughters >0 )
109 {
111 // copy daughters name
112 for (G4int index=0; index<numberOfDaughters; ++index)
113 {
114 daughters_name[index] = new G4String(*right.daughters_name[index]);
115 }
116 }
117
118 // particle table
120
122
125}
126
127// --------------------------------------------------------------------
129{
130 if (this != &right)
131 {
134 rbranch = right.rbranch;
135 rangeMass = right.rangeMass;
137 // copy parent name
138 delete parent_name;
139 parent_name = new G4String(*right.parent_name);
140
141 // clear daughters_name array
143
144 // recreate array
146 if ( numberOfDaughters >0 )
147 {
149 // copy daughters name
150 for (G4int index=0; index<numberOfDaughters; ++index)
151 {
152 daughters_name[index] = new G4String(*right.daughters_name[index]);
153 }
154 }
155 }
156
157 // particle table
159
162
163 return *this;
164}
165
166// --------------------------------------------------------------------
168{
170 delete parent_name;
171 parent_name = nullptr;
172 delete [] G4MT_daughters_mass;
173 G4MT_daughters_mass = nullptr;
174 delete [] G4MT_daughters_width;
175 G4MT_daughters_width = nullptr;
178}
179
180// --------------------------------------------------------------------
182{
184 if ( daughters_name != nullptr)
185 {
186 if (numberOfDaughters>0)
187 {
188#ifdef G4VERBOSE
189 if (verboseLevel>1)
190 {
191 G4cout << "G4VDecayChannel::ClearDaughtersName() "
192 << " for " << *parent_name << G4endl;
193 }
194#endif
195 for (G4int index=0; index<numberOfDaughters; ++index)
196 {
197 delete daughters_name[index];
198 }
199 }
200 delete [] daughters_name;
201 daughters_name = nullptr;
202 }
203
204 delete [] G4MT_daughters;
205 delete [] G4MT_daughters_mass;
206 delete [] G4MT_daughters_width;
207 G4MT_daughters_width = nullptr;
208 G4MT_daughters = nullptr;
209 G4MT_daughters_mass = nullptr;
210
212}
213
214// --------------------------------------------------------------------
216{
217 if (size >0)
218 {
219 // remove old contents
221 // cleate array
222 daughters_name = new G4String*[size];
223 for (G4int index=0; index<size; ++index)
224 {
225 daughters_name[index] = nullptr;
226 }
227 numberOfDaughters = size;
228 }
229}
230
231// --------------------------------------------------------------------
232void G4VDecayChannel::SetDaughter(G4int anIndex, const G4String& particle_name)
233{
234 // check numberOfDaughters is positive
235 if (numberOfDaughters<=0)
236 {
237#ifdef G4VERBOSE
238 if (verboseLevel>0)
239 {
240 G4cout << "G4VDecayChannel::SetDaughter() - "
241 << "Number of daughters is not defined" << G4endl;
242 }
243#endif
244 return;
245 }
246
247 // An analysis of this code, shows that this method is called
248 // only in the constructor of derived classes.
249 // The general idea of this method is probably to support
250 // the possibility to re-define daughters on the fly, however
251 // this design is extremely problematic for MT mode, we thus
252 // require (as practically happens) that the method is called only
253 // at construction, i.e. when G4MT_daughters == 0
254 // moreover this method can be called only after SetNumberOfDaugthers()
255 // has been called (see previous if), in such a case daughters_name != nullptr
256 //
257 if ( daughters_name == nullptr )
258 {
259 G4Exception("G4VDecayChannel::SetDaughter()", "PART112", FatalException,
260 "Trying to add a daughter without specifying number of secondaries!");
261 return;
262 }
263 if ( G4MT_daughters != nullptr )
264 {
265 G4Exception("G4VDecayChannel::SetDaughter()", "PART111", FatalException,
266 "Trying to modify a daughter of a decay channel, \
267 but decay channel already has daughters.");
268 return;
269 }
270
271 // check an index
272 if ( (anIndex<0) || (anIndex>=numberOfDaughters) )
273 {
274#ifdef G4VERBOSE
275 if (verboseLevel>0)
276 {
277 G4cout << "G4VDecayChannel::SetDaughter() - "
278 << "index out of range " << anIndex << G4endl;
279 }
280#endif
281 }
282 else
283 {
284 // fill the name
285 daughters_name[anIndex] = new G4String(particle_name);
286#ifdef G4VERBOSE
287 if (verboseLevel>1)
288 {
289 G4cout << "G4VDecayChannel::SetDaughter[" << anIndex <<"] :";
290 G4cout << daughters_name[anIndex] << ":"
291 << *daughters_name[anIndex]<<G4endl;
292 }
293#endif
294 }
295}
296
297// --------------------------------------------------------------------
299SetDaughter(G4int anIndex, const G4ParticleDefinition* parent_type)
300{
301 if (parent_type != nullptr)
302 SetDaughter(anIndex, parent_type->GetParticleName());
303}
304
305// --------------------------------------------------------------------
306void G4VDecayChannel::FillDaughters()
307{
309
310 // Double check, check again if another thread has already filled this, in
311 // case do not need to do anything
312 if ( G4MT_daughters != nullptr ) return;
313
314 G4int index;
315
316#ifdef G4VERBOSE
317 if (verboseLevel>1) G4cout << "G4VDecayChannel::FillDaughters()" << G4endl;
318#endif
319 if (G4MT_daughters != nullptr)
320 {
321 delete [] G4MT_daughters;
322 G4MT_daughters = nullptr;
323 }
324
325 // parent mass
327 G4double parentmass = G4MT_parent->GetPDGMass();
328
329 //
330 G4double sumofdaughtermass = 0.0;
331 G4double sumofdaughterwidthsq = 0.0;
332
333 if ((numberOfDaughters <=0) || (daughters_name == nullptr) )
334 {
335#ifdef G4VERBOSE
336 if (verboseLevel>0)
337 {
338 G4cout << "G4VDecayChannel::FillDaughters() - "
339 << "[ " << G4MT_parent->GetParticleName() << " ]"
340 << "numberOfDaughters is not defined yet";
341 }
342#endif
343 G4MT_daughters = nullptr;
344 G4Exception("G4VDecayChannel::FillDaughters()",
345 "PART011", FatalException,
346 "Cannot fill daughters: numberOfDaughters is not defined yet");
347 }
348
349 // create and set the array of pointers to daughter particles
351 if (G4MT_daughters_mass != nullptr) delete [] G4MT_daughters_mass;
352 if (G4MT_daughters_width != nullptr) delete [] G4MT_daughters_width;
355 // loop over all daughters
356 for (index=0; index<numberOfDaughters; ++index)
357 {
358 if (daughters_name[index] == nullptr)
359 {
360 // daughter name is not defined
361#ifdef G4VERBOSE
362 if (verboseLevel>0)
363 {
364 G4cout << "G4VDecayChannel::FillDaughters() - "
365 << "[ " << G4MT_parent->GetParticleName() << " ]"
366 << index << "-th daughter is not defined yet" << G4endl;
367 }
368#endif
369 G4MT_daughters[index] = nullptr;
370 G4Exception("G4VDecayChannel::FillDaughters()",
371 "PART011", FatalException,
372 "Cannot fill daughters: name of daughter is not defined yet");
373 }
374 // search daughter particles in the particle table
376 if (G4MT_daughters[index] == nullptr )
377 {
378 // cannot find the daughter particle
379#ifdef G4VERBOSE
380 if (verboseLevel>0)
381 {
382 G4cout << "G4VDecayChannel::FillDaughters() - "
383 << "[ " << G4MT_parent->GetParticleName() << " ]"
384 << index << ":" << *daughters_name[index]
385 << " is not defined !!" << G4endl;
386 G4cout << " The BR of this decay mode is set to zero." << G4endl;
387 }
388#endif
389 SetBR(0.0);
390 return;
391 }
392#ifdef G4VERBOSE
393 if (verboseLevel>1)
394 {
395 G4cout << index << ":" << *daughters_name[index];
396 G4cout << ":" << G4MT_daughters[index] << G4endl;
397 }
398#endif
400 G4double d_width = G4MT_daughters[index]->GetPDGWidth();
401 G4MT_daughters_width[index] = d_width;
402 sumofdaughtermass += G4MT_daughters[index]->GetPDGMass();
403 sumofdaughterwidthsq += d_width*d_width;
404 } // end loop over all daughters
405
406 // check sum of daghter mass
407 G4double widthMass = std::sqrt(G4MT_parent->GetPDGWidth()
409 +sumofdaughterwidthsq);
410 if ( (G4MT_parent->GetParticleType() != "nucleus")
411 && (numberOfDaughters !=1)
412 && (sumofdaughtermass > parentmass + rangeMass*widthMass) )
413 {
414 // !!! illegal mass !!!
415#ifdef G4VERBOSE
416 if (GetVerboseLevel()>0)
417 {
418 G4cout << "G4VDecayChannel::FillDaughters() - "
419 << "[ " << G4MT_parent->GetParticleName() << " ]"
420 << " Energy/Momentum conserevation breaks " << G4endl;
421 if (GetVerboseLevel()>1)
422 {
423 G4cout << " parent:" << *parent_name
424 << " mass:" << parentmass/GeV << "[GeV/c/c]" << G4endl;
425 for (index=0; index<numberOfDaughters; ++index)
426 {
427 G4cout << " daughter " << index << ":" << *daughters_name[index]
428 << " mass:" << G4MT_daughters[index]->GetPDGMass()/GeV
429 << "[GeV/c/c]" << G4endl;
430 }
431 }
432 }
433#endif
434 }
435}
436
437// --------------------------------------------------------------------
438void G4VDecayChannel::FillParent()
439{
440 G4AutoLock lock(&parentMutex);
441
442 // Double check, check again if another thread has already filled this, in
443 // case do not need to do anything
444 if ( G4MT_parent != nullptr ) return;
445
446 if (parent_name == nullptr)
447 {
448 // parent name is not defined
449#ifdef G4VERBOSE
450 if (verboseLevel>0)
451 {
452 G4cout << "G4VDecayChannel::FillParent() - "
453 << "parent name is not defined !!" << G4endl;
454 }
455#endif
456 G4MT_parent = nullptr;
457 G4Exception("G4VDecayChannel::FillParent()",
458 "PART012", FatalException,
459 "Cannot fill parent: parent name is not defined yet");
460 return;
461 }
462
463 // search parent particle in the particle table
465 if (G4MT_parent == nullptr)
466 {
467 // parent particle does not exist
468#ifdef G4VERBOSE
469 if (verboseLevel>0)
470 {
471 G4cout << "G4VDecayChannel::FillParent() - "
472 << *parent_name << " does not exist !!" << G4endl;
473 }
474#endif
475 G4Exception("G4VDecayChannel::FillParent()",
476 "PART012", FatalException,
477 "Cannot fill parent: parent does not exist");
478 return;
479 }
481}
482
483// --------------------------------------------------------------------
485{
486 if (parent_type != nullptr) SetParent(parent_type->GetParticleName());
487}
488
489// --------------------------------------------------------------------
491{
492 // determine angular momentum
493
494 // fill pointers to daughter particles if not yet set
496
497 const G4int PiSpin = G4MT_parent->GetPDGiSpin();
498 const G4int PParity = G4MT_parent->GetPDGiParity();
499 if (2==numberOfDaughters) // up to now we can only handle two particle decays
500 {
501 const G4int D1iSpin = G4MT_daughters[0]->GetPDGiSpin();
502 const G4int D1Parity = G4MT_daughters[0]->GetPDGiParity();
503 const G4int D2iSpin = G4MT_daughters[1]->GetPDGiSpin();
504 const G4int D2Parity = G4MT_daughters[1]->GetPDGiParity();
505 const G4int MiniSpin = std::abs (D1iSpin - D2iSpin);
506 const G4int MaxiSpin = D1iSpin + D2iSpin;
507 const G4int lMax = (PiSpin+D1iSpin+D2iSpin)/2; // l is always int
508 G4int lMin;
509#ifdef G4VERBOSE
510 if (verboseLevel>1)
511 {
512 G4cout << "iSpin: " << PiSpin << " -> " << D1iSpin
513 << " + " << D2iSpin << G4endl;
514 G4cout << "2*jmin, 2*jmax, lmax " << MiniSpin << " " << MaxiSpin
515 << " " << lMax << G4endl;
516 }
517#endif
518 for (G4int j=MiniSpin; j<=MaxiSpin; j+=2) // loop over all possible
519 { // spin couplings
520 lMin = std::abs(PiSpin-j)/2;
521#ifdef G4VERBOSE
522 if (verboseLevel>1)
523 G4cout << "-> checking 2*j=" << j << G4endl;
524#endif
525 for (G4int l=lMin; l<=lMax; ++l)
526 {
527#ifdef G4VERBOSE
528 if (verboseLevel>1)
529 G4cout << " checking l=" << l << G4endl;
530#endif
531 if (l%2==0)
532 {
533 if (PParity == D1Parity*D2Parity) // check parity for this l
534 return l;
535 }
536 else
537 {
538 if (PParity == -1*D1Parity*D2Parity) // check parity for this l
539 return l;
540 }
541 }
542 }
543 }
544 else
545 {
546 G4Exception("G4VDecayChannel::GetAngularMomentum()",
547 "PART111", JustWarning,
548 "Sorry, can't handle 3 particle decays (up to now)");
549 return 0;
550 }
551 G4Exception ("G4VDecayChannel::GetAngularMomentum()",
552 "PART111", JustWarning,
553 "Can't find angular momentum for this decay");
554 return 0;
555}
556
557// --------------------------------------------------------------------
559{
560 G4cout << " BR: " << rbranch << " [" << kinematics_name << "]";
561 G4cout << " : " ;
562 for (G4int index=0; index<numberOfDaughters; ++index)
563 {
564 if(daughters_name[index] != nullptr)
565 {
566 G4cout << " " << *(daughters_name[index]);
567 }
568 else
569 {
570 G4cout << " not defined ";
571 }
572 }
573 G4cout << G4endl;
574}
575
576// --------------------------------------------------------------------
577const G4String& G4VDecayChannel::GetNoName() const
578{
579 return noName;
580}
581
582// --------------------------------------------------------------------
584DynamicalMass(G4double massPDG, G4double width, G4double maxDev ) const
585{
586 if (width<=0.0) return massPDG;
587 if (maxDev >rangeMass) maxDev = rangeMass;
588 if (maxDev <=-1.*rangeMass) return massPDG; // cannot calculate
589
590 G4double x = G4UniformRand()*(maxDev+rangeMass) - rangeMass;
592 const std::size_t MAX_LOOP=10000;
593 for (std::size_t loop_counter=0; loop_counter <MAX_LOOP; ++loop_counter)
594 {
595 if ( y * (width*width*x*x + massPDG*massPDG*width*width)
596 <= massPDG*massPDG*width*width ) break;
597 x = G4UniformRand()*(maxDev+rangeMass) - rangeMass;
598 y = G4UniformRand();
599 }
600 G4double mass = massPDG + x*width;
601 return mass;
602}
603
604// --------------------------------------------------------------------
606{
607 G4double sumOfDaughterMassMin = 0.0;
610 // skip one body decay
611 if (numberOfDaughters==1) return true;
612
613 for (G4int index=0; index<numberOfDaughters; ++index)
614 {
615 sumOfDaughterMassMin +=
617 }
618 return (parentMass >= sumOfDaughterMassMin);
619}
620
621// --------------------------------------------------------------------
623{
624 rbranch = value;
625 if (rbranch <0. ) rbranch = 0.0;
626 else if (rbranch >1.0 ) rbranch = 1.0;
627}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4MUTEXDESTROY(mutex)
Definition: G4Threading.hh:90
#define G4MUTEXINIT(mutex)
Definition: G4Threading.hh:87
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
const G4String & GetParticleType() const
G4double GetPDGWidth() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4ParticleDefinition ** G4MT_daughters
virtual ~G4VDecayChannel()
G4String * parent_name
void SetBR(G4double value)
G4String ** daughters_name
G4double DynamicalMass(G4double massPDG, G4double width, G4double maxDev=1.0) const
G4int GetVerboseLevel() const
G4double G4MT_parent_mass
void SetNumberOfDaughters(G4int value)
static const G4String noName
G4ParticleDefinition * G4MT_parent
virtual G4bool IsOKWithParentMass(G4double parentMass)
void CheckAndFillDaughters()
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
G4double * G4MT_daughters_mass
G4VDecayChannel & operator=(const G4VDecayChannel &)
G4int GetAngularMomentum()
G4ParticleTable * particletable
G4ThreeVector parent_polarization
G4String kinematics_name
G4double * G4MT_daughters_width
void SetParent(const G4ParticleDefinition *particle_type)