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
00014 #include<algorithm>
00015
00016
00017
00018 StackingAction::StackingAction(const edm::ParameterSet & p) {
00019
00020 trackNeutrino = p.getParameter<bool>("TrackNeutrino");
00021 killHeavy = p.getParameter<bool>("KillHeavy");
00022 kmaxIon = p.getParameter<double>("IonThreshold")*MeV;
00023 kmaxProton = p.getParameter<double>("ProtonThreshold")*MeV;
00024 kmaxNeutron = p.getParameter<double>("NeutronThreshold")*MeV;
00025 killDeltaRay = p.getParameter<bool>("KillDeltaRay");
00026 maxTrackTime = p.getParameter<double>("MaxTrackTime")*ns;
00027 maxTrackTimes = p.getParameter<std::vector<double> >("MaxTrackTimes");
00028 maxTimeNames = p.getParameter<std::vector<std::string> >("MaxTimeNames");
00029 savePDandCinTracker = p.getUntrackedParameter<bool>("SavePrimaryDecayProductsAndConversionsInTracker",false);
00030 savePDandCinCalo = p.getUntrackedParameter<bool>("SavePrimaryDecayProductsAndConversionsInCalo",false);
00031 savePDandCinMuon = p.getUntrackedParameter<bool>("SavePrimaryDecayProductsAndConversionsInMuon",false);
00032 saveFirstSecondary = p.getUntrackedParameter<bool>("SaveFirstLevelSecondary",false);
00033
00034 edm::LogInfo("SimG4CoreApplication") << "StackingAction initiated with"
00035 << " flag for saving decay products in "
00036 << " Tracker: " << savePDandCinTracker
00037 << " in Calo: " << savePDandCinCalo
00038 << " in Muon: " << savePDandCinMuon
00039 << "\n saveFirstSecondary"
00040 << ": " << saveFirstSecondary
00041 << " Flag for tracking neutrino: "
00042 << trackNeutrino << " Killing Flag "
00043 << killHeavy << " protons below "
00044 << kmaxProton <<" MeV, neutrons below "
00045 << kmaxNeutron << " MeV and ions"
00046 << " below " << kmaxIon << " MeV\n"
00047 << " kill tracks with "
00048 << "time larger than " << maxTrackTime
00049 << " ns and kill Delta Ray flag set to "
00050 << killDeltaRay;
00051 for (unsigned int i=0; i<maxTrackTimes.size(); i++) {
00052 maxTrackTimes[i] *= ns;
00053 edm::LogInfo("SimG4CoreApplication") << "StackingAction::MaxTrackTime for "
00054 << maxTimeNames[i] << " is "
00055 << maxTrackTimes[i];
00056 }
00057 initPointer();
00058 }
00059
00060 StackingAction::~StackingAction() {}
00061
00062 G4ClassificationOfNewTrack StackingAction::ClassifyNewTrack(const G4Track * aTrack) {
00063
00064
00065 G4ClassificationOfNewTrack classification = fUrgent;
00066 int flag = 0;
00067
00068 NewTrackAction newTA;
00069 if (aTrack->GetCreatorProcess()==0 || aTrack->GetParentID()==0) {
00070 newTA.primary(aTrack);
00071 } else if (aTrack->GetTouchable() == 0) {
00072 edm::LogError("SimG4CoreApplication")
00073 << "StackingAction: no touchable for track " << aTrack->GetTrackID()
00074 << " from " << aTrack->GetParentID()
00075 << " with PDG code " << aTrack->GetDefinition()->GetParticleName();
00076 classification = fKill;
00077 } else {
00078 const G4Track * mother = CurrentG4Track::track();
00079 if ((savePDandCinTracker && isThisVolume(aTrack->GetTouchable(),tracker))||
00080 (savePDandCinCalo && isThisVolume(aTrack->GetTouchable(),calo)) ||
00081 (savePDandCinMuon && isThisVolume(aTrack->GetTouchable(),muon)))
00082 flag = isItPrimaryDecayProductOrConversion(aTrack, *mother);
00083 if (saveFirstSecondary) flag = isItFromPrimary(*mother, flag);
00084 newTA.secondary(aTrack, *mother, flag);
00085
00086 if (aTrack->GetTrackStatus() == fStopAndKill) classification = fKill;
00087 if (killHeavy) {
00088 int pdg = aTrack->GetDefinition()->GetPDGEncoding();
00089 double ke = aTrack->GetKineticEnergy()/MeV;
00090 if (((pdg/1000000000 == 1) && (((pdg/10000)%100) > 0) &&
00091 (((pdg/10)%100) > 0) && (ke<kmaxIon)) ||
00092 ((pdg == 2212) && (ke < kmaxProton)) ||
00093 ((pdg == 2112) && (ke < kmaxNeutron))) classification = fKill;
00094 }
00095 if (!trackNeutrino) {
00096 int pdg = std::abs(aTrack->GetDefinition()->GetPDGEncoding());
00097 if (pdg == 12 || pdg == 14 || pdg == 16 || pdg == 18)
00098 classification = fKill;
00099 }
00100 if (isItLongLived(aTrack)) classification = fKill;
00101 if (killDeltaRay) {
00102 if (aTrack->GetCreatorProcess()->GetProcessType() == fElectromagnetic &&
00103 aTrack->GetCreatorProcess()->GetProcessSubType() == fIonisation)
00104 classification = fKill;
00105 }
00106 #ifdef DebugLog
00107 LogDebug("SimG4CoreApplication") << "StackingAction:Classify Track "
00108 << aTrack->GetTrackID() << " Parent "
00109 << aTrack->GetParentID() << " Type "
00110 << aTrack->GetDefinition()->GetParticleName()
00111 << " K.E. " << aTrack->GetKineticEnergy()/MeV
00112 << " MeV from process/subprocess "
00113 << aTrack->GetCreatorProcess()->GetProcessType() << "|"
00114 << aTrack->GetCreatorProcess()->GetProcessSubType()
00115 << " as " << classification << " Flag " << flag;
00116 #endif
00117 }
00118 return classification;
00119 }
00120
00121 void StackingAction::NewStage() {}
00122
00123 void StackingAction::PrepareNewEvent() {}
00124
00125 void StackingAction::initPointer() {
00126
00127 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
00128 if (lvs) {
00129 std::vector<G4LogicalVolume*>::const_iterator lvcite;
00130 for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
00131 if (savePDandCinTracker) {
00132 if (strcmp("Tracker",(*lvcite)->GetName().c_str()) == 0) tracker.push_back(*lvcite);
00133 if (strcmp("BEAM",(*lvcite)->GetName().substr(0,4).c_str()) == 0) tracker.push_back(*lvcite);
00134 }
00135 if (savePDandCinCalo) {
00136 if (strcmp("CALO",(*lvcite)->GetName().c_str()) == 0) calo.push_back(*lvcite);
00137 if (strcmp("VCAL",(*lvcite)->GetName().c_str()) == 0) calo.push_back(*lvcite);
00138 }
00139 if (savePDandCinMuon) {
00140 if (strcmp("MUON",(*lvcite)->GetName().c_str()) == 0) muon.push_back(*lvcite);
00141 }
00142 }
00143 edm::LogInfo("SimG4CoreApplication") << "# of LV for Tracker "
00144 << tracker.size() << " for Calo "
00145 << calo.size() << " for Muon "
00146 << muon.size();
00147 for (unsigned int i=0; i<tracker.size(); ++i)
00148 edm::LogInfo("SimG4CoreApplication") << "Tracker vol " << i << " name "
00149 << tracker[i]->GetName();
00150 for (unsigned int i=0; i<calo.size(); ++i)
00151 edm::LogInfo("SimG4CoreApplication") << "Calorimeter vol " <<i <<" name "
00152 << calo[i]->GetName();
00153 for (unsigned int i=0; i<muon.size(); ++i)
00154 edm::LogInfo("SimG4CoreApplication") << "Muon vol " << i << " name "
00155 << muon[i]->GetName();
00156 }
00157
00158 const G4RegionStore * rs = G4RegionStore::GetInstance();
00159 unsigned int num = maxTimeNames.size();
00160 if (num > 0) {
00161 std::vector<double> tofs;
00162 if (rs) {
00163 std::vector<G4Region*>::const_iterator rcite;
00164 for (rcite = rs->begin(); rcite != rs->end(); rcite++) {
00165 for (unsigned int i=0; i<num; i++) {
00166 if ((*rcite)->GetName() == (G4String)(maxTimeNames[i])) {
00167 maxTimeRegions.push_back(*rcite);
00168 tofs.push_back(maxTrackTimes[i]);
00169 break;
00170 }
00171 }
00172 if (tofs.size() == num) break;
00173 }
00174 }
00175 for (unsigned int i=0; i<tofs.size(); i++) {
00176 maxTrackTimes[i] = tofs[i];
00177 G4String name = "Unknown";
00178 if (maxTimeRegions[i]) name = maxTimeRegions[i]->GetName();
00179 edm::LogInfo("SimG4CoreApplication") << name << " with pointer "
00180 << maxTimeRegions[i]<<" KE cut off "
00181 << maxTrackTimes[i];
00182 }
00183 }
00184
00185 }
00186
00187 bool StackingAction::isThisVolume(const G4VTouchable* touch,
00188 std::vector<G4LogicalVolume*> & lvs) const {
00189
00190 bool flag = false;
00191 if (lvs.size() > 0 && touch !=0) {
00192 int level = ((touch->GetHistoryDepth())+1);
00193 if (level >= 3) {
00194 unsigned int ii = (unsigned int)(level - 3);
00195 flag = (std::count(lvs.begin(),lvs.end(),(touch->GetVolume(ii)->GetLogicalVolume())) != 0);
00196 }
00197 }
00198 return flag;
00199 }
00200
00201 int StackingAction::isItPrimaryDecayProductOrConversion(const G4Track * aTrack,
00202 const G4Track & mother) const {
00203
00204 int flag = 0;
00205 TrackInformationExtractor extractor;
00206 const TrackInformation & motherInfo(extractor(mother));
00207
00208 if (motherInfo.isPrimary()) {
00209 if (aTrack->GetCreatorProcess()->GetProcessType() == fDecay) flag = 1;
00210 else if (aTrack->GetCreatorProcess()->GetProcessType() == fElectromagnetic &&
00211 aTrack->GetCreatorProcess()->GetProcessSubType() == fGammaConversion) flag = 2;
00212 }
00213 return flag;
00214 }
00215
00216 int StackingAction::isItFromPrimary(const G4Track & mother, int flagIn) const {
00217
00218 int flag = flagIn;
00219 if (flag != 1) {
00220 TrackInformationExtractor extractor;
00221 const TrackInformation & motherInfo(extractor(mother));
00222 if (motherInfo.isPrimary()) flag = 3;
00223 }
00224 return flag;
00225 }
00226
00227 bool StackingAction::isItLongLived(const G4Track * aTrack) const {
00228
00229 bool flag = false;
00230 double time = (aTrack->GetGlobalTime())/nanosecond;
00231 double tofM = maxTrackTime;
00232 if (maxTimeRegions.size() > 0) {
00233 G4Region* reg = aTrack->GetTouchable()->GetVolume(0)->GetLogicalVolume()->GetRegion();
00234 for (unsigned int i=0; i<maxTimeRegions.size(); i++) {
00235 if (reg == maxTimeRegions[i]) {
00236 tofM = maxTrackTimes[i];
00237 break;
00238 }
00239 }
00240 }
00241 if (time > tofM) flag = true;
00242 return flag;
00243 }