00001 #include "SimG4CMS/Muon/interface/MuonSensitiveDetector.h"
00002 #include "SimG4CMS/Muon/interface/MuonSlaveSD.h"
00003 #include "SimG4CMS/Muon//interface/MuonEndcapFrameRotation.h"
00004 #include "SimG4CMS/Muon/interface/MuonGemFrameRotation.h"
00005 #include "SimG4CMS/Muon/interface/MuonRpcFrameRotation.h"
00006 #include "Geometry/MuonNumbering/interface/MuonSubDetector.h"
00007
00008 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00009
00010 #include "SimG4CMS/Muon/interface/SimHitPrinter.h"
00011 #include "SimDataFormats/TrackingHit/interface/UpdatablePSimHit.h"
00012
00013 #include "SimG4CMS/Muon/interface/MuonG4Numbering.h"
00014 #include "Geometry/MuonNumbering/interface/MuonBaseNumber.h"
00015 #include "Geometry/MuonNumbering/interface/MuonSimHitNumberingScheme.h"
00016
00017 #include "SimG4Core/Notification/interface/TrackInformation.h"
00018 #include "SimG4Core/Notification/interface/G4TrackToParticleID.h"
00019 #include "SimG4Core/Physics/interface/G4ProcessTypeEnumerator.h"
00020
00021 #include "G4SDManager.hh"
00022 #include "G4VProcess.hh"
00023 #include "G4EventManager.hh"
00024
00025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00026
00027 #include <iostream>
00028
00029 MuonSensitiveDetector::MuonSensitiveDetector(std::string name,
00030 const DDCompactView & cpv,
00031 SensitiveDetectorCatalog & clg,
00032 edm::ParameterSet const & p,
00033 const SimTrackManager* manager)
00034 : SensitiveTkDetector(name, cpv, clg, p),
00035 thePV(0), theHit(0), theDetUnitId(0), theTrackID(0), theManager(manager)
00036 {
00037 edm::ParameterSet m_MuonSD = p.getParameter<edm::ParameterSet>("MuonSD");
00038 STenergyPersistentCut = m_MuonSD.getParameter<double>("EnergyThresholdForPersistency");
00039 STallMuonsPersistent = m_MuonSD.getParameter<bool>("AllMuonsPersistent");
00040 printHits = m_MuonSD.getParameter<bool>("PrintHits");
00041
00042
00043
00044
00045
00046 LogDebug("MuonSimDebug") << "create MuonSubDetector "<<name<<std::endl;
00047
00048 detector = new MuonSubDetector(name);
00049
00050 LogDebug("MuonSimDebug") << "create MuonFrameRotation"<<std::endl;
00051
00052 if (detector->isEndcap()) {
00053
00054 theRotation=new MuonEndcapFrameRotation();
00055 } else if (detector->isRpc()) {
00056
00057 theRotation=new MuonRpcFrameRotation( cpv );
00058 } else if (detector->isGem()) {
00059
00060 theRotation=new MuonGemFrameRotation( cpv );
00061 } else {
00062 theRotation = 0;
00063 }
00064 LogDebug("MuonSimDebug") << "create MuonSlaveSD"<<std::endl;
00065 slaveMuon = new MuonSlaveSD(detector,theManager);
00066 LogDebug("MuonSimDebug") << "create MuonSimHitNumberingScheme"<<std::endl;
00067 numbering = new MuonSimHitNumberingScheme(detector, cpv);
00068 g4numbering = new MuonG4Numbering(cpv);
00069
00070
00071
00072
00073
00074 std::vector<std::string> lvNames = clg.logicalNames(name);
00075 this->Register();
00076 for (std::vector<std::string>::iterator it = lvNames.begin(); it != lvNames.end(); it++){
00077 LogDebug("MuonSimDebug") << name << " MuonSensitiveDetector:: attaching SD to LV " << *it << std::endl;
00078 this->AssignSD(*it);
00079 }
00080
00081 if (printHits) {
00082 thePrinter = new SimHitPrinter("HitPositionOSCAR.dat");
00083 }
00084
00085
00086 LogDebug("MuonSimDebug") << " EnergyThresholdForPersistency " << STenergyPersistentCut << " AllMuonsPersistent " << STallMuonsPersistent << std::endl;
00087
00088 theG4ProcessTypeEnumerator = new G4ProcessTypeEnumerator;
00089 myG4TrackToParticleID = new G4TrackToParticleID;
00090
00091 }
00092
00093
00094 MuonSensitiveDetector::~MuonSensitiveDetector() {
00095 delete g4numbering;
00096 delete numbering;
00097 delete slaveMuon;
00098 delete theRotation;
00099 delete detector;
00100
00101 delete theG4ProcessTypeEnumerator;
00102
00103 delete myG4TrackToParticleID;
00104 }
00105
00106 void MuonSensitiveDetector::update(const BeginOfEvent * i){
00107 clearHits();
00108
00109
00110 thePV = 0;
00111 theDetUnitId = 0;
00112 theTrackID = 0;
00113
00114 }
00115
00116 void MuonSensitiveDetector::update(const ::EndOfEvent * ev)
00117 {
00118
00119 }
00120
00121
00122 void MuonSensitiveDetector::clearHits()
00123 {
00124 LogDebug("MuonSimDebug") << "MuonSensitiveDetector::clearHits"<<std::endl;
00125 slaveMuon->Initialize();
00126 }
00127
00128 bool MuonSensitiveDetector::ProcessHits(G4Step * aStep, G4TouchableHistory * ROhist)
00129 {
00130 LogDebug("MuonSimDebug") <<" MuonSensitiveDetector::ProcessHits "<<InitialStepPosition(aStep,WorldCoordinates)<<std::endl;
00131
00132
00133
00134 if (aStep->GetTotalEnergyDeposit()>0.){
00135
00136 if( aStep->GetTrack()->GetDynamicParticle()->GetCharge() != 0 ){
00137
00138 if (newHit(aStep)) {
00139 saveHit();
00140 createHit(aStep);
00141 } else {
00142 updateHit(aStep);
00143 }
00144 return true;
00145 } else {
00146 storeVolumeAndTrack(aStep);
00147 return false;
00148 }
00149 }
00150 return false;
00151 }
00152
00153 uint32_t MuonSensitiveDetector::setDetUnitId(G4Step * aStep)
00154 {
00155
00156 MuonBaseNumber num = g4numbering->PhysicalVolumeToBaseNumber(aStep);
00157 return numbering->baseNumberToUnitNumber(num);
00158 }
00159
00160
00161 Local3DPoint MuonSensitiveDetector::toOrcaRef(Local3DPoint in ,G4Step * s){
00162 if (theRotation !=0 ) {
00163 return theRotation->transformPoint(in,s);
00164 }
00165 return (in);
00166 }
00167
00168 Local3DPoint MuonSensitiveDetector::toOrcaUnits(Local3DPoint in){
00169 return Local3DPoint(in.x()/cm,in.y()/cm,in.z()/cm);
00170 }
00171
00172 Global3DPoint MuonSensitiveDetector::toOrcaUnits(Global3DPoint in){
00173 return Global3DPoint(in.x()/cm,in.y()/cm,in.z()/cm);
00174 }
00175
00176 void MuonSensitiveDetector::storeVolumeAndTrack(G4Step * aStep) {
00177 G4VPhysicalVolume* pv = aStep->GetPreStepPoint()->GetPhysicalVolume();
00178 G4Track * t = aStep->GetTrack();
00179 thePV=pv;
00180 theTrackID=t->GetTrackID();
00181 }
00182
00183 bool MuonSensitiveDetector::newHit(G4Step * aStep){
00184
00185 G4VPhysicalVolume* pv = aStep->GetPreStepPoint()->GetPhysicalVolume();
00186 G4Track * t = aStep->GetTrack();
00187 uint32_t currentUnitId=setDetUnitId(aStep);
00188 unsigned int currentTrackID=t->GetTrackID();
00189
00190 bool changed=((pv!=thePV) ||
00191 (currentUnitId!=theDetUnitId) ||
00192 (currentTrackID!=theTrackID));
00193 return changed;
00194 }
00195
00196 void MuonSensitiveDetector::createHit(G4Step * aStep){
00197
00198 G4Track * theTrack = aStep->GetTrack();
00199
00200 Local3DPoint theEntryPoint;
00201 Local3DPoint theExitPoint;
00202
00203 if (detector->isBarrel()) {
00204 theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPositionVsParent(aStep,1),aStep));
00205 theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,1),aStep));
00206 } else if (detector->isEndcap()) {
00207
00208 theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPosition(aStep,LocalCoordinates),aStep));
00209 theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
00210 float zentry = theEntryPoint.z();
00211 float zexit = theExitPoint.z();
00212 Local3DPoint tempEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPositionVsParent(aStep,4),aStep));
00213 Local3DPoint tempExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,4),aStep));
00214
00215 theEntryPoint = Local3DPoint( tempEntryPoint.x(), tempEntryPoint.y(), zentry );
00216 theExitPoint = Local3DPoint( tempExitPoint.x(), tempExitPoint.y(), zexit );
00217 } else {
00218 theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPosition(aStep,LocalCoordinates),aStep));
00219 theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
00220 }
00221
00222 float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
00223 float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
00224 float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
00225
00226 int theParticleType = myG4TrackToParticleID->particleID(theTrack);
00227 G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
00228 G4ThreeVector lmd = ((const G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()
00229 ->GetTopTransform().TransformAxis(gmd);
00230 Local3DPoint lnmd = toOrcaRef(ConvertToLocal3DPoint(lmd),aStep);
00231 float theThetaAtEntry = lnmd.theta();
00232 float thePhiAtEntry = lnmd.phi();
00233
00234 storeVolumeAndTrack( aStep );
00235 theDetUnitId = setDetUnitId(aStep);
00236
00237 Global3DPoint theGlobalPos;
00238 if (printHits) {
00239 Local3DPoint theGlobalHelp = InitialStepPosition(aStep,WorldCoordinates);
00240 theGlobalEntry = toOrcaUnits(Global3DPoint (theGlobalHelp.x(),theGlobalHelp.y(),theGlobalHelp.z()));
00241
00242 G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
00243 const G4TouchableHistory * theTouchable=(const G4TouchableHistory *)
00244 (preStepPoint->GetTouchable());
00245 theGlobalHelp=ConvertToLocal3DPoint(theTouchable->GetTranslation());
00246 theGlobalPos = toOrcaUnits(Global3DPoint (theGlobalHelp.x(),theGlobalHelp.y(),theGlobalHelp.z()));
00247
00248 }
00249
00250 LogDebug("MuonSimDebug") << "MuonSensitiveDetector::createHit UpdatablePSimHit"<<std::endl;
00251
00252 theHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
00253 theEnergyLoss,theParticleType,theDetUnitId,
00254 theTrackID,theThetaAtEntry,thePhiAtEntry,
00255 theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
00256
00257
00258 LogDebug("MuonSimDebug") <<"=== NEW ==================> ELOSS = "<<theEnergyLoss<<" "
00259 <<thePV->GetLogicalVolume()->GetName()<<std::endl;
00260 const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep();
00261 const G4VProcess* p2 = aStep->GetPreStepPoint()->GetProcessDefinedStep();
00262 if (p)
00263 LogDebug("MuonSimDebug") <<" POST PROCESS = "<<p->GetProcessName()<<std::endl;
00264 if (p2)
00265 LogDebug("MuonSimDebug") <<" PRE PROCESS = "<<p2->GetProcessName()<<std::endl;
00266 LogDebug("MuonSimDebug") << "newhit theta " << theThetaAtEntry<<std::endl;
00267 LogDebug("MuonSimDebug") << "newhit phi " << thePhiAtEntry<<std::endl;
00268 LogDebug("MuonSimDebug") << "newhit pabs " << thePabs<<std::endl;
00269 LogDebug("MuonSimDebug") << "newhit tof " << theTof<<std::endl;
00270 LogDebug("MuonSimDebug") << "newhit track " << theTrackID<<std::endl;
00271 LogDebug("MuonSimDebug") << "newhit entry " << theEntryPoint<<std::endl;
00272 LogDebug("MuonSimDebug") << "newhit exit " << theExitPoint<<std::endl;
00273 LogDebug("MuonSimDebug") << "newhit eloss " << theEnergyLoss << std::endl;
00274 LogDebug("MuonSimDebug") << "newhit detid " << theDetUnitId<<std::endl;
00275 LogDebug("MuonSimDebug") << "newhit delta " << (theExitPoint-theEntryPoint)<<std::endl;
00276 LogDebug("MuonSimDebug") << "newhit deltu " << (theExitPoint-theEntryPoint).unit();
00277 LogDebug("MuonSimDebug") << " " << (theExitPoint-theEntryPoint).mag()<<std::endl;
00278 LogDebug("MuonSimDebug") << "newhit glob " << theGlobalEntry<<std::endl;
00279 LogDebug("MuonSimDebug") << "newhit dpos " << theGlobalPos<<std::endl;
00280 LogDebug("MuonSimDebug") << "newhit drot " << std::endl;
00281
00282
00283
00284
00285
00286
00287 int thePID = theTrack->GetDefinition()->GetPDGEncoding();
00288 LogDebug("MuonSimDebug") << " checking simtrack " << thePID << " " << thePabs << " STenergyPersistentCut " << STenergyPersistentCut << std::endl;
00289
00290 if( thePabs*GeV > STenergyPersistentCut
00291 || ( abs(thePID) == 13 && STallMuonsPersistent ) ){
00292 TrackInformation* info = getOrCreateTrackInformation(theTrack);
00293 LogDebug("MuonSimDebug") <<" track leaving hit in muons made selected for persistency"<<std::endl;
00294
00295 info->storeTrack(true);
00296 }
00297
00298 }
00299
00300 void MuonSensitiveDetector::updateHit(G4Step * aStep){
00301
00302
00303
00304
00305 Local3DPoint theExitPoint;
00306
00307 if (detector->isBarrel()) {
00308 theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,1),aStep));
00309 } else if (detector->isEndcap()) {
00310
00311 theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
00312 float zexit = theExitPoint.z();
00313 Local3DPoint tempExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,4),aStep));
00314 theExitPoint = Local3DPoint( tempExitPoint.x(), tempExitPoint.y(), zexit );
00315 } else {
00316 theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
00317 }
00318
00319 float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
00320
00321 if( theHit == 0 ){
00322 std::cerr << "!!ERRROR in MuonSensitiveDetector::updateHit. It is called when there is no hit " << std::endl;
00323 }
00324
00325 theHit->updateExitPoint(theExitPoint);
00326 theHit->addEnergyLoss(theEnergyLoss);
00327
00328 LogDebug("MuonSimDebug") <<"=== UPDATE ===============> ELOSS = "<<theEnergyLoss<<" "
00329 <<thePV->GetLogicalVolume()->GetName()<<std::endl;
00330 const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep();
00331 const G4VProcess* p2 = aStep->GetPreStepPoint()->GetProcessDefinedStep();
00332 if (p)
00333 LogDebug("MuonSimDebug") <<" POST PROCESS = "<<p->GetProcessName()<<std::endl;
00334 if (p2)
00335 LogDebug("MuonSimDebug") <<" PRE PROCESS = "<<p2->GetProcessName()<<std::endl;
00336 LogDebug("MuonSimDebug") << "updhit exit " << theExitPoint<<std::endl;
00337 LogDebug("MuonSimDebug") << "updhit eloss " << theHit->energyLoss() <<std::endl;
00338 LogDebug("MuonSimDebug") << "updhit detid " << theDetUnitId<<std::endl;
00339 LogDebug("MuonSimDebug") << "updhit delta " << (theExitPoint-theHit->entryPoint())<<std::endl;
00340 LogDebug("MuonSimDebug") << "updhit deltu " << (theExitPoint-theHit->entryPoint()).unit();
00341 LogDebug("MuonSimDebug") << " " << (theExitPoint-theHit->entryPoint()).mag()<<std::endl;
00342
00343 }
00344
00345 void MuonSensitiveDetector::saveHit(){
00346
00347 if (theHit) {
00348 if (printHits) {
00349 thePrinter->startNewSimHit(detector->name());
00350 thePrinter->printId(theHit->detUnitId());
00351
00352
00353
00354 thePrinter->printLocal(theHit->entryPoint(),theHit->exitPoint());
00355 thePrinter->printGlobal(theGlobalEntry);
00356 }
00357 slaveMuon->processHits(*theHit);
00358
00359
00360 delete theHit;
00361 theHit = 0;
00362 }
00363
00364 }
00365
00366 TrackInformation* MuonSensitiveDetector::getOrCreateTrackInformation( const G4Track* gTrack)
00367 {
00368 G4VUserTrackInformation* temp = gTrack->GetUserInformation();
00369 if (temp == 0){
00370 std::cerr <<" ERROR: no G4VUserTrackInformation available"<<std::endl;
00371 abort();
00372 }else{
00373 TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
00374 if (info ==0){
00375 std::cerr <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation"<<std::endl;
00376 abort();
00377 }
00378 return info;
00379 }
00380 }
00381
00382 void MuonSensitiveDetector::EndOfEvent(G4HCofThisEvent*)
00383 {
00384
00385
00386 saveHit();
00387 }
00388
00389
00390 void MuonSensitiveDetector::fillHits(edm::PSimHitContainer& c, std::string n){
00391
00392
00393
00394
00395 if (slaveMuon->name() == n) c=slaveMuon->hits();
00396
00397 }
00398
00399 std::vector<std::string> MuonSensitiveDetector::getNames(){
00400 std::vector<std::string> temp;
00401 temp.push_back(slaveMuon->name());
00402 return temp;
00403 }
00404
00405 Local3DPoint MuonSensitiveDetector::InitialStepPositionVsParent(G4Step * currentStep, G4int levelsUp) {
00406
00407 G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
00408 G4ThreeVector globalCoordinates = preStepPoint->GetPosition();
00409
00410 const G4TouchableHistory * theTouchable=(const G4TouchableHistory *)
00411 (preStepPoint->GetTouchable());
00412
00413 G4int depth = theTouchable->GetHistory()->GetDepth();
00414 G4ThreeVector localCoordinates = theTouchable->GetHistory()
00415 ->GetTransform(depth-levelsUp).TransformPoint(globalCoordinates);
00416
00417 return ConvertToLocal3DPoint(localCoordinates);
00418 }
00419
00420 Local3DPoint MuonSensitiveDetector::FinalStepPositionVsParent(G4Step * currentStep, G4int levelsUp) {
00421
00422 G4StepPoint * postStepPoint = currentStep->GetPostStepPoint();
00423 G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
00424 G4ThreeVector globalCoordinates = postStepPoint->GetPosition();
00425
00426 const G4TouchableHistory * theTouchable = (const G4TouchableHistory *)
00427 (preStepPoint->GetTouchable());
00428
00429 G4int depth = theTouchable->GetHistory()->GetDepth();
00430 G4ThreeVector localCoordinates = theTouchable->GetHistory()
00431 ->GetTransform(depth-levelsUp).TransformPoint(globalCoordinates);
00432
00433 return ConvertToLocal3DPoint(localCoordinates);
00434 }