CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/SimG4Core/CheckSecondary/src/TreatSecondary.cc

Go to the documentation of this file.
00001 #include "SimG4Core/CheckSecondary/interface/TreatSecondary.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 #include "G4Step.hh"
00006 #include "G4Track.hh"
00007 #include "G4VProcess.hh"
00008 #include "G4HCofThisEvent.hh"
00009 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00010 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00011 
00012 #include <cmath>
00013 #include <iostream>
00014 #include <iomanip>
00015 
00016 TreatSecondary::TreatSecondary(const edm::ParameterSet &p): typeEnumerator(0) {
00017 
00018   verbosity     = p.getUntrackedParameter<int>("Verbosity",0);
00019   killAfter     = p.getUntrackedParameter<int>("KillAfter", -1);
00020   minDeltaE     = p.getUntrackedParameter<double>("MinimumDeltaE", 10.0)*MeV;
00021   
00022   edm::LogInfo("CheckSecondary") << "Instantiate CheckSecondary with Flag"
00023                                  << " for Killing track after "<< killAfter 
00024                                  << " hadronic interactions\nDefine inelastic"
00025                                  << " if > 1 seondary or change in KE > "
00026                                  << minDeltaE << " MeV\n";
00027 
00028   typeEnumerator = new G4ProcessTypeEnumerator();
00029 } 
00030    
00031 TreatSecondary::~TreatSecondary() {
00032   if (typeEnumerator) delete typeEnumerator;
00033 }
00034 
00035 void TreatSecondary::initTrack(const G4Track * thTk) {
00036 
00037   step   = 0;
00038   nsecL  = 0;
00039   nHad   = 0;
00040   eTrack = thTk->GetKineticEnergy()/MeV;
00041   LogDebug("CheckSecondary") << "TreatSecondary::initTrack:Track: " 
00042                              << thTk->GetTrackID() << " Type: " 
00043                              << thTk->GetDefinition()->GetParticleName() 
00044                              << " KE "    << thTk->GetKineticEnergy()/GeV 
00045                              << " GeV p " << thTk->GetMomentum().mag()/GeV
00046                              << " GeV daughter of particle " 
00047                              << thTk->GetParentID();
00048 }
00049 
00050 std::vector<math::XYZTLorentzVector> TreatSecondary::tracks(const G4Step*aStep,
00051                                                             std::string & name,
00052                                                             int & procid,
00053                                                             bool & hadrInt,
00054                                                             double & deltaE,
00055                                                             std::vector<int> & charges) {
00056 
00057   step++;
00058   procid  = -1;
00059   name    = "Unknown";
00060   hadrInt = false;
00061   deltaE  = 0;
00062   std::vector<math::XYZTLorentzVector> secondaries;
00063   charges.clear();
00064 
00065   if (aStep != NULL) {
00066     G4TrackVector* tkV  = const_cast<G4TrackVector*>(aStep->GetSecondary());
00067     G4Track*       thTk = aStep->GetTrack();
00068     const G4StepPoint* preStepPoint  = aStep->GetPreStepPoint();
00069     const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
00070     double eTrackNew = thTk->GetKineticEnergy()/MeV;
00071     deltaE = eTrack-eTrackNew;
00072     eTrack = eTrackNew;
00073     if (tkV != 0) {
00074       int nsc  = (*tkV).size();
00075       const G4VProcess*  proc = 0;
00076       if (postStepPoint) proc = postStepPoint->GetProcessDefinedStep();
00077       procid = typeEnumerator->processIdLong(proc);
00078       G4ProcessType type   = fNotDefined;
00079       if (proc) {
00080         type = proc->GetProcessType();
00081         name = proc->GetProcessName();
00082       }
00083       int           sec   = nsc - nsecL;
00084       LogDebug("CheckSecondary") << sec << " secondaries in step " 
00085                                  << thTk->GetCurrentStepNumber() 
00086                                  << " of track "  << thTk->GetTrackID() 
00087                                  << " from " << name << " of type " 
00088                                  << type << " ID " << procid << " (" 
00089                                  << typeEnumerator->processG4Name(procid) 
00090                                  << ")"; 
00091 
00092       G4TrackStatus state = thTk->GetTrackStatus();
00093       if (state == fAlive || state == fStopButAlive) sec++;
00094 
00095       if (type == fHadronic || type == fPhotolepton_hadron || type == fDecay) {
00096         if (deltaE > minDeltaE || sec > 1) {
00097           hadrInt = true;
00098           nHad++;
00099         }
00100         LogDebug("CheckSecondary") << "Hadronic Interaction " << nHad
00101                                    << " of Type " << type << " with "
00102                                    << sec << " secondaries from process "
00103                                    << proc->GetProcessName() << " Delta E "
00104                                    << deltaE << " Flag " << hadrInt;
00105       }
00106       if (hadrInt) {
00107         math::XYZTLorentzVector secondary;
00108         if (state == fAlive || state == fStopButAlive) {
00109           G4ThreeVector pp    = postStepPoint->GetMomentum();
00110           double        ee    = postStepPoint->GetTotalEnergy();
00111           secondary = math::XYZTLorentzVector(pp.x(),pp.y(),pp.z(),ee);
00112           secondaries.push_back(secondary);
00113           int           charge = (int)(postStepPoint->GetCharge());
00114           charges.push_back(charge);
00115         }
00116         for (int i=nsecL; i<nsc; i++) {
00117           G4Track*      tk = (*tkV)[i];
00118           G4ThreeVector pp = tk->GetMomentum();
00119           double        ee = tk->GetTotalEnergy();
00120           secondary = math::XYZTLorentzVector(pp.x(),pp.y(),pp.z(),ee);
00121           secondaries.push_back(secondary);
00122           int           charge = (int)(tk->GetDefinition()->GetPDGCharge());
00123           charges.push_back(charge);
00124         }
00125       }
00126       
00127       if (killAfter >= 0 && nHad >= killAfter)
00128         thTk->SetTrackStatus(fStopAndKill);
00129 
00130       if (verbosity > 0) {
00131         sec = 0;
00132         if (state == fAlive || state == fStopButAlive) {
00133           sec++;
00134           LogDebug("CheckSecondary") << "Secondary: " << sec << " ID " 
00135                                      << thTk->GetTrackID() << " Status " 
00136                                      << thTk->GetTrackStatus() << " Particle " 
00137                                      <<thTk->GetDefinition()->GetParticleName()
00138                                      << " Position "
00139                                      << postStepPoint->GetPosition() << " KE " 
00140                                      << postStepPoint->GetKineticEnergy() 
00141                                      << " Time " 
00142                                      << postStepPoint->GetGlobalTime();
00143         }
00144         for (int i=nsecL; i<nsc; i++) {
00145           sec++;
00146           G4Track* tk = (*tkV)[i];
00147           LogDebug("CheckSecondary") << "Secondary: " << sec << " ID " 
00148                                      << tk->GetTrackID() << " Status " 
00149                                      << tk->GetTrackStatus() << " Particle " 
00150                                      << tk->GetDefinition()->GetParticleName()
00151                                      << " Position "  << tk->GetPosition()
00152                                      << " KE " << tk->GetKineticEnergy() 
00153                                      << " Time " << tk->GetGlobalTime();
00154         }
00155       }
00156       nsecL  = nsc;
00157     }
00158 
00159     if (verbosity > 1)
00160       LogDebug("CheckSecondary") << "Track: " << thTk->GetTrackID() 
00161                                  << " Status " << thTk->GetTrackStatus()
00162                                  << " Particle " 
00163                                  << thTk->GetDefinition()->GetParticleName()
00164                                  << " at " << preStepPoint->GetPosition()
00165                                  << " Step: " << step << " KE " 
00166                                  << thTk->GetKineticEnergy()/GeV << " GeV; p " 
00167                                  << thTk->GetMomentum().mag()/GeV 
00168                                  << " GeV/c; Step Length "
00169                                  << aStep->GetStepLength()<< " Energy Deposit "
00170                                  << aStep->GetTotalEnergyDeposit()/MeV
00171                                  << " MeV; Interaction " << hadrInt;
00172   }
00173   return secondaries;
00174 }