00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015 #include <cmath>
00016 #include <iostream>
00017 #include <iomanip>
00018
00019
00020 #include "SimG4CMS/HcalTestBeam/interface/HcalTB06Analysis.h"
00021
00022 #include "SimG4Core/Notification/interface/BeginOfRun.h"
00023 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
00024 #include "SimG4Core/Notification/interface/EndOfEvent.h"
00025 #include "SimG4Core/Watcher/interface/SimWatcherFactory.h"
00026
00027
00028 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00029 #include "SimG4CMS/Calo/interface/CaloG4Hit.h"
00030 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
00031 #include "DataFormats/Math/interface/Point3D.h"
00032
00033 #include "FWCore/Framework/interface/Event.h"
00034 #include "FWCore/Framework/interface/EventSetup.h"
00035 #include "FWCore/Framework/interface/ESHandle.h"
00036 #include "FWCore/Framework/interface/MakerMacros.h"
00037 #include "FWCore/PluginManager/interface/ModuleDef.h"
00038 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00039
00040 #include "G4SDManager.hh"
00041 #include "G4VProcess.hh"
00042 #include "G4HCofThisEvent.hh"
00043 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00044 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00045
00046
00047
00048
00049
00050 HcalTB06Analysis::HcalTB06Analysis(const edm::ParameterSet &p): histo(0) {
00051
00052 edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("HcalTB06Analysis");
00053 names = m_Anal.getParameter<std::vector<std::string> >("Names");
00054 beamOffset =-m_Anal.getParameter<double>("BeamPosition")*cm;
00055 double fMinEta = m_Anal.getParameter<double>("MinEta");
00056 double fMaxEta = m_Anal.getParameter<double>("MaxEta");
00057 double fMinPhi = m_Anal.getParameter<double>("MinPhi");
00058 double fMaxPhi = m_Anal.getParameter<double>("MaxPhi");
00059 double beamEta = (fMaxEta+fMinEta)/2.;
00060 double beamPhi = (fMaxPhi+fMinPhi)/2.;
00061 double beamThet= 2*atan(exp(-beamEta));
00062 if (beamPhi < 0) beamPhi += twopi;
00063 iceta = (int)(beamEta/0.087) + 1;
00064 icphi = (int)(fabs(beamPhi)/0.087) + 5;
00065 if (icphi > 72) icphi -= 73;
00066
00067 produces<PHcalTB06Info>();
00068
00069 beamline_RM = new G4RotationMatrix;
00070 beamline_RM->rotateZ(-beamPhi);
00071 beamline_RM->rotateY(-beamThet);
00072
00073 edm::LogInfo("HcalTBSim") << "HcalTB06:: Initialised as observer of BeginOf"
00074 << "Job/BeginOfRun/BeginOfEvent/G4Step/EndOfEvent"
00075 << " with Parameter values:\n \tbeamOffset = "
00076 << beamOffset << "\ticeta = " << iceta
00077 << "\ticphi = " << icphi << "\n\tbeamline_RM = "
00078 << *beamline_RM;
00079
00080 init();
00081
00082 histo = new HcalTB06Histo(m_Anal);
00083 }
00084
00085 HcalTB06Analysis::~HcalTB06Analysis() {
00086
00087 edm::LogInfo("HcalTBSim") << "\n --------> Total number of selected entries"
00088 << " : " << count << "\nPointers:: Histo " <<histo;
00089 if (histo) {
00090 delete histo;
00091 histo = 0;
00092 }
00093 }
00094
00095
00096
00097
00098
00099 void HcalTB06Analysis::produce(edm::Event& e, const edm::EventSetup&) {
00100
00101 std::auto_ptr<PHcalTB06Info> product(new PHcalTB06Info);
00102 fillEvent(*product);
00103 e.put(product);
00104 }
00105
00106 void HcalTB06Analysis::init() {
00107
00108
00109 count = 0;
00110 evNum = 0;
00111 clear();
00112 }
00113
00114 void HcalTB06Analysis::update(const BeginOfRun * run) {
00115
00116 int irun = (*run)()->GetRunID();
00117 edm::LogInfo("HcalTBSim") <<" =====> Begin of Run = " << irun;
00118
00119 }
00120
00121 void HcalTB06Analysis::update(const BeginOfEvent * evt) {
00122
00123 evNum = (*evt) ()->GetEventID ();
00124 clear();
00125 edm::LogInfo("HcalTBSim") << "HcalTB06Analysis: =====> Begin of event = "
00126 << evNum;
00127 }
00128
00129 void HcalTB06Analysis::update(const G4Step * aStep) {
00130
00131 if (aStep != NULL) {
00132
00133 G4ThreeVector thePreStepPoint = aStep->GetPreStepPoint()->GetPosition();
00134 G4ThreeVector thePostStepPoint;
00135
00136
00137 G4Track* aTrack = aStep->GetTrack();
00138 int trackID = aTrack->GetTrackID();
00139 int parentID = aTrack->GetParentID();
00140 G4ThreeVector position = aTrack->GetPosition();
00141 G4ThreeVector momentum = aTrack->GetMomentum();
00142 G4String partType = aTrack->GetDefinition()->GetParticleType();
00143 G4String partSubType = aTrack->GetDefinition()->GetParticleSubType();
00144 int partPDGEncoding = aTrack->GetDefinition()->GetPDGEncoding();
00145 #ifdef ddebug
00146 bool isPDGStable = aTrack->GetDefinition()->GetPDGStable();
00147 #endif
00148 double pDGlifetime = aTrack->GetDefinition()->GetPDGLifeTime();
00149 double gammaFactor = aStep->GetPreStepPoint()->GetGamma();
00150
00151 if (!pvFound) {
00152 double stepDeltaEnergy = aStep->GetDeltaEnergy ();
00153 double kinEnergy = aTrack->GetKineticEnergy ();
00154
00155
00156 if (trackID == 1 && parentID == 0 &&
00157 ((kinEnergy == 0.) || (fabs (stepDeltaEnergy / kinEnergy) > 0.1))) {
00158 pvType = -1;
00159 if (kinEnergy == 0.) {
00160 pvType = 0;
00161 } else {
00162 if (fabs (stepDeltaEnergy / kinEnergy) > 0.1) pvType = 1;
00163 }
00164 pvFound = true;
00165 pvPosition = position;
00166 pvMomentum = momentum;
00167
00168 pvUVW = (*beamline_RM)*(pvPosition);
00169
00170
00171 G4String thePostPVname = "NoName";
00172 G4StepPoint * thePostPoint = aStep->GetPostStepPoint ();
00173 if (thePostPoint) {
00174 thePostStepPoint = thePostPoint->GetPosition();
00175 G4VPhysicalVolume * thePostPV = thePostPoint->GetPhysicalVolume ();
00176 if (thePostPV) thePostPVname = thePostPV->GetName ();
00177 }
00178 #ifdef ddebug
00179 LogDebug("HcalTBSim") << "HcalTB06Analysis:: V1 found at: "
00180 << thePostStepPoint << " G4VPhysicalVolume: "
00181 << thePostPVname;
00182 #endif
00183 LogDebug("HcalTBSim") << "HcalTB06Analysis::fill_v1Pos: Primary Track "
00184 << "momentum: " << pvMomentum << " psoition "
00185 << pvPosition << " u/v/w " << pvUVW;
00186 }
00187 } else {
00188
00189 if ((trackID != 1 && parentID == 1 &&
00190 (aTrack->GetCurrentStepNumber () == 1) &&
00191 (thePreStepPoint == pvPosition)) ||
00192 (trackID == 1 && thePreStepPoint == pvPosition)) {
00193 #ifdef ddebug
00194 LogDebug("HcalTBSim") << "HcalTB06Analysis::A secondary... PDG:"
00195 << partPDGEncoding << " TrackID:" << trackID
00196 << " ParentID:" << parentID << " stable: "
00197 << isPDGStable << " Tau: " << pDGlifetime
00198 << " cTauGamma="
00199 << c_light*pDGlifetime*gammaFactor*1000.
00200 << "um" << " GammaFactor: " << gammaFactor;
00201 #endif
00202 secTrackID.push_back(trackID);
00203 secPartID.push_back(partPDGEncoding);
00204 secMomentum.push_back(momentum);
00205 secEkin.push_back(aTrack->GetKineticEnergy());
00206
00207
00208 double ctaugamma_um = c_light * pDGlifetime * gammaFactor * 1000.;
00209 if ((ctaugamma_um>0.) && (ctaugamma_um<100.)) {
00210 shortLivedSecondaries.push_back(trackID);
00211 } else {
00212
00213 }
00214 }
00215
00216
00217 if (aTrack->GetCurrentStepNumber() == 1) {
00218 if (shortLivedSecondaries.size() > 0) {
00219 int pid = parentID;
00220 std::vector<int>::iterator pos1= shortLivedSecondaries.begin();
00221 std::vector<int>::iterator pos2 = shortLivedSecondaries.end();
00222 std::vector<int>::iterator pos;
00223 for (pos = pos1; pos != pos2; pos++) {
00224 if (*pos == pid) {
00225
00226 #ifdef ddebug
00227 LogDebug("HcalTBSim") << "HcalTB06Analysis:: A tertiary... PDG:"
00228 << partPDGEncoding << " TrackID:" <<trackID
00229 << " ParentID:" << parentID << " stable: "
00230 << isPDGStable << " Tau: " << pDGlifetime
00231 << " cTauGamma="
00232 << c_light*pDGlifetime*gammaFactor*1000.
00233 << "um GammaFactor: " << gammaFactor;
00234 #endif
00235 }
00236 }
00237 }
00238 }
00239 }
00240 }
00241 }
00242
00243 void HcalTB06Analysis::update(const EndOfEvent * evt) {
00244
00245 count++;
00246
00247
00248 LogDebug("HcalTBSim") << "HcalTB06Analysis::Fill event "
00249 << (*evt)()->GetEventID();
00250 fillBuffer (evt);
00251
00252
00253 LogDebug("HcalTBSim") << "HcalTB06Analysis::Final analysis";
00254 finalAnalysis();
00255
00256 int iEvt = (*evt)()->GetEventID();
00257 if (iEvt < 10)
00258 edm::LogInfo("HcalTBSim") << "HcalTB06Analysis:: Event " << iEvt;
00259 else if ((iEvt < 100) && (iEvt%10 == 0))
00260 edm::LogInfo("HcalTBSim") << "HcalTB06Analysis:: Event " << iEvt;
00261 else if ((iEvt < 1000) && (iEvt%100 == 0))
00262 edm::LogInfo("HcalTBSim") << "HcalTB06Analysis:: Event " << iEvt;
00263 else if ((iEvt < 10000) && (iEvt%1000 == 0))
00264 edm::LogInfo("HcalTBSim") << "HcalTB06Analysis:: Event " << iEvt;
00265 }
00266
00267 void HcalTB06Analysis::fillBuffer(const EndOfEvent * evt) {
00268
00269 std::vector<CaloHit> hhits;
00270 int idHC, j;
00271 CaloG4HitCollection* theHC;
00272 std::map<int,float,std::less<int> > primaries;
00273 double etot1=0, etot2=0;
00274
00275
00276 G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
00277 std::string sdName = names[0];
00278 idHC = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
00279 theHC = (CaloG4HitCollection*) allHC->GetHC(idHC);
00280 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Hit Collection for " << sdName
00281 << " of ID " << idHC << " is obtained at " << theHC;
00282
00283 if (idHC >= 0 && theHC > 0) {
00284 hhits.reserve(theHC->entries());
00285 for (j = 0; j < theHC->entries(); j++) {
00286 CaloG4Hit* aHit = (*theHC)[j];
00287 double e = aHit->getEnergyDeposit()/GeV;
00288 double time = aHit->getTimeSlice();
00289 math::XYZPoint pos = aHit->getEntry();
00290 unsigned int id = aHit->getUnitID();
00291 double theta = pos.theta();
00292 double eta = -log(tan(theta * 0.5));
00293 double phi = pos.phi();
00294 CaloHit hit(2,1,e,eta,phi,time,id);
00295 hhits.push_back(hit);
00296 primaries[aHit->getTrackID()]+= e;
00297 etot1 += e;
00298 #ifdef ddebug
00299 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Hcal Hit i/p " << j
00300 << " ID 0x" << std::hex << id << std::dec
00301 << " time " << std::setw(6) << time << " theta "
00302 << std::setw(8) << theta << " eta " << std::setw(8)
00303 << eta << " phi " << std::setw(8) << phi << " e "
00304 << std::setw(8) << e;
00305 #endif
00306 }
00307 }
00308
00309
00310 std::vector<CaloHit>::iterator itr;
00311 int nHit = hhits.size();
00312 std::vector<CaloHit*> hits(nHit);
00313 for (j = 0, itr = hhits.begin(); itr != hhits.end(); j++, itr++) {
00314 hits[j] = &hhits[j];
00315 }
00316 sort(hits.begin(),hits.end(),CaloHitIdMore());
00317 std::vector<CaloHit*>::iterator k1, k2;
00318 int nhit = 0;
00319 for (k1 = hits.begin(); k1 != hits.end(); k1++) {
00320 int det = (**k1).det();
00321 int layer = (**k1).layer();
00322 double ehit = (**k1).e();
00323 double eta = (**k1).eta();
00324 double phi = (**k1).phi();
00325 double jitter = (**k1).t();
00326 uint32_t unitID = (**k1).id();
00327 int jump = 0;
00328 for (k2 = k1+1; k2 != hits.end() && fabs(jitter-(**k2).t())<1 &&
00329 unitID==(**k2).id(); k2++) {
00330 ehit += (**k2).e();
00331 jump++;
00332 }
00333 nhit++;
00334 CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
00335 hcalHitCache.push_back(hit);
00336 etot2 += ehit;
00337 k1 += jump;
00338 #ifdef ddebug
00339 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Hcal Hit store " << nhit
00340 << " ID 0x" << std::hex << unitID << std::dec
00341 << " time " << std::setw(6) << jitter << " eta "
00342 << std::setw(8) << eta << " phi " << std::setw(8)
00343 << phi << " e " << std::setw(8) << ehit;
00344 #endif
00345 }
00346 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Stores " << nhit << " HCal hits"
00347 << " from " << nHit << " input hits E(Hcal) " << etot1
00348 << " " << etot2;
00349
00350
00351 std::vector<CaloHit> ehits;
00352 sdName= names[1];
00353 idHC = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
00354 theHC = (CaloG4HitCollection*) allHC->GetHC(idHC);
00355 etot1 = etot2 = 0;
00356 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Hit Collection for " << sdName
00357 << " of ID " << idHC << " is obtained at " << theHC;
00358 if (idHC >= 0 && theHC > 0) {
00359 ehits.reserve(theHC->entries());
00360 for (j = 0; j < theHC->entries(); j++) {
00361 CaloG4Hit* aHit = (*theHC)[j];
00362 double e = aHit->getEnergyDeposit()/GeV;
00363 double time = aHit->getTimeSlice();
00364 math::XYZPoint pos = aHit->getEntry();
00365 unsigned int id = aHit->getUnitID();
00366 double theta = pos.theta();
00367 double eta = -log(tan(theta * 0.5));
00368 double phi = pos.phi();
00369 if (e < 0 || e > 100000.) e = 0;
00370 CaloHit hit(1,0,e,eta,phi,time,id);
00371 ehits.push_back(hit);
00372 primaries[aHit->getTrackID()]+= e;
00373 etot1 += e;
00374 #ifdef ddebug
00375 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Ecal Hit i/p " << j
00376 << " ID 0x" << std::hex << id << std::dec
00377 << " time " << std::setw(6) << time << " theta "
00378 << std::setw(8) << theta << " eta " <<std::setw(8)
00379 << eta << " phi " << std::setw(8) << phi << " e "
00380 << std::setw(8) << e;
00381 #endif
00382 }
00383 }
00384
00385
00386 nHit = ehits.size();
00387 std::vector<CaloHit*> hite(nHit);
00388 for (j = 0, itr = ehits.begin(); itr != ehits.end(); j++, itr++) {
00389 hite[j] = &ehits[j];
00390 }
00391 sort(hite.begin(),hite.end(),CaloHitIdMore());
00392 nhit = 0;
00393 for (k1 = hite.begin(); k1 != hite.end(); k1++) {
00394 int det = (**k1).det();
00395 int layer = (**k1).layer();
00396 double ehit = (**k1).e();
00397 double eta = (**k1).eta();
00398 double phi = (**k1).phi();
00399 double jitter = (**k1).t();
00400 uint32_t unitID = (**k1).id();
00401 int jump = 0;
00402 for (k2 = k1+1; k2 != hite.end() && fabs(jitter-(**k2).t())<1 &&
00403 unitID==(**k2).id(); k2++) {
00404 ehit += (**k2).e();
00405 jump++;
00406 }
00407 nhit++;
00408 CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
00409 ecalHitCache.push_back(hit);
00410 etot2 += ehit;
00411 k1 += jump;
00412 #ifdef ddebug
00413 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Ecal Hit store " << nhit
00414 << " ID 0x" << std::hex << unitID << std::dec
00415 << " time " << std::setw(6) << jitter << " eta "
00416 << std::setw(8) << eta << " phi " << std::setw(8)
00417 << phi << " e " << std::setw(8) << ehit;
00418 #endif
00419 }
00420 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Stores " << nhit << " ECal hits"
00421 << " from " << nHit << " input hits E(Ecal) " << etot1
00422 << " " << etot2;
00423
00424
00425 nPrimary = (int)(primaries.size());
00426 int trackID = 0;
00427 G4PrimaryParticle* thePrim=0;
00428 int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
00429 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Event has " << nvertex
00430 << " verteices";
00431 if (nvertex<=0)
00432 edm::LogInfo("HcalTBSim") << "HcalTB06Analysis::EndOfEvent ERROR: no "
00433 << "vertex found for event " << evNum;
00434
00435 for (int i = 0 ; i<nvertex; i++) {
00436 G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
00437 if (avertex == 0) {
00438 edm::LogInfo("HcalTBSim") << "HcalTB06Analysis::EndOfEvent ERR: pointer "
00439 << "to vertex = 0 for event " << evNum;
00440 } else {
00441 LogDebug("HcalTBSim") << "HcalTB06Analysis::Vertex number :" << i << " "
00442 << avertex->GetPosition();
00443 int npart = avertex->GetNumberOfParticle();
00444 if (npart == 0)
00445 edm::LogWarning("HcalTBSim") << "HcalTB06Analysis::End Of Event ERR: "
00446 << "no primary!";
00447 if (thePrim==0) thePrim=avertex->GetPrimary(trackID);
00448 }
00449 }
00450
00451 if (thePrim != 0) {
00452 double px = thePrim->GetPx();
00453 double py = thePrim->GetPy();
00454 double pz = thePrim->GetPz();
00455 double p = std::sqrt(pow(px,2.)+pow(py,2.)+pow(pz,2.));
00456 pInit = p/GeV;
00457 if (p==0)
00458 edm::LogWarning("HcalTBSim") << "HcalTB06Analysis:: EndOfEvent ERR: "
00459 << "primary has p=0 ";
00460 else {
00461 double costheta = pz/p;
00462 double theta = acos(std::min(std::max(costheta,-1.),1.));
00463 etaInit = -log(tan(theta/2));
00464 if (px != 0 || py != 0) phiInit = atan2(py,px);
00465 }
00466 particleType = thePrim->GetPDGcode();
00467 } else
00468 edm::LogWarning("HcalTBSim") << "HcalTB06Analysis::EndOfEvent ERR: could "
00469 << "not find primary";
00470
00471 }
00472
00473 void HcalTB06Analysis::finalAnalysis() {
00474
00475
00476 histo->fillPrimary(pInit, etaInit, phiInit);
00477
00478
00479 eecals = ehcals = 0.;
00480 for (unsigned int i=0; i<hcalHitCache.size(); i++) {
00481 ehcals += hcalHitCache[i].e();
00482 }
00483 for (unsigned int i=0; i<ecalHitCache.size(); i++) {
00484 eecals += ecalHitCache[i].e();
00485 }
00486 etots = eecals + ehcals;
00487 LogDebug("HcalTBSim") << "HcalTB06Analysis:: Energy deposit at Sim Level "
00488 << "(Total) " << etots << " (ECal) " << eecals
00489 << " (HCal) " << ehcals;
00490 histo->fillEdep(etots, eecals, ehcals);
00491 }
00492
00493
00494 void HcalTB06Analysis::fillEvent (PHcalTB06Info& product) {
00495
00496
00497 product.setPrimary(nPrimary, particleType, pInit, etaInit, phiInit);
00498
00499
00500 product.setEdep(etots, eecals, ehcals);
00501
00502
00503 for (unsigned int i=0; i<hcalHitCache.size(); i++)
00504 product.saveHit(hcalHitCache[i].id(), hcalHitCache[i].eta(),
00505 hcalHitCache[i].phi(), hcalHitCache[i].e(),
00506 hcalHitCache[i].t());
00507 for (unsigned int i=0; i<ecalHitCache.size(); i++)
00508 product.saveHit(ecalHitCache[i].id(), ecalHitCache[i].eta(),
00509 ecalHitCache[i].phi(), ecalHitCache[i].e(),
00510 ecalHitCache[i].t());
00511
00512
00513 product.setVtxPrim(evNum, pvType, pvPosition.x(), pvPosition.y(),
00514 pvPosition.z(), pvUVW.x(), pvUVW.y(), pvUVW.z(),
00515 pvMomentum.x(), pvMomentum.y(), pvMomentum.z());
00516 for (unsigned int i = 0; i < secTrackID.size(); i++) {
00517 product.setVtxSec(secTrackID[i], secPartID[i], secMomentum[i].x(),
00518 secMomentum[i].y(), secMomentum[i].z(), secEkin[i]);
00519 }
00520 }
00521
00522 void HcalTB06Analysis::clear() {
00523
00524 pvFound = false;
00525 pvType =-2;
00526 pvPosition = G4ThreeVector();
00527 pvMomentum = G4ThreeVector();
00528 pvUVW = G4ThreeVector();
00529 secTrackID.clear();
00530 secPartID.clear();
00531 secMomentum.clear();
00532 secEkin.clear();
00533 shortLivedSecondaries.clear();
00534
00535 ecalHitCache.erase(ecalHitCache.begin(), ecalHitCache.end());
00536 hcalHitCache.erase(hcalHitCache.begin(), hcalHitCache.end());
00537 nPrimary = particleType = 0;
00538 pInit = etaInit = phiInit = 0;
00539 }
00540
00541 DEFINE_SIMWATCHER (HcalTB06Analysis);