Go to the documentation of this file.00001 #include "SimG4CMS/Forward/interface/PLTSensitiveDetector.h"
00002
00003 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00004 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00005
00006 #include "SimG4Core/Notification/interface/TrackInformation.h"
00007 #include "SimG4Core/Notification/interface/G4TrackToParticleID.h"
00008 #include "SimG4Core/Physics/interface/G4ProcessTypeEnumerator.h"
00009
00010 #include "SimDataFormats/TrackingHit/interface/UpdatablePSimHit.h"
00011 #include "SimDataFormats/SimHitMaker/interface/TrackingSlaveSD.h"
00012
00013 #include "SimG4Core/Notification/interface/TrackInformation.h"
00014
00015 #include "FWCore/Framework/interface/ESTransientHandle.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 #include "FWCore/Framework/interface/EventSetup.h"
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019
00020 #include "G4Track.hh"
00021 #include "G4SDManager.hh"
00022 #include "G4VProcess.hh"
00023 #include "G4EventManager.hh"
00024 #include "G4Event.hh"
00025 #include "G4VProcess.hh"
00026
00027 #include <string>
00028 #include <vector>
00029 #include <iostream>
00030
00031
00032 PLTSensitiveDetector::PLTSensitiveDetector(std::string name,
00033 const DDCompactView & cpv,
00034 SensitiveDetectorCatalog & clg,
00035 edm::ParameterSet const & p,
00036 const SimTrackManager* manager) :
00037 SensitiveTkDetector(name, cpv, clg, p), myName(name), mySimHit(0),
00038 oldVolume(0), lastId(0), lastTrack(0), eventno(0) {
00039
00040 edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("PltSD");
00041 energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*GeV;
00042 energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*GeV;
00043
00044 edm::LogInfo("PLTSensitiveDetector") <<"Criteria for Saving Tracker SimTracks: \n "
00045 <<" History: "<<energyHistoryCut<< " MeV ; Persistency: "<< energyCut<<" MeV\n"
00046 <<" Constructing a PLTSensitiveDetector with ";
00047
00048 slave = new TrackingSlaveSD(name);
00049
00050
00051 std::vector<std::string> lvNames = clg.logicalNames(name);
00052 this->Register();
00053 for (std::vector<std::string>::iterator it = lvNames.begin();
00054 it != lvNames.end(); it++) {
00055 edm::LogInfo("PLTSensitiveDetector")<< name << " attaching LV " << *it;
00056 this->AssignSD(*it);
00057 }
00058
00059 theG4ProcessTypeEnumerator = new G4ProcessTypeEnumerator;
00060 myG4TrackToParticleID = new G4TrackToParticleID;
00061 }
00062
00063 PLTSensitiveDetector::~PLTSensitiveDetector() {
00064 delete slave;
00065 delete theG4ProcessTypeEnumerator;
00066 delete myG4TrackToParticleID;
00067 }
00068
00069 bool PLTSensitiveDetector::ProcessHits(G4Step * aStep, G4TouchableHistory *) {
00070
00071 LogDebug("PLTSensitiveDetector") << " Entering a new Step "
00072 << aStep->GetTotalEnergyDeposit() << " "
00073 << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
00074
00075 G4Track * gTrack = aStep->GetTrack();
00076 if ((unsigned int)(gTrack->GetTrackID()) != lastTrack) {
00077
00078 if (gTrack->GetKineticEnergy() > energyCut){
00079 TrackInformation* info = getOrCreateTrackInformation(gTrack);
00080 info->storeTrack(true);
00081 }
00082 if (gTrack->GetKineticEnergy() > energyHistoryCut){
00083 TrackInformation* info = getOrCreateTrackInformation(gTrack);
00084 info->putInHistory();
00085 }
00086 }
00087
00088 if (aStep->GetTotalEnergyDeposit()>0.) {
00089 if (newHit(aStep) == true) {
00090 sendHit();
00091 createHit(aStep);
00092 } else {
00093 updateHit(aStep);
00094 }
00095 return true;
00096 }
00097 return false;
00098 }
00099
00100 uint32_t PLTSensitiveDetector::setDetUnitId(G4Step * ) {
00101
00102 unsigned int detId = 0;
00103
00104 LogDebug("PLTSensitiveDetector")<< " DetID = "<<detId;
00105
00106 return detId;
00107 }
00108
00109 void PLTSensitiveDetector::EndOfEvent(G4HCofThisEvent *) {
00110
00111 LogDebug("PLTSensitiveDetector")<< " Saving the last hit in a ROU " << myName;
00112
00113 if (mySimHit == 0) return;
00114 sendHit();
00115 }
00116
00117 void PLTSensitiveDetector::fillHits(edm::PSimHitContainer& c, std::string n){
00118 if (slave->name() == n) c=slave->hits();
00119 }
00120
00121 void PLTSensitiveDetector::sendHit() {
00122 if (mySimHit == 0) return;
00123 LogDebug("PLTSensitiveDetector") << " Storing PSimHit: " << pname << " " << mySimHit->detUnitId()
00124 << " " << mySimHit->trackId() << " " << mySimHit->energyLoss()
00125 << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
00126
00127 slave->processHits(*mySimHit);
00128
00129
00130 delete mySimHit;
00131 mySimHit = 0;
00132 lastTrack = 0;
00133 lastId = 0;
00134 }
00135
00136 void PLTSensitiveDetector::updateHit(G4Step * aStep) {
00137
00138 Local3DPoint theExitPoint = SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates);
00139 float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
00140 mySimHit->setExitPoint(theExitPoint);
00141 LogDebug("PLTSensitiveDetector")<< " Before : " << mySimHit->energyLoss();
00142 mySimHit->addEnergyLoss(theEnergyLoss);
00143 globalExitPoint = SensitiveDetector::FinalStepPosition(aStep,WorldCoordinates);
00144
00145 LogDebug("PLTSensitiveDetector") << " Updating: new exitpoint " << pname << " "
00146 << theExitPoint << " new energy loss " << theEnergyLoss
00147 << "\n Updated PSimHit: " << mySimHit->detUnitId()
00148 << " " << mySimHit->trackId()
00149 << " " << mySimHit->energyLoss() << " "
00150 << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
00151 }
00152
00153 bool PLTSensitiveDetector::newHit(G4Step * aStep) {
00154
00155 G4Track * theTrack = aStep->GetTrack();
00156 uint32_t theDetUnitId = setDetUnitId(aStep);
00157 unsigned int theTrackID = theTrack->GetTrackID();
00158
00159 LogDebug("PLTSensitiveDetector") << " OLD (d,t) = (" << lastId << "," << lastTrack
00160 << "), new = (" << theDetUnitId << "," << theTrackID << ") return "
00161 << ((theTrackID == lastTrack) && (lastId == theDetUnitId));
00162 if ((mySimHit != 0) && (theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep))
00163 return false;
00164 return true;
00165 }
00166
00167 bool PLTSensitiveDetector::closeHit(G4Step * aStep) {
00168
00169 if (mySimHit == 0) return false;
00170 const float tolerance = 0.05 * mm;
00171
00172 Local3DPoint theEntryPoint = SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates);
00173 LogDebug("PLTSensitiveDetector")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag();
00174
00175 if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true;
00176 return false;
00177 }
00178
00179 void PLTSensitiveDetector::createHit(G4Step * aStep) {
00180
00181 if (mySimHit != 0) {
00182 delete mySimHit;
00183 mySimHit=0;
00184 }
00185
00186 G4Track * theTrack = aStep->GetTrack();
00187 G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
00188
00189 Local3DPoint theEntryPoint = SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates);
00190 Local3DPoint theExitPoint = SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates);
00191
00192 float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
00193 float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
00194 float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
00195 int theParticleType = myG4TrackToParticleID->particleID(theTrack);
00196 uint32_t theDetUnitId = setDetUnitId(aStep);
00197
00198 globalEntryPoint = SensitiveDetector::InitialStepPosition(aStep,WorldCoordinates);
00199 globalExitPoint = SensitiveDetector::FinalStepPosition(aStep,WorldCoordinates);
00200 pname = theTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
00201
00202 unsigned int theTrackID = theTrack->GetTrackID();
00203
00204 G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
00205
00206 G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd);
00207 Local3DPoint lnmd = ConvertToLocal3DPoint(lmd);
00208 float theThetaAtEntry = lnmd.theta();
00209 float thePhiAtEntry = lnmd.phi();
00210
00211 mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
00212 theEnergyLoss,theParticleType,theDetUnitId,
00213 theTrackID,theThetaAtEntry,thePhiAtEntry,
00214 theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
00215
00216 LogDebug("PLTSensitiveDetector") << " Created PSimHit: " << pname << " "
00217 << mySimHit->detUnitId() << " " << mySimHit->trackId()
00218 << " " << mySimHit->energyLoss() << " "
00219 << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
00220 lastId = theDetUnitId;
00221 lastTrack = theTrackID;
00222 oldVolume = v;
00223 }
00224
00225 void PLTSensitiveDetector::update(const BeginOfJob * ) { }
00226
00227 void PLTSensitiveDetector::update(const BeginOfEvent * i) {
00228
00229 clearHits();
00230 eventno = (*i)()->GetEventID();
00231 mySimHit = 0;
00232 }
00233
00234 void PLTSensitiveDetector::update(const BeginOfTrack *bot) {
00235
00236 const G4Track* gTrack = (*bot)();
00237 pname = gTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
00238 }
00239
00240 void PLTSensitiveDetector::clearHits() {
00241 slave->Initialize();
00242 }
00243
00244 TrackInformation* PLTSensitiveDetector::getOrCreateTrackInformation( const G4Track* gTrack) {
00245 G4VUserTrackInformation* temp = gTrack->GetUserInformation();
00246 if (temp == 0){
00247 edm::LogError("PLTSensitiveDetector") <<" ERROR: no G4VUserTrackInformation available";
00248 abort();
00249 }else{
00250 TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
00251 if (info == 0){
00252 edm::LogError("PLTSensitiveDetector") <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation";
00253 abort();
00254 }
00255 return info;
00256 }
00257 }