CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/SimG4CMS/FP420/plugins/FP420SD.cc

Go to the documentation of this file.
00001 
00002 // File: FP420SD.cc
00003 // Date: 02.2006
00004 // Description: Sensitive Detector class for FP420
00005 // Modifications: 
00007  
00008 //#include "Geometry/Vector/interface/LocalPoint.h"
00009 //#include "Geometry/Vector/interface/LocalVector.h"
00010 
00011 #include "SimG4Core/SensitiveDetector/interface/FrameRotation.h"
00012 #include "SimG4Core/Notification/interface/TrackInformation.h"
00013 #include "SimG4Core/Notification/interface/G4TrackToParticleID.h"
00014 #include "SimG4Core/Physics/interface/G4ProcessTypeEnumerator.h"
00015 
00016 #include "SimDataFormats/SimHitMaker/interface/TrackingSlaveSD.h"
00017 #include "SimDataFormats/TrackingHit/interface/UpdatablePSimHit.h"
00018 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00019 
00020 #include "SimG4CMS/FP420/interface/FP420SD.h"
00021 #include "SimG4CMS/FP420/interface/FP420G4Hit.h"
00022 #include "SimG4CMS/FP420/interface/FP420G4HitCollection.h"
00023 #include "SimG4CMS/FP420/interface/FP420NumberingScheme.h"
00024 
00025 #include "FWCore/Framework/interface/ESHandle.h"
00026 #include "FWCore/Framework/interface/EventSetup.h"
00027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00028 
00029 #include "SimG4Core/Notification/interface/TrackInformation.h"
00030 
00031 #include "G4Track.hh"
00032 #include "G4SDManager.hh"
00033 #include "G4VProcess.hh"
00034 #include "G4EventManager.hh"
00035 #include "G4Step.hh"
00036 #include "G4ParticleTable.hh"
00037 
00038 #include <string>
00039 #include <vector>
00040 #include <iostream>
00041 
00042 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00043 
00044 #ifdef DUMPPROCESSES
00045 #include "G4VProcess.hh"
00046 #endif
00047 
00048 using std::cout;
00049 using std::endl;
00050 using std::vector;
00051 using std::string;
00052 
00053 //#define debug
00054 //-------------------------------------------------------------------
00055 FP420SD::FP420SD(std::string name, const DDCompactView & cpv,
00056                  SensitiveDetectorCatalog & clg, 
00057                  edm::ParameterSet const & p, const SimTrackManager* manager) :
00058   SensitiveTkDetector(name, cpv, clg, p), numberingScheme(0), name(name),
00059   hcID(-1), theHC(0), theManager(manager), currentHit(0), theTrack(0), 
00060   currentPV(0), unitID(0),  previousUnitID(0), preStepPoint(0), 
00061   postStepPoint(0), eventno(0){
00062 //-------------------------------------------------------------------
00063 /*
00064 FP420SD::FP420SD(G4String name, const DDCompactView & cpv,
00065                  edm::ParameterSet const & p, const SimTrackManager* manager) :
00066   CaloSD(name, cpv, p, manager), numberingScheme(0), name(name),hcID(-1),
00067   theHC(0), currentHit(0), theTrack(0), currentPV(0), 
00068   unitID(0),  previousUnitID(0), preStepPoint(0), postStepPoint(0), eventno(0){
00069 */
00070 //-------------------------------------------------------------------
00071 
00072 
00073     
00074     //Add FP420 Sentitive Detector Name
00075     collectionName.insert(name);
00076     
00077     
00078     //Parameters
00079       edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("FP420SD");
00080        int verbn = m_p.getUntrackedParameter<int>("Verbosity");
00081     //int verbn = 1;
00082     
00083     SetVerboseLevel(verbn);
00084     LogDebug("FP420Sim") 
00085       << "*******************************************************\n"
00086       << "*                                                     *\n"
00087       << "* Constructing a FP420SD  with name " << name << "\n"
00088       << "*                                                     *\n"
00089       << "*******************************************************";
00090     
00091     
00092     slave  = new TrackingSlaveSD(name);
00093     
00094     //
00095     // attach detectors (LogicalVolumes)
00096     //
00097     std::vector<std::string> lvNames = clg.logicalNames(name);
00098 
00099     this->Register();
00100 
00101     for (std::vector<std::string>::iterator it=lvNames.begin();  
00102          it !=lvNames.end(); it++) {
00103       this->AssignSD(*it);
00104       edm::LogInfo("FP420Sim") << "FP420SD : Assigns SD to LV " << (*it);
00105     }
00106     
00107     if      (name == "FP420SI") {
00108       if (verbn > 0) {
00109       edm::LogInfo("FP420Sim") << "name = FP420SI and  new FP420NumberingSchem";
00110     std::cout << "name = FP420SI and  new FP420NumberingSchem"<< std::endl;
00111       }
00112       numberingScheme = new FP420NumberingScheme() ;
00113     } else {
00114       edm::LogWarning("FP420Sim") << "FP420SD: ReadoutName not supported\n";
00115     std::cout << "FP420SD: ReadoutName not supported"<< std::endl;
00116     }
00117     
00118     edm::LogInfo("FP420Sim") << "FP420SD: Instantiation completed";
00119     std::cout << "FP420SD: Instantiation completed"<< std::endl;
00120   }
00121 
00122 
00123 
00124 
00125 FP420SD::~FP420SD() { 
00126   //AZ:
00127   if (slave) delete slave; 
00128 
00129   if (numberingScheme)
00130     delete numberingScheme;
00131 
00132 }
00133 
00134 double FP420SD::getEnergyDeposit(G4Step* aStep) {
00135   return aStep->GetTotalEnergyDeposit();
00136 }
00137 
00138 void FP420SD::Initialize(G4HCofThisEvent * HCE) { 
00139 #ifdef debug
00140     LogDebug("FP420Sim") << "FP420SD : Initialize called for " << name << std::endl;
00141 #endif
00142 
00143   theHC = new FP420G4HitCollection(name, collectionName[0]);
00144   if (hcID<0) 
00145     hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
00146   HCE->AddHitsCollection(hcID, theHC);
00147 
00148   tsID   = -2;
00149   //  primID = -2;
00150   primID = 0;
00151 
00153 }
00154 
00155 
00156 bool FP420SD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
00157 
00158   if (aStep == NULL) {
00159     return true;
00160   } else {
00161     GetStepInfo(aStep);
00162   //   LogDebug("FP420Sim") << edeposit <<std::endl;
00163 
00164     //AZ
00165 #ifdef debug
00166   LogDebug("FP420Sim") << "FP420SD :  number of hits = " << theHC->entries() << std::endl;
00167 #endif
00168 
00169     if (HitExists() == false && edeposit>0. &&  theHC->entries()< 200 ){ 
00170       CreateNewHit();
00171     return true;
00172     }
00173   }
00174     return true;
00175 } 
00176 
00177 void FP420SD::GetStepInfo(G4Step* aStep) {
00178   
00179   preStepPoint = aStep->GetPreStepPoint(); 
00180   postStepPoint= aStep->GetPostStepPoint(); 
00181   theTrack     = aStep->GetTrack();   
00182   hitPoint     = preStepPoint->GetPosition();   
00183   currentPV    = preStepPoint->GetPhysicalVolume();
00184   hitPointExit = postStepPoint->GetPosition();  
00185 
00186   hitPointLocal = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
00187   hitPointLocalExit = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPointExit);
00188 
00189 
00190   G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
00191   if (particleCode == emPDG ||
00192       particleCode == epPDG ||
00193       particleCode == gammaPDG ) {
00194     edepositEM  = getEnergyDeposit(aStep); edepositHAD = 0.;
00195   } else {
00196     edepositEM  = 0.; edepositHAD = getEnergyDeposit(aStep);
00197   }
00198   edeposit = aStep->GetTotalEnergyDeposit();
00199   tSlice    = (postStepPoint->GetGlobalTime() )/nanosecond;
00200   tSliceID  = (int) tSlice;
00201   unitID    = setDetUnitId(aStep);
00202 #ifdef debug
00203   LogDebug("FP420Sim") << "unitID=" << unitID <<std::endl;
00204 #endif
00205   primaryID    = theTrack->GetTrackID();
00206   //  Position     = hitPoint;
00207   Pabs         = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
00208   //Tof          = 1400. + aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
00209   Tof          = aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
00210   Eloss        = aStep->GetTotalEnergyDeposit()/GeV;
00211   ParticleType = theTrack->GetDefinition()->GetPDGEncoding();      
00212   ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta()/deg;
00213   PhiAtEntry   = aStep->GetPreStepPoint()->GetPosition().phi()/deg;
00214 
00215   ParentId = theTrack->GetParentID();
00216   Vx = theTrack->GetVertexPosition().x();
00217   Vy = theTrack->GetVertexPosition().y();
00218   Vz = theTrack->GetVertexPosition().z();
00219   X  = hitPoint.x();
00220   Y  = hitPoint.y();
00221   Z  = hitPoint.z();
00222 }
00223 
00224 uint32_t FP420SD::setDetUnitId(G4Step * aStep) { 
00225 
00226   return (numberingScheme == 0 ? 0 : numberingScheme->getUnitID(aStep));
00227 }
00228 
00229 
00230 G4bool FP420SD::HitExists() {
00231   if (primaryID<1) {
00232     edm::LogWarning("FP420Sim") << "***** FP420SD error: primaryID = " 
00233                                   << primaryID
00234                                   << " maybe detector name changed";
00235   }
00236 
00237   // Update if in the same detector, time-slice and for same track   
00238   //  if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
00239   if (tSliceID == tsID && unitID==previousUnitID) {
00240   //AZ:
00241     UpdateHit();
00242     return true;
00243   }
00244   // Reset entry point for new primary
00245   if (primaryID != primID)
00246     ResetForNewPrimary();
00247    
00248   //look in the HitContainer whether a hit with the same primID, unitID,
00249   //tSliceID already exists:
00250    
00251   G4bool found = false;
00252 
00253   //    LogDebug("FP420Sim") << "FP420SD: HCollection=  " << theHC->entries()    <<std::endl;
00254   
00255   for (int j=0; j<theHC->entries()&&!found; j++) {
00256     FP420G4Hit* aPreviousHit = (*theHC)[j];
00257     if (aPreviousHit->getTrackID()     == primaryID &&
00258         aPreviousHit->getTimeSliceID() == tSliceID  &&
00259         aPreviousHit->getUnitID()      == unitID       ) {
00260   //AZ:
00261       currentHit = aPreviousHit;
00262       found      = true;
00263     }
00264   }          
00265 
00266   if (found) {
00267   //AZ:
00268     UpdateHit();
00269     return true;
00270   } else {
00271     return false;
00272   }    
00273 }
00274 
00275 
00276 void FP420SD::ResetForNewPrimary() {
00277   
00278   entrancePoint  = SetToLocal(hitPoint);
00279   exitPoint      = SetToLocalExit(hitPointExit);
00280   incidentEnergy = preStepPoint->GetKineticEnergy();
00281 
00282 }
00283 
00284 
00285 void FP420SD::StoreHit(FP420G4Hit* hit){
00286 
00287   //  if (primID<0) return;
00288   if (hit == 0 ) {
00289     edm::LogWarning("FP420Sim") << "FP420SD: hit to be stored is NULL !!";
00290     return;
00291   }
00292 
00293   theHC->insert( hit );
00294 }
00295 
00296 
00297 void FP420SD::CreateNewHit() {
00298 
00299 #ifdef debug
00300   //       << " MVid = " << currentPV->GetMother()->GetCopyNo()
00301   LogDebug("FP420Sim") << "FP420SD CreateNewHit for"
00302        << " PV "     << currentPV->GetName()
00303        << " PVid = " << currentPV->GetCopyNo()
00304        << " Unit "   << unitID <<std::endl;
00305   LogDebug("FP420Sim") << " primary "    << primaryID
00306        << " time slice " << tSliceID 
00307        << " For Track  " << theTrack->GetTrackID()
00308        << " which is a " << theTrack->GetDefinition()->GetParticleName();
00309            
00310   if (theTrack->GetTrackID()==1) {
00311     LogDebug("FP420Sim") << " of energy "     << theTrack->GetTotalEnergy();
00312   } else {
00313     LogDebug("FP420Sim") << " daughter of part. " << theTrack->GetParentID();
00314   }
00315 
00316   LogDebug("FP420Sim")  << " and created by " ;
00317   if (theTrack->GetCreatorProcess()!=NULL)
00318     LogDebug("FP420Sim") << theTrack->GetCreatorProcess()->GetProcessName() ;
00319   else 
00320     LogDebug("FP420Sim") << "NO process";
00321   LogDebug("FP420Sim") << std::endl;
00322 #endif          
00323     
00324 
00325   currentHit = new FP420G4Hit;
00326   currentHit->setTrackID(primaryID);
00327   currentHit->setTimeSlice(tSlice);
00328   currentHit->setUnitID(unitID);
00329   currentHit->setIncidentEnergy(incidentEnergy);
00330 
00331   currentHit->setPabs(Pabs);
00332   currentHit->setTof(Tof);
00333   currentHit->setEnergyLoss(Eloss);
00334   currentHit->setParticleType(ParticleType);
00335   currentHit->setThetaAtEntry(ThetaAtEntry);
00336   currentHit->setPhiAtEntry(PhiAtEntry);
00337 
00338 // currentHit->setEntry(entrancePoint);
00339   currentHit->setEntry(hitPoint);
00340 
00341   currentHit->setEntryLocalP(hitPointLocal);
00342   currentHit->setExitLocalP(hitPointLocalExit);
00343 
00344   currentHit->setParentId(ParentId);
00345   currentHit->setVx(Vx);
00346   currentHit->setVy(Vy);
00347   currentHit->setVz(Vz);
00348 
00349   currentHit->setX(X);
00350   currentHit->setY(Y);
00351   currentHit->setZ(Z);
00352   //AZ:12.10.2007
00353   //  UpdateHit();
00354   // buffer for next steps:
00355   tsID           = tSliceID;
00356   primID         = primaryID;
00357   previousUnitID = unitID;
00358   
00359   StoreHit(currentHit);
00360 }        
00361  
00362 
00363 void FP420SD::UpdateHit() {
00364   //
00365   if (Eloss > 0.) {
00366   currentHit->addEnergyDeposit(edepositEM,edepositHAD);
00367 
00368 #ifdef debug
00369     LogDebug("FP420Sim") << "updateHit: add eloss " << Eloss <<std::endl;
00370     LogDebug("FP420Sim") << "CurrentHit=" << currentHit
00371          << ", PostStepPoint=" << postStepPoint->GetPosition() << std::endl;
00372 #endif
00373     //AZ
00374     //   currentHit->setEnergyLoss(Eloss);
00375         currentHit->addEnergyLoss(Eloss);
00376   }  
00377 
00378 
00379 
00380   // buffer for next steps:
00381   tsID           = tSliceID;
00382   primID         = primaryID;
00383   previousUnitID = unitID;
00384 }
00385 
00386 
00387 G4ThreeVector FP420SD::SetToLocal(G4ThreeVector global){
00388 
00389   const G4VTouchable* touch= preStepPoint->GetTouchable();
00390   theEntryPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
00391   return theEntryPoint;  
00392 }
00393      
00394 
00395 G4ThreeVector FP420SD::SetToLocalExit(G4ThreeVector globalPoint){
00396 
00397   const G4VTouchable* touch= postStepPoint->GetTouchable();
00398   theExitPoint = touch->GetHistory()->GetTopTransform().TransformPoint(globalPoint);
00399   return theExitPoint;  
00400 }
00401      
00402 
00403 void FP420SD::EndOfEvent(G4HCofThisEvent* ) {
00404   
00405   // here we loop over transient hits and make them persistent
00406   
00407   //  if(theHC->entries() > 100){
00408   //    LogDebug("FP420Sim") << "FP420SD: warning!!! Number of hits exceed 100 and =" << theHC->entries() << "\n";
00409   //  }
00410   //  for (int j=0; j<theHC->entries() && j<100; j++) {
00411   int nhitsHPS240=0;
00412   int nhitsFP420=0;
00413   for (int j=0; j<theHC->entries(); j++) {
00414     FP420G4Hit* aHit = (*theHC)[j];
00415     if((fabs(aHit->getTof())> 780. && fabs(aHit->getTof())< 840.)) ++nhitsHPS240;
00416     if((fabs(aHit->getTof())>1380. && fabs(aHit->getTof())<1450.)) ++nhitsFP420;
00417     //    if(fabs(aHit->getTof()) < 1700.) {
00418     if((fabs(aHit->getTof())>780. && fabs(aHit->getTof())<840. && nhitsHPS240<200.) || ( fabs(aHit->getTof())>1380. && fabs(aHit->getTof())<1450. && nhitsFP420<200.)) {
00419       
00420       
00421 #ifdef ddebug
00422       //    LogDebug("FP420SD") << " FP420Hit " << j << " " << *aHit << std::endl;
00423       LogDebug("FP420Sim") << "hit number" << j << "unit ID = "<<aHit->getUnitID()<< "\n";
00424       LogDebug("FP420Sim") << "entry z " << aHit->getEntry().z()<< "\n";
00425       LogDebug("FP420Sim") << "entr theta " << aHit->getThetaAtEntry()<< "\n";
00426 #endif
00427       
00428       //    Local3DPoint locExitPoint(0,0,0);
00429       //  Local3DPoint locEntryPoint(aHit->getEntry().x(),
00430       //         aHit->getEntry().y(),
00431       //         aHit->getEntry().z());
00432       Local3DPoint locExitPoint(aHit->getExitLocalP().x(),
00433                                 aHit->getExitLocalP().y(),
00434                                 aHit->getExitLocalP().z());
00435       Local3DPoint locEntryPoint(aHit->getEntryLocalP().x(),
00436                                  aHit->getEntryLocalP().y(),
00437                                  aHit->getEntryLocalP().z());
00438       // implicit conversion (slicing) to PSimHit!!!
00439       slave->processHits(PSimHit(locEntryPoint,locExitPoint,//entryPoint(), exitPoint()  Local3DPoint
00440                                  aHit->getPabs(),// pabs()  float
00441                                  aHit->getTof(), // tof() float
00442                                  aHit->getEnergyLoss(), // energyLoss() float
00443                                  aHit->getParticleType(),// particleType()   int
00444                                  aHit->getUnitID(), // detUnitId() unsigned int 
00445                                  aHit->getTrackID(),// trackId() unsigned int 
00446                                  aHit->getThetaAtEntry(),//  thetaAtEntry()   float
00447                                  aHit->getPhiAtEntry())); //  phiAtEntry()   float  
00448       
00449       //PSimHit( const Local3DPoint& entry, const Local3DPoint& exit, 
00450       //           float pabs, float tof, float eloss, int particleType,
00451       //           unsigned int detId, unsigned int trackId,
00452       //           float theta, float phi, unsigned short processType=0) :
00453         
00454       //  LocalVector direction = hit.exitPoint() - hit.entryPoint(); 
00455       //hit.energyLoss
00456       
00457       
00458       /*      
00459               aHit->getEM(),               -
00460               aHit->getHadr(),             -
00461               aHit->getIncidentEnergy(),   -
00462               aHit->getTimeSlice(),       -
00463               aHit->getEntry(),     -
00464               aHit->getParentId(),     
00465               aHit->getEntryLocalP(),  -
00466               aHit->getExitLocalP(),   -
00467               aHit->getX(),    -
00468               aHit->getY(),   -
00469               aHit->getZ(),   -
00470               aHit->getVx(),  -
00471               aHit->getVy(),  -
00472               aHit->getVz()));  -
00473       */
00474     }//if Tof<1600. if nhits<100
00475   } // for loop on hits
00476   
00477   Summarize();
00478 }
00479 
00480 
00481 void FP420SD::Summarize() {
00482 }
00483 
00484 
00485 void FP420SD::clear() {
00486 } 
00487 
00488 
00489 void FP420SD::DrawAll() {
00490 } 
00491 
00492 
00493 void FP420SD::PrintAll() {
00494   LogDebug("FP420Sim") << "FP420SD: Collection " << theHC->GetName() << "\n";
00495   theHC->PrintAllHits();
00496 } 
00497 
00498 
00499 //void FP420SD::SetNumberingScheme(FP420NumberingScheme* scheme){
00500 //
00501 //  if (numberingScheme)
00502 //    delete numberingScheme;
00503 //  numberingScheme = scheme;
00504 //
00505 //}
00506 
00507 void FP420SD::fillHits(edm::PSimHitContainer& c, std::string n) {
00508   if (slave->name() == n) {
00509     c=slave->hits();
00510     std::cout << "FP420SD: fillHits to PSimHitContainer for name= " <<  name << std::endl;
00511   }
00512 }
00513 
00514 void FP420SD::update (const BeginOfEvent * i) {
00515   LogDebug("ForwardSim") << " Dispatched BeginOfEvent for " << GetName()
00516                        << " !" ;
00517    clearHits();
00518    eventno = (*i)()->GetEventID();
00519 }
00520 
00521 void FP420SD::update(const BeginOfRun *) {
00522 
00523   G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00524   G4String particleName;
00525   emPDG = theParticleTable->FindParticle(particleName="e-")->GetPDGEncoding();
00526   epPDG = theParticleTable->FindParticle(particleName="e+")->GetPDGEncoding();
00527   gammaPDG = theParticleTable->FindParticle(particleName="gamma")->GetPDGEncoding();
00528 
00529 }
00530 
00531 void FP420SD::update (const ::EndOfEvent*) {
00532 }
00533 
00534 void FP420SD::clearHits(){
00535   //AZ:
00536     slave->Initialize();
00537 }
00538 
00539 std::vector<std::string> FP420SD::getNames(){
00540   std::vector<std::string> temp;
00541   temp.push_back(slave->name());
00542   return temp;
00543 }
00544