CMS 3D CMS Logo

TrackingVerboseAction.cc
Go to the documentation of this file.
1 // File: TrackingVerboseAction.cc
3 // Creation: P.Arce 09/01
4 // Modifications: porting to CMSSW by M. Stavrianakou 22/03/06
5 // Description:
7 
9 
16 
17 #include "G4Track.hh"
18 #include "G4Event.hh"
19 #include "G4ios.hh"
20 #include "G4TrackingManager.hh"
21 #include "G4EventManager.hh"
22 #include "G4VSteppingVerbose.hh"
23 #include "G4UnitsTable.hh"
24 
25 #include <algorithm>
26 
27 using namespace CLHEP;
28 
30  : theTrackingManager(nullptr), fVerbose(nullptr) {
31  fLarge = int(1E10);
32  fDEBUG = p.getUntrackedParameter<bool>("DEBUG", false);
33  fHighEtPhotons = p.getUntrackedParameter<bool>("CheckForHighEtPhotons", false);
34  fG4Verbose = p.getUntrackedParameter<bool>("G4Verbose", false);
35 
36  //----- Set which events are verbose
37  fTVEventMin = p.getUntrackedParameter<int>("EventMin", 0);
38  fTVEventMax = p.getUntrackedParameter<int>("EventMax", fLarge);
39  fTVEventStep = p.getUntrackedParameter<int>("EventStep", 1);
40 
41  //----- Set which tracks of those events are verbose
42  fTVTrackMin = p.getUntrackedParameter<int>("TrackMin", 0);
43  fTVTrackMax = p.getUntrackedParameter<int>("TrackMax", fLarge);
44  fTVTrackStep = p.getUntrackedParameter<int>("TrackStep", 1);
45 
46  //----- Set the verbosity level
47  fVerboseLevel = p.getUntrackedParameter<int>("VerboseLevel", 1);
48  fPdgIds = p.getUntrackedParameter<std::vector<int> >("PDGids");
49  if (fDEBUG) {
50  G4cout << "TV: fTVTrackMin " << fTVTrackMin << " fTVTrackMax " << fTVTrackMax << " fTVTrackStep " << fTVTrackStep
51  << " fTVEventMin " << fTVEventMin << " fTVEventMax " << fTVEventMax << " fTVEventStep " << fTVEventStep
52  << " fVerboseLevel " << fVerboseLevel << " fG4Verbose " << fG4Verbose << " PDGIds " << fPdgIds.size()
53  << G4endl;
54  for (unsigned int ii = 0; ii < fPdgIds.size(); ++ii)
55  G4cout << "TV: PDGId[" << ii << "] = " << fPdgIds[ii] << G4endl;
56  }
57 
58  //----- Set verbosity off to start
59  fTrackingVerboseON = false;
60  fTkVerbThisEventON = false;
61 
62  G4cout << " TrackingVerbose constructed " << G4endl;
63 }
64 
66 
68  /*
69  TrackingAction * ta =
70  dynamic_cast<TrackingAction*>(G4EventManager::GetEventManager()->GetUserTrackingAction());
71  theTrackingManager = ta->getTrackManager();
72  */
73  fVerbose = G4VSteppingVerbose::GetInstance();
74  if (fDEBUG)
75  G4cout << " TV: Get the Tracking Manager: " << theTrackingManager << " and the SteppingVerbose: " << fVerbose
76  << G4endl;
77 }
78 
80  if (evt == nullptr)
81  return;
82  const G4Event* anEvent = (*evt)();
83  if (anEvent == nullptr)
84  return;
85 
86  //----------- Set /tracking/verbose for this event
87  int eventNo = anEvent->GetEventID();
88  if (fDEBUG)
89  G4cout << "TV: trackID: NEW EVENT " << eventNo << G4endl;
90 
91  fTkVerbThisEventON = false;
92  //----- Check if event is in the selected range
93  if (eventNo >= fTVEventMin && eventNo <= fTVEventMax) {
94  if ((eventNo - fTVEventMin) % fTVEventStep == 0)
95  fTkVerbThisEventON = true;
96  }
97 
98  if (fDEBUG)
99  G4cout << " TV: fTkVerbThisEventON " << fTkVerbThisEventON << " fTrackingVerboseON " << fTrackingVerboseON
100  << " fTVEventMin " << fTVEventMin << " fTVEventMax " << fTVEventMax << G4endl;
101  //----- check if verbosity has to be changed
103  if (fTVTrackMin == 0 && fTVTrackMax == fLarge && fTVTrackStep != 1) {
105  fTrackingVerboseON = true;
106  if (fDEBUG)
107  G4cout << "TV: VERBOSEet1 " << eventNo << G4endl;
108  }
109  } else if ((!fTkVerbThisEventON) && (fTrackingVerboseON)) {
111  fTrackingVerboseON = false;
112  if (fDEBUG)
113  G4cout << "TV: VERBOSEet0 " << eventNo << G4endl;
114  }
115 }
116 
118  const G4Track* aTrack = (*trk)();
119 
120  //----- High ET photon printout
121  //---------- Set /tracking/verbose
122  //----- track is verbose only if event is verbose
123  double tkP = aTrack->GetMomentum().mag();
124  double tkPx = aTrack->GetMomentum().x();
125  double tkPy = aTrack->GetMomentum().y();
126  double tkPz = aTrack->GetMomentum().z();
127 
128  double tvtx = aTrack->GetVertexPosition().x();
129  double tvty = aTrack->GetVertexPosition().y();
130  double tvtz = aTrack->GetVertexPosition().z();
131 
132  double g4t_phi = atan2(tkPy, tkPx);
133 
134  double drpart = sqrt(tkPx * tkPx + tkPy * tkPy);
135 
136  double mythetapart = acos(tkPz / sqrt(drpart * drpart + tkPz * tkPz));
137 
138  double g4t_eta = -log(tan(mythetapart / 2.));
139  G4int MytrackNo = aTrack->GetTrackID();
140 
141  if (fHighEtPhotons) {
142  if (aTrack->GetDefinition()->GetParticleName() == "gamma" && aTrack->GetParentID() != 0) {
143  if ((tkPx * tkPx + tkPy * tkPy + tkPz * tkPz) > 1000.0 * 1000.0 &&
144  aTrack->GetCreatorProcess()->GetProcessName() == "LCapture") {
145  G4cout << "MY NEW GAMMA " << G4endl;
146  G4cout << "**********************************************************************" << G4endl;
147  G4cout << "MY NEW TRACK ID = " << MytrackNo << "(" << aTrack->GetDefinition()->GetParticleName() << ")"
148  << " PARENT =" << aTrack->GetParentID() << G4endl;
149  G4cout << "Primary particle: " << aTrack->GetDynamicParticle()->GetPrimaryParticle() << G4endl;
150  G4cout << "Process type: " << aTrack->GetCreatorProcess()->GetProcessType()
151  << " Process name: " << aTrack->GetCreatorProcess()->GetProcessName() << G4endl;
152  G4cout << "ToT E = " << aTrack->GetTotalEnergy() << " KineE = " << aTrack->GetKineticEnergy()
153  << " Tot P = " << tkP << " Pt = " << sqrt(tkPx * tkPx + tkPy * tkPy) << " VTX=(" << tvtx << "," << tvty
154  << "," << tvtz << ")" << G4endl;
155  if (aTrack->GetKineticEnergy() > 1. * GeV && aTrack->GetCreatorProcess()->GetProcessName() != "LCapture")
156  G4cout << " KineE > 1 GeV !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << G4endl;
157  const G4VTouchable* touchable = aTrack->GetTouchable();
158  if (touchable != nullptr && touchable->GetVolume() != nullptr &&
159  touchable->GetVolume()->GetLogicalVolume() != nullptr) {
160  G4Material* material = touchable->GetVolume()->GetLogicalVolume()->GetMaterial();
161  G4cout << "G4LCapture Gamma E(GeV) " << aTrack->GetTotalEnergy() / GeV << " " << material->GetName() << " "
162  << touchable->GetVolume()->GetName() << G4endl;
163  G4cout << "G4LCapture Gamma position(m): " << aTrack->GetPosition() / m << G4endl;
164  G4cout << "G4LCapture created Gamma direction " << aTrack->GetMomentumDirection() << G4endl;
165  G4cout << "G4LCapture gamma (eta,phi) = "
166  << "(" << g4t_eta << "," << g4t_phi << ")" << G4endl;
167  }
168  aTrack->GetUserInformation()->Print();
169  G4cout << "**********************************************************************" << G4endl;
170  }
171  }
172 
173  if (aTrack->GetDefinition()->GetParticleName() == "gamma") {
174  const G4VProcess* proc = aTrack->GetCreatorProcess();
175  double Tgamma = aTrack->GetKineticEnergy();
176  std::string ProcName;
177  const std::string nullstr("Null_prc");
178  if (proc)
179  ProcName = proc->GetProcessName();
180  else
181  ProcName = nullstr;
182  if (Tgamma > 2.5 * GeV) { //&& ProcName!="Decay" && ProcName!="eBrem")
183  std::string volumeName("_Unknown_Vol_");
184  std::string materialName("_Unknown_Mat_");
185  G4Material* material = nullptr;
186  G4VPhysicalVolume* pvolume = nullptr;
187  G4LogicalVolume* lvolume = nullptr;
188  const G4VTouchable* touchable = aTrack->GetTouchable();
189  if (touchable)
190  pvolume = touchable->GetVolume();
191  if (pvolume) {
192  volumeName = pvolume->GetName();
193  lvolume = pvolume->GetLogicalVolume();
194  }
195  if (lvolume)
196  material = lvolume->GetMaterial();
197  if (material)
198  materialName = material->GetName();
199  G4cout << "**** ALL photons > 2.5 GeV ****" << G4endl;
200  G4cout << ProcName << "**** ALL photons: gamma E(GeV) " << aTrack->GetTotalEnergy() / GeV << " "
201  << materialName << " " << volumeName << G4endl;
202  G4cout << ProcName << "**** ALL photons: gamma position(m): " << aTrack->GetPosition() / m << G4endl;
203  G4cout << ProcName << "**** ALL photons: gamma direction " << aTrack->GetMomentumDirection() << G4endl;
204  G4cout << "**********************************************************************" << G4endl;
205  }
206  }
207  }
208 
209  //---------- Set /tracking/verbose
210  //----- track is verbose only if event is verbose
211  if (fTkVerbThisEventON) {
212  bool trackingVerboseThisTrack = checkTrackingVerbose(aTrack);
213 
214  //----- Set the /tracking/verbose for this track
215  if ((trackingVerboseThisTrack) && (!fTrackingVerboseON)) {
217  fTrackingVerboseON = true;
218  if (fDEBUG)
219  G4cout << "TV: VERBOSEtt1 " << aTrack->GetTrackID() << G4endl;
220  printTrackInfo(aTrack);
221  } else if ((!trackingVerboseThisTrack) && (fTrackingVerboseON)) {
223  fTrackingVerboseON = false;
224  if (fDEBUG)
225  G4cout << "TV: VERBOSEtt0 " << aTrack->GetTrackID() << G4endl;
226  }
227  }
228 }
229 
231  const G4Track* aTrack = (*trk)();
232  if (fTkVerbThisEventON) {
233  bool trackingVerboseThisTrack = checkTrackingVerbose(aTrack);
234  if ((trackingVerboseThisTrack) && (fTrackingVerboseON) && (fTVTrackMax < fLarge || fTVTrackStep != 1)) {
236  fTrackingVerboseON = false;
237  if (fDEBUG)
238  G4cout << "TV: VERBOSEtt0 " << aTrack->GetTrackID() << G4endl;
239  }
240  }
241 }
242 
243 void TrackingVerboseAction::update(const G4Step* fStep) {
244  if ((fG4Verbose) && (fTrackingVerboseON)) {
245  G4Track* fTrack = fStep->GetTrack();
246  G4cout << std::setw(5) << fTrack->GetCurrentStepNumber() << " " << std::setw(8)
247  << G4BestUnit(fTrack->GetPosition().x(), "Length") << " " << std::setw(8)
248  << G4BestUnit(fTrack->GetPosition().y(), "Length") << " " << std::setw(8)
249  << G4BestUnit(fTrack->GetPosition().z(), "Length") << " " << std::setw(9)
250  << G4BestUnit(fTrack->GetKineticEnergy(), "Energy") << " " << std::setw(8)
251  << G4BestUnit(fStep->GetTotalEnergyDeposit(), "Energy") << " " << std::setw(8)
252  << G4BestUnit(fStep->GetStepLength(), "Length") << " " << std::setw(9)
253  << G4BestUnit(fTrack->GetTrackLength(), "Length") << " " << std::setw(9)
254  << G4BestUnit(fTrack->GetGlobalTime(), "Time") << " ";
255 
256  // Put cut comment here
257  if (fTrack->GetNextVolume() != nullptr) {
258  G4cout << std::setw(11) << fTrack->GetNextVolume()->GetName() << " ";
259  } else {
260  G4cout << std::setw(11) << "OutOfWorld"
261  << " ";
262  }
263  if (fStep->GetPostStepPoint()->GetProcessDefinedStep() != nullptr) {
264  G4cout << fStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
265  } else {
266  G4cout << "User Limit";
267  }
268  G4cout << G4endl;
269  }
270 }
271 
273  if (fDEBUG)
274  G4cout << " setting verbose level " << verblev << G4endl;
275  if (theTrackingManager != nullptr)
276  theTrackingManager->SetVerboseLevel(verblev);
277 }
278 
279 bool TrackingVerboseAction::checkTrackingVerbose(const G4Track* aTrack) {
280  int trackNo = aTrack->GetTrackID();
281  bool trackingVerboseThisTrack = false;
282  //----- Check if track is in the selected range
283  if (trackNo >= fTVTrackMin && trackNo <= fTVTrackMax) {
284  if ((trackNo - fTVTrackMin) % fTVTrackStep == 0)
285  trackingVerboseThisTrack = true;
286  }
287  if (trackingVerboseThisTrack && (!fPdgIds.empty())) {
288  int pdgId = aTrack->GetDefinition()->GetPDGEncoding();
289  if (std::count(fPdgIds.begin(), fPdgIds.end(), pdgId) == 0)
290  trackingVerboseThisTrack = false;
291  }
292  return trackingVerboseThisTrack;
293 }
294 
295 void TrackingVerboseAction::printTrackInfo(const G4Track* aTrack) {
296  G4cout << G4endl << "*******************************************************"
297  << "**************************************************" << G4endl << "* G4Track Information: "
298  << " Particle = " << aTrack->GetDefinition()->GetParticleName() << ","
299  << " Track ID = " << aTrack->GetTrackID() << ","
300  << " Parent ID = " << aTrack->GetParentID() << G4endl
301  << "*******************************************************"
302  << "**************************************************" << G4endl << G4endl;
303  if (fVerbose)
304  fVerbose->TrackingStarted();
305 }
std::vector< int > fPdgIds
bool checkTrackingVerbose(const G4Track *)
void printTrackInfo(const G4Track *)
G4TrackingManager * theTrackingManager
T sqrt(T t)
Definition: SSEVec.h:23
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
void update(const BeginOfRun *) override
This routine will be called when the appropriate signal arrives.
ii
Definition: cuy.py:589
G4VSteppingVerbose * fVerbose
void setTrackingVerbose(int verblev)
TrackingVerboseAction(edm::ParameterSet const &p)