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/HcalTB04Analysis.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/ECalSD.h"
00030 #include "SimG4CMS/Calo/interface/HCalSD.h"
00031 #include "SimG4CMS/Calo/interface/CaloG4Hit.h"
00032 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
00033 #include "SimG4CMS/Calo/interface/HcalTestNumberingScheme.h"
00034 #include "SimG4CMS/HcalTestBeam/interface/HcalTBNumberingScheme.h"
00035 #include "SimG4CMS/HcalTestBeam/interface/HcalTB04XtalNumberingScheme.h"
00036 #include "DataFormats/Math/interface/Point3D.h"
00037
00038 #include "FWCore/Framework/interface/Event.h"
00039 #include "FWCore/Framework/interface/EventSetup.h"
00040 #include "FWCore/Framework/interface/ESHandle.h"
00041 #include "FWCore/Framework/interface/MakerMacros.h"
00042 #include "FWCore/PluginManager/interface/ModuleDef.h"
00043 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00044
00045 #include "FWCore/ServiceRegistry/interface/Service.h"
00046 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00047 #include "CLHEP/Random/RandGaussQ.h"
00048 #include "FWCore/Utilities/interface/Exception.h"
00049
00050 #include "G4SDManager.hh"
00051 #include "G4VProcess.hh"
00052 #include "G4HCofThisEvent.hh"
00053 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00054 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00055
00056
00057
00058
00059
00060 HcalTB04Analysis::HcalTB04Analysis(const edm::ParameterSet &p): myQie(0),
00061 histo(0) {
00062
00063 edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("HcalTB04Analysis");
00064 hcalOnly = m_Anal.getParameter<bool>("HcalOnly");
00065 mode = m_Anal.getParameter<int>("Mode");
00066 type = m_Anal.getParameter<int>("Type");
00067 ecalNoise = m_Anal.getParameter<double>("EcalNoise");
00068 scaleHB0 = m_Anal.getParameter<double>("ScaleHB0");
00069 scaleHB16 = m_Anal.getParameter<double>("ScaleHB16");
00070 scaleHO = m_Anal.getParameter<double>("ScaleHO");
00071 scaleHE0 = m_Anal.getParameter<double>("ScaleHE0");
00072 names = m_Anal.getParameter<std::vector<std::string> >("Names");
00073 beamOffset =-m_Anal.getParameter<double>("BeamPosition")*cm;
00074 double fMinEta = m_Anal.getParameter<double>("MinEta");
00075 double fMaxEta = m_Anal.getParameter<double>("MaxEta");
00076 double fMinPhi = m_Anal.getParameter<double>("MinPhi");
00077 double fMaxPhi = m_Anal.getParameter<double>("MaxPhi");
00078 double beamEta = (fMaxEta+fMinEta)/2.;
00079 double beamPhi = (fMaxPhi+fMinPhi)/2.;
00080 double beamThet= 2*atan(exp(-beamEta));
00081 if (beamPhi < 0) beamPhi += twopi;
00082 iceta = (int)(beamEta/0.087) + 1;
00083 icphi = (int)(fabs(beamPhi)/0.087) + 5;
00084 if (icphi > 72) icphi -= 73;
00085
00086 produces<PHcalTB04Info>();
00087
00088 beamline_RM = new G4RotationMatrix;
00089 beamline_RM->rotateZ(-beamPhi);
00090 beamline_RM->rotateY(-beamThet);
00091
00092 edm::LogInfo("HcalTBSim") << "HcalTB04:: Initialised as observer of BeginOf"
00093 << "Job/BeginOfRun/BeginOfEvent/G4Step/EndOfEvent"
00094 << " with Parameter values:\n \thcalOnly = "
00095 << hcalOnly << "\tecalNoise = " << ecalNoise
00096 << "\n\tMode = " << mode << " (0: HB2 Standard; "
00097 << "1:HB2 Segmented)" << "\tType = " << type
00098 << " (0: HB; 1 HE; 2 HB+HE)\n\tbeamOffset = "
00099 << beamOffset << "\ticeta = " << iceta
00100 << "\ticphi = " << icphi << "\n\tbeamline_RM = "
00101 << *beamline_RM;
00102
00103 init();
00104
00105 myQie = new HcalQie(p);
00106 histo = new HcalTB04Histo(m_Anal);
00107 }
00108
00109 HcalTB04Analysis::~HcalTB04Analysis() {
00110
00111 edm::LogInfo("HcalTBSim") << "\n --------> Total number of selected entries"
00112 << " : " << count << "\nPointers:: QIE " << myQie
00113 << " Histo " << histo;
00114 if (myQie) {
00115 delete myQie;
00116 myQie = 0;
00117 }
00118 if (histo) {
00119 delete histo;
00120 histo = 0;
00121 }
00122 }
00123
00124
00125
00126
00127
00128 void HcalTB04Analysis::produce(edm::Event& e, const edm::EventSetup&) {
00129
00130 std::auto_ptr<PHcalTB04Info> product(new PHcalTB04Info);
00131 fillEvent(*product);
00132 e.put(product);
00133 }
00134
00135 void HcalTB04Analysis::init() {
00136
00137 idTower = HcalTBNumberingScheme::getUnitIDs(type, mode);
00138 nTower = idTower.size();
00139 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis:: Save information from "
00140 << nTower << " HCal towers";
00141 idHcal.reserve(nTower);
00142 for (int i=0; i<nTower; i++) {
00143 int id = unitID(idTower[i]);
00144 idHcal.push_back(id);
00145 LogDebug("HcalTBSim") << "\tTower[" << i << "] Original " << std::hex
00146 << idTower[i] << " Stored " << idHcal[i] << std::dec;
00147 }
00148
00149 if (!hcalOnly) {
00150 int det = 10;
00151 uint32_t id1;
00152 nCrystal = 0;
00153 for (int lay=1; lay<8; lay++) {
00154 for (int icr=1; icr<8; icr++) {
00155 id1 = HcalTestNumbering::packHcalIndex(det,0,1,icr,lay,1);
00156 int id = unitID(id1);
00157 idEcal.push_back(id1);
00158 idXtal.push_back(id);
00159 nCrystal++;
00160 }
00161 }
00162 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis:: Save information from "
00163 << nCrystal << " ECal Crystals";
00164 for (int i=0; i<nCrystal; i++) {
00165 LogDebug("HcalTBSim") << "\tCrystal[" << i << "] Original " << std::hex
00166 << idEcal[i] << " Stored " << idXtal[i] <<std::dec;
00167 }
00168 }
00169
00170 eseta.reserve(5);
00171 eqeta.reserve(5);
00172 esphi.reserve(3);
00173 eqphi.reserve(3);
00174 eslay.reserve(20);
00175 eqlay.reserve(20);
00176 for (int i=0; i<5; i++) {
00177 eseta.push_back(0.);
00178 eqeta.push_back(0.);
00179 }
00180 for (int i=0; i<3; i++) {
00181 esphi.push_back(0.);
00182 eqphi.push_back(0.);
00183 }
00184 for (int i=0; i<20; i++) {
00185 eslay.push_back(0.);
00186 eqlay.push_back(0.);
00187 }
00188
00189
00190 count = 0;
00191 evNum = 0;
00192 clear();
00193 }
00194
00195 void HcalTB04Analysis::update(const BeginOfRun * run) {
00196
00197 int irun = (*run)()->GetRunID();
00198 edm::LogInfo("HcalTBSim") <<" =====> Begin of Run = " << irun;
00199
00200 G4SDManager* sd = G4SDManager::GetSDMpointerIfExist();
00201 if (sd != 0) {
00202 std::string sdname = names[0];
00203 G4VSensitiveDetector* aSD = sd->FindSensitiveDetector(sdname);
00204 if (aSD==0) {
00205 edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::beginOfRun: No SD"
00206 << " with name " << sdname << " in this "
00207 << "Setup";
00208 } else {
00209 HCalSD* theCaloSD = dynamic_cast<HCalSD*>(aSD);
00210 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis::beginOfRun: Finds SD "
00211 << "with name " << theCaloSD->GetName()
00212 << " in this Setup";
00213 HcalNumberingScheme* org = new HcalTestNumberingScheme(false);
00214 theCaloSD->setNumberingScheme(org);
00215 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis::beginOfRun: set a "
00216 << "new numbering scheme";
00217 }
00218 if (!hcalOnly) {
00219 sdname = names[1];
00220 aSD = sd->FindSensitiveDetector(sdname);
00221 if (aSD==0) {
00222 edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::beginOfRun: No SD"
00223 << " with name " << sdname << " in this "
00224 << "Setup";
00225 } else {
00226 ECalSD* theCaloSD = dynamic_cast<ECalSD*>(aSD);
00227 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis::beginOfRun: Finds SD "
00228 << "with name " << theCaloSD->GetName()
00229 << " in this Setup";
00230 EcalNumberingScheme* org = new HcalTB04XtalNumberingScheme();
00231 theCaloSD->setNumberingScheme(org);
00232 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis::beginOfRun: set a "
00233 << "new numbering scheme";
00234 }
00235 }
00236 } else {
00237 edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::beginOfRun: Could "
00238 << "not get SD Manager!";
00239 }
00240
00241 }
00242
00243 void HcalTB04Analysis::update(const BeginOfEvent * evt) {
00244
00245 evNum = (*evt) ()->GetEventID ();
00246 clear();
00247 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis: =====> Begin of event = "
00248 << evNum;
00249 }
00250
00251 void HcalTB04Analysis::update(const G4Step * aStep) {
00252
00253 if (aStep != NULL) {
00254
00255 G4ThreeVector thePreStepPoint = aStep->GetPreStepPoint()->GetPosition();
00256 G4ThreeVector thePostStepPoint;
00257
00258
00259 G4Track* aTrack = aStep->GetTrack();
00260 int trackID = aTrack->GetTrackID();
00261 int parentID = aTrack->GetParentID();
00262 G4ThreeVector position = aTrack->GetPosition();
00263 G4ThreeVector momentum = aTrack->GetMomentum();
00264 G4String partType = aTrack->GetDefinition()->GetParticleType();
00265 G4String partSubType = aTrack->GetDefinition()->GetParticleSubType();
00266 int partPDGEncoding = aTrack->GetDefinition()->GetPDGEncoding();
00267 #ifdef ddebug
00268 bool isPDGStable = aTrack->GetDefinition()->GetPDGStable();
00269 #endif
00270 double pDGlifetime = aTrack->GetDefinition()->GetPDGLifeTime();
00271 double gammaFactor = aStep->GetPreStepPoint()->GetGamma();
00272
00273 if (!pvFound) {
00274 double stepDeltaEnergy = aStep->GetDeltaEnergy ();
00275 double kinEnergy = aTrack->GetKineticEnergy ();
00276
00277
00278 if (trackID == 1 && parentID == 0 &&
00279 ((kinEnergy == 0.) || (fabs (stepDeltaEnergy / kinEnergy) > 0.1))) {
00280 pvType = -1;
00281 if (kinEnergy == 0.) {
00282 pvType = 0;
00283 } else {
00284 if (fabs (stepDeltaEnergy / kinEnergy) > 0.1) pvType = 1;
00285 }
00286 pvFound = true;
00287 pvPosition = position;
00288 pvMomentum = momentum;
00289
00290 pvUVW = (*beamline_RM)*(pvPosition);
00291
00292
00293 G4String thePostPVname = "NoName";
00294 G4StepPoint * thePostPoint = aStep->GetPostStepPoint ();
00295 if (thePostPoint) {
00296 thePostStepPoint = thePostPoint->GetPosition();
00297 G4VPhysicalVolume * thePostPV = thePostPoint->GetPhysicalVolume ();
00298 if (thePostPV) thePostPVname = thePostPV->GetName ();
00299 }
00300 #ifdef ddebug
00301 LogDebug("HcalTBSim") << "HcalTB04Analysis:: V1 found at: "
00302 << thePostStepPoint << " G4VPhysicalVolume: "
00303 << thePostPVname;
00304 #endif
00305 LogDebug("HcalTBSim") << "HcalTB04Analysis::fill_v1Pos: Primary Track "
00306 << "momentum: " << pvMomentum << " psoition "
00307 << pvPosition << " u/v/w " << pvUVW;
00308 }
00309 } else {
00310
00311 if ((trackID != 1 && parentID == 1 &&
00312 (aTrack->GetCurrentStepNumber () == 1) &&
00313 (thePreStepPoint == pvPosition)) ||
00314 (trackID == 1 && thePreStepPoint == pvPosition)) {
00315 #ifdef ddebug
00316 LogDebug("HcalTBSim") << "HcalTB04Analysis::A secondary... PDG:"
00317 << partPDGEncoding << " TrackID:" << trackID
00318 << " ParentID:" << parentID << " stable: "
00319 << isPDGStable << " Tau: " << pDGlifetime
00320 << " cTauGamma="
00321 << c_light*pDGlifetime*gammaFactor*1000.
00322 << "um" << " GammaFactor: " << gammaFactor;
00323 #endif
00324 secTrackID.push_back(trackID);
00325 secPartID.push_back(partPDGEncoding);
00326 secMomentum.push_back(momentum);
00327 secEkin.push_back(aTrack->GetKineticEnergy());
00328
00329
00330 double ctaugamma_um = c_light * pDGlifetime * gammaFactor * 1000.;
00331 if ((ctaugamma_um>0.) && (ctaugamma_um<100.)) {
00332 shortLivedSecondaries.push_back(trackID);
00333 } else {
00334
00335 }
00336 }
00337
00338
00339 if (aTrack->GetCurrentStepNumber() == 1) {
00340 if (shortLivedSecondaries.size() > 0) {
00341 int pid = parentID;
00342 std::vector<int>::iterator pos1= shortLivedSecondaries.begin();
00343 std::vector<int>::iterator pos2 = shortLivedSecondaries.end();
00344 std::vector<int>::iterator pos;
00345 for (pos = pos1; pos != pos2; pos++) {
00346 if (*pos == pid) {
00347
00348 #ifdef ddebug
00349 LogDebug("HcalTBSim") << "HcalTB04Analysis:: A tertiary... PDG:"
00350 << partPDGEncoding << " TrackID:" <<trackID
00351 << " ParentID:" << parentID << " stable: "
00352 << isPDGStable << " Tau: " << pDGlifetime
00353 << " cTauGamma="
00354 << c_light*pDGlifetime*gammaFactor*1000.
00355 << "um GammaFactor: " << gammaFactor;
00356 #endif
00357 }
00358 }
00359 }
00360 }
00361 }
00362 }
00363 }
00364
00365 void HcalTB04Analysis::update(const EndOfEvent * evt) {
00366
00367 count++;
00368
00369
00370 LogDebug("HcalTBSim") << "HcalTB04Analysis::Fill event "
00371 << (*evt)()->GetEventID();
00372 fillBuffer (evt);
00373
00374
00375 LogDebug("HcalTBSim") << "HcalTB04Analysis::Do QIE analysis with "
00376 << hcalHitCache.size() << " hits";
00377 qieAnalysis();
00378
00379
00380 if (!hcalOnly) {
00381 LogDebug("HcalTBSim") << "HcalTB04Analysis::Do Xtal analysis with "
00382 << ecalHitCache.size() << " hits";
00383 xtalAnalysis();
00384 }
00385
00386
00387 LogDebug("HcalTBSim") << "HcalTB04Analysis::Final analysis";
00388 finalAnalysis();
00389
00390 int iEvt = (*evt)()->GetEventID();
00391 if (iEvt < 10)
00392 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis:: Event " << iEvt;
00393 else if ((iEvt < 100) && (iEvt%10 == 0))
00394 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis:: Event " << iEvt;
00395 else if ((iEvt < 1000) && (iEvt%100 == 0))
00396 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis:: Event " << iEvt;
00397 else if ((iEvt < 10000) && (iEvt%1000 == 0))
00398 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis:: Event " << iEvt;
00399 }
00400
00401 void HcalTB04Analysis::fillBuffer(const EndOfEvent * evt) {
00402
00403 std::vector<CaloHit> hhits, hhitl;
00404 int idHC, j;
00405 CaloG4HitCollection* theHC;
00406 std::map<int,float,std::less<int> > primaries;
00407 double etot1=0, etot2=0;
00408
00409
00410 G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
00411 std::string sdName = names[0];
00412 idHC = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
00413 theHC = (CaloG4HitCollection*) allHC->GetHC(idHC);
00414 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Hit Collection for " << sdName
00415 << " of ID " << idHC << " is obtained at " << theHC;
00416
00417 if (idHC >= 0 && theHC > 0) {
00418 hhits.reserve(theHC->entries());
00419 hhitl.reserve(theHC->entries());
00420 for (j = 0; j < theHC->entries(); j++) {
00421 CaloG4Hit* aHit = (*theHC)[j];
00422 double e = aHit->getEnergyDeposit()/GeV;
00423 double time = aHit->getTimeSlice();
00424 math::XYZPoint pos = aHit->getEntry();
00425 unsigned int id = aHit->getUnitID();
00426 double theta = pos.theta();
00427 double eta = -log(tan(theta * 0.5));
00428 double phi = pos.phi();
00429 int det, z, group, ieta, iphi, layer;
00430 HcalTestNumbering::unpackHcalIndex(id,det,z,group,ieta,iphi,layer);
00431 double jitter = time-timeOfFlight(det,layer,eta);
00432 if (jitter<0) jitter = 0;
00433 if (e < 0 || e > 1.) e = 0;
00434 double escl = e * scale(det,layer);
00435 unsigned int idx= HcalTBNumberingScheme::getUnitID(id,mode);
00436 CaloHit hit(det,layer,escl,eta,phi,jitter,idx);
00437 hhits.push_back(hit);
00438 CaloHit hitl(det,layer,escl,eta,phi,jitter,id);
00439 hhitl.push_back(hitl);
00440 primaries[aHit->getTrackID()]+= e;
00441 etot1 += escl;
00442 #ifdef ddebug
00443 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Hcal Hit i/p " << j
00444 << " ID 0x" << std::hex << id << " 0x" << idx
00445 << std::dec << " time " << std::setw(6) << time
00446 << " " << std::setw(6) << jitter << " theta "
00447 << std::setw(8) << theta << " eta " << std::setw(8)
00448 << eta << " phi " << std::setw(8) << phi << " e "
00449 << std::setw(8) << e << " " << std::setw(8) <<escl;
00450 #endif
00451 }
00452 }
00453
00454
00455 std::vector<CaloHit>::iterator itr;
00456 int nHit = hhits.size();
00457 std::vector<CaloHit*> hits(nHit);
00458 for (j = 0, itr = hhits.begin(); itr != hhits.end(); j++, itr++) {
00459 hits[j] = &hhits[j];
00460 }
00461 sort(hits.begin(),hits.end(),CaloHitIdMore());
00462 std::vector<CaloHit*>::iterator k1, k2;
00463 int nhit = 0;
00464 for (k1 = hits.begin(); k1 != hits.end(); k1++) {
00465 int det = (**k1).det();
00466 int layer = (**k1).layer();
00467 double ehit = (**k1).e();
00468 double eta = (**k1).eta();
00469 double phi = (**k1).phi();
00470 double jitter = (**k1).t();
00471 uint32_t unitID = (**k1).id();
00472 int jump = 0;
00473 for (k2 = k1+1; k2 != hits.end() && fabs(jitter-(**k2).t())<1 &&
00474 unitID==(**k2).id(); k2++) {
00475 ehit += (**k2).e();
00476 jump++;
00477 }
00478 nhit++;
00479 CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
00480 hcalHitCache.push_back(hit);
00481 etot2 += ehit;
00482 k1 += jump;
00483 #ifdef ddebug
00484 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Hcal Hit store " << nhit
00485 << " ID 0x" << std::hex << unitID << std::dec
00486 << " time " << std::setw(6) << jitter << " eta "
00487 << std::setw(8) << eta << " phi " << std::setw(8)
00488 << phi << " e " << std::setw(8) << ehit;
00489 #endif
00490 }
00491 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Stores " << nhit << " HCal hits"
00492 << " from " << nHit << " input hits E(Hcal) " << etot1
00493 << " " << etot2;
00494
00495
00496 for (j = 0, itr = hhitl.begin(); itr != hhitl.end(); j++, itr++) {
00497 hits[j] = &hhitl[j];
00498 }
00499 sort(hits.begin(),hits.end(),CaloHitIdMore());
00500 int nhitl = 0;
00501 double etotl = 0;
00502 for (k1 = hits.begin(); k1 != hits.end(); k1++) {
00503 int det = (**k1).det();
00504 int layer = (**k1).layer();
00505 double ehit = (**k1).e();
00506 double eta = (**k1).eta();
00507 double phi = (**k1).phi();
00508 double jitter = (**k1).t();
00509 uint32_t unitID = (**k1).id();
00510 int jump = 0;
00511 for (k2 = k1+1; k2 != hits.end() && fabs(jitter-(**k2).t())<1 &&
00512 unitID==(**k2).id(); k2++) {
00513 ehit += (**k2).e();
00514 jump++;
00515 }
00516 nhitl++;
00517 CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
00518 hcalHitLayer.push_back(hit);
00519 etotl += ehit;
00520 k1 += jump;
00521 #ifdef ddebug
00522 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Hcal Hit store " << nhitl
00523 << " ID 0x" << std::hex << unitID << std::dec
00524 << " time " << std::setw(6) << jitter << " eta "
00525 << std::setw(8) << eta << " phi " << std::setw(8)
00526 << phi << " e " << std::setw(8) << ehit;
00527 #endif
00528 }
00529 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Stores " << nhitl << " HCal "
00530 << "hits from " << nHit << " input hits E(Hcal) "
00531 << etot1 << " " << etotl;
00532
00533
00534 std::vector<CaloHit> ehits;
00535 sdName= names[1];
00536 idHC = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
00537 theHC = (CaloG4HitCollection*) allHC->GetHC(idHC);
00538 etot1 = etot2 = 0;
00539 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Hit Collection for " << sdName
00540 << " of ID " << idHC << " is obtained at " << theHC;
00541 if (idHC >= 0 && theHC > 0) {
00542 ehits.reserve(theHC->entries());
00543 for (j = 0; j < theHC->entries(); j++) {
00544 CaloG4Hit* aHit = (*theHC)[j];
00545 double e = aHit->getEnergyDeposit()/GeV;
00546 double time = aHit->getTimeSlice();
00547 math::XYZPoint pos = aHit->getEntry();
00548 unsigned int id = aHit->getUnitID();
00549 double theta = pos.theta();
00550 double eta = -log(tan(theta * 0.5));
00551 double phi = pos.phi();
00552 if (e < 0 || e > 100000.) e = 0;
00553 int det, z, group, ieta, iphi, layer;
00554 HcalTestNumbering::unpackHcalIndex(id,det,z,group,ieta,iphi,layer);
00555 CaloHit hit(det,0,e,eta,phi,time,id);
00556 ehits.push_back(hit);
00557 primaries[aHit->getTrackID()]+= e;
00558 etot1 += e;
00559 #ifdef ddebug
00560 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Ecal Hit i/p " << j
00561 << " ID 0x" << std::hex << id << std::dec
00562 << " time " << std::setw(6) << time << " theta "
00563 << std::setw(8) << theta << " eta " <<std::setw(8)
00564 << eta << " phi " << std::setw(8) << phi << " e "
00565 << std::setw(8) << e;
00566 #endif
00567 }
00568 }
00569
00570
00571 nHit = ehits.size();
00572 std::vector<CaloHit*> hite(nHit);
00573 for (j = 0, itr = ehits.begin(); itr != ehits.end(); j++, itr++) {
00574 hite[j] = &ehits[j];
00575 }
00576 sort(hite.begin(),hite.end(),CaloHitIdMore());
00577 nhit = 0;
00578 for (k1 = hite.begin(); k1 != hite.end(); k1++) {
00579 int det = (**k1).det();
00580 int layer = (**k1).layer();
00581 double ehit = (**k1).e();
00582 double eta = (**k1).eta();
00583 double phi = (**k1).phi();
00584 double jitter = (**k1).t();
00585 uint32_t unitID = (**k1).id();
00586 int jump = 0;
00587 for (k2 = k1+1; k2 != hite.end() && fabs(jitter-(**k2).t())<1 &&
00588 unitID==(**k2).id(); k2++) {
00589 ehit += (**k2).e();
00590 jump++;
00591 }
00592 nhit++;
00593 CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
00594 ecalHitCache.push_back(hit);
00595 etot2 += ehit;
00596 k1 += jump;
00597 #ifdef ddebug
00598 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Ecal Hit store " << nhit
00599 << " ID 0x" << std::hex << unitID << std::dec
00600 << " time " << std::setw(6) << jitter << " eta "
00601 << std::setw(8) << eta << " phi " << std::setw(8)
00602 << phi << " e " << std::setw(8) << ehit;
00603 #endif
00604 }
00605 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Stores " << nhit << " ECal hits"
00606 << " from " << nHit << " input hits E(Ecal) " << etot1
00607 << " " << etot2;
00608
00609
00610 nPrimary = (int)(primaries.size());
00611 int trackID = 0;
00612 G4PrimaryParticle* thePrim=0;
00613 int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
00614 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Event has " << nvertex
00615 << " verteices";
00616 if (nvertex<=0)
00617 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis::EndOfEvent ERROR: no "
00618 << "vertex found for event " << evNum;
00619
00620 for (int i = 0 ; i<nvertex; i++) {
00621 G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
00622 if (avertex == 0) {
00623 edm::LogInfo("HcalTBSim") << "HcalTB04Analysis::EndOfEvent ERR: pointer "
00624 << "to vertex = 0 for event " << evNum;
00625 } else {
00626 LogDebug("HcalTBSim") << "HcalTB04Analysis::Vertex number :" << i << " "
00627 << avertex->GetPosition();
00628 int npart = avertex->GetNumberOfParticle();
00629 if (npart == 0)
00630 edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::End Of Event ERR: "
00631 << "no primary!";
00632 if (thePrim==0) thePrim=avertex->GetPrimary(trackID);
00633 }
00634 }
00635
00636 if (thePrim != 0) {
00637 double px = thePrim->GetPx();
00638 double py = thePrim->GetPy();
00639 double pz = thePrim->GetPz();
00640 double p = std::sqrt(pow(px,2.)+pow(py,2.)+pow(pz,2.));
00641 pInit = p/GeV;
00642 if (p==0)
00643 edm::LogWarning("HcalTBSim") << "HcalTB04Analysis:: EndOfEvent ERR: "
00644 << "primary has p=0 ";
00645 else {
00646 double costheta = pz/p;
00647 double theta = acos(std::min(std::max(costheta,-1.),1.));
00648 etaInit = -log(tan(theta/2));
00649 if (px != 0 || py != 0) phiInit = atan2(py,px);
00650 }
00651 particleType = thePrim->GetPDGcode();
00652 } else
00653 edm::LogWarning("HcalTBSim") << "HcalTB04Analysis::EndOfEvent ERR: could "
00654 << "not find primary";
00655
00656 }
00657
00658 void HcalTB04Analysis::qieAnalysis() {
00659
00660 int hittot = hcalHitCache.size();
00661 if (hittot<=0) hittot = 1;
00662 std::vector<CaloHit> hits(hittot);
00663 std::vector<int> todo(nTower,0);
00664
00665 LogDebug("HcalTBSim") << "HcalTB04Analysis::qieAnalysis: Size "
00666 << hits.size() << " " << todo.size() << " "
00667 << idTower.size() << " " << esimh.size() << " "
00668 << eqie.size();
00669
00670 for (unsigned int k1 = 0; k1 < hcalHitCache.size(); k1++) {
00671 CaloHit hit = hcalHitCache[k1];
00672 uint32_t id = hit.id();
00673 int nhit = 0;
00674 double esim = hit.e();
00675 hits[nhit] = hit;
00676 for (unsigned int k2 = k1+1; k2 < hcalHitCache.size(); k2++) {
00677 hit = hcalHitCache[k2];
00678 if (hit.id() == id) {
00679 nhit++;
00680 hits[nhit] = hit;
00681 esim += hit.e();
00682 }
00683 }
00684 k1 += nhit;
00685 nhit++;
00686 std::vector<int> cd = myQie->getCode(nhit,hits);
00687 double eq = myQie->getEnergy(cd);
00688 LogDebug("HcalTBSim") << "HcalTB04Analysis:: ID 0x" << std::hex << id
00689 << std::dec << " registers " << esim << " energy "
00690 << "from " << nhit << " hits starting with hit # "
00691 << k1 << " energy with noise " << eq;
00692 for (int k2 = 0; k2 < nTower; k2++) {
00693 if (id == idTower[k2]) {
00694 todo[k2] = 1;
00695 esimh[k2] = esim;
00696 eqie[k2] = eq;
00697 }
00698 }
00699 }
00700
00701
00702 for (int k2 = 0; k2 < nTower; k2++) {
00703 if (todo[k2] == 0) {
00704 std::vector<int> cd = myQie->getCode(0,hits);
00705 double eq = myQie->getEnergy(cd);
00706 esimh[k2] = 0;
00707 eqie[k2] = eq;
00708 #ifdef ddebug
00709 LogDebug("HcalTBSim") << "HcalTB04Analysis:: ID 0x" << std::hex
00710 << idTower[k2] << std::dec << " registers "
00711 << esimh[k2] << " energy from hits and energy "
00712 << "after QIE analysis " << eqie[k2];
00713 #endif
00714 }
00715 }
00716 }
00717
00718 void HcalTB04Analysis::xtalAnalysis() {
00719
00720 edm::Service<edm::RandomNumberGenerator> rng;
00721 if ( ! rng.isAvailable()) {
00722 throw cms::Exception("Configuration")
00723 << "HcalTB04Analysis requires the RandomNumberGeneratorService\n"
00724 << "which is not present in the configuration file. "
00725 << "You must add the service\n in the configuration file or "
00726 << "remove the modules that require it.";
00727 }
00728 CLHEP::RandGaussQ randGauss(rng->getEngine());
00729
00730
00731 std::vector<int> iok(nCrystal,0);
00732 LogDebug("HcalTBSim") << "HcalTB04Analysis::xtalAnalysis: Size " <<iok.size()
00733 << " " << idEcal.size() << " " << esime.size() << " "
00734 << enois.size();
00735 for (unsigned int k1 = 0; k1 < ecalHitCache.size(); k1++) {
00736 uint32_t id = ecalHitCache[k1].id();
00737 int nhit = 0;
00738 double esim = ecalHitCache[k1].e();
00739 for (unsigned int k2 = k1+1; k2 < ecalHitCache.size(); k2++) {
00740 if (ecalHitCache[k2].id() == id) {
00741 nhit++;
00742 esim += ecalHitCache[k2].e();
00743 }
00744 }
00745 k1 += nhit;
00746 nhit++;
00747 double eq = esim + randGauss.fire(0., ecalNoise);
00748 #ifdef ddebug
00749 LogDebug("HcalTBSim") << "HcalTB04Analysis:: ID 0x" << std::hex << id
00750 << std::dec << " registers " << esim << " energy "
00751 << "from " << nhit << " hits starting with hit # "
00752 << k1 << " energy with noise " << eq;
00753 #endif
00754 for (int k2 = 0; k2 < nCrystal; k2++) {
00755 if (id == idEcal[k2]) {
00756 iok[k2] = 1;
00757 esime[k2] = esim;
00758 enois[k2] = eq;
00759 }
00760 }
00761 }
00762
00763
00764 for (int k2 = 0; k2 < nCrystal; k2++) {
00765 if (iok[k2] == 0) {
00766 esime[k2] = 0;
00767 enois[k2] = randGauss.fire(0., ecalNoise);
00768 #ifdef ddebug
00769 LogDebug("HcalTBSim") << "HcalTB04Analysis:: ID 0x" << std::hex
00770 << idEcal[k2] << std::dec << " registers "
00771 << esime[k2] << " energy from hits and energy from"
00772 << " noise " << enois[k2];
00773 #endif
00774 }
00775 }
00776 }
00777
00778 void HcalTB04Analysis::finalAnalysis() {
00779
00780
00781 histo->fillPrimary(pInit, etaInit, phiInit);
00782
00783
00784 eecals = ehcals = eecalq = ehcalq = 0.;
00785 for (int i=0; i<nTower; i++) {
00786 ehcals += esimh[i];
00787 ehcalq += eqie[i];
00788 }
00789 for (int i=0; i<nCrystal; i++) {
00790 eecals += esime[i];
00791 eecalq += enois[i];
00792 }
00793 etots = eecals + ehcals;
00794 etotq = eecalq + ehcalq;
00795 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Energy deposit at Sim Level "
00796 << "(Total) " << etots << " (ECal) " << eecals
00797 << " (HCal) " << ehcals << "\nHcalTB04Analysis:: "
00798 << "Energy deposit at Qie Level (Total) " << etotq
00799 << " (ECal) " << eecalq << " (HCal) " << ehcalq;
00800 histo->fillEdep(etots, eecals, ehcals, etotq, eecalq, ehcalq);
00801
00802
00803 for (int i=0; i<5; i++) {
00804 eseta[i] = 0.;
00805 eqeta[i] = 0.;
00806 }
00807 for (int i=0; i<3; i++) {
00808 esphi[i] = 0.;
00809 eqphi[i] = 0.;
00810 }
00811 double e1=0, e2=0;
00812 unsigned int id;
00813 for (int i=0; i<nTower; i++) {
00814 int det, z, group, ieta, iphi, layer;
00815 id = idTower[i];
00816 HcalTestNumbering::unpackHcalIndex(id,det,z,group,ieta,iphi,layer);
00817 iphi -= (icphi - 1);
00818 if (icphi > 4) {
00819 if (ieta == 0) ieta = 2;
00820 else ieta =-1;
00821 } else {
00822 ieta = ieta - iceta + 2;
00823 }
00824 if (iphi >= 0 && iphi < 3 && ieta >= 0 && ieta < 5) {
00825 eseta[ieta] += esimh[i];
00826 esphi[iphi] += esimh[i];
00827 e1 += esimh[i];
00828 eqeta[ieta] += eqie[i];
00829 eqphi[iphi] += eqie[i];
00830 e2 += eqie[i];
00831 }
00832 }
00833 for (int i=0; i<3; i++) {
00834 if (e1>0) esphi[i] /= e1;
00835 if (e2>0) eqphi[i] /= e2;
00836 }
00837 for (int i=0; i<5; i++) {
00838 if (e1>0) eseta[i] /= e1;
00839 if (e2>0) eqeta[i] /= e2;
00840 }
00841 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Energy fraction along Eta and"
00842 << " Phi (Sim/Qie)";
00843 for (int i=0; i<5; i++)
00844 LogDebug("HcalTBSim") << "HcalTB04Analysis:: [" << i << "] Eta Sim = "
00845 << eseta[i] << " Qie = " << eqeta[i] << " Phi Sim = "
00846 << esphi[i] << " Qie = " << eqphi[i];
00847 histo->fillTrnsProf(eseta,eqeta,esphi,eqphi);
00848
00849
00850 for (int i=0; i<20; i++) {
00851 eslay[i] = 0.;
00852 eqlay[i] = 0.;
00853 }
00854 e1=0; e2=0;
00855 for (int i=0; i<nTower; i++) {
00856 int det, z, group, ieta, iphi, layer;
00857 id = idTower[i];
00858 HcalTestNumbering::unpackHcalIndex(id,det,z,group,ieta,iphi,layer);
00859 iphi -= (icphi - 1);
00860 layer -= 1;
00861 if (iphi >= 0 && iphi < 3 && layer >= 0 && layer < 20) {
00862 eslay[layer] += esimh[i];
00863 e1 += esimh[i];
00864 eqlay[layer] += eqie[i];
00865 e2 += eqie[i];
00866 }
00867 }
00868 for (int i=0; i<20; i++) {
00869 if (e1>0) eslay[i] /= e1;
00870 if (e2>0) eqlay[i] /= e2;
00871 }
00872 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Energy fraction along Layer";
00873 for (int i=0; i<20; i++)
00874 LogDebug("HcalTBSim") << "HcalTB04Analysis:: [" << i << "] Sim = "
00875 << eslay[i] << " Qie = " << eqlay[i];
00876 histo->fillLongProf(eslay, eqlay);
00877 }
00878
00879
00880 void HcalTB04Analysis::fillEvent (PHcalTB04Info& product) {
00881
00882
00883 product.setIDs(idHcal, idXtal);
00884
00885
00886 product.setPrimary(nPrimary, particleType, pInit, etaInit, phiInit);
00887
00888
00889 product.setEdepHcal(esimh, eqie);
00890 product.setEdepHcal(esime, enois);
00891
00892
00893 product.setEdep(etots, eecals, ehcals, etotq, eecalq, ehcalq);
00894
00895
00896 product.setTrnsProf(eseta,eqeta,esphi,eqphi);
00897
00898
00899 product.setLongProf(eslay, eqlay);
00900
00901
00902 int i, nhit=0;
00903 std::vector<CaloHit>::iterator itr;
00904 for (i=0, itr=ecalHitCache.begin(); itr!=ecalHitCache.end(); i++,itr++) {
00905 uint32_t id = itr->id();
00906 int det, z, group, ieta, iphi, lay;
00907 HcalTestNumbering::unpackHcalIndex(id,det,z,group,ieta,iphi,lay);
00908 product.saveHit(det, lay, ieta, iphi, itr->e(), itr->t());
00909 nhit++;
00910 #ifdef debug
00911 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Save Hit " << std::setw(3)
00912 << i+1 << " ID 0x" << std::hex << group << std::dec
00913 << " " << std::setw(2) << det << " " << std::setw(2)
00914 << lay << " " << std::setw(1) << z << " "
00915 << std::setw(3) << ieta << " " << std::setw(3) <<iphi
00916 << " T " << std::setw(6) << itr->t() << " E "
00917 << std::setw(6) << itr->e();
00918 #endif
00919 }
00920 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Saves " << nhit
00921 << " hits from Crystals";
00922 int hit = nhit;
00923 nhit = 0;
00924
00925 for (i=hit, itr=hcalHitCache.begin(); itr!=hcalHitCache.end(); i++,itr++) {
00926 uint32_t id = itr->id();
00927 int det, z, group, ieta, iphi, lay;
00928 HcalTestNumbering::unpackHcalIndex(id,det,z,group,ieta,iphi,lay);
00929 product.saveHit(det, lay, ieta, iphi, itr->e(), itr->t());
00930 nhit++;
00931 #ifdef debug
00932 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Save Hit " << std::setw(3)
00933 << i+1 << " ID 0x" << std::hex << group << std::dec
00934 << " " << std::setw(2) << det << " " << std::setw(2)
00935 << lay << " " << std::setw(1) << z << " "
00936 << std::setw(3) << ieta << " " << std::setw(3) <<iphi
00937 << " T " << std::setw(6) << itr->t() << " E "
00938 << std::setw(6) << itr->e();
00939 #endif
00940 }
00941 LogDebug("HcalTBSim") << "HcalTB04Analysis:: Saves " << nhit
00942 << " hits from HCal";
00943
00944
00945 product.setVtxPrim(evNum, pvType, pvPosition.x(), pvPosition.y(),
00946 pvPosition.z(), pvUVW.x(), pvUVW.y(), pvUVW.z(),
00947 pvMomentum.x(), pvMomentum.y(), pvMomentum.z());
00948 for (unsigned int i = 0; i < secTrackID.size(); i++) {
00949 product.setVtxSec(secTrackID[i], secPartID[i], secMomentum[i].x(),
00950 secMomentum[i].y(), secMomentum[i].z(), secEkin[i]);
00951 }
00952 }
00953
00954 void HcalTB04Analysis::clear(){
00955 pvFound = false;
00956 pvType =-2;
00957 pvPosition = G4ThreeVector();
00958 pvMomentum = G4ThreeVector();
00959 pvUVW = G4ThreeVector();
00960 secTrackID.clear();
00961 secPartID.clear();
00962 secMomentum.clear();
00963 secEkin.clear();
00964 shortLivedSecondaries.clear();
00965
00966 ecalHitCache.erase(ecalHitCache.begin(), ecalHitCache.end());
00967 hcalHitCache.erase(hcalHitCache.begin(), hcalHitCache.end());
00968 hcalHitLayer.erase(hcalHitLayer.begin(), hcalHitLayer.end());
00969 nPrimary = particleType = 0;
00970 pInit = etaInit = phiInit = 0;
00971
00972 esimh.clear();
00973 eqie.clear();
00974 esimh.reserve(nTower);
00975 eqie.reserve(nTower);
00976 for (int i=0; i<nTower; i++) {
00977 esimh.push_back(0.);
00978 eqie.push_back(0.);
00979 }
00980 esime.clear();
00981 enois.clear();
00982 esime.reserve(nCrystal);
00983 enois.reserve(nCrystal);
00984 for (int i=0; i<nCrystal; i++) {
00985 esime.push_back(0.);
00986 enois.push_back(0.);
00987 }
00988 }
00989
00990 int HcalTB04Analysis::unitID(uint32_t id) {
00991
00992 int det, z, group, ieta, iphi, lay;
00993 HcalTestNumbering::unpackHcalIndex(id,det,z,group,ieta,iphi,lay);
00994 group = (det&15)<<20;
00995 group += ((lay-1)&31)<<15;
00996 group += (z&1)<<14;
00997 group += (ieta&127)<<7;
00998 group += (iphi&127);
00999 return group;
01000 }
01001
01002 double HcalTB04Analysis::scale(int det, int layer) {
01003
01004 double tmp = 1.;
01005 if (det == static_cast<int>(HcalBarrel)) {
01006 if (layer == 1) tmp = scaleHB0;
01007 else if (layer == 17) tmp = scaleHB16;
01008 else if (layer > 17) tmp = scaleHO;
01009 } else {
01010 if (layer <= 2) tmp = scaleHE0;
01011 }
01012 return tmp;
01013 }
01014
01015 double HcalTB04Analysis::timeOfFlight(int det, int layer, double eta) {
01016
01017 double theta = 2.0*atan(exp(-eta));
01018 double dist = beamOffset;
01019 if (det == static_cast<int>(HcalBarrel)) {
01020 const double rLay[19] = {
01021 1836.0, 1902.0, 1962.0, 2022.0, 2082.0, 2142.0, 2202.0, 2262.0, 2322.0,
01022 2382.0, 2448.0, 2514.0, 2580.0, 2646.0, 2712.0, 2776.0, 2862.5, 3847.0,
01023 4052.0};
01024 if (layer>0 && layer<=19) dist += rLay[layer-1]*mm/sin(theta);
01025 } else {
01026 const double zLay[19] = {
01027 4034.0, 4032.0, 4123.0, 4210.0, 4297.0, 4384.0, 4471.0, 4558.0, 4645.0,
01028 4732.0, 4819.0, 4906.0, 4993.0, 5080.0, 5167.0, 5254.0, 5341.0, 5428.0,
01029 5515.0};
01030 if (layer>0 && layer<=19) dist += zLay[layer-1]*mm/cos(theta);
01031 }
01032
01033 double tmp = dist/c_light/ns;
01034 #ifdef ddebug
01035 LogDebug("HcalTBSim") << "HcalTB04Analysis::timeOfFlight " << tmp
01036 << " for det/lay " << det << " " << layer
01037 << " eta/theta " << eta << " " << theta/deg
01038 << " dist " << dist;
01039 #endif
01040 return tmp;
01041 }
01042
01043 DEFINE_SIMWATCHER (HcalTB04Analysis);