CMS 3D CMS Logo

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