]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4PhysicsConstructorSpecialCuts.cxx
added comment lines separating methods
[u/mrichter/AliRoot.git] / TGeant4 / TG4PhysicsConstructorSpecialCuts.cxx
1 // $Id$
2 // Category: physics
3 //
4 // See the class description in the header file.
5
6 #include "TG4PhysicsConstructorSpecialCuts.h"
7 #include "TG4G3PhysicsManager.h"
8 #include "TG4SpecialCuts.h"
9
10 #include <G4ParticleDefinition.hh>
11 #include <G4ProcessManager.hh>
12 #include <G4VProcess.hh>
13
14
15 //_____________________________________________________________________________
16 TG4PhysicsConstructorSpecialCuts::TG4PhysicsConstructorSpecialCuts(
17                                      const G4String& name)
18   : G4VPhysicsConstructor(name)
19 {
20 //
21   SetVerboseLevel(1);
22 }
23
24 //_____________________________________________________________________________
25 TG4PhysicsConstructorSpecialCuts::~TG4PhysicsConstructorSpecialCuts() {
26 //
27 }
28
29 // protected methods
30
31 //_____________________________________________________________________________
32 void TG4PhysicsConstructorSpecialCuts::ConstructParticle()
33 {
34 // The particles are constructed in the 
35 // TG4ModularPhysicsList.
36 // ---
37 }
38
39 //_____________________________________________________________________________
40 void TG4PhysicsConstructorSpecialCuts::ConstructProcess()
41 {
42 // Adds TG4SpecialCuts "process" that activates
43 // the kinetic energy cuts defined in 
44 // the vector of cuts (PhysicsManager::fCutVector) or in TG4Limits.
45 // ---
46
47   TG4G3PhysicsManager* g3PhysicsManager 
48     = TG4G3PhysicsManager::Instance();
49
50   if (g3PhysicsManager->IsSpecialCuts())
51   {
52     TG4G3CutVector* cutVector
53       = g3PhysicsManager->GetCutVector(); 
54     TG4boolVector* isCutVector 
55       = g3PhysicsManager->GetIsCutVector(); 
56
57     theParticleIterator->reset();
58     while ((*theParticleIterator)())
59     {
60       G4ParticleDefinition* particle = theParticleIterator->value();
61       TG4G3ParticleWSP particleWSP 
62         = g3PhysicsManager->GetG3ParticleWSP(particle);
63       G4String name;
64       g3PhysicsManager->GetG3ParticleWSPName(particleWSP, name);
65       
66       // uncomment this to see all particles "WSP"
67       //G4cout << "Iterating particle: " 
68       //       << particle->GetParticleName() << " " << particleWSP << " "
69       //       << name << G4endl;
70
71       // special process is created in case
72       // cutVector (vector of kinetic energy cuts) is set
73       // or the special cut is set by TG4Limits
74       if ((particleWSP !=kNofParticlesWSP) && 
75           ((*isCutVector)[particleWSP])) {
76         // check if process already exists
77         G4String processName = "specialCutFor" + name;
78         G4VProcess* process = g3PhysicsManager->FindProcess(processName);
79         if (!process) {
80           process = new TG4SpecialCuts(particleWSP, cutVector, processName);
81         }  
82         //particle->GetProcessManager()->AddProcess(process, 0, -1, 1);
83         particle->GetProcessManager()->AddDiscreteProcess(process);
84       }
85     }
86
87     if (verboseLevel>0) {
88       G4cout << "TG4PhysicsList::ConstructSpecialCuts: " << G4endl;
89       if (cutVector)
90         G4cout << "   Global kinetic energy cuts are set." << G4endl;
91       G4cout << "   Special cuts process is defined for: " << G4endl 
92              << "   ";
93       for (G4int i=0; i<kAny; i++) {
94         G4String name;
95         g3PhysicsManager->GetG3ParticleWSPName(i, name);
96         if ((*isCutVector)[i]) G4cout << name << " ";
97       }  
98       G4cout << G4endl;
99     }  
100   }
101 }
102
103