Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/eventgenerator/pythia/py8decayer/src/Py8DecayerPhysics.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 /// \file eventgenerator/pythia/pythia8decayer/src/Py8DecayerPhysics.cc
 28 /// \brief Implementation of the Py8DecayerPhysics class
 29 ///
 30 /// \author J. Yarba; FNAL
 31 
 32 #include "Py8DecayerPhysics.hh"
 33 
 34 #include "Py8Decayer.hh"
 35 
 36 #include "G4Decay.hh"
 37 #include "G4DecayTable.hh"
 38 #include "G4ParticleDefinition.hh"
 39 #include "G4ProcessManager.hh"
 40 
 41 // factory
 42 //
 43 #include "G4PhysicsConstructorFactory.hh"
 44 //
 45 // register it with contructor factory
 46 //
 47 G4_DECLARE_PHYSCONSTR_FACTORY(Py8DecayerPhysics);
 48 
 49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 50 
 51 Py8DecayerPhysics::Py8DecayerPhysics(G4int) : G4VPhysicsConstructor("Py8DecayerPhysics") {}
 52 
 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 54 
 55 Py8DecayerPhysics::~Py8DecayerPhysics() {}
 56 
 57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 58 
 59 void Py8DecayerPhysics::ConstructParticle()
 60 {
 61   // Nothing needs to be done here
 62 }
 63 
 64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 65 
 66 void Py8DecayerPhysics::ConstructProcess()
 67 {
 68   // Adding external decayer to G4Decay process (per each thread).
 69   // G4Decay will use the external decayer if G4Decay process is
 70   // assigned to an unstable particle and that particle does not
 71   // have its decay table.
 72 
 73   // Loop over all particles instantiated and remove already-assigned
 74   // decay table for tau's and B+/- so that they will decay through
 75   // the external decayer (Pythia8).
 76 
 77   // NOTE: The extDecayer will be deleted in G4Decay destructor
 78 
 79   Py8Decayer* extDecayer = new Py8Decayer();
 80   G4bool setOnce = true;
 81 
 82   auto particleIterator = GetParticleIterator();
 83   particleIterator->reset();
 84   while ((*particleIterator)()) {
 85     G4ParticleDefinition* particle = particleIterator->value();
 86 
 87     // remove native/existing decay table for
 88     // a)tau's
 89     // b) B+/-
 90     // so that G4Decay will use the external decayer
 91     if (std::abs(particle->GetPDGEncoding()) == 15 || std::abs(particle->GetPDGEncoding()) == 521) {
 92       if (particle->GetDecayTable()) {
 93         delete particle->GetDecayTable();
 94         particle->SetDecayTable(nullptr);
 95         /*
 96                   if ( verboseLevel > 1 ) {
 97                      G4cout << "Use ext decayer for: "
 98                         <<  particleIterator->value()->GetParticleName()
 99                         << G4endl;
100                   }
101         */
102       }
103     }
104 
105     if (setOnce)
106     // One G4Decay object is shared by all unstable particles (per thread).
107     // Thus, we set the external decayer only once.
108     {
109       G4ProcessManager* pmanager = particle->GetProcessManager();
110       G4ProcessVector* processVector = pmanager->GetProcessList();
111       for (size_t i = 0; i < processVector->length(); ++i) {
112         G4Decay* decay = dynamic_cast<G4Decay*>((*processVector)[i]);
113         if (decay) {
114           decay->SetExtDecayer(extDecayer);
115           setOnce = false;
116         }
117       }
118     }
119   }
120 
121   return;
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125