00001 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00002 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00003
00004 #include "SimG4Core/SensitiveDetector/interface/FrameRotation.h"
00005 #include "SimG4Core/Notification/interface/TrackInformation.h"
00006 #include "SimG4Core/Notification/interface/G4TrackToParticleID.h"
00007 #include "SimG4Core/Physics/interface/G4ProcessTypeEnumerator.h"
00008
00009 #include "SimDataFormats/TrackingHit/interface/UpdatablePSimHit.h"
00010 #include "SimDataFormats/SimHitMaker/interface/TrackingSlaveSD.h"
00011
00012 #include "SimG4CMS/Tracker/interface/TkAccumulatingSensitiveDetector.h"
00013 #include "SimG4CMS/Tracker/interface/FakeFrameRotation.h"
00014 #include "SimG4CMS/Tracker/interface/StandardFrameRotation.h"
00015 #include "SimG4CMS/Tracker/interface/TrackerFrameRotation.h"
00016 #include "SimG4CMS/Tracker/interface/TkSimHitPrinter.h"
00017 #include "SimG4CMS/Tracker/interface/TrackerG4SimHitNumberingScheme.h"
00018
00019 #include "FWCore/Framework/interface/ESTransientHandle.h"
00020 #include "FWCore/Framework/interface/ESHandle.h"
00021 #include "FWCore/Framework/interface/EventSetup.h"
00022
00023 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00024 #include "Geometry/TrackerNumberingBuilder/interface/GeometricDet.h"
00025
00026 #include "SimG4Core/Notification/interface/TrackInformation.h"
00027
00028 #include "G4Track.hh"
00029 #include "G4SDManager.hh"
00030 #include "G4VProcess.hh"
00031 #include "G4EventManager.hh"
00032 #include "G4Event.hh"
00033
00034 #include <string>
00035 #include <vector>
00036 #include <iostream>
00037
00038 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00039
00040
00041
00042
00043 #ifdef DUMPPROCESSES
00044 #include "G4VProcess.hh"
00045 #endif
00046
00047 using std::cout;
00048 using std::endl;
00049 using std::vector;
00050 using std::string;
00051
00052 static
00053 TrackerG4SimHitNumberingScheme&
00054 numberingScheme(const DDCompactView& cpv, const GeometricDet& det)
00055 {
00056 static TrackerG4SimHitNumberingScheme s_scheme(cpv, det);
00057 return s_scheme;
00058 }
00059
00060
00061 TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector(string name,
00062 const DDCompactView & cpv,
00063 SensitiveDetectorCatalog & clg,
00064 edm::ParameterSet const & p,
00065 const SimTrackManager* manager) :
00066 SensitiveTkDetector(name, cpv, clg, p), myName(name), myRotation(0), mySimHit(0),theManager(manager),
00067 oldVolume(0), lastId(0), lastTrack(0), eventno(0) ,rTracker(1200.*mm),zTracker(3000.*mm),
00068 numberingScheme_(0)
00069 {
00070
00071 edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("TrackerSD");
00072 allowZeroEnergyLoss = m_TrackerSD.getParameter<bool>("ZeroEnergyLoss");
00073 neverAccumulate = m_TrackerSD.getParameter<bool>("NeverAccumulate");
00074 printHits = m_TrackerSD.getParameter<bool>("PrintHits");
00075 theSigma = m_TrackerSD.getParameter<double>("ElectronicSigmaInNanoSeconds");
00076 energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*GeV;
00077 energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*GeV;
00078
00079 edm::LogInfo("TrackerSimInfo") <<"Criteria for Saving Tracker SimTracks: ";
00080 edm::LogInfo("TrackerSimInfo")<<" History: "<<energyHistoryCut<< " MeV ; Persistency: "<< energyCut<<" MeV ";
00081 edm::LogInfo("TrackerSimInfo")<<" Constructing a TkAccumulatingSensitiveDetector with ";
00082
00083 #ifndef FAKEFRAMEROTATION
00084
00085 if (name.find("TrackerHits") != string::npos)
00086 {
00087 edm::LogInfo("TrackerSimInfo")<<" TkAccumulatingSensitiveDetector: using TrackerFrameRotation for "<<myName;
00088 myRotation = new TrackerFrameRotation;
00089 }
00090
00091 if (myRotation == 0)
00092 {
00093 edm::LogInfo("TrackerSimInfo")<<" TkAccumulatingSensitiveDetector: using StandardFrameRotation for "<<myName;
00094 myRotation = new StandardFrameRotation;
00095 }
00096 #else
00097 myRotation = new FakeFrameRotation;
00098 edm::LogWarning("TrackerSimInfo")<< " WARNING - Using FakeFrameRotation in TkAccumulatingSensitiveDetector;";
00099 #endif
00100
00101 slaveLowTof = new TrackingSlaveSD(name+"LowTof");
00102 slaveHighTof = new TrackingSlaveSD(name+"HighTof");
00103
00104
00105 vector<string> lvNames = clg.logicalNames(name);
00106 this->Register();
00107 for (vector<string>::iterator it = lvNames.begin(); it != lvNames.end(); it++)
00108 {
00109 edm::LogInfo("TrackerSimInfo")<< name << " attaching LV " << *it;
00110 this->AssignSD(*it);
00111 }
00112
00113 theG4ProcessTypeEnumerator = new G4ProcessTypeEnumerator;
00114 myG4TrackToParticleID = new G4TrackToParticleID;
00115 }
00116
00117 TkAccumulatingSensitiveDetector::~TkAccumulatingSensitiveDetector()
00118 {
00119 delete slaveLowTof;
00120 delete slaveHighTof;
00121 delete theG4ProcessTypeEnumerator;
00122 delete myG4TrackToParticleID;
00123 }
00124
00125
00126 bool TkAccumulatingSensitiveDetector::ProcessHits(G4Step * aStep, G4TouchableHistory * ROhist)
00127 {
00128
00129 LogDebug("TrackerSimDebug")<< " Entering a new Step " << aStep->GetTotalEnergyDeposit() << " "
00130 << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
00131
00132
00133 if (aStep->GetTotalEnergyDeposit()>0. || allowZeroEnergyLoss == true)
00134 {
00135 if (newHit(aStep) == true)
00136 {
00137 sendHit();
00138 createHit(aStep);
00139 }
00140 else
00141 {
00142 updateHit(aStep);
00143 }
00144 return true;
00145 }
00146 return false;
00147 }
00148
00149 uint32_t TkAccumulatingSensitiveDetector::setDetUnitId(G4Step * s)
00150 {
00151 unsigned int detId = numberingScheme_->g4ToNumberingScheme(s->GetPreStepPoint()->GetTouchable());
00152
00153 LogDebug("TrackerSimDebug")<< " DetID = "<<detId;
00154
00155 return detId;
00156 }
00157
00158
00159 int TkAccumulatingSensitiveDetector::tofBin(float tof)
00160 {
00161
00162 float threshold = 3. * theSigma;
00163 if (tof < threshold) return 1;
00164 return 2;
00165 }
00166
00167 Local3DPoint TkAccumulatingSensitiveDetector::toOrcaRef(Local3DPoint in ,G4VPhysicalVolume * v)
00168 {
00169 if (myRotation !=0) return myRotation->transformPoint(in,v);
00170 return (in);
00171 }
00172
00173
00174 void TkAccumulatingSensitiveDetector::update(const BeginOfTrack *bot){
00175 const G4Track* gTrack = (*bot)();
00176 #ifdef DUMPPROCESSES
00177 edm::LogInfo("TrackerSimInfo") <<" -> process creator pointer "<<gTrack->GetCreatorProcess();
00178 if (gTrack->GetCreatorProcess())
00179 edm::LogInfo("TrackerSimInfo")<<" -> PROCESS CREATOR : "<<gTrack->GetCreatorProcess()->GetProcessName();
00180 #endif
00181
00182
00183
00184
00185
00186 const G4ThreeVector pos = gTrack->GetPosition();
00187
00188
00189
00190
00191
00192
00193
00194 if (pos.perp() < rTracker &&std::fabs(pos.z()) < zTracker){
00195
00196
00197
00198 LogDebug("TrackerSimDebug")<<" INSIDE TRACKER";
00199
00200 if (gTrack->GetKineticEnergy() > energyCut){
00201 TrackInformation* info = getOrCreateTrackInformation(gTrack);
00202
00203
00204
00205
00206 info->storeTrack(true);
00207 }
00208
00209
00210
00211 if (gTrack->GetKineticEnergy() > energyHistoryCut){
00212 TrackInformation* info = getOrCreateTrackInformation(gTrack);
00213 info->putInHistory();
00214
00215
00216
00217 }
00218
00219 }
00220 }
00221
00222
00223
00224 void TkAccumulatingSensitiveDetector::sendHit()
00225 {
00226 if (mySimHit == 0) return;
00227 LogDebug("TrackerSimDebug")<< " Storing PSimHit: " << pname << " " << mySimHit->detUnitId()
00228 << " " << mySimHit->trackId() << " " << mySimHit->energyLoss()
00229 << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
00230 if (printHits)
00231 {
00232 TkSimHitPrinter thePrinter("TkHitPositionOSCAR.dat");
00233 thePrinter.startNewSimHit(myName,oldVolume->GetLogicalVolume()->GetName(),
00234 mySimHit->detUnitId(),mySimHit->trackId(),eventno);
00235 thePrinter.printLocal(mySimHit->entryPoint(),mySimHit->exitPoint());
00236 thePrinter.printGlobal(globalEntryPoint,globalExitPoint);
00237 thePrinter.printHitData(mySimHit->energyLoss(),mySimHit->timeOfFlight());
00238 thePrinter.printMomentumOfTrack(mySimHit->pabs(),pname,
00239 thePrinter.getPropagationSign(globalEntryPoint,globalExitPoint));
00240 thePrinter.printGlobalMomentum(px,py,pz);
00241 }
00242
00243 if (tofBin(mySimHit->timeOfFlight()) == 1)
00244 slaveLowTof->processHits(*mySimHit);
00245 else
00246 slaveHighTof->processHits(*mySimHit);
00247
00248
00249 delete mySimHit;
00250 mySimHit = 0;
00251 lastTrack = 0;
00252 lastId = 0;
00253 }
00254
00255 void TkAccumulatingSensitiveDetector::createHit(G4Step * aStep)
00256 {
00257 if (mySimHit != 0)
00258 {
00259 delete mySimHit;
00260 mySimHit=0;
00261 }
00262
00263 G4Track * theTrack = aStep->GetTrack();
00264
00265 G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
00266 Local3DPoint theEntryPoint;
00267 Local3DPoint theExitPoint =
00268 toOrcaRef(SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates),v);
00269
00270
00271
00272
00273
00274 if(0.0 == theTrack->GetDefinition()->GetPDGCharge()) {
00275 theEntryPoint = theExitPoint;
00276 } else {
00277 theEntryPoint = toOrcaRef(SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates),v);
00278 }
00279
00280
00281
00282
00283 checkExitPoint(theExitPoint);
00284 float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
00285 float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
00286 float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
00287 int theParticleType = myG4TrackToParticleID->particleID(theTrack);
00288 uint32_t theDetUnitId = setDetUnitId(aStep);
00289
00290
00291
00292
00293
00294 globalEntryPoint = SensitiveDetector::InitialStepPosition(aStep,WorldCoordinates);
00295 globalExitPoint = SensitiveDetector::FinalStepPosition(aStep,WorldCoordinates);
00296
00297 pname = theTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
00298
00299 if (theDetUnitId == 0)
00300 {
00301 edm::LogError("TrackerSimInfo") << " Error: theDetUnitId is not valid.";
00302 abort();
00303 }
00304 unsigned int theTrackID = theTrack->GetTrackID();
00305
00306
00307
00308
00309 unsigned int theTrackIDInsideTheSimHit=theTrackID;
00310
00311
00312 G4VUserTrackInformation * info = theTrack->GetUserInformation();
00313 if (info == 0) edm::LogError("TrackerSimInfo")<< " Error: no UserInformation available ";
00314 else
00315 {
00316 TrackInformation * temp = dynamic_cast<TrackInformation* >(info);
00317 if (temp ==0) edm::LogError("TrackerSimInfo")<< " Error:G4VUserTrackInformation is not a TrackInformation.";
00318 if (temp->storeTrack() == false)
00319 {
00320
00321 LogDebug("TrackerSimDebug")<< " TkAccumulatingSensitiveDetector:createHit(): setting the TrackID from "
00322 << theTrackIDInsideTheSimHit;
00323 theTrackIDInsideTheSimHit = theTrack->GetParentID();
00324 LogDebug("TrackerSimDebug")<< " to the mother one " << theTrackIDInsideTheSimHit << " " << theEnergyLoss;
00325 }
00326 else
00327 {
00328 LogDebug("TrackerSimDebug")<< " TkAccumulatingSensitiveDetector:createHit(): leaving the current TrackID "
00329 << theTrackIDInsideTheSimHit;
00330 }
00331 }
00332
00333 px = aStep->GetPreStepPoint()->GetMomentum().x()/GeV;
00334 py = aStep->GetPreStepPoint()->GetMomentum().y()/GeV;
00335 pz = aStep->GetPreStepPoint()->GetMomentum().z()/GeV;
00336
00337 G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
00338
00339 G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()
00340 ->GetTopTransform().TransformAxis(gmd);
00341 Local3DPoint lnmd = toOrcaRef(ConvertToLocal3DPoint(lmd),v);
00342 float theThetaAtEntry = lnmd.theta();
00343 float thePhiAtEntry = lnmd.phi();
00344
00345 mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
00346 theEnergyLoss,theParticleType,theDetUnitId,
00347 theTrackIDInsideTheSimHit,theThetaAtEntry,thePhiAtEntry,
00348 theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
00349 LogDebug("TrackerSimDebug")<< " Created PSimHit: " << pname << " " << mySimHit->detUnitId() << " " << mySimHit->trackId()
00350 << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint()
00351 << " " << mySimHit->exitPoint();
00352 lastId = theDetUnitId;
00353 lastTrack = theTrackID;
00354 oldVolume = v;
00355 }
00356
00357 void TkAccumulatingSensitiveDetector::updateHit(G4Step * aStep)
00358 {
00359 G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
00360 Local3DPoint theExitPoint = toOrcaRef(SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates),v);
00361
00362
00363
00364 checkExitPoint(theExitPoint);
00365 float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
00366 mySimHit->setExitPoint(theExitPoint);
00367 LogDebug("TrackerSimDebug")<< " Before : " << mySimHit->energyLoss();
00368 mySimHit->addEnergyLoss(theEnergyLoss);
00369 globalExitPoint = SensitiveDetector::FinalStepPosition(aStep,WorldCoordinates);
00370
00371 LogDebug("TrackerSimDebug")<< " Updating: new exitpoint " << pname << " " << theExitPoint
00372 << " new energy loss " << theEnergyLoss;
00373 LogDebug("TrackerSimDebug")<< " Updated PSimHit: " << mySimHit->detUnitId() << " " << mySimHit->trackId()
00374 << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint()
00375 << " " << mySimHit->exitPoint();
00376 }
00377
00378 bool TkAccumulatingSensitiveDetector::newHit(G4Step * aStep)
00379 {
00380 if (neverAccumulate == true) return true;
00381 G4Track * theTrack = aStep->GetTrack();
00382
00383
00384 if(0.0 == theTrack->GetDefinition()->GetPDGCharge()) return true;
00385
00386 uint32_t theDetUnitId = setDetUnitId(aStep);
00387 unsigned int theTrackID = theTrack->GetTrackID();
00388
00389 LogDebug("TrackerSimDebug")<< " OLD (d,t) = (" << lastId << "," << lastTrack
00390 << "), new = (" << theDetUnitId << "," << theTrackID << ") return "
00391 << ((theTrackID == lastTrack) && (lastId == theDetUnitId));
00392 if ((mySimHit != 0) && (theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep))
00393 return false;
00394 return true;
00395 }
00396
00397 bool TkAccumulatingSensitiveDetector::closeHit(G4Step * aStep)
00398 {
00399 if (mySimHit == 0) return false;
00400 const float tolerance = 0.05 * mm;
00401
00402 G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
00403 Local3DPoint theEntryPoint = toOrcaRef(SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates),v);
00404 LogDebug("TrackerSimDebug")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag();
00405
00406 if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true;
00407 return false;
00408 }
00409
00410 void TkAccumulatingSensitiveDetector::EndOfEvent(G4HCofThisEvent *)
00411 {
00412
00413 LogDebug("TrackerSimDebug")<< " Saving the last hit in a ROU " << myName;
00414
00415 if (mySimHit == 0) return;
00416 sendHit();
00417 }
00418
00419 void TkAccumulatingSensitiveDetector::update(const BeginOfEvent * i)
00420 {
00421 clearHits();
00422 eventno = (*i)()->GetEventID();
00423 mySimHit = 0;
00424 }
00425
00426 void TkAccumulatingSensitiveDetector::update(const BeginOfJob * i)
00427 {
00428 edm::ESHandle<GeometricDet> pDD;
00429 const edm::EventSetup* es = (*i)();
00430 es->get<IdealGeometryRecord>().get( pDD );
00431
00432 edm::ESTransientHandle<DDCompactView> pView;
00433 es->get<IdealGeometryRecord>().get(pView);
00434
00435 numberingScheme_=&(numberingScheme(*pView,*pDD));
00436 }
00437
00438 void TkAccumulatingSensitiveDetector::clearHits()
00439 {
00440 slaveLowTof->Initialize();
00441 slaveHighTof->Initialize();
00442 }
00443
00444 void TkAccumulatingSensitiveDetector::checkExitPoint(Local3DPoint p)
00445 {
00446 double z = p.z();
00447 if (std::abs(z)<0.3*mm) return;
00448 bool sendExc = false;
00449
00450 edm::LogWarning("TrackerSimInfo")<< " ************ Hit outside the detector ; Local z " << z
00451 << "; skipping event = " << sendExc;
00452 if (sendExc == true)
00453 {
00454 G4EventManager::GetEventManager()->AbortCurrentEvent();
00455 G4EventManager::GetEventManager()->GetNonconstCurrentEvent()->SetEventAborted();
00456 }
00457 }
00458
00459 TrackInformation* TkAccumulatingSensitiveDetector::getOrCreateTrackInformation( const G4Track* gTrack){
00460 G4VUserTrackInformation* temp = gTrack->GetUserInformation();
00461 if (temp == 0){
00462 edm::LogError("TrackerSimInfo") <<" ERROR: no G4VUserTrackInformation available";
00463 abort();
00464 }else{
00465 TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
00466 if (info ==0){
00467 edm::LogError("TrackerSimInfo")<<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation";
00468 abort();
00469 }
00470 return info;
00471 }
00472 }
00473
00474 void TkAccumulatingSensitiveDetector::fillHits(edm::PSimHitContainer& c, std::string n){
00475
00476
00477
00478
00479 if (slaveLowTof->name() == n) c=slaveLowTof->hits();
00480 if (slaveHighTof->name() == n) c=slaveHighTof->hits();
00481
00482 }
00483
00484 std::vector<string> TkAccumulatingSensitiveDetector::getNames(){
00485 std::vector<string> temp;
00486 temp.push_back(slaveLowTof->name());
00487 temp.push_back(slaveHighTof->name());
00488 return temp;
00489 }