00001
00002
00003
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
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
00032 collectionName.insert(name);
00033
00034
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
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
00172 if (currentID == previousID) {
00173 updateHit(currentHit);
00174 } else {
00175 posGlobal = aSpot->GetEnergySpot()->GetPosition();
00176
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
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
00221
00222 cleanHitCollection();
00223
00224 #ifdef DebugLog
00225 edm::LogInfo("CaloSim") << "CaloSD: EndofEvent entered with " << theHC->entries()
00226 << " entries";
00227 #endif
00228
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 double time = (aStep->GetPostStepPoint()->GetGlobalTime())/nanosecond;
00253 unsigned int unitID= setDetUnitId(aStep);
00254 uint16_t depth = getDepth(aStep);
00255 int primaryID = getTrackID(theTrack);
00256
00257 bool flag = (unitID > 0);
00258 if (flag) {
00259 currentID.setID(unitID, time, primaryID, depth);
00260 #ifdef DebugLog
00261 G4TouchableHistory* touch =(G4TouchableHistory*)(theTrack->GetTouchable());
00262 LogDebug("CaloSim") << "CaloSD:: GetStepInfo for"
00263 << " PV " << touch->GetVolume(0)->GetName()
00264 << " PVid = " << touch->GetReplicaNumber(0)
00265 << " MVid = " << touch->GetReplicaNumber(1)
00266 << " Unit " << currentID.unitID()
00267 << " Edeposit = " << edepositEM << " " << edepositHAD;
00268 } else {
00269 G4TouchableHistory* touch =(G4TouchableHistory*)(theTrack->GetTouchable());
00270 LogDebug("CaloSim") << "CaloSD:: GetStepInfo for"
00271 << " PV " << touch->GetVolume(0)->GetName()
00272 << " PVid = " << touch->GetReplicaNumber(0)
00273 << " MVid = " << touch->GetReplicaNumber(1)
00274 << " Unit " << std::hex << unitID << std::dec
00275 << " Edeposit = " << edepositEM << " " << edepositHAD;
00276 #endif
00277 }
00278
00279 G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
00280 if (particleCode == emPDG ||
00281 particleCode == epPDG ||
00282 particleCode == gammaPDG ) {
00283 edepositEM = getEnergyDeposit(aStep);
00284 edepositHAD = 0.;
00285 } else {
00286 edepositEM = 0.;
00287 edepositHAD = getEnergyDeposit(aStep);
00288 }
00289
00290 return flag;
00291 }
00292
00293 G4ThreeVector CaloSD::setToLocal(G4ThreeVector global, const G4VTouchable* touch) {
00294
00295 G4ThreeVector localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
00296
00297 return localPoint;
00298 }
00299
00300 G4ThreeVector CaloSD::setToGlobal(G4ThreeVector local, const G4VTouchable* touch) {
00301
00302 G4ThreeVector globalPoint = touch->GetHistory()->GetTopTransform().Inverse().TransformPoint(local);
00303
00304 return globalPoint;
00305 }
00306
00307 G4bool CaloSD::hitExists() {
00308 #ifdef DebugLog
00309 if (currentID.trackID()<1)
00310 edm::LogWarning("CaloSim") << "***** CaloSD error: primaryID = "
00311 << currentID.trackID()
00312 << " maybe detector name changed";
00313 #endif
00314
00315 if (currentID == previousID) {
00316 updateHit(currentHit);
00317 return true;
00318 }
00319
00320
00321 posGlobal = preStepPoint->GetPosition();
00322 if (currentID.trackID() != previousID.trackID())
00323 resetForNewPrimary(preStepPoint->GetPosition(), preStepPoint->GetKineticEnergy());
00324
00325 return checkHit();
00326 }
00327
00328 G4bool CaloSD::checkHit() {
00329
00330 bool found = false;
00331 if (useMap) {
00332 std::map<CaloHitID,CaloG4Hit*>::const_iterator it = hitMap.find(currentID);
00333 if (it != hitMap.end()) {
00334 currentHit = it->second;
00335 found = true;
00336 }
00337 } else {
00338 if (checkHits <= 0) return false;
00339 int minhit= (theHC->entries()>checkHits ? theHC->entries()-checkHits : 0);
00340 int maxhit= theHC->entries()-1;
00341
00342 for (int j=maxhit; j>minhit&&!found; --j) {
00343 if ((*theHC)[j]->getID() == currentID) {
00344 currentHit = (*theHC)[j];
00345 found = true;
00346 }
00347 }
00348 }
00349
00350 if (found) {
00351 updateHit(currentHit);
00352 return true;
00353 } else {
00354 return false;
00355 }
00356 }
00357
00358 int CaloSD::getNumberOfHits() { return theHC->entries(); }
00359
00360 CaloG4Hit* CaloSD::createNewHit() {
00361 #ifdef DebugLog
00362 LogDebug("CaloSim") << "CaloSD::CreateNewHit for"
00363 << " Unit " << currentID.unitID()
00364 << " " << currentID.depth()
00365 << " Edeposit = " << edepositEM << " " << edepositHAD;
00366 LogDebug("CaloSim") << " primary " << currentID.trackID()
00367 << " time slice " << currentID.timeSliceID()
00368 << " For Track " << theTrack->GetTrackID()
00369 << " which is a " <<theTrack->GetDefinition()->GetParticleName()
00370 << " of energy " << theTrack->GetKineticEnergy()/GeV
00371 << " " << theTrack->GetMomentum().mag()/GeV
00372 << " daughter of part. " << theTrack->GetParentID()
00373 << " and created by " ;
00374
00375 if (theTrack->GetCreatorProcess()!=NULL)
00376 LogDebug("CaloSim") << theTrack->GetCreatorProcess()->GetProcessName() ;
00377 else
00378 LogDebug("CaloSim") << "NO process";
00379 #endif
00380
00381 CaloG4Hit* aHit;
00382 if (reusehit.size() > 0) {
00383 aHit = reusehit[0];
00384 aHit->setEM(0.);
00385 aHit->setHadr(0.);
00386 reusehit.erase(reusehit.begin());
00387 } else {
00388 aHit = new CaloG4Hit;
00389 }
00390
00391 aHit->setID(currentID);
00392 aHit->setEntry(entrancePoint.x(),entrancePoint.y(),entrancePoint.z());
00393 aHit->setEntryLocal(entranceLocal.x(),entranceLocal.y(),entranceLocal.z());
00394 aHit->setPosition(posGlobal.x(),posGlobal.y(),posGlobal.z());
00395 aHit->setIncidentEnergy(incidentEnergy);
00396 updateHit(aHit);
00397
00398 storeHit(aHit);
00399 double etrack = 0;
00400 if (currentID.trackID() == primIDSaved) {
00401 } else if (currentID.trackID() == theTrack->GetTrackID()) {
00402 etrack= theTrack->GetKineticEnergy();
00403
00404
00405 if (etrack >= energyCut || forceSave) {
00406 TrackInformation* trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
00407 trkInfo->storeTrack(true);
00408 trkInfo->putInHistory();
00409
00410 #ifdef DebugLog
00411 LogDebug("CaloSim") << "CaloSD: set save the track " << currentID.trackID()
00412 << " with Hit";
00413 #endif
00414 }
00415 } else {
00416 TrackWithHistory * trkh = tkMap[currentID.trackID()];
00417 #ifdef DebugLog
00418 LogDebug("CaloSim") << "CaloSD : TrackwithHistory pointer for "
00419 << currentID.trackID() << " is " << trkh;
00420 #endif
00421 if (trkh != NULL) {
00422 etrack = sqrt(trkh->momentum().Mag2());
00423 if (etrack >= energyCut) {
00424 trkh->save();
00425 #ifdef DebugLog
00426 LogDebug("CaloSim") << "CaloSD: set save the track "
00427 << currentID.trackID() << " with Hit";
00428 #endif
00429 }
00430 }
00431 }
00432 primIDSaved = currentID.trackID();
00433 if (useMap) totalHits++;
00434 return aHit;
00435 }
00436
00437 void CaloSD::updateHit(CaloG4Hit* aHit) {
00438 if (edepositEM+edepositHAD != 0) {
00439 aHit->addEnergyDeposit(edepositEM,edepositHAD);
00440 #ifdef DebugLog
00441 LogDebug("CaloSim") << "CaloSD: Add energy deposit in " << currentID
00442 << " em " << edepositEM/MeV << " hadronic "
00443 << edepositHAD/MeV << " MeV";
00444 #endif
00445 }
00446
00447
00448 previousID = currentID;
00449 }
00450
00451 void CaloSD::resetForNewPrimary(G4ThreeVector point, double energy) {
00452 entrancePoint = point;
00453 entranceLocal = setToLocal(entrancePoint, preStepPoint->GetTouchable());
00454 incidentEnergy = energy;
00455 #ifdef DebugLog
00456 LogDebug("CaloSim") << "CaloSD: Incident energy " << incidentEnergy/GeV
00457 << " GeV and" << " entrance point " << entrancePoint
00458 << " (Global) " << entranceLocal << " (Local)";
00459 #endif
00460 }
00461
00462 double CaloSD::getAttenuation(G4Step* aStep, double birk1, double birk2, double birk3) {
00463 double weight = 1.;
00464 double charge = aStep->GetPreStepPoint()->GetCharge();
00465
00466 if (charge != 0. && aStep->GetStepLength() > 0) {
00467 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
00468 double density = mat->GetDensity();
00469 double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
00470 double rkb = birk1/density;
00471 double c = birk2*rkb*rkb;
00472 if (std::abs(charge) >= 2.) rkb /= birk3;
00473 weight = 1./(1.+rkb*dedx+c*dedx*dedx);
00474 #ifdef DebugLog
00475 LogDebug("CaloSim") << "CaloSD::getAttenuation in " << mat->GetName()
00476 << " Charge " << charge << " dE/dx " << dedx
00477 << " Birk Const " << rkb << ", " << c << " Weight = "
00478 << weight << " dE " << aStep->GetTotalEnergyDeposit();
00479 #endif
00480 }
00481 return weight;
00482 }
00483
00484 void CaloSD::update(const BeginOfRun *) {
00485 G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00486 G4String particleName;
00487 emPDG = theParticleTable->FindParticle(particleName="e-")->GetPDGEncoding();
00488 epPDG = theParticleTable->FindParticle(particleName="e+")->GetPDGEncoding();
00489 gammaPDG = theParticleTable->FindParticle(particleName="gamma")->GetPDGEncoding();
00490 #ifdef DebugLog
00491 LogDebug("CaloSim") << "CaloSD: Particle code for e- = " << emPDG
00492 << " for e+ = " << epPDG << " for gamma = " << gammaPDG;
00493 #endif
00494 initRun();
00495 runInit = true;
00496 }
00497
00498 void CaloSD::update(const BeginOfEvent *) {
00499 #ifdef DebugLog
00500 LogDebug("CaloSim") << "CaloSD: Dispatched BeginOfEvent for " << GetName()
00501 << " !" ;
00502 #endif
00503 clearHits();
00504 }
00505
00506 void CaloSD::update(const EndOfTrack * trk) {
00507 int id = (*trk)()->GetTrackID();
00508 TrackInformation *trkI =(TrackInformation *)((*trk)()->GetUserInformation());
00509 int lastTrackID = -1;
00510 if (trkI) lastTrackID = trkI->getIDonCaloSurface();
00511 if (id == lastTrackID) {
00512 const TrackContainer * trksForThisEvent = m_trackManager->trackContainer();
00513 if (trksForThisEvent != NULL) {
00514 int it = (int)(trksForThisEvent->size()) - 1;
00515 if (it >= 0) {
00516 TrackWithHistory * trkH = (*trksForThisEvent)[it];
00517 if (trkH->trackID() == (unsigned int)(id)) tkMap[id] = trkH;
00518 #ifdef DebugLog
00519 LogDebug("CaloSim") << "CaloSD: get track " << it << " from "
00520 << "Container of size " << trksForThisEvent->size()
00521 << " with ID " << trkH->trackID();
00522 } else {
00523 LogDebug("CaloSim") << "CaloSD: get track " << it << " from "
00524 << "Container of size " << trksForThisEvent->size()
00525 << " with no ID";
00526 #endif
00527 }
00528 }
00529 }
00530 }
00531
00532 void CaloSD::update(const ::EndOfEvent * ) {
00533 int count = 0, wrong = 0;
00534 bool ok;
00535
00536 slave->ReserveMemory(theHC->entries());
00537
00538 for (int i=0; i<theHC->entries(); ++i) {
00539 ok = saveHit((*theHC)[i]);
00540 ++count;
00541 if (!ok) ++wrong;
00542 }
00543
00544 edm::LogInfo("CaloSim") << "CaloSD: " << GetName() << " store " << count
00545 << " hits recorded with " << wrong
00546 << " track IDs not given properly and "
00547 << totalHits-count << " hits not passing cuts";
00548 summarize();
00549
00550 tkMap.erase (tkMap.begin(), tkMap.end());
00551 }
00552
00553 void CaloSD::clearHits() {
00554 if (useMap) hitMap.erase (hitMap.begin(), hitMap.end());
00555 for (unsigned int i = 0; i<reusehit.size(); ++i) delete reusehit[i];
00556 std::vector<CaloG4Hit*>().swap(reusehit);
00557 cleanIndex = 0;
00558 previousID.reset();
00559 primIDSaved = -99;
00560 #ifdef DebugLog
00561 LogDebug("CaloSim") << "CaloSD: Clears hit vector for " << GetName() << " " << slave;
00562 #endif
00563 slave->Initialize();
00564 #ifdef DebugLog
00565 LogDebug("CaloSim") << "CaloSD: Initialises slave SD for " << GetName();
00566 #endif
00567 }
00568
00569 void CaloSD::initRun() {}
00570
00571 int CaloSD::getTrackID(G4Track* aTrack) {
00572 int primaryID = 0;
00573 forceSave = false;
00574 TrackInformation* trkInfo=(TrackInformation *)(aTrack->GetUserInformation());
00575 if (trkInfo) {
00576 primaryID = trkInfo->getIDonCaloSurface();
00577 #ifdef DebugLog
00578 LogDebug("CaloSim") << "CaloSD: hit update from track Id on Calo Surface "
00579 << trkInfo->getIDonCaloSurface();
00580 #endif
00581 } else {
00582 primaryID = aTrack->GetTrackID();
00583 #ifdef DebugLog
00584 edm::LogWarning("CaloSim") << "CaloSD: Problem with primaryID **** set by "
00585 << "force to TkID **** " << primaryID << " in "
00586 << preStepPoint->GetTouchable()->GetVolume(0)->GetName();
00587 #endif
00588 }
00589 return primaryID;
00590 }
00591
00592 uint16_t CaloSD::getDepth(G4Step*) { return 0; }
00593
00594 bool CaloSD::filterHit(CaloG4Hit* hit, double time) {
00595 double emin(eminHit);
00596 if (hit->getDepth() > 0) emin = eminHitD;
00597 #ifdef DebugLog
00598 LogDebug("CaloSim") << "Depth " << hit->getDepth() << " Emin = " << emin << " ("
00599 << eminHit << ", " << eminHitD << ")";
00600 #endif
00601 return ((time <= tmaxHit) && (hit->getEnergyDeposit() > emin));
00602 }
00603
00604 double CaloSD::getResponseWt(G4Track* aTrack) {
00605 if (meanResponse) {
00606 TrackInformation * trkInfo = (TrackInformation *)(aTrack->GetUserInformation());
00607 return meanResponse->getWeight(trkInfo->genParticlePID(), trkInfo->genParticleP());
00608 } else {
00609 return 1;
00610 }
00611 }
00612
00613 void CaloSD::storeHit(CaloG4Hit* hit) {
00614 if (previousID.trackID()<0) return;
00615 if (hit == 0) {
00616 edm::LogWarning("CaloSim") << "CaloSD: hit to be stored is NULL !!";
00617 return;
00618 }
00619
00620 theHC->insert(hit);
00621 if (useMap) hitMap.insert(std::pair<CaloHitID,CaloG4Hit*>(previousID,hit));
00622 }
00623
00624 bool CaloSD::saveHit(CaloG4Hit* aHit) {
00625 int tkID;
00626 bool ok = true;
00627 if (m_trackManager) {
00628 tkID = m_trackManager->giveMotherNeeded(aHit->getTrackID());
00629 if (tkID == 0) {
00630 if (m_trackManager->trackExists(aHit->getTrackID())) tkID = (aHit->getTrackID());
00631 else {
00632 ok = false;
00633 }
00634 }
00635 } else {
00636 tkID = aHit->getTrackID();
00637 ok = false;
00638 }
00639
00640 #ifdef DebugLog
00641 LogDebug("CaloSim") << "CalosD: Track ID " << aHit->getTrackID()
00642 << " changed to " << tkID << " by SimTrackManager"
00643 << " Status " << ok;
00644 #endif
00645 double time = aHit->getTimeSlice();
00646 if (corrTOFBeam) time += correctT;
00647 slave->processHits(aHit->getUnitID(), aHit->getEM()/GeV,
00648 aHit->getHadr()/GeV, time, tkID, aHit->getDepth());
00649 #ifdef DebugLog
00650 LogDebug("CaloSim") << "CaloSD: Store Hit at " << std::hex
00651 << aHit->getUnitID() << std::dec << " "
00652 << aHit->getDepth() << " due to " << tkID
00653 << " in time " << time << " of energy "
00654 << aHit->getEM()/GeV << " GeV (EM) and "
00655 << aHit->getHadr()/GeV << " GeV (Hadr)";
00656 #endif
00657 return ok;
00658 }
00659
00660 void CaloSD::summarize() {}
00661
00662 void CaloSD::update(const BeginOfTrack * trk) {
00663 int primary = -1;
00664 TrackInformation * trkInfo = (TrackInformation *)((*trk)()->GetUserInformation());
00665 if ( trkInfo->isPrimary() ) primary = (*trk)()->GetTrackID();
00666
00667 #ifdef DebugLog
00668 LogDebug("CaloSim") << "New track: isPrimary " << trkInfo->isPrimary()
00669 << " primary ID = " << primary
00670 << " primary ancestor ID " << primAncestor;
00671 #endif
00672
00673
00674
00675 if (primary > 0 && primary != primAncestor) {
00676 primAncestor = primary;
00677
00678
00679
00680 if (theHC->entries()>0) cleanHitCollection();
00681
00682 }
00683 }
00684
00685 void CaloSD::cleanHitCollection() {
00686 std::vector<CaloG4Hit*>* theCollection = theHC->GetVector();
00687
00688 #ifdef DebugLog
00689 LogDebug("CaloSim") << "CaloSD: collection before merging, size = " << theHC->entries();
00690 #endif
00691
00692 selIndex.reserve(theHC->entries()-cleanIndex);
00693 if ( reusehit.size() == 0 ) reusehit.reserve(theHC->entries()-cleanIndex);
00694
00695
00696 if ( !useMap ) {
00697 hitvec.swap(*theCollection);
00698 sort((hitvec.begin()+cleanIndex), hitvec.end(), CaloG4HitLess());
00699 #ifdef DebugLog
00700 LogDebug("CaloSim") << "CaloSD::cleanHitCollection: sort hits in buffer "
00701 << "starting from element = " << cleanIndex;
00702 for (unsigned int i = 0; i<hitvec.size(); ++i)
00703 LogDebug("CaloSim")<<i<<" "<<*hitvec[i];
00704 #endif
00705 unsigned int i, j;
00706 CaloG4HitEqual equal;
00707 for (i=cleanIndex; i<hitvec.size(); ++i) {
00708 selIndex.push_back(i-cleanIndex);
00709 int jump = 0;
00710 for (j = i+1; j <hitvec.size() && equal(hitvec[i], hitvec[j]); ++j) {
00711 ++jump;
00712
00713 (*hitvec[i]).addEnergyDeposit(*hitvec[j]);
00714 (*hitvec[j]).setEM(0.);
00715 (*hitvec[j]).setHadr(0.);
00716 reusehit.push_back(hitvec[j]);
00717 }
00718 i+=jump;
00719 }
00720 #ifdef DebugLog
00721 LogDebug("CaloSim") << "CaloSD: cleanHitCollection merge the hits in buffer ";
00722 for (unsigned int i = 0; i<hitvec.size(); ++i)
00723 LogDebug("CaloSim")<<i<<" "<<*hitvec[i];
00724 #endif
00725 for ( unsigned int i = cleanIndex; i < cleanIndex+selIndex.size(); ++i ) {
00726 hitvec[i] = hitvec[selIndex[i-cleanIndex]+cleanIndex];
00727 }
00728 hitvec.resize(cleanIndex+selIndex.size());
00729 #ifdef DebugLog
00730 LogDebug("CaloSim") << "CaloSD::cleanHitCollection: remove the merged hits in buffer,"
00731 << " new size = " << hitvec.size();
00732 for (unsigned int i = 0; i<hitvec.size(); ++i)
00733 LogDebug("CaloSim")<<i<<" "<<*hitvec[i];
00734 #endif
00735 hitvec.swap(*theCollection);
00736 std::vector<CaloG4Hit*>().swap(hitvec);
00737 selIndex.clear();
00738 totalHits = theHC->entries();
00739 }
00740
00741 #ifdef DebugLog
00742 LogDebug("CaloSim") << "CaloSD: collection after merging, size = " << theHC->entries();
00743 #endif
00744
00745 int addhit = 0;
00746
00747 #ifdef DebugLog
00748 LogDebug("CaloSim") << "CaloSD: Size of reusehit after merge = " << reusehit.size();
00749 LogDebug("CaloSim") << "CaloSD: Starting hit selection from index = " << cleanIndex;
00750 #endif
00751
00752 selIndex.reserve(theCollection->size()-cleanIndex);
00753 for (unsigned int i = cleanIndex; i<theCollection->size(); ++i) {
00754 CaloG4Hit* aHit((*theCollection)[i]);
00755
00756
00757
00758 double time = aHit->getTimeSlice();
00759 if (corrTOFBeam) time += correctT;
00760 if (!filterHit(aHit,time)) {
00761 #ifdef DebugLog
00762 LogDebug("CaloSim") << "CaloSD: dropped CaloG4Hit " << " " << *aHit;
00763 #endif
00764
00765
00766
00767 reusehit.push_back((*theCollection)[i]);
00768 ++addhit;
00769 } else {
00770 selIndex.push_back(i-cleanIndex);
00771 }
00772 }
00773
00774 #ifdef DebugLog
00775 LogDebug("CaloSim") << "CaloSD: Size of reusehit after selection = " << reusehit.size()
00776 << " Number of added hit = " << addhit;
00777 #endif
00778 if (useMap) {
00779 if ( addhit>0 ) {
00780 int offset = reusehit.size()-addhit;
00781 for (int ii = addhit-1; ii>=0; --ii) {
00782 CaloHitID theID = reusehit[offset+ii]->getID();
00783 hitMap.erase(theID);
00784 }
00785 }
00786 }
00787 for (unsigned int j = 0; j<selIndex.size(); ++j) {
00788 (*theCollection)[cleanIndex+j] = (*theCollection)[cleanIndex+selIndex[j]];
00789 }
00790
00791 theCollection->resize(cleanIndex+selIndex.size());
00792 std::vector<unsigned int>().swap(selIndex);
00793
00794 #ifdef DebugLog
00795 LogDebug("CaloSim") << "CaloSD: hit collection after selection, size = "
00796 << theHC->entries();
00797 theHC->PrintAllHits();
00798 #endif
00799
00800 cleanIndex = theHC->entries();
00801 }