00001 #include "SimG4Core/Notification/interface/BeginOfJob.h"
00002 #include "SimG4Core/Notification/interface/BeginOfRun.h"
00003 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
00004 #include "SimG4Core/Notification/interface/EndOfEvent.h"
00005
00006
00007 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00008 #include "SimG4CMS/Calo/interface/HcalTestAnalysis.h"
00009 #include "SimG4CMS/Calo/interface/HCalSD.h"
00010 #include "SimG4CMS/Calo/interface/CaloG4Hit.h"
00011 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
00012 #include "DataFormats/Math/interface/Point3D.h"
00013
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00017 #include "DetectorDescription/Core/interface/DDCompactView.h"
00018
00019 #include "G4SDManager.hh"
00020 #include "G4Step.hh"
00021 #include "G4Track.hh"
00022 #include "G4VProcess.hh"
00023 #include "G4HCofThisEvent.hh"
00024 #include "CLHEP/Units/SystemOfUnits.h"
00025 #include "CLHEP/Units/PhysicalConstants.h"
00026
00027 #include <cmath>
00028 #include <iostream>
00029 #include <iomanip>
00030
00031 HcalTestAnalysis::HcalTestAnalysis(const edm::ParameterSet &p):
00032 myqie(0), addTower(3), tuplesManager(0), tuples(0), numberingFromDDD(0), org(0) {
00033
00034 edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("HcalTestAnalysis");
00035 eta0 = m_Anal.getParameter<double>("Eta0");
00036 phi0 = m_Anal.getParameter<double>("Phi0");
00037 int laygroup = m_Anal.getParameter<int>("LayerGrouping");
00038 centralTower = m_Anal.getParameter<int>("CentralTower");
00039 names = m_Anal.getParameter<std::vector<std::string> >("Names");
00040 fileName = m_Anal.getParameter<std::string>("FileName");
00041
00042 edm::LogInfo("HcalSim") << "HcalTestAnalysis:: Initialised as observer of "
00043 << "begin/end events and of G4step";
00044
00045 count = 0;
00046 group_ = layerGrouping(laygroup);
00047 nGroup = 0;
00048 for (unsigned int i=0; i<group_.size(); i++)
00049 if (group_[i]>nGroup) nGroup = group_[i];
00050 tower_ = towersToAdd(centralTower, addTower);
00051 nTower = tower_.size()/2;
00052
00053 edm::LogInfo("HcalSim") << "HcalTestAnalysis:: initialised for " << nGroup
00054 << " Longitudinal groups and " << nTower
00055 << " towers";
00056
00057
00058 myqie = new HcalQie(p);
00059 }
00060
00061 HcalTestAnalysis::~HcalTestAnalysis() {
00062 edm::LogInfo("HcalSim") << "HcalTestAnalysis: --------> Total number of "
00063 << "selected entries : " << count;
00064 edm::LogInfo("HcalSim") << "HcalTestAnalysis: Pointers:: HcalQie " << myqie
00065 << ", HistoClass " << tuples << ", Numbering Scheme "
00066 << org << " and FromDDD " << numberingFromDDD;
00067 if (myqie) {
00068 edm::LogInfo("HcalSim") << "HcalTestAnalysis: Delete HcalQie";
00069 delete myqie;
00070 }
00071 if (numberingFromDDD) {
00072 edm::LogInfo("HcalSim") << "HcalTestAnalysis: Delete HcalNumberingFromDDD";
00073 delete numberingFromDDD;
00074 }
00075 }
00076
00077 std::vector<int> HcalTestAnalysis::layerGrouping(int group) {
00078
00079 std::vector<int> temp(19);
00080 if (group <= 1) {
00081 int grp[19] = {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 2, 2};
00082 for (int i=0; i<19; i++)
00083 temp[i] = grp[i];
00084 } else if (group == 2) {
00085 int grp[19] = {1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19};
00086 for (int i=0; i<19; i++)
00087 temp[i] = grp[i];
00088 } else if (group == 3) {
00089 int grp[19] = {1, 1, 2, 2, 2, 3, 3, 3, 4, 4, 4, 5, 5, 5, 6, 6, 6, 7, 7};
00090 for (int i=0; i<19; i++)
00091 temp[i] = grp[i];
00092 } else if (group == 4) {
00093 int grp[19] = {1, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 5, 6, 6, 6, 6, 6, 7, 7};
00094 for (int i=0; i<19; i++)
00095 temp[i] = grp[i];
00096 } else {
00097 int grp[19] = {1, 1, 2, 3, 3, 4, 4, 5, 5, 5, 6, 6, 6, 6, 6, 6, 6, 7, 7};
00098 for (int i=0; i<19; i++)
00099 temp[i] = grp[i];
00100 }
00101
00102 edm::LogInfo("HcalSim") << "HcalTestAnalysis:: Layer Grouping ";
00103 for (int i=0; i<19; i++)
00104 edm::LogInfo("HcalSim") << "HcalTestAnalysis: Group[" << i << "] = "
00105 << temp[i];
00106 return temp;
00107 }
00108
00109 std::vector<int> HcalTestAnalysis::towersToAdd(int centre, int nadd) {
00110
00111 int etac = (centre/100)%100;
00112 int phic = (centre%100);
00113 int etamin, etamax, phimin, phimax;
00114 if (etac>0) {
00115 etamin = etac-nadd;
00116 etamax = etac+nadd;
00117 } else {
00118 etamin = etac;
00119 etamax = etac;
00120 }
00121 if (phic>0) {
00122 phimin = phic-nadd;
00123 phimax = phic+nadd;
00124 } else {
00125 phimin = phic;
00126 phimax = phic;
00127 }
00128
00129 int nbuf, kount=0;
00130 nbuf = (etamax-etamin+1)*(phimax-phimin+1);
00131 std::vector<int> temp(2*nbuf);
00132 for (int eta=etamin; eta<=etamax; eta++) {
00133 for (int phi=phimin; phi<=phimax; phi++) {
00134 temp[kount] = (eta*100 + phi);
00135 temp[kount+nbuf] = std::max(abs(eta-etac),abs(phi-phic));
00136 kount++;
00137 }
00138 }
00139
00140 edm::LogInfo("HcalSim") << "HcalTestAnalysis:: Towers to be considered for"
00141 << " Central " << centre << " and " << nadd
00142 << " on either side";
00143 for (int i=0; i<nbuf; i++)
00144 edm::LogInfo("HcalSim") << "HcalTestAnalysis: Tower[" << std::setw(3) << i
00145 << "] " << temp[i] << " " << temp[nbuf+i];
00146 return temp;
00147 }
00148
00149
00150 void HcalTestAnalysis::update(const BeginOfJob * job) {
00151
00152
00153 edm::ESHandle<DDCompactView> pDD;
00154 (*job)()->get<IdealGeometryRecord>().get(pDD);
00155 edm::LogInfo("HcalSim") << "HcalTestAnalysis:: Initialise "
00156 << "HcalNumberingFromDDD for " << names[0];
00157 numberingFromDDD = new HcalNumberingFromDDD(names[0], (*pDD));
00158
00159
00160 tuplesManager.reset(new HcalTestHistoManager(fileName));
00161
00162
00163 org = new HcalTestNumberingScheme(false);
00164
00165 }
00166
00167
00168 void HcalTestAnalysis::update(const BeginOfRun * run) {
00169
00170 int irun = (*run)()->GetRunID();
00171 edm::LogInfo("HcalSim") << "HcalTestAnalysis:: Begin of Run = " << irun;
00172
00173 bool loop = true, eta = true, phi = true;
00174 int etac = (centralTower/100)%100;
00175 if (etac == 0) {
00176 etac = 1;
00177 eta = false;
00178 }
00179 int phic = (centralTower%100);
00180 if (phic == 0) {
00181 phic = 1;
00182 phi = false;
00183 }
00184 int idet = static_cast<int>(HcalBarrel);
00185 while (loop) {
00186 HcalCellType::HcalCell tmp = numberingFromDDD->cell(idet,1,1,etac,phic);
00187 if (tmp.ok) {
00188 if (eta) eta0 = tmp.eta;
00189 if (phi) phi0 = tmp.phi;
00190 loop = false;
00191 } else if (idet == static_cast<int>(HcalBarrel)) {
00192 idet = static_cast<int>(HcalEndcap);
00193 } else if (idet == static_cast<int>(HcalEndcap)) {
00194 idet = static_cast<int>(HcalForward);
00195 } else {
00196 loop = false;
00197 }
00198 }
00199
00200 edm::LogInfo("HcalSim") << "HcalTestAnalysis:: Central Tower "
00201 << centralTower << " corresponds to eta0 = " << eta0
00202 << " phi0 = " << phi0;
00203
00204 std::string sdname = names[0];
00205 G4SDManager* sd = G4SDManager::GetSDMpointerIfExist();
00206 if (sd != 0) {
00207 G4VSensitiveDetector* aSD = sd->FindSensitiveDetector(sdname);
00208 if (aSD==0) {
00209 edm::LogWarning("HcalSim") << "HcalTestAnalysis::beginOfRun: No SD with "
00210 << "name " << sdname << " in this Setup";
00211 } else {
00212 HCalSD* theCaloSD = dynamic_cast<HCalSD*>(aSD);
00213 edm::LogInfo("HcalSim") << "HcalTestAnalysis::beginOfRun: Finds SD with "
00214 << "name " << theCaloSD->GetName()
00215 << " in this Setup";
00216 if (org) {
00217 theCaloSD->setNumberingScheme(org);
00218 edm::LogInfo("HcalSim") << "HcalTestAnalysis::beginOfRun: set a new "
00219 << "numbering scheme";
00220 }
00221 }
00222 } else {
00223 edm::LogWarning("HcalSim") << "HcalTestAnalysis::beginOfRun: Could not get"
00224 << " SD Manager!";
00225 }
00226
00227 }
00228
00229
00230 void HcalTestAnalysis::update(const BeginOfEvent * evt) {
00231
00232
00233 tuples = new HcalTestHistoClass();
00234
00235 tuples->setCounters();
00236
00237 int i = 0;
00238 edepEB = edepEE = edepHB = edepHE = edepHO = 0.;
00239 for (i = 0; i < 20; i++) edepl[i] = 0.;
00240 for (i = 0; i < 20; i++) mudist[i] = -1.;
00241
00242 int iev = (*evt)()->GetEventID();
00243 LogDebug("HcalSim") <<"HcalTestAnalysis: Begin of event = " << iev;
00244 }
00245
00246
00247 void HcalTestAnalysis::update(const G4Step * aStep) {
00248
00249 if (aStep != NULL) {
00250 G4VPhysicalVolume* curPV = aStep->GetPreStepPoint()->GetPhysicalVolume();
00251 G4String name = curPV->GetName();
00252 name.assign(name,0,3);
00253 double edeposit = aStep->GetTotalEnergyDeposit();
00254 int layer=-1;
00255 if (name == "EBR") {
00256 edepEB += edeposit;
00257 } else if (name == "EFR") {
00258 edepEE += edeposit;
00259 } else if (name == "HBS") {
00260 layer = (curPV->GetCopyNo()/10)%100;
00261 if (layer >= 0 && layer < 17) {
00262 edepHB += edeposit;
00263 } else {
00264 edm::LogWarning("HcalSim") << "HcalTestAnalysis::Error in HB "
00265 << curPV->GetName() << curPV->GetCopyNo();
00266 layer = -1;
00267 }
00268 } else if (name == "HES") {
00269 layer = (curPV->GetCopyNo()/10)%100;
00270 if (layer >= 0 && layer < 19) {
00271 edepHE += edeposit;
00272 } else {
00273 edm::LogWarning("HcalSim") << "HcalTestAnalysis::Error in HE "
00274 << curPV->GetName() << curPV->GetCopyNo();
00275 layer = -1;
00276 }
00277 } else if (name == "HTS") {
00278 layer = (curPV->GetCopyNo()/10)%100;
00279 if (layer >= 17 && layer < 20) {
00280 edepHO += edeposit;
00281 } else {
00282 edm::LogWarning("HcalSim") << "HcalTestAnalysis::Error in HO "
00283 << curPV->GetName() << curPV->GetCopyNo();
00284 layer = -1;
00285 }
00286 }
00287 if (layer >= 0 && layer < 20) {
00288 edepl[layer] += edeposit;
00289
00290
00291 G4String part = aStep->GetTrack()->GetDefinition()->GetParticleName();
00292 if ((part == "mu-" || part == "mu+") && mudist[layer] < 0) {
00293 math::XYZPoint pos(aStep->GetPreStepPoint()->GetPosition().x(),
00294 aStep->GetPreStepPoint()->GetPosition().y(),
00295 aStep->GetPreStepPoint()->GetPosition().z());
00296 double theta = pos.theta();
00297 double eta = -log(tan(theta * 0.5));
00298 double phi = pos.phi();
00299 double dist = sqrt ((eta-eta0)*(eta-eta0) + (phi-phi0)*(phi-phi0));
00300 mudist[layer] = dist*std::sqrt(pos.perp2());
00301 }
00302 }
00303
00304 if (layer >= 0 && layer < 20) {
00305 LogDebug("HcalSim") << "HcalTestAnalysis:: G4Step: " << name << " Layer "
00306 << std::setw(3) << layer << " Edep " << std::setw(6)
00307 << edeposit/MeV << " MeV";
00308 }
00309 } else {
00310 edm::LogInfo("HcalSim") << "HcalTestAnalysis:: G4Step: Null Step";
00311 }
00312 }
00313
00314
00315 void HcalTestAnalysis::update(const EndOfEvent * evt) {
00316
00317 count++;
00318
00319 fill(evt);
00320 LogDebug("HcalSim") << "HcalTestAnalysis:: --- after Fill";
00321
00322
00323 qieAnalysis();
00324 LogDebug("HcalSim") << "HcalTestAnalysis:: --- after QieAnalysis";
00325
00326
00327 layerAnalysis();
00328 LogDebug("HcalSim") << "HcalTestAnalysis:: --- after LayerAnalysis";
00329
00330
00331 tuplesManager->fillTree(tuples);
00332 tuples = 0;
00333 LogDebug("HcalSim") << "HcalTestAnalysis:: --- after fillTree";
00334
00335 }
00336
00337
00338 void HcalTestAnalysis::fill(const EndOfEvent * evt) {
00339
00340 LogDebug("HcalSim") << "HcalTestAnalysis: Fill event "
00341 << (*evt)()->GetEventID();
00342
00343
00344 G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
00345
00346 int nhc = 0, neb = 0, nef = 0, j = 0;
00347 caloHitCache.erase (caloHitCache.begin(), caloHitCache.end());
00348
00349
00350 int HCHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names[0]);
00351 CaloG4HitCollection* theHCHC = (CaloG4HitCollection*) allHC->GetHC(HCHCid);
00352 LogDebug("HcalSim") << "HcalTestAnalysis :: Hit Collection for " << names[0]
00353 << " of ID " << HCHCid << " is obtained at " << theHCHC;
00354 if (HCHCid >= 0 && theHCHC > 0) {
00355 for (j = 0; j < theHCHC->entries(); j++) {
00356
00357 CaloG4Hit* aHit = (*theHCHC)[j];
00358
00359 double e = aHit->getEnergyDeposit()/GeV;
00360 double time = aHit->getTimeSlice();
00361
00362 math::XYZPoint pos = aHit->getPosition();
00363 double theta = pos.theta();
00364 double eta = -log(tan(theta * 0.5));
00365 double phi = pos.phi();
00366
00367 uint32_t unitID = aHit->getUnitID();
00368 int subdet, zside, layer, etaIndex, phiIndex, lay;
00369 org->unpackHcalIndex(unitID,subdet,zside,layer,etaIndex,phiIndex,lay);
00370 double jitter = time-timeOfFlight(subdet,lay,eta);
00371 if (jitter<0) jitter = 0;
00372 CaloHit hit(subdet,lay,e,eta,phi,jitter,unitID);
00373 caloHitCache.push_back(hit);
00374 nhc++;
00375
00376 std::string det = "HB";
00377 if (subdet == static_cast<int>(HcalForward)) {
00378 det = "HF";
00379 } else if (subdet == static_cast<int>(HcalEndcap)) {
00380 if (etaIndex <= 20) {
00381 det = "HES";
00382 } else {
00383 det = "HED";
00384 }
00385 }
00386
00387 LogDebug("HcalSim") << "HcalTest: " << det << " layer " << std::setw(2)
00388 << layer << " time " << std::setw(6) << time
00389 << " theta " << std::setw(8) << theta << " eta "
00390 << std::setw(8) << eta << " phi " << std::setw(8)
00391 << phi << " e " << std::setw(8) << e;
00392 }
00393 }
00394 LogDebug("HcalSim") << "HcalTestAnalysis::HCAL hits : " << nhc << "\n";
00395
00396
00397 int EBHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names[1]);
00398 CaloG4HitCollection* theEBHC = (CaloG4HitCollection*) allHC->GetHC(EBHCid);
00399 LogDebug("HcalSim") << "HcalTestAnalysis :: Hit Collection for " << names[1]
00400 << " of ID " << EBHCid << " is obtained at " << theEBHC;
00401 if (EBHCid >= 0 && theEBHC > 0) {
00402 for (j = 0; j < theEBHC->entries(); j++) {
00403
00404 CaloG4Hit* aHit = (*theEBHC)[j];
00405
00406 double e = aHit->getEnergyDeposit()/GeV;
00407 double time = aHit->getTimeSlice();
00408 std::string det = "EB";
00409
00410 math::XYZPoint pos = aHit->getPosition();
00411 double theta = pos.theta();
00412 double eta = -log(tan(theta/2.));
00413 double phi = pos.phi();
00414
00415 HcalNumberingFromDDD::HcalID id = numberingFromDDD->unitID(eta,phi,1,1);
00416 uint32_t unitID = org->getUnitID(id);
00417 int subdet, zside, layer, ieta, iphi, lay;
00418 org->unpackHcalIndex(unitID,subdet,zside,layer,ieta,iphi,lay);
00419 subdet = 10;
00420 layer = 0;
00421 unitID = org->packHcalIndex(subdet,zside,layer,ieta,iphi,lay);
00422 CaloHit hit(subdet,lay,e,eta,phi,time,unitID);
00423 caloHitCache.push_back(hit);
00424 neb++;
00425 LogDebug("HcalSim") << "HcalTest: " << det << " layer " << std::setw(2)
00426 << layer << " time " << std::setw(6) << time
00427 << " theta " << std::setw(8) << theta << " eta "
00428 << std::setw(8) << eta << " phi " << std::setw(8)
00429 << phi << " e " << std::setw(8) << e;
00430 }
00431 }
00432 LogDebug("HcalSim") << "HcalTestAnalysis::EB hits : " << neb << "\n";
00433
00434
00435 int EEHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names[2]);
00436 CaloG4HitCollection* theEEHC = (CaloG4HitCollection*) allHC->GetHC(EEHCid);
00437 LogDebug("HcalSim") << "HcalTestAnalysis :: Hit Collection for " << names[2]
00438 << " of ID " << EEHCid << " is obtained at " << theEEHC;
00439 if (EEHCid >= 0 && theEEHC > 0) {
00440 for (j = 0; j < theEEHC->entries(); j++) {
00441
00442 CaloG4Hit* aHit = (*theEEHC)[j];
00443
00444 double e = aHit->getEnergyDeposit()/GeV;
00445 double time = aHit->getTimeSlice();
00446 std::string det = "EE";
00447
00448 math::XYZPoint pos = aHit->getPosition();
00449 double theta = pos.theta();
00450 double eta = -log(tan(theta/2.));
00451 double phi = pos.phi();
00452
00453 HcalNumberingFromDDD::HcalID id = numberingFromDDD->unitID(eta,phi,1,1);
00454 uint32_t unitID = org->getUnitID(id);
00455 int subdet, zside, layer, ieta, iphi, lay;
00456 org->unpackHcalIndex(unitID,subdet,zside,layer,ieta,iphi,lay);
00457 subdet = 11;
00458 layer = 0;
00459 unitID = org->packHcalIndex(subdet,zside,layer,ieta,iphi,lay);
00460 CaloHit hit(subdet,lay,e,eta,phi,time,unitID);
00461 caloHitCache.push_back(hit);
00462 nef++;
00463 LogDebug("HcalSim") << "HcalTest: " << det << " layer " << std::setw(2)
00464 << layer << " time " << std::setw(6) << time
00465 << " theta " << std::setw(8) << theta << " eta "
00466 << std::setw(8) << eta << " phi " << std::setw(8)
00467 << phi << " e " << std::setw(8) << e;
00468 }
00469 }
00470 LogDebug("HcalSim") << "HcalTestAnalysis::EE hits : " << nef << "\n";
00471 }
00472
00473
00474 void HcalTestAnalysis::qieAnalysis() {
00475
00476
00477 int hittot = caloHitCache.size();
00478 tuples->fillHits(caloHitCache);
00479
00480
00481 HcalNumberingFromDDD::HcalID id = numberingFromDDD->unitID(eta0,phi0,1,1);
00482 uint32_t unitID = org->getUnitID(id);
00483 int subdet, zside, layer, ieta, iphi, lay;
00484 org->unpackHcalIndex(unitID,subdet,zside,layer,ieta,iphi,lay);
00485 int laymax = 0;
00486 std::string det = "Unknown";
00487 if (subdet == static_cast<int>(HcalBarrel)) {
00488 laymax = 4; det = "HB";
00489 } else if (subdet == static_cast<int>(HcalEndcap)) {
00490 laymax = 2; det = "HES";
00491 }
00492 LogDebug("HcalSim") << "HcalTestAnalysis::Qie: " << det << " Eta " << ieta
00493 << " Phi " << iphi << " Laymax " << laymax << " Hits "
00494 << hittot;
00495
00496 if (laymax>0 && hittot>0) {
00497 std::vector<CaloHit> hits(hittot);
00498 std::vector<double> eqielay(80,0.0), esimlay(80,0.0), esimtot(4,0.0);
00499 std::vector<double> eqietow(200,0.0), esimtow(200,0.0), eqietot(4,0.0);
00500 int etac = (centralTower/100)%100;
00501 int phic = (centralTower%100);
00502
00503 for (int layr=0; layr<nGroup; layr++) {
00504 int layx, layy=20;
00505 for (int i=0; i<20; i++)
00506 if (group_[i] == layr+1 && i < layy) layy = i+1;
00507 if (subdet == static_cast<int>(HcalBarrel)) {
00508 if (layy < 2) layx = 0;
00509 else if (layy < 17) layx = 1;
00510 else if (layy == 17) layx = 2;
00511 else layx = 3;
00512 } else {
00513 if (layy < 2) layx = 0;
00514 else layx = 1;
00515 }
00516
00517 for (int it=0; it<nTower; it++) {
00518 int nhit = 0;
00519 double esim = 0;
00520 for (int k1 = 0; k1 < hittot; k1++) {
00521 CaloHit hit = caloHitCache[k1];
00522 int subdetc = hit.det();
00523 int layer = hit.layer();
00524 int group = 0;
00525 if (layer > 0 && layer < 20) group = group_[layer];
00526 if (subdetc == subdet && group == layr+1) {
00527 int zsidec, ietac, iphic, idx;
00528 unitID = hit.id();
00529 org->unpackHcalIndex(unitID,subdetc,zsidec,layer,ietac,iphic,lay);
00530 if (etac > 0 && phic > 0) {
00531 idx = ietac*100 + iphic;
00532 } else if (etac > 0) {
00533 idx = ietac*100;
00534 } else if (phic > 0) {
00535 idx = iphic;
00536 } else {
00537 idx = 0;
00538 }
00539 if (zsidec==zside && idx==tower_[it]) {
00540 hits[nhit] = hit;
00541 LogDebug("HcalSim") << "HcalTest: Hit " << nhit << " " << hit;
00542 nhit++;
00543 esim += hit.e();
00544 }
00545 }
00546 }
00547
00548 std::vector<int> cd = myqie->getCode(nhit,hits);
00549 double eqie = myqie->getEnergy(cd);
00550
00551 LogDebug("HcalSim") << "HcalTestAnalysis::Qie: Energy in layer "
00552 << layr << " Sim " << esim << " After QIE " <<eqie;
00553 for (int i=0; i<4; i++) {
00554 if (tower_[nTower+it] <= i) {
00555 esimtot[i] += esim;
00556 eqietot[i] += eqie;
00557 esimlay[20*i+layr] += esim;
00558 eqielay[20*i+layr] += eqie;
00559 esimtow[50*i+it] += esim;
00560 eqietow[50*i+it] += eqie;
00561 }
00562 }
00563 }
00564 }
00565 LogDebug("HcalSim") << "HcalTestAnalysis::Qie: Total energy " << esimtot[3]
00566 << " (SimHit) " << eqietot[3] << " (After QIE)";
00567
00568 std::vector<double> latphi(10);
00569 int nt = 2*addTower + 1;
00570 for (int it=0; it<nt; it++)
00571 latphi[it] = it-addTower;
00572 for (int i=0; i<4; i++) {
00573 double scals=1, scalq=1;
00574 std::vector<double> latfs(10,0.), latfq(10,0.), longs(20), longq(20);
00575 if (esimtot[i]>0) scals = 1./esimtot[i];
00576 if (eqietot[i]>0) scalq = 1./eqietot[i];
00577 for (int it=0; it<nTower; it++) {
00578 int phib = it%nt;
00579 latfs[phib] += scals*esimtow[50*i+it];
00580 latfq[phib] += scalq*eqietow[50*i+it];
00581 }
00582 for (int layr=0; layr<=nGroup; layr++) {
00583 longs[layr] = scals*esimlay[20*i+layr];
00584 longq[layr] = scalq*eqielay[20*i+layr];
00585 }
00586 tuples->fillQie(i,esimtot[i],eqietot[i],nGroup,longs,longq,
00587 nt,latphi,latfs,latfq);
00588 }
00589 }
00590 }
00591
00592
00593 void HcalTestAnalysis::layerAnalysis(){
00594
00595 int i = 0;
00596 LogDebug("HcalSim") << "\n ===>>> HcalTestAnalysis: Energy deposit in MeV "
00597 << "\n at EB : " << std::setw(6) << edepEB/MeV
00598 << "\n at EE : " << std::setw(6) << edepEE/MeV
00599 << "\n at HB : " << std::setw(6) << edepHB/MeV
00600 << "\n at HE : " << std::setw(6) << edepHE/MeV
00601 << "\n at HO : " << std::setw(6) << edepHO/MeV
00602 << "\n ---- HcalTestAnalysis: Energy deposit in Layers";
00603 for (i = 0; i < 20; i++)
00604 LogDebug("HcalSim") << " Layer " << std::setw(2) << i << " E "
00605 << std::setw(8) << edepl[i]/MeV << " MeV";
00606
00607 tuples->fillLayers(edepl, edepHO, edepHB+edepHE, mudist);
00608 }
00609
00610
00611
00612 double HcalTestAnalysis::timeOfFlight(int det, int layer, double eta) {
00613
00614 double theta = 2.0*atan(exp(-eta));
00615 double dist = 0.;
00616 if (det == static_cast<int>(HcalBarrel)) {
00617 const double rLay[19] = {
00618 1836.0, 1902.0, 1962.0, 2022.0, 2082.0, 2142.0, 2202.0, 2262.0, 2322.0,
00619 2382.0, 2448.0, 2514.0, 2580.0, 2646.0, 2712.0, 2776.0, 2862.5, 3847.0,
00620 4052.0};
00621 if (layer>0 && layer<20) dist += rLay[layer-1]*mm/sin(theta);
00622 } else {
00623 const double zLay[19] = {
00624 4034.0, 4032.0, 4123.0, 4210.0, 4297.0, 4384.0, 4471.0, 4558.0, 4645.0,
00625 4732.0, 4819.0, 4906.0, 4993.0, 5080.0, 5167.0, 5254.0, 5341.0, 5428.0,
00626 5515.0};
00627 if (layer>0 && layer<20) dist += zLay[layer-1]*mm/cos(theta);
00628 }
00629 double tmp = dist/c_light/ns;
00630 LogDebug("HcalSim") << "HcalTestAnalysis::timeOfFlight " << tmp
00631 << " for det/lay " << det << " " << layer
00632 << " eta/theta " << eta << " " << theta/deg << " dist "
00633 << dist;
00634 return tmp;
00635 }