CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/SimG4Core/Application/src/StackingAction.cc

Go to the documentation of this file.
00001 #include "SimG4Core/Application/interface/StackingAction.h"
00002 #include "SimG4Core/Notification/interface/CurrentG4Track.h"
00003 #include "SimG4Core/Notification/interface/NewTrackAction.h"
00004 #include "SimG4Core/Notification/interface/TrackInformation.h"
00005 #include "SimG4Core/Notification/interface/TrackInformationExtractor.h"
00006 
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 #include "G4VProcess.hh"
00010 #include "G4EmProcessSubType.hh"
00011 #include "G4LogicalVolumeStore.hh"
00012 #include "G4RegionStore.hh"
00013 #include "Randomize.hh"
00014 
00015 #include<algorithm>
00016 
00017 //#define DebugLog
00018 
00019 StackingAction::StackingAction(const edm::ParameterSet & p) {
00020 
00021   trackNeutrino  = p.getParameter<bool>("TrackNeutrino");
00022   killHeavy      = p.getParameter<bool>("KillHeavy");
00023   kmaxIon        = p.getParameter<double>("IonThreshold")*MeV;
00024   kmaxProton     = p.getParameter<double>("ProtonThreshold")*MeV;
00025   kmaxNeutron    = p.getParameter<double>("NeutronThreshold")*MeV;
00026   killDeltaRay   = p.getParameter<bool>("KillDeltaRay");
00027   maxTrackTime   = p.getParameter<double>("MaxTrackTime")*ns;
00028   maxTrackTimes  = p.getParameter<std::vector<double> >("MaxTrackTimes");
00029   maxTimeNames   = p.getParameter<std::vector<std::string> >("MaxTimeNames");
00030   savePDandCinTracker = p.getUntrackedParameter<bool>("SavePrimaryDecayProductsAndConversionsInTracker",false);
00031   savePDandCinCalo    = p.getUntrackedParameter<bool>("SavePrimaryDecayProductsAndConversionsInCalo",false);
00032   savePDandCinMuon    = p.getUntrackedParameter<bool>("SavePrimaryDecayProductsAndConversionsInMuon",false);
00033   saveFirstSecondary  = p.getUntrackedParameter<bool>("SaveFirstLevelSecondary",false);
00034   killInCalo = false;
00035   killInCaloEfH = false;
00036 
00037   // Russian Roulette
00038   nRusRoEcal = p.getParameter<double>("RusRoEcalNeutron");
00039   nRusRoHcal = p.getParameter<double>("RusRoHcalNeutron");
00040   nRusRoEcalLim = p.getParameter<double>("RusRoEcalNeutronLimit");
00041   nRusRoHcalLim = p.getParameter<double>("RusRoHcalNeutronLimit");
00042   pRusRoEcal = p.getParameter<double>("RusRoEcalProton");
00043   pRusRoHcal = p.getParameter<double>("RusRoHcalProton");
00044   pRusRoEcalLim = p.getParameter<double>("RusRoEcalProtonLimit");
00045   pRusRoHcalLim = p.getParameter<double>("RusRoHcalProtonLimit");
00046   regionEcal = 0;
00047   regionHcal = 0;
00048 
00049   if ( p.exists("TestKillingOptions") ) {
00050 
00051     killInCalo = (p.getParameter<edm::ParameterSet>("TestKillingOptions")).getParameter<bool>("KillInCalo");
00052     killInCaloEfH = (p.getParameter<edm::ParameterSet>("TestKillingOptions")).getParameter<bool>("KillInCaloEfH");
00053     edm::LogWarning("SimG4CoreApplication") << " *** Activating special test killing options in StackingAction \n"
00054                                             << " *** Kill secondaries in Calorimetetrs volume = " << killInCalo << "\n"
00055                                             << " *** Kill electromagnetic secondaries from hadrons in Calorimeters volume = " << killInCaloEfH;
00056 
00057   }
00058 
00059   edm::LogInfo("SimG4CoreApplication") << "StackingAction initiated with"
00060                                        << " flag for saving decay products in "
00061                                        << " Tracker: " << savePDandCinTracker
00062                                        << " in Calo: " << savePDandCinCalo
00063                                        << " in Muon: " << savePDandCinMuon
00064                                        << "\n               saveFirstSecondary"
00065                                        << ": " << saveFirstSecondary
00066                                        << " Flag for tracking neutrino: "
00067                                        << trackNeutrino << " Killing Flag "
00068                                        << killHeavy << " protons below " 
00069                                        << kmaxProton <<" MeV, neutrons below "
00070                                        << kmaxNeutron << " MeV and ions"
00071                                        << " below " << kmaxIon << " MeV\n"
00072                                        << "               kill tracks with "
00073                                        << "time larger than " << maxTrackTime
00074                                        << " ns and kill Delta Ray flag set to "
00075                                        << killDeltaRay;
00076   for (unsigned int i=0; i<maxTrackTimes.size(); i++) {
00077     maxTrackTimes[i] *= ns;
00078     edm::LogInfo("SimG4CoreApplication") << "StackingAction::MaxTrackTime for "
00079                                          << maxTimeNames[i] << " is " 
00080                                          << maxTrackTimes[i];
00081   }
00082   initPointer();
00083 }
00084 
00085 StackingAction::~StackingAction() {}
00086 
00087 G4ClassificationOfNewTrack StackingAction::ClassifyNewTrack(const G4Track * aTrack) {
00088 
00089   // G4 interface part
00090   G4ClassificationOfNewTrack classification = fUrgent;
00091   int flag = 0;
00092 
00093   NewTrackAction newTA;
00094   if (aTrack->GetCreatorProcess()==0 || aTrack->GetParentID()==0) {
00095     /*
00096     std::cout << "StackingAction: primary weight= " 
00097               << aTrack->GetWeight() << " "
00098               << aTrack->GetDefinition()->GetParticleName()
00099               << " " << aTrack->GetKineticEnergy()
00100               << " Id=" << aTrack->GetTrackID()
00101               << "  trackInfo " << aTrack->GetUserInformation() 
00102               << std::endl;
00103     */
00104     newTA.primary(aTrack);
00105     if (!trackNeutrino) {
00106       int    pdg = std::abs(aTrack->GetDefinition()->GetPDGEncoding());
00107       if (pdg == 12 || pdg == 14 || pdg == 16 || pdg == 18) {
00108         classification = fKill;
00109       }
00110     }
00111   } else if (aTrack->GetTouchable() == 0) {
00112     edm::LogError("SimG4CoreApplication")
00113       << "StackingAction: no touchable for track " << aTrack->GetTrackID()
00114       << " from " << aTrack->GetParentID()
00115       << " with PDG code " << aTrack->GetDefinition()->GetParticleName();
00116     classification = fKill;
00117   } else {
00118     const G4Track * mother = CurrentG4Track::track();
00119     if ((savePDandCinTracker && isThisVolume(aTrack->GetTouchable(),tracker))||
00120         (savePDandCinCalo && isThisVolume(aTrack->GetTouchable(),calo)) ||
00121         (savePDandCinMuon && isThisVolume(aTrack->GetTouchable(),muon)))
00122       flag = isItPrimaryDecayProductOrConversion(aTrack, *mother);
00123     if (saveFirstSecondary) flag = isItFromPrimary(*mother, flag);
00124     newTA.secondary(aTrack, *mother, flag);
00125 
00126     if (aTrack->GetTrackStatus() == fStopAndKill) classification = fKill;
00127     if (killHeavy && classification != fKill) {
00128       int    pdg = aTrack->GetDefinition()->GetPDGEncoding();
00129       double ke  = aTrack->GetKineticEnergy()/MeV;
00130       if (((pdg/1000000000 == 1) && (((pdg/10000)%100) > 0) && 
00131            (((pdg/10)%100) > 0) && (ke<kmaxIon)) || 
00132           ((pdg == 2212) && (ke < kmaxProton)) ||
00133           ((pdg == 2112) && (ke < kmaxNeutron))) classification = fKill;
00134     }
00135     if (!trackNeutrino  && classification != fKill) {
00136       int    pdg = std::abs(aTrack->GetDefinition()->GetPDGEncoding());
00137       if (pdg == 12 || pdg == 14 || pdg == 16 || pdg == 18) 
00138         classification = fKill;
00139     }
00140     if (isItLongLived(aTrack)) classification = fKill;
00141     if (killDeltaRay) {
00142       if (aTrack->GetCreatorProcess()->GetProcessType() == fElectromagnetic &&
00143           aTrack->GetCreatorProcess()->GetProcessSubType() == fIonisation)
00144         classification = fKill;
00145     }
00146     if (killInCalo && classification != fKill) {
00147       if ( isThisVolume(aTrack->GetTouchable(),calo)) { 
00148         classification = fKill; 
00149       }
00150     }
00151     if (killInCaloEfH && classification != fKill) {
00152       int    pdg = aTrack->GetDefinition()->GetPDGEncoding();
00153       int    pdgMother = mother->GetDefinition()->GetPDGEncoding();
00154       if ( (pdg == 22 || std::abs(pdg) == 11) && (std::abs(pdgMother) < 11 || std::abs(pdgMother) > 17) && pdgMother != 22  ) {
00155         if ( isThisVolume(aTrack->GetTouchable(),calo)) { 
00156           classification = fKill; 
00157         }
00158       }
00159     }
00160     if((regionEcal || regionHcal) && classification != fKill) {
00161       int    pdg = aTrack->GetDefinition()->GetPDGEncoding();
00162       if(2112 == pdg || 2212 == pdg) {
00163         G4Region* reg = aTrack->GetVolume()->GetLogicalVolume()->GetRegion();
00164         double prob = 1.0;
00165         double elim = 0.0;
00166         if(reg == regionEcal) {
00167           if(2112 == pdg) {
00168             prob = nRusRoEcal;
00169             elim = nRusRoEcalLim;
00170           } else {
00171             prob = pRusRoEcal;
00172             elim = pRusRoEcalLim;
00173           }
00174         } else if(reg == regionHcal) {
00175           if(2112 == pdg) {
00176             prob = nRusRoHcal;
00177             elim = nRusRoHcalLim;
00178           } else {
00179             prob = pRusRoHcal;
00180             elim = pRusRoHcalLim;
00181           }
00182         }
00183         if(prob < 1.0 && aTrack->GetKineticEnergy() < elim) {
00184           if(G4UniformRand() < prob) {
00185             const_cast<G4Track*>(aTrack)->SetWeight(aTrack->GetWeight()/prob);
00186           } else {
00187             classification = fKill;
00188           }
00189         }
00190       }  
00191     }
00192   /*
00193   double wt2 = aTrack->GetWeight();
00194   if(wt2 != 1.0) { 
00195     G4Region* reg = aTrack->GetVolume()->GetLogicalVolume()->GetRegion();
00196     std::cout << "StackingAction: weight= " << wt2 << " "
00197               << aTrack->GetDefinition()->GetParticleName()
00198               << " " << aTrack->GetKineticEnergy()
00199               << " Id=" << aTrack->GetTrackID()
00200               << " IdP=" << aTrack->GetParentID();
00201     const G4VProcess* pr = aTrack->GetCreatorProcess();
00202     if(pr) std::cout << " from  " << pr->GetProcessName();
00203     if(reg) std::cout << " in  " << reg->GetName()
00204                       << "  trackInfo " << aTrack->GetUserInformation(); 
00205     std::cout << std::endl;
00206   }
00207   */
00208 #ifdef DebugLog
00209     LogDebug("SimG4CoreApplication") << "StackingAction:Classify Track "
00210                                      << aTrack->GetTrackID() << " Parent " 
00211                                      << aTrack->GetParentID() << " Type "
00212                                      << aTrack->GetDefinition()->GetParticleName() 
00213                                      << " K.E. " << aTrack->GetKineticEnergy()/MeV
00214                                      << " MeV from process/subprocess " 
00215                                      << aTrack->GetCreatorProcess()->GetProcessType() << "|"
00216                                      << aTrack->GetCreatorProcess()->GetProcessSubType()
00217                                      << " as " << classification << " Flag " << flag;
00218 #endif
00219   }
00220   return classification;
00221 }
00222 
00223 void StackingAction::NewStage() {}
00224 
00225 void StackingAction::PrepareNewEvent() {}
00226 
00227 void StackingAction::initPointer() {
00228 
00229   const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
00230   if (lvs) {
00231     std::vector<G4LogicalVolume*>::const_iterator lvcite;
00232     for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
00233       if (savePDandCinTracker) {
00234         if (strcmp("Tracker",(*lvcite)->GetName().c_str()) == 0) tracker.push_back(*lvcite);
00235         if (strcmp("BEAM",(*lvcite)->GetName().substr(0,4).c_str()) == 0) tracker.push_back(*lvcite);
00236       }
00237       if (savePDandCinCalo || killInCalo || killInCaloEfH ) {
00238         if (strcmp("CALO",(*lvcite)->GetName().c_str()) == 0) calo.push_back(*lvcite);
00239         if (strcmp("VCAL",(*lvcite)->GetName().c_str()) == 0) calo.push_back(*lvcite);
00240       }
00241       if (savePDandCinMuon) {
00242         if (strcmp("MUON",(*lvcite)->GetName().c_str()) == 0) muon.push_back(*lvcite);
00243       }
00244     }
00245     edm::LogInfo("SimG4CoreApplication") << "# of LV for Tracker " 
00246                                          << tracker.size() << " for Calo "
00247                                          << calo.size() << " for Muon "
00248                                          << muon.size();
00249     for (unsigned int i=0; i<tracker.size(); ++i)
00250       edm::LogInfo("SimG4CoreApplication") << "Tracker vol " << i << " name "
00251                                            << tracker[i]->GetName();
00252     for (unsigned int i=0; i<calo.size(); ++i)
00253       edm::LogInfo("SimG4CoreApplication") << "Calorimeter vol " <<i <<" name "
00254                                            << calo[i]->GetName();
00255     for (unsigned int i=0; i<muon.size(); ++i)
00256       edm::LogInfo("SimG4CoreApplication") << "Muon vol " << i << " name "
00257                                            << muon[i]->GetName();
00258   }
00259 
00260   const G4RegionStore * rs = G4RegionStore::GetInstance();
00261   unsigned int num = maxTimeNames.size();
00262   if (num > 0) {
00263     std::vector<double> tofs;
00264     if (rs) {
00265       std::vector<G4Region*>::const_iterator rcite;
00266       for (rcite = rs->begin(); rcite != rs->end(); rcite++) {
00267         if ((nRusRoEcal < 1.0 || pRusRoEcal < 1.0) && 
00268             (*rcite)->GetName() == "EcalRegion") {  regionEcal = (*rcite); }
00269         if ((nRusRoHcal < 1.0 || pRusRoHcal < 1.0) && 
00270             (*rcite)->GetName() == "HcalRegion") {  regionHcal = (*rcite); }
00271         for (unsigned int i=0; i<num; i++) {
00272           if ((*rcite)->GetName() == (G4String)(maxTimeNames[i])) {
00273             maxTimeRegions.push_back(*rcite);
00274             tofs.push_back(maxTrackTimes[i]);
00275             break;
00276           }
00277         }
00278         if (tofs.size() == num) break;
00279       }
00280     }
00281     for (unsigned int i=0; i<tofs.size(); i++) {
00282       maxTrackTimes[i] = tofs[i];
00283       G4String name = "Unknown";
00284       if (maxTimeRegions[i]) name = maxTimeRegions[i]->GetName();
00285       edm::LogInfo("SimG4CoreApplication") << name << " with pointer " 
00286                                            << maxTimeRegions[i]<<" KE cut off "
00287                                            << maxTrackTimes[i];
00288     }
00289   }
00290 }
00291 
00292 bool StackingAction::isThisVolume(const G4VTouchable* touch, 
00293                                   std::vector<G4LogicalVolume*> & lvs) const {
00294 
00295   bool flag = false;
00296   if (lvs.size() > 0 && touch !=0) {
00297     int level = ((touch->GetHistoryDepth())+1);
00298     if (level >= 3) {
00299       unsigned int  ii = (unsigned int)(level - 3);
00300       flag    = (std::count(lvs.begin(),lvs.end(),(touch->GetVolume(ii)->GetLogicalVolume())) != 0);
00301     }
00302   }
00303   return flag;
00304 }
00305 
00306 int StackingAction::isItPrimaryDecayProductOrConversion(const G4Track * aTrack,
00307                                                         const G4Track & mother) const {
00308 
00309   int flag = 0;
00310   TrackInformationExtractor extractor;
00311   const TrackInformation & motherInfo(extractor(mother));
00312   // Check whether mother is a primary
00313   if (motherInfo.isPrimary()) {
00314     if (aTrack->GetCreatorProcess()->GetProcessType() == fDecay) flag = 1;
00315     else if (aTrack->GetCreatorProcess()->GetProcessType() == fElectromagnetic &&
00316              aTrack->GetCreatorProcess()->GetProcessSubType() == fGammaConversion) flag = 2;
00317   }
00318   return flag;
00319 }
00320 
00321 int StackingAction::isItFromPrimary(const G4Track & mother, int flagIn) const {
00322 
00323   int flag = flagIn;
00324   if (flag != 1) {
00325     TrackInformationExtractor extractor;
00326     const TrackInformation & motherInfo(extractor(mother));
00327     if (motherInfo.isPrimary()) flag = 3;
00328   }
00329   return flag;
00330 }
00331 
00332 bool StackingAction::isItLongLived(const G4Track * aTrack) const {
00333 
00334   bool   flag = false;
00335   double time = (aTrack->GetGlobalTime())/nanosecond;
00336   double tofM = maxTrackTime;
00337   if (maxTimeRegions.size() > 0) {
00338     G4Region* reg = aTrack->GetTouchable()->GetVolume(0)->GetLogicalVolume()->GetRegion();
00339     for (unsigned int i=0; i<maxTimeRegions.size(); i++) {
00340       if (reg == maxTimeRegions[i]) {
00341         tofM = maxTrackTimes[i];
00342         break;
00343       }
00344     }
00345   }
00346   if (time > tofM) flag = true;
00347   return flag;
00348 }
00349 
00350