CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/SimG4CMS/Calo/src/CaloSD.cc

Go to the documentation of this file.
00001 
00002 // File: CaloSD.cc
00003 // Description: Sensitive Detector class for calorimeters
00005 
00006 #include "SimG4CMS/Calo/interface/CaloSD.h"
00007 #include "SimDataFormats/SimHitMaker/interface/CaloSlaveSD.h"
00008 #include "SimG4Core/Notification/interface/TrackInformation.h"
00009 #include "SimG4Core/Application/interface/EventAction.h"
00010 
00011 #include "G4EventManager.hh"
00012 #include "G4SDManager.hh"
00013 #include "G4Step.hh"
00014 #include "G4Track.hh"
00015 #include "G4VProcess.hh"
00016 #include "G4GFlashSpot.hh"
00017 #include "G4ParticleTable.hh"
00018 
00019 //#define DebugLog
00020 
00021 CaloSD::CaloSD(G4String name, const DDCompactView & cpv,
00022         SensitiveDetectorCatalog & clg, 
00023         edm::ParameterSet const & p, const SimTrackManager* manager,
00024         int tSlice, bool ignoreTkID) : 
00025   SensitiveCaloDetector(name, cpv, clg, p),
00026   G4VGFlashSensitiveDetector(), theTrack(0), preStepPoint(0), eminHit(0), 
00027   eminHitD(0), m_trackManager(manager), currentHit(0), runInit(false),
00028   timeSlice(tSlice), ignoreTrackID(ignoreTkID), hcID(-1), theHC(0), 
00029   meanResponse(0) {
00030   //Add Hcal Sentitive Detector Names
00031 
00032   collectionName.insert(name);
00033 
00034   //Parameters
00035   edm::ParameterSet m_CaloSD = p.getParameter<edm::ParameterSet>("CaloSD");
00036   energyCut    = m_CaloSD.getParameter<double>("EminTrack")*GeV;
00037   tmaxHit      = m_CaloSD.getParameter<double>("TmaxHit")*ns;
00038   std::vector<double> eminHits = m_CaloSD.getParameter<std::vector<double> >("EminHits");
00039   std::vector<double> tmaxHits = m_CaloSD.getParameter<std::vector<double> >("TmaxHits");
00040   std::vector<std::string> hcn = m_CaloSD.getParameter<std::vector<std::string> >("HCNames");
00041   std::vector<int>   useResMap = m_CaloSD.getParameter<std::vector<int> >("UseResponseTables");
00042   std::vector<double> eminHitX = m_CaloSD.getParameter<std::vector<double> >("EminHitsDepth");
00043   suppressHeavy= m_CaloSD.getParameter<bool>("SuppressHeavy");
00044   kmaxIon      = m_CaloSD.getParameter<double>("IonThreshold")*MeV;
00045   kmaxProton   = m_CaloSD.getParameter<double>("ProtonThreshold")*MeV;
00046   kmaxNeutron  = m_CaloSD.getParameter<double>("NeutronThreshold")*MeV;
00047   checkHits    = m_CaloSD.getUntrackedParameter<int>("CheckHits", 25);
00048   useMap       = m_CaloSD.getUntrackedParameter<bool>("UseMap", true);
00049   int verbn    = m_CaloSD.getUntrackedParameter<int>("Verbosity", 0);
00050   corrTOFBeam  = m_CaloSD.getParameter<bool>("CorrectTOFBeam");
00051   double beamZ = m_CaloSD.getParameter<double>("BeamPosition")*cm;
00052   correctT     = beamZ/c_light/nanosecond;
00053 
00054   SetVerboseLevel(verbn);
00055   for (unsigned int k=0; k<hcn.size(); ++k) {
00056     if (name == (G4String)(hcn[k])) {
00057       if (k < eminHits.size()) eminHit = eminHits[k]*MeV;
00058       if (k < eminHitX.size()) eminHitD= eminHitX[k]*MeV;
00059       if (k < tmaxHits.size()) tmaxHit = tmaxHits[k]*ns;
00060       if (k < useResMap.size() && useResMap[k] > 0) meanResponse = new CaloMeanResponse(p);
00061       break;
00062     }
00063   }
00064 #ifdef DebugLog
00065   LogDebug("CaloSim") << "***************************************************" 
00066                       << "\n"
00067                       << "*                                                 *" 
00068                       << "\n"
00069                       << "* Constructing a CaloSD  with name " << GetName()
00070                       << "\n"
00071                       << "*                                                 *" 
00072                       << "\n"
00073                       << "***************************************************";
00074 #endif
00075   slave      = new CaloSlaveSD(name);
00076   currentID  = CaloHitID(timeSlice, ignoreTrackID);
00077   previousID = CaloHitID(timeSlice, ignoreTrackID);
00078   
00079   primAncestor = 0;
00080   cleanIndex = 0;
00081   totalHits = 0;
00082   forceSave = false;
00083 
00084   //
00085   // Now attach the right detectors (LogicalVolumes) to me
00086   //
00087   std::vector<std::string> lvNames = clg.logicalNames(name);
00088   this->Register();
00089   for (std::vector<std::string>::iterator it=lvNames.begin(); it !=lvNames.end(); ++it) {
00090     this->AssignSD(*it);
00091 #ifdef DebugLog
00092     LogDebug("CaloSim") << "CaloSD : Assigns SD to LV " << (*it);
00093 #endif
00094   }
00095 
00096   edm::LogInfo("CaloSim") << "CaloSD: Minimum energy of track for saving it " 
00097                           << energyCut/GeV  << " GeV" << "\n"
00098                           << "        Use of HitID Map " << useMap << "\n"
00099                           << "        Check last " << checkHits 
00100                           << " before saving the hit\n" 
00101                           << "        Correct TOF globally by " << correctT
00102                           << " ns (Flag =" << corrTOFBeam << ")\n"
00103                           << "        Save hits recorded before " << tmaxHit
00104                           << " ns and if energy is above " << eminHit/MeV
00105                           << " MeV (for depth 0) or " << eminHitD/MeV
00106                           << " MeV (for nonzero depths); Time Slice Unit " 
00107                           << timeSlice << " Ignore TrackID Flag " << ignoreTrackID;
00108 }
00109 
00110 CaloSD::~CaloSD() { 
00111   if (slave)           delete slave; 
00112   if (theHC)           delete theHC;
00113   if (meanResponse)    delete meanResponse;
00114 }
00115 
00116 bool CaloSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
00117   
00118   NaNTrap( aStep ) ;
00119   
00120   if (aStep == NULL) {
00121     return true;
00122   } else {
00123     if (getStepInfo(aStep)) {
00124       if (hitExists() == false && edepositEM+edepositHAD>0.) 
00125         currentHit = createNewHit();
00126     }
00127   }
00128   return true;
00129 } 
00130 
00131 bool CaloSD::ProcessHits(G4GFlashSpot* aSpot, G4TouchableHistory*) { 
00132 
00133   if (aSpot != NULL) {   
00134     theTrack = const_cast<G4Track *>(aSpot->GetOriginatorTrack()->GetPrimaryTrack());
00135     G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
00136     
00137     if (particleCode == emPDG ||
00138         particleCode == epPDG ||
00139         particleCode == gammaPDG ) {
00140       edepositEM  = aSpot->GetEnergySpot()->GetEnergy();
00141       edepositHAD = 0.;
00142     } else {
00143       edepositEM  = 0.;
00144       edepositHAD = 0.;
00145     }
00146  
00147     if (edepositEM>0.) {
00148       G4Step *      fFakeStep          = new G4Step();
00149       preStepPoint                     = fFakeStep->GetPreStepPoint();
00150       G4StepPoint * fFakePostStepPoint = fFakeStep->GetPostStepPoint();
00151       preStepPoint->SetPosition(aSpot->GetPosition());
00152       fFakePostStepPoint->SetPosition(aSpot->GetPosition());
00153       
00154       G4TouchableHandle fTouchableHandle   = aSpot->GetTouchableHandle();
00155       preStepPoint->SetTouchableHandle(fTouchableHandle);
00156       fFakeStep->SetTotalEnergyDeposit(aSpot->GetEnergySpot()->GetEnergy());
00157       
00158       double       time   = 0;
00159       unsigned int unitID = setDetUnitId(fFakeStep);
00160       int          primaryID = getTrackID(theTrack);
00161       uint16_t     depth = getDepth(fFakeStep);
00162 
00163       if (unitID > 0) {
00164         currentID.setID(unitID, time, primaryID, depth);
00165 #ifdef DebugLog
00166         LogDebug("CaloSim") << "CaloSD:: GetSpotInfo for"
00167                             << " Unit 0x" << std::hex << currentID.unitID() 
00168                             << std::dec << " Edeposit = " << edepositEM << " " 
00169                             << edepositHAD;
00170 #endif
00171         // Update if in the same detector, time-slice and for same track   
00172         if (currentID == previousID) {
00173           updateHit(currentHit);
00174         } else {
00175           posGlobal = aSpot->GetEnergySpot()->GetPosition();
00176           // Reset entry point for new primary
00177           if (currentID.trackID() != previousID.trackID()) {
00178             entrancePoint  = aSpot->GetPosition();
00179             entranceLocal  = aSpot->GetTouchableHandle()->GetHistory()->
00180                                       GetTopTransform().TransformPoint(entrancePoint);
00181             incidentEnergy = theTrack->GetKineticEnergy();
00182 #ifdef DebugLog
00183             LogDebug("CaloSim") << "CaloSD: Incident energy " 
00184                                 << incidentEnergy/GeV << " GeV and" 
00185                                 << " entrance point " << entrancePoint 
00186                                 << " (Global) " << entranceLocal << " (Local)";
00187 #endif
00188           }
00189 
00190           if (checkHit() == false) currentHit = createNewHit();
00191         }
00192       }
00193       delete  fFakeStep;
00194     }
00195     return true;
00196   } 
00197   return false;
00198 }                                   
00199 
00200 double CaloSD::getEnergyDeposit(G4Step* aStep) {
00201   return aStep->GetTotalEnergyDeposit();
00202 }
00203 
00204 void CaloSD::Initialize(G4HCofThisEvent * HCE) { 
00205   totalHits = 0;
00206   
00207 #ifdef DebugLog
00208   edm::LogInfo("CaloSim") << "CaloSD : Initialize called for " << GetName(); 
00209 #endif
00210   
00211   //This initialization is performed at the beginning of an event
00212   //------------------------------------------------------------
00213   theHC = new CaloG4HitCollection(GetName(), collectionName[0]);
00214   
00215   if (hcID<0) hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
00216   HCE->AddHitsCollection(hcID, theHC);
00217 }
00218 
00219 void CaloSD::EndOfEvent(G4HCofThisEvent* ) {
00220   // clean the hits for the last tracks
00221   
00222   cleanHitCollection();
00223   
00224 #ifdef DebugLog
00225   edm::LogInfo("CaloSim") << "CaloSD: EndofEvent entered with " << theHC->entries()
00226                           << " entries";
00227 #endif
00228   //  TimeMe("CaloSD:sortAndMergeHits",false);
00229 }
00230 
00231 void CaloSD::clear() {} 
00232 
00233 void CaloSD::DrawAll() {} 
00234 
00235 void CaloSD::PrintAll() {
00236 #ifdef DebugLog
00237   edm::LogInfo("CaloSim") << "CaloSD: Collection " << theHC->GetName();
00238 #endif
00239   theHC->PrintAllHits();
00240 } 
00241 
00242 void CaloSD::fillHits(edm::PCaloHitContainer& c, std::string n) {
00243   if (slave->name() == n) c=slave->hits();
00244   slave->Clean();
00245 }
00246 
00247 bool CaloSD::getStepInfo(G4Step* aStep) {  
00248 
00249   preStepPoint = aStep->GetPreStepPoint(); 
00250   theTrack     = aStep->GetTrack();   
00251   
00252   G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
00253   if (particleCode == emPDG ||
00254       particleCode == epPDG ||
00255       particleCode == gammaPDG ) {
00256     edepositEM  = getEnergyDeposit(aStep);
00257     edepositHAD = 0.;
00258   } else {
00259     edepositEM  = 0.;
00260     edepositHAD = getEnergyDeposit(aStep);
00261   }
00262   
00263   double       time  = (aStep->GetPostStepPoint()->GetGlobalTime())/nanosecond;
00264   unsigned int unitID= setDetUnitId(aStep);
00265   uint16_t     depth = getDepth(aStep);
00266   int          primaryID = getTrackID(theTrack);
00267   
00268   bool flag = (unitID > 0);
00269   if (flag) {
00270     currentID.setID(unitID, time, primaryID, depth);
00271 #ifdef DebugLog
00272     G4TouchableHistory* touch =(G4TouchableHistory*)(theTrack->GetTouchable());
00273     LogDebug("CaloSim") << "CaloSD:: GetStepInfo for"
00274                         << " PV "     << touch->GetVolume(0)->GetName()
00275                         << " PVid = " << touch->GetReplicaNumber(0)
00276                         << " MVid = " << touch->GetReplicaNumber(1)
00277                         << " Unit   " << currentID.unitID() 
00278                         << " Edeposit = " << edepositEM << " " << edepositHAD;
00279   } else {
00280     G4TouchableHistory* touch =(G4TouchableHistory*)(theTrack->GetTouchable());
00281     LogDebug("CaloSim") << "CaloSD:: GetStepInfo for"
00282                         << " PV "     << touch->GetVolume(0)->GetName()
00283                         << " PVid = " << touch->GetReplicaNumber(0)
00284                         << " MVid = " << touch->GetReplicaNumber(1)
00285                         << " Unit   " << std::hex << unitID << std::dec 
00286                         << " Edeposit = " << edepositEM << " " << edepositHAD;
00287 #endif
00288   }
00289   return flag;
00290 }
00291 
00292 G4ThreeVector CaloSD::setToLocal(G4ThreeVector global, const G4VTouchable* touch) {
00293 
00294   G4ThreeVector localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
00295   
00296   return localPoint;  
00297 }
00298 
00299 G4bool CaloSD::hitExists() {
00300 #ifdef DebugLog
00301   if (currentID.trackID()<1)
00302     edm::LogWarning("CaloSim") << "***** CaloSD error: primaryID = " 
00303                                << currentID.trackID()
00304                                << " maybe detector name changed";
00305 #endif  
00306   // Update if in the same detector, time-slice and for same track   
00307   if (currentID == previousID) {
00308     updateHit(currentHit);
00309     return true;
00310   }
00311   
00312   // Reset entry point for new primary
00313   posGlobal = preStepPoint->GetPosition();
00314   if (currentID.trackID() != previousID.trackID()) 
00315     resetForNewPrimary(preStepPoint->GetPosition(), preStepPoint->GetKineticEnergy());
00316   
00317   return checkHit();
00318 }
00319 
00320 G4bool CaloSD::checkHit() {  
00321   //look in the HitContainer whether a hit with the same ID already exists:
00322   bool       found = false;
00323   if (useMap) {
00324     std::map<CaloHitID,CaloG4Hit*>::const_iterator it = hitMap.find(currentID);
00325     if (it != hitMap.end()) {
00326       currentHit = it->second;
00327       found      = true;
00328     }
00329   } else {
00330     if (checkHits <= 0) return false;
00331     int  minhit= (theHC->entries()>checkHits ? theHC->entries()-checkHits : 0);
00332     int  maxhit= theHC->entries()-1;
00333     
00334     for (int j=maxhit; j>minhit&&!found; --j) {
00335       if ((*theHC)[j]->getID() == currentID) {
00336         currentHit = (*theHC)[j];
00337         found      = true;
00338       }
00339     }          
00340   }
00341   
00342   if (found) {
00343     updateHit(currentHit);
00344     return true;
00345   } else {
00346     return false;
00347   }
00348 }
00349 
00350 int CaloSD::getNumberOfHits() { return theHC->entries(); }
00351 
00352 CaloG4Hit* CaloSD::createNewHit() {
00353 #ifdef DebugLog
00354   LogDebug("CaloSim") << "CaloSD::CreateNewHit for"
00355                       << " Unit " << currentID.unitID() 
00356                       << " " << currentID.depth()
00357                       << " Edeposit = " << edepositEM << " " << edepositHAD;
00358   LogDebug("CaloSim") << " primary "    << currentID.trackID()
00359                       << " time slice " << currentID.timeSliceID()
00360                       << " For Track  " << theTrack->GetTrackID()
00361                       << " which is a " <<theTrack->GetDefinition()->GetParticleName()
00362                       << " of energy "  << theTrack->GetKineticEnergy()/GeV
00363                       << " " << theTrack->GetMomentum().mag()/GeV
00364                       << " daughter of part. " << theTrack->GetParentID()
00365                       << " and created by " ;
00366   
00367   if (theTrack->GetCreatorProcess()!=NULL)
00368     LogDebug("CaloSim") << theTrack->GetCreatorProcess()->GetProcessName() ;
00369   else 
00370     LogDebug("CaloSim") << "NO process";
00371 #endif  
00372   
00373   CaloG4Hit* aHit;
00374   if (reusehit.size() > 0) {
00375     aHit = reusehit[0];
00376     aHit->setEM(0.);
00377     aHit->setHadr(0.);
00378     reusehit.erase(reusehit.begin());
00379   } else {
00380     aHit = new CaloG4Hit;
00381   }
00382   
00383   aHit->setID(currentID);
00384   aHit->setEntry(entrancePoint.x(),entrancePoint.y(),entrancePoint.z());
00385   aHit->setEntryLocal(entranceLocal.x(),entranceLocal.y(),entranceLocal.z());
00386   aHit->setPosition(posGlobal.x(),posGlobal.y(),posGlobal.z());
00387   aHit->setIncidentEnergy(incidentEnergy);
00388   updateHit(aHit);
00389   
00390   storeHit(aHit);
00391   double etrack = 0;
00392   if (currentID.trackID() == primIDSaved) { // The track is saved; nothing to be done
00393   } else if (currentID.trackID() == theTrack->GetTrackID()) {
00394     etrack= theTrack->GetKineticEnergy();
00395     //edm::LogInfo("CaloSim") << "CaloSD: set save the track " << currentID.trackID()
00396     //      << " etrack " << etrack << " eCut " << energyCut << " flag " << forceSave;
00397     if (etrack >= energyCut || forceSave) {
00398       TrackInformation* trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
00399       trkInfo->storeTrack(true);
00400       trkInfo->putInHistory();
00401       //      trkInfo->setAncestor();
00402 #ifdef DebugLog
00403       LogDebug("CaloSim") << "CaloSD: set save the track " << currentID.trackID()
00404                           << " with Hit";
00405 #endif
00406     }
00407   } else {
00408     TrackWithHistory * trkh = tkMap[currentID.trackID()];
00409 #ifdef DebugLog
00410     LogDebug("CaloSim") << "CaloSD : TrackwithHistory pointer for " 
00411                         << currentID.trackID() << " is " << trkh;
00412 #endif
00413     if (trkh != NULL) {
00414       etrack = sqrt(trkh->momentum().Mag2());
00415       if (etrack >= energyCut) {
00416         trkh->save();
00417 #ifdef DebugLog
00418         LogDebug("CaloSim") << "CaloSD: set save the track " 
00419                             << currentID.trackID() << " with Hit";
00420 #endif
00421       }
00422     }
00423   }
00424   primIDSaved = currentID.trackID();
00425   if (useMap) totalHits++;
00426   return aHit;
00427 }  
00428 
00429 void CaloSD::updateHit(CaloG4Hit* aHit) {
00430   if (edepositEM+edepositHAD != 0) {
00431     aHit->addEnergyDeposit(edepositEM,edepositHAD);
00432 #ifdef DebugLog
00433     LogDebug("CaloSim") << "CaloSD: Add energy deposit in " << currentID 
00434                         << " em " << edepositEM/MeV << " hadronic " 
00435                         << edepositHAD/MeV << " MeV"; 
00436 #endif
00437   }
00438 
00439   // buffer for next steps:
00440   previousID = currentID;
00441 }
00442 
00443 void CaloSD::resetForNewPrimary(G4ThreeVector point, double energy) { 
00444   entrancePoint  = point;
00445   entranceLocal  = setToLocal(entrancePoint, preStepPoint->GetTouchable());
00446   incidentEnergy = energy;
00447 #ifdef DebugLog
00448   LogDebug("CaloSim") << "CaloSD: Incident energy " << incidentEnergy/GeV 
00449                       << " GeV and" << " entrance point " << entrancePoint 
00450                       << " (Global) " << entranceLocal << " (Local)";
00451 #endif
00452 }
00453 
00454 double CaloSD::getAttenuation(G4Step* aStep, double birk1, double birk2, double birk3) {
00455   double weight = 1.;
00456   double charge = aStep->GetPreStepPoint()->GetCharge();
00457 
00458   if (charge != 0. && aStep->GetStepLength() > 0) {
00459     G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
00460     double density = mat->GetDensity();
00461     double dedx    = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
00462     double rkb     = birk1/density;
00463     double c       = birk2*rkb*rkb;
00464     if (std::abs(charge) >= 2.) rkb /= birk3; // based on alpha particle data
00465     weight = 1./(1.+rkb*dedx+c*dedx*dedx);
00466 #ifdef DebugLog
00467     LogDebug("CaloSim") << "CaloSD::getAttenuation in " << mat->GetName() 
00468                         << " Charge " << charge << " dE/dx " << dedx 
00469                         << " Birk Const " << rkb << ", " << c << " Weight = " 
00470                         << weight << " dE " << aStep->GetTotalEnergyDeposit();
00471 #endif
00472   }
00473   return weight;
00474 }
00475 
00476 void CaloSD::update(const BeginOfRun *) {
00477   G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00478   G4String particleName;
00479   emPDG = theParticleTable->FindParticle(particleName="e-")->GetPDGEncoding();
00480   epPDG = theParticleTable->FindParticle(particleName="e+")->GetPDGEncoding();
00481   gammaPDG = theParticleTable->FindParticle(particleName="gamma")->GetPDGEncoding();
00482 #ifdef DebugLog
00483   LogDebug("CaloSim") << "CaloSD: Particle code for e- = " << emPDG
00484                       << " for e+ = " << epPDG << " for gamma = " << gammaPDG;
00485 #endif
00486   initRun();
00487   runInit = true;
00488 } 
00489 
00490 void CaloSD::update(const BeginOfEvent *) {
00491 #ifdef DebugLog
00492   LogDebug("CaloSim")  << "CaloSD: Dispatched BeginOfEvent for " << GetName() 
00493                        << " !" ;
00494 #endif
00495   clearHits();
00496 }
00497 
00498 void CaloSD::update(const EndOfTrack * trk) {
00499   int id = (*trk)()->GetTrackID();
00500   TrackInformation *trkI =(TrackInformation *)((*trk)()->GetUserInformation());
00501   int lastTrackID = -1;
00502   if (trkI) lastTrackID = trkI->getIDonCaloSurface();
00503   if (id == lastTrackID) {
00504     const TrackContainer * trksForThisEvent = m_trackManager->trackContainer();
00505     if (trksForThisEvent != NULL) {
00506       int it = (int)(trksForThisEvent->size()) - 1;
00507       if (it >= 0) {
00508         TrackWithHistory * trkH = (*trksForThisEvent)[it];
00509         if (trkH->trackID() == (unsigned int)(id)) tkMap[id] = trkH;
00510 #ifdef DebugLog
00511         LogDebug("CaloSim") << "CaloSD: get track " << it << " from "
00512                             << "Container of size " << trksForThisEvent->size()
00513                             << " with ID " << trkH->trackID();
00514       } else {
00515         LogDebug("CaloSim") << "CaloSD: get track " << it << " from "
00516                             << "Container of size " << trksForThisEvent->size()
00517                             << " with no ID";
00518 #endif
00519       }
00520     }
00521   }
00522 }
00523 
00524 void CaloSD::update(const ::EndOfEvent * ) {
00525   int count = 0, wrong = 0;
00526   bool ok;
00527   
00528   slave->ReserveMemory(theHC->entries());
00529 
00530   for (int i=0; i<theHC->entries(); ++i) {
00531     ok = saveHit((*theHC)[i]);
00532     ++count;
00533     if (!ok)  ++wrong;
00534   }
00535   
00536   edm::LogInfo("CaloSim") << "CaloSD: " << GetName() << " store " << count
00537                           << " hits recorded with " << wrong 
00538                           << " track IDs not given properly and "
00539                           << totalHits-count << " hits not passing cuts";
00540   summarize();
00541 
00542   tkMap.erase (tkMap.begin(), tkMap.end());
00543 }
00544 
00545 void CaloSD::clearHits() {  
00546   if (useMap) hitMap.erase (hitMap.begin(), hitMap.end());
00547   for (unsigned int i = 0; i<reusehit.size(); ++i) delete reusehit[i];
00548   std::vector<CaloG4Hit*>().swap(reusehit);
00549   cleanIndex  = 0;
00550   previousID.reset();
00551   primIDSaved = -99;
00552 #ifdef DebugLog
00553   LogDebug("CaloSim") << "CaloSD: Clears hit vector for " << GetName() << " " << slave;
00554 #endif
00555   slave->Initialize();
00556 #ifdef DebugLog
00557   LogDebug("CaloSim") << "CaloSD: Initialises slave SD for " << GetName();
00558 #endif
00559 }
00560 
00561 void CaloSD::initRun() {}
00562 
00563 int CaloSD::getTrackID(G4Track* aTrack) {
00564   int primaryID = 0;
00565   forceSave = false;
00566   TrackInformation* trkInfo=(TrackInformation *)(aTrack->GetUserInformation());
00567   if (trkInfo) {
00568     primaryID = trkInfo->getIDonCaloSurface(); 
00569 #ifdef DebugLog
00570     LogDebug("CaloSim") << "CaloSD: hit update from track Id on Calo Surface " 
00571                         << trkInfo->getIDonCaloSurface();
00572 #endif   
00573   } else {
00574     primaryID = aTrack->GetTrackID();
00575 #ifdef DebugLog
00576     edm::LogWarning("CaloSim") << "CaloSD: Problem with primaryID **** set by "
00577                                << "force to TkID **** " << primaryID << " in "
00578                                << preStepPoint->GetTouchable()->GetVolume(0)->GetName();
00579 #endif
00580   }
00581   return primaryID;
00582 }
00583 
00584 uint16_t CaloSD::getDepth(G4Step*) { return 0; }
00585 
00586 bool CaloSD::filterHit(CaloG4Hit* hit, double time) {
00587   double emin(eminHit);
00588   if (hit->getDepth() > 0) emin = eminHitD;
00589 #ifdef DebugLog
00590   LogDebug("CaloSim") << "Depth " << hit->getDepth() << " Emin = " << emin << " ("
00591                       << eminHit << ", " << eminHitD << ")";
00592 #endif   
00593   return ((time <= tmaxHit) && (hit->getEnergyDeposit() > emin));
00594 }
00595 
00596 double CaloSD::getResponseWt(G4Track* aTrack) {
00597   if (meanResponse) {
00598     TrackInformation * trkInfo = (TrackInformation *)(aTrack->GetUserInformation());
00599     return meanResponse->getWeight(trkInfo->genParticlePID(), trkInfo->genParticleP());
00600   } else {
00601     return 1;
00602   }
00603 }
00604 
00605 void CaloSD::storeHit(CaloG4Hit* hit) {
00606   if (previousID.trackID()<0) return;
00607   if (hit == 0) {
00608     edm::LogWarning("CaloSim") << "CaloSD: hit to be stored is NULL !!";
00609     return;
00610   }
00611   
00612   theHC->insert(hit);
00613   if (useMap) hitMap.insert(std::pair<CaloHitID,CaloG4Hit*>(previousID,hit));
00614 }
00615 
00616 bool CaloSD::saveHit(CaloG4Hit* aHit) {  
00617   int tkID;
00618   bool ok   = true;
00619   if (m_trackManager) {
00620     tkID = m_trackManager->giveMotherNeeded(aHit->getTrackID());
00621     if (tkID == 0) {
00622       if (m_trackManager->trackExists(aHit->getTrackID())) tkID = (aHit->getTrackID());
00623       else {
00624         ok = false;
00625       }
00626     }
00627   } else {
00628     tkID = aHit->getTrackID();
00629     ok = false;
00630   }
00631   //  edm::LogInfo("CaloSim") << "CalosD: Track ID " << aHit->getTrackID() << " changed to " << tkID << " by SimTrackManager" << " Status " << ok;
00632 #ifdef DebugLog
00633   LogDebug("CaloSim") << "CalosD: Track ID " << aHit->getTrackID() 
00634                       << " changed to " << tkID << " by SimTrackManager"
00635                       << " Status " << ok;
00636 #endif
00637   double time = aHit->getTimeSlice();
00638   if (corrTOFBeam) time += correctT;
00639   slave->processHits(aHit->getUnitID(), aHit->getEM()/GeV, 
00640                      aHit->getHadr()/GeV, time, tkID, aHit->getDepth());
00641 #ifdef DebugLog
00642   LogDebug("CaloSim") << "CaloSD: Store Hit at " << std::hex 
00643                       << aHit->getUnitID() << std::dec << " " 
00644                       << aHit->getDepth() << " due to " << tkID 
00645                       << " in time " << time << " of energy " 
00646                       << aHit->getEM()/GeV << " GeV (EM) and " 
00647                       << aHit->getHadr()/GeV << " GeV (Hadr)";
00648 #endif
00649   return ok;
00650 }
00651 
00652 void CaloSD::summarize() {}
00653 
00654 void CaloSD::update(const BeginOfTrack * trk) {
00655   int primary = -1;
00656   TrackInformation * trkInfo = (TrackInformation *)((*trk)()->GetUserInformation());
00657   if ( trkInfo->isPrimary() ) primary = (*trk)()->GetTrackID();
00658   
00659 #ifdef DebugLog
00660   LogDebug("CaloSim") << "New track: isPrimary " << trkInfo->isPrimary() 
00661                       << " primary ID = " << primary 
00662                       << " primary ancestor ID " << primAncestor;
00663 #endif
00664   
00665   // update the information if a different primary track ID 
00666   
00667   if (primary > 0 && primary != primAncestor) {
00668     primAncestor = primary;
00669     
00670     // clean the hits information
00671     
00672     if (theHC->entries()>0) cleanHitCollection();
00673     
00674   }
00675 }
00676 
00677 void CaloSD::cleanHitCollection() {
00678   std::vector<CaloG4Hit*>* theCollection = theHC->GetVector();
00679 
00680 #ifdef DebugLog
00681   LogDebug("CaloSim") << "CaloSD: collection before merging, size = " << theHC->entries();
00682 #endif
00683   
00684   selIndex.reserve(theHC->entries()-cleanIndex);
00685   if ( reusehit.size() == 0 ) reusehit.reserve(theHC->entries()-cleanIndex); 
00686 
00687   // if no map used, merge before hits to have the save situation as a map
00688   if ( !useMap ) {
00689     hitvec.swap(*theCollection);
00690     sort((hitvec.begin()+cleanIndex), hitvec.end(), CaloG4HitLess());
00691 #ifdef DebugLog
00692     LogDebug("CaloSim") << "CaloSD::cleanHitCollection: sort hits in buffer "
00693                         << "starting from element = " << cleanIndex;
00694     for (unsigned int i = 0; i<hitvec.size(); ++i) 
00695       LogDebug("CaloSim")<<i<<" "<<*hitvec[i];
00696 #endif
00697     unsigned int i, j;
00698     CaloG4HitEqual equal;
00699     for (i=cleanIndex; i<hitvec.size(); ++i) {
00700       selIndex.push_back(i-cleanIndex);
00701       int jump = 0;
00702       for (j = i+1; j <hitvec.size() && equal(hitvec[i], hitvec[j]); ++j) {
00703         ++jump;
00704         // merge j to i
00705         (*hitvec[i]).addEnergyDeposit(*hitvec[j]);
00706         (*hitvec[j]).setEM(0.);
00707         (*hitvec[j]).setHadr(0.);
00708         reusehit.push_back(hitvec[j]);
00709       }
00710       i+=jump;
00711     }
00712 #ifdef DebugLog
00713     LogDebug("CaloSim") << "CaloSD: cleanHitCollection merge the hits in buffer ";
00714     for (unsigned int i = 0; i<hitvec.size(); ++i) 
00715       LogDebug("CaloSim")<<i<<" "<<*hitvec[i];
00716 #endif
00717     for ( unsigned int i = cleanIndex; i < cleanIndex+selIndex.size(); ++i ) {
00718       hitvec[i] = hitvec[selIndex[i-cleanIndex]+cleanIndex];
00719     }
00720     hitvec.resize(cleanIndex+selIndex.size());
00721 #ifdef DebugLog
00722     LogDebug("CaloSim") << "CaloSD::cleanHitCollection: remove the merged hits in buffer,"
00723                         << " new size = " << hitvec.size();
00724     for (unsigned int i = 0; i<hitvec.size(); ++i) 
00725       LogDebug("CaloSim")<<i<<" "<<*hitvec[i];
00726 #endif
00727     hitvec.swap(*theCollection);
00728     std::vector<CaloG4Hit*>().swap(hitvec);
00729     selIndex.clear();
00730     totalHits = theHC->entries();
00731   }
00732 
00733 #ifdef DebugLog
00734   LogDebug("CaloSim") << "CaloSD: collection after merging, size = " << theHC->entries();
00735 #endif
00736 
00737   int addhit = 0;
00738 
00739 #ifdef DebugLog
00740   LogDebug("CaloSim") << "CaloSD: Size of reusehit after merge = " << reusehit.size();
00741   LogDebug("CaloSim") << "CaloSD: Starting hit selection from index = " << cleanIndex;
00742 #endif
00743   
00744   selIndex.reserve(theCollection->size()-cleanIndex);
00745   for (unsigned int i = cleanIndex; i<theCollection->size(); ++i) {   
00746     CaloG4Hit* aHit((*theCollection)[i]);
00747     
00748     // selection
00749     
00750     double time = aHit->getTimeSlice();
00751     if (corrTOFBeam) time += correctT;
00752     if (!filterHit(aHit,time)) {
00753 #ifdef DebugLog
00754       LogDebug("CaloSim") << "CaloSD: dropped CaloG4Hit " << " " << *aHit; 
00755 #endif
00756       
00757       // create the list of hits to be reused
00758       
00759       reusehit.push_back((*theCollection)[i]);
00760       ++addhit;
00761     } else {
00762       selIndex.push_back(i-cleanIndex);
00763     }
00764   }
00765 
00766 #ifdef DebugLog
00767   LogDebug("CaloSim") << "CaloSD: Size of reusehit after selection = " << reusehit.size()  
00768                       << " Number of added hit = " << addhit;
00769 #endif
00770   if (useMap) {
00771     if ( addhit>0 ) {
00772       int offset = reusehit.size()-addhit;
00773       for (int ii = addhit-1; ii>=0; --ii) {
00774         CaloHitID theID = reusehit[offset+ii]->getID();
00775         hitMap.erase(theID);
00776       }
00777     }
00778   }
00779   for (unsigned int j = 0; j<selIndex.size(); ++j) {
00780     (*theCollection)[cleanIndex+j] = (*theCollection)[cleanIndex+selIndex[j]];
00781   }
00782 
00783   theCollection->resize(cleanIndex+selIndex.size());
00784   std::vector<unsigned int>().swap(selIndex);
00785 
00786 #ifdef DebugLog
00787   LogDebug("CaloSim") << "CaloSD: hit collection after selection, size = "
00788                       << theHC->entries();
00789   theHC->PrintAllHits();
00790 #endif
00791     
00792   cleanIndex = theHC->entries();
00793 }