File indexing completed on 2023-10-25 10:06:34
0001
0002
0003
0004
0005 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
0006 #include "SimG4Core/Notification/interface/BeginOfJob.h"
0007 #include "SimG4Core/Notification/interface/BeginOfRun.h"
0008 #include "SimG4Core/Notification/interface/EndOfEvent.h"
0009 #include "SimG4Core/Notification/interface/Observer.h"
0010 #include "SimG4Core/Watcher/interface/SimProducer.h"
0011
0012
0013 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0014 #include "DataFormats/Math/interface/GeantUnits.h"
0015 #include "DataFormats/Math/interface/Point3D.h"
0016 #include "SimDataFormats/CaloHit/interface/CaloHit.h"
0017 #include "SimDataFormats/ValidationFormats/interface/PValidationFormats.h"
0018
0019 #include "SimG4CMS/Calo/interface/CaloG4Hit.h"
0020 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
0021 #include "SimG4CMS/Calo/interface/HCalSD.h"
0022 #include "SimG4CMS/Calo/interface/HcalTestNumberingScheme.h"
0023
0024 #include "FWCore/Framework/interface/Event.h"
0025 #include "FWCore/Framework/interface/EventSetup.h"
0026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0028
0029 #include "Geometry/HcalCommonData/interface/HcalDDDSimConstants.h"
0030 #include "Geometry/HcalCommonData/interface/HcalNumberingFromDDD.h"
0031 #include "Geometry/Records/interface/HcalSimNumberingRecord.h"
0032
0033 #include "Validation/HcalHits/interface/SimG4HcalHitCluster.h"
0034 #include "Validation/HcalHits/interface/SimG4HcalHitJetFinder.h"
0035
0036 #include "G4HCofThisEvent.hh"
0037 #include "G4SDManager.hh"
0038 #include "G4Step.hh"
0039 #include "G4Track.hh"
0040 #include "G4VProcess.hh"
0041
0042 #include <iostream>
0043 #include <memory>
0044 #include <string>
0045 #include <vector>
0046
0047 using namespace geant_units::operators;
0048
0049 class SimG4HcalValidation : public SimProducer,
0050 public Observer<const BeginOfRun *>,
0051 public Observer<const BeginOfEvent *>,
0052 public Observer<const EndOfEvent *>,
0053 public Observer<const G4Step *> {
0054 public:
0055 SimG4HcalValidation(const edm::ParameterSet &p);
0056 SimG4HcalValidation(const SimG4HcalValidation &) = delete;
0057 const SimG4HcalValidation &operator=(const SimG4HcalValidation &) = delete;
0058 ~SimG4HcalValidation() override;
0059
0060 void registerConsumes(edm::ConsumesCollector) override;
0061 void produce(edm::Event &, const edm::EventSetup &) override;
0062 void beginRun(edm::EventSetup const &) override;
0063
0064 private:
0065 void init();
0066
0067
0068 void update(const BeginOfRun *run) override;
0069 void update(const BeginOfEvent *evt) override;
0070 void update(const G4Step *step) override;
0071 void update(const EndOfEvent *evt) override;
0072
0073
0074 void fill(const EndOfEvent *ev);
0075 void layerAnalysis(PHcalValidInfoLayer &);
0076 void nxNAnalysis(PHcalValidInfoNxN &);
0077 void jetAnalysis(PHcalValidInfoJets &);
0078 void fetchHits(PHcalValidInfoLayer &);
0079 void clear();
0080 void collectEnergyRdir(const double, const double);
0081 double getHcalScale(std::string, int) const;
0082
0083 private:
0084 edm::ESGetToken<HcalDDDSimConstants, HcalSimNumberingRecord> ddconsToken_;
0085
0086
0087 std::unique_ptr<SimG4HcalHitJetFinder> jetf;
0088
0089
0090 std::unique_ptr<HcalNumberingFromDDD> numberingFromDDD;
0091
0092
0093 std::unique_ptr<HcalTestNumberingScheme> org;
0094
0095
0096 std::vector<CaloHit> hitcache;
0097
0098
0099 std::vector<float> scaleHB;
0100 std::vector<float> scaleHE;
0101 std::vector<float> scaleHF;
0102
0103
0104 std::vector<std::string> names;
0105 double coneSize, ehitThreshold, hhitThreshold;
0106 float timeLowlim, timeUplim, eta0, phi0, jetThreshold;
0107 bool applySampling, hcalOnly;
0108 int infolevel;
0109 std::string labelLayer, labelNxN, labelJets;
0110
0111
0112 std::vector<double> dEta;
0113 std::vector<double> dPhi;
0114
0115
0116 unsigned int count;
0117 double edepEB, edepEE, edepHB, edepHE, edepHO;
0118 double edepd[5], edepl[20];
0119 double een, hen, hoen;
0120 double vhitec, vhithc, enEcal, enHcal;
0121 };
0122
0123 SimG4HcalValidation::SimG4HcalValidation(const edm::ParameterSet &p) {
0124 edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("SimG4HcalValidation");
0125 infolevel = m_Anal.getParameter<int>("InfoLevel");
0126 hcalOnly = m_Anal.getParameter<bool>("HcalClusterOnly");
0127 applySampling = m_Anal.getParameter<bool>("HcalSampling");
0128 coneSize = m_Anal.getParameter<double>("ConeSize");
0129 ehitThreshold = m_Anal.getParameter<double>("EcalHitThreshold");
0130 hhitThreshold = m_Anal.getParameter<double>("HcalHitThreshold");
0131 timeLowlim = m_Anal.getParameter<double>("TimeLowLimit");
0132 timeUplim = m_Anal.getParameter<double>("TimeUpLimit");
0133 jetThreshold = m_Anal.getParameter<double>("JetThreshold");
0134 eta0 = m_Anal.getParameter<double>("Eta0");
0135 phi0 = m_Anal.getParameter<double>("Phi0");
0136 names = m_Anal.getParameter<std::vector<std::string>>("Names");
0137 labelLayer = m_Anal.getUntrackedParameter<std::string>("LabelLayerInfo", "HcalInfoLayer");
0138 labelNxN = m_Anal.getUntrackedParameter<std::string>("LabelNxNInfo", "HcalInfoNxN");
0139 labelJets = m_Anal.getUntrackedParameter<std::string>("LabelJetsInfo", "HcalInfoJets");
0140
0141 produces<PHcalValidInfoLayer>(labelLayer);
0142 if (infolevel > 0)
0143 produces<PHcalValidInfoNxN>(labelNxN);
0144 if (infolevel > 1)
0145 produces<PHcalValidInfoJets>(labelJets);
0146
0147 edm::LogVerbatim("ValidHcal") << "HcalTestAnalysis:: Initialized as observer of begin/end events and "
0148 << "of G4step with Parameter values: \n\tInfoLevel = " << infolevel
0149 << "\n\thcalOnly = " << hcalOnly << "\n\tapplySampling = " << applySampling
0150 << "\n\tconeSize = " << coneSize << "\n\tehitThreshold = " << ehitThreshold
0151 << "\n\thhitThreshold = " << hhitThreshold << "\n\tttimeLowlim = " << timeLowlim
0152 << "\n\tttimeUplim = " << timeUplim << "\n\teta0 = " << eta0
0153 << "\n\tphi0 = " << phi0 << "\nLabels (Layer): " << labelLayer
0154 << " (NxN): " << labelNxN << " (Jets): " << labelJets;
0155
0156 init();
0157 }
0158
0159 SimG4HcalValidation::~SimG4HcalValidation() {
0160 edm::LogVerbatim("ValidHcal") << "\n --------> Total number of selected entries"
0161 << " : " << count;
0162 }
0163
0164 void SimG4HcalValidation::registerConsumes(edm::ConsumesCollector cc) {
0165 ddconsToken_ = cc.esConsumes<HcalDDDSimConstants, HcalSimNumberingRecord, edm::Transition::BeginRun>();
0166 edm::LogVerbatim("ValidHcal") << "SimG4HcalValidation::Initialize ESGetToken for HcalDDDSimConstants";
0167 }
0168
0169 void SimG4HcalValidation::produce(edm::Event &e, const edm::EventSetup &) {
0170 std::unique_ptr<PHcalValidInfoLayer> productLayer(new PHcalValidInfoLayer);
0171 layerAnalysis(*productLayer);
0172 e.put(std::move(productLayer), labelLayer);
0173
0174 if (infolevel > 0) {
0175 std::unique_ptr<PHcalValidInfoNxN> productNxN(new PHcalValidInfoNxN);
0176 nxNAnalysis(*productNxN);
0177 e.put(std::move(productNxN), labelNxN);
0178 }
0179
0180 if (infolevel > 1) {
0181 std::unique_ptr<PHcalValidInfoJets> productJets(new PHcalValidInfoJets);
0182 jetAnalysis(*productJets);
0183 e.put(std::move(productJets), labelJets);
0184 }
0185 }
0186
0187
0188 void SimG4HcalValidation::init() {
0189 float sHB[4] = {117., 117., 178., 217.};
0190 float sHE[3] = {178., 178., 178.};
0191 float sHF[3] = {2.84, 2.09, 0.};
0192
0193 float deta[4] = {0.0435, 0.1305, 0.2175, 0.3045};
0194 float dphi[4] = {0.0436, 0.1309, 0.2182, 0.3054};
0195
0196 int i = 0;
0197 for (i = 0; i < 4; i++) {
0198 scaleHB.push_back(sHB[i]);
0199 }
0200 for (i = 0; i < 3; i++) {
0201 scaleHE.push_back(sHE[i]);
0202 }
0203 for (i = 0; i < 3; i++) {
0204 scaleHF.push_back(sHF[i]);
0205 }
0206
0207
0208 for (i = 0; i < 4; i++) {
0209 dEta.push_back(deta[i]);
0210 dPhi.push_back(dphi[i]);
0211 }
0212
0213
0214 jetf = std::make_unique<SimG4HcalHitJetFinder>(coneSize);
0215
0216
0217 count = 0;
0218 }
0219
0220 void SimG4HcalValidation::beginRun(edm::EventSetup const &es) {
0221
0222 const HcalDDDSimConstants *hcons = &es.getData(ddconsToken_);
0223 edm::LogVerbatim("ValidHcal") << "HcalTestAnalysis:: Initialise HcalNumberingFromDDD";
0224 numberingFromDDD = std::make_unique<HcalNumberingFromDDD>(hcons);
0225
0226
0227 org = std::make_unique<HcalTestNumberingScheme>(false);
0228 }
0229
0230 void SimG4HcalValidation::update(const BeginOfRun *run) {
0231 int irun = (*run)()->GetRunID();
0232
0233 edm::LogVerbatim("ValidHcal") << " =====> Begin of Run = " << irun;
0234
0235 std::string sdname = names[0];
0236 G4SDManager *sd = G4SDManager::GetSDMpointerIfExist();
0237 if (sd != nullptr) {
0238 G4VSensitiveDetector *aSD = sd->FindSensitiveDetector(sdname);
0239 if (aSD == nullptr) {
0240 edm::LogWarning("ValidHcal") << "SimG4HcalValidation::beginOfRun: No SD"
0241 << " with name " << sdname << " in this "
0242 << "Setup";
0243 } else {
0244 HCalSD *theCaloSD = dynamic_cast<HCalSD *>(aSD);
0245 edm::LogVerbatim("ValidHcal") << "SimG4HcalValidation::beginOfRun: Finds SD with name " << theCaloSD->GetName()
0246 << " in this Setup";
0247 if (org.get()) {
0248 theCaloSD->setNumberingScheme(org.get());
0249 edm::LogVerbatim("ValidHcal") << "SimG4HcalValidation::beginOfRun: set a new numbering scheme";
0250 }
0251 }
0252 } else {
0253 edm::LogWarning("ValidHcal") << "SimG4HcalValidation::beginOfRun: Could "
0254 << "not get SD Manager!";
0255 }
0256 }
0257
0258
0259 void SimG4HcalValidation::update(const BeginOfEvent *evt) {
0260 int i = 0;
0261 edepEB = edepEE = edepHB = edepHE = edepHO = 0.;
0262 for (i = 0; i < 5; i++)
0263 edepd[i] = 0.;
0264 for (i = 0; i < 20; i++)
0265 edepl[i] = 0.;
0266 vhitec = vhithc = enEcal = enHcal = 0;
0267
0268 clear();
0269
0270 int iev = (*evt)()->GetEventID();
0271 LogDebug("ValidHcal") << "SimG4HcalValidation: =====> Begin of event = " << iev;
0272 }
0273
0274
0275 void SimG4HcalValidation::update(const G4Step *aStep) {
0276 if (aStep != nullptr) {
0277 G4VPhysicalVolume *curPV = aStep->GetPreStepPoint()->GetPhysicalVolume();
0278 G4String name = curPV->GetName();
0279 name.assign(name, 0, 3);
0280 double edeposit = aStep->GetTotalEnergyDeposit();
0281 int layer = -1, depth = -1;
0282 if (name == "EBR") {
0283 depth = 0;
0284 edepEB += edeposit;
0285 } else if (name == "EFR") {
0286 depth = 0;
0287 edepEE += edeposit;
0288 } else if (name == "HBS") {
0289 layer = (curPV->GetCopyNo() / 10) % 100;
0290 depth = (curPV->GetCopyNo()) % 10 + 1;
0291 if (depth > 0 && depth < 4 && layer >= 0 && layer < 17) {
0292 edepHB += edeposit;
0293 } else {
0294 edm::LogWarning("ValidHcal") << "SimG4HcalValidation:Error " << curPV->GetName() << curPV->GetCopyNo();
0295 depth = -1;
0296 layer = -1;
0297 }
0298 } else if (name == "HES") {
0299 layer = (curPV->GetCopyNo() / 10) % 100;
0300 depth = (curPV->GetCopyNo()) % 10 + 1;
0301 if (depth > 0 && depth < 3 && layer >= 0 && layer < 19) {
0302 edepHE += edeposit;
0303 } else {
0304 edm::LogWarning("ValidHcal") << "SimG4HcalValidation:Error " << curPV->GetName() << curPV->GetCopyNo();
0305 depth = -1;
0306 layer = -1;
0307 }
0308 } else if (name == "HTS") {
0309 layer = (curPV->GetCopyNo() / 10) % 100;
0310 depth = (curPV->GetCopyNo()) % 10 + 1;
0311 if (depth > 3 && depth < 5 && layer >= 17 && layer < 20) {
0312 edepHO += edeposit;
0313 } else {
0314 edm::LogWarning("ValidHcal") << "SimG4HcalValidation:Error " << curPV->GetName() << curPV->GetCopyNo();
0315 depth = -1;
0316 layer = -1;
0317 }
0318 }
0319 if (depth >= 0 && depth < 5)
0320 edepd[depth] += edeposit;
0321 if (layer >= 0 && layer < 20)
0322 edepl[layer] += edeposit;
0323
0324 if (layer >= 0 && layer < 20) {
0325 LogDebug("ValidHcal") << "SimG4HcalValidation:: G4Step: " << name << " Layer " << std::setw(3) << layer
0326 << " Depth " << std::setw(2) << depth << " Edep " << std::setw(6) << edeposit / MeV
0327 << " MeV";
0328 }
0329 }
0330 }
0331
0332
0333 void SimG4HcalValidation::update(const EndOfEvent *evt) {
0334 count++;
0335
0336
0337 fill(evt);
0338 LogDebug("ValidHcal") << "SimG4HcalValidation:: --- after Fill ";
0339 }
0340
0341
0342 void SimG4HcalValidation::fill(const EndOfEvent *evt) {
0343 LogDebug("ValidHcal") << "SimG4HcalValidation:Fill event " << (*evt)()->GetEventID();
0344
0345
0346 G4HCofThisEvent *allHC = (*evt)()->GetHCofThisEvent();
0347
0348 int nhc = 0, j = 0;
0349
0350
0351 int HCHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names[0]);
0352 CaloG4HitCollection *theHCHC = (CaloG4HitCollection *)allHC->GetHC(HCHCid);
0353 LogDebug("ValidHcal") << "SimG4HcalValidation :: Hit Collection for " << names[0] << " of ID " << HCHCid
0354 << " is obtained at " << theHCHC;
0355 if (HCHCid >= 0 && theHCHC != nullptr) {
0356 int nhits = theHCHC->entries();
0357 for (j = 0; j < nhits; j++) {
0358 CaloG4Hit *aHit = (*theHCHC)[j];
0359
0360 double e = aHit->getEnergyDeposit() / GeV;
0361 double time = aHit->getTimeSlice();
0362
0363 math::XYZPoint pos = aHit->getPosition();
0364 double theta = pos.theta();
0365 double eta = -log(tan(theta * 0.5));
0366 double phi = pos.phi();
0367
0368 uint32_t unitID = aHit->getUnitID();
0369 int subdet, zside, layer, etaIndex, phiIndex, lay;
0370 org->unpackHcalIndex(unitID, subdet, zside, layer, etaIndex, phiIndex, lay);
0371
0372
0373 layer--;
0374 std::string det = "HB";
0375 if (subdet == static_cast<int>(HcalForward)) {
0376 det = "HF";
0377 uint16_t depth = aHit->getDepth();
0378 if (depth != 0) {
0379 layer += 2;
0380 lay += 2;
0381 }
0382 if (layer == 1)
0383 vhithc += e;
0384 else if (layer == 0)
0385 vhitec += e;
0386 else
0387 edm::LogVerbatim("ValidHcal") << "SimG4HcalValidation::HitPMT " << subdet << " " << (2 * zside - 1) * etaIndex
0388 << " " << phiIndex << " " << layer + 1 << " R " << pos.rho() << " Phi "
0389 << convertRadToDeg(phi) << " Edep " << e << " Time " << time;
0390 } else if (subdet == static_cast<int>(HcalEndcap)) {
0391 if (etaIndex <= 20) {
0392 det = "HES";
0393 } else {
0394 det = "HED";
0395 }
0396 }
0397 LogDebug("ValidHcal") << "SimG4HcalValidation::debugFill Hcal " << det << " layer " << std::setw(2) << layer
0398 << " lay " << std::setw(2) << lay << " time " << std::setw(6) << time << " theta "
0399 << std::setw(8) << theta << " eta " << std::setw(8) << eta << " phi " << std::setw(8) << phi
0400 << " e " << std::setw(8) << e << " ID 0x" << std::hex << unitID << " ID dec " << std::dec
0401 << (int)unitID;
0402
0403
0404 if (applySampling)
0405 e *= getHcalScale(det, layer);
0406
0407
0408 if (time >= timeLowlim && time <= timeUplim && e > hhitThreshold) {
0409 enHcal += e;
0410 CaloHit ahit(subdet, lay, e, eta, phi, time, unitID);
0411 hitcache.push_back(ahit);
0412 ++nhc;
0413 }
0414 }
0415 }
0416 LogDebug("ValidHcal") << "SimG4HcalValidation:: HCAL hits : " << nhc;
0417
0418 if (!hcalOnly) {
0419 int ndets = names.size();
0420 if (ndets > 3)
0421 ndets = 3;
0422 for (int idty = 1; idty < ndets; idty++) {
0423 std::string det = "EB";
0424 if (idty == 2)
0425 det = "EE";
0426 else if (idty > 2)
0427 det = "ES";
0428
0429 int nec = 0;
0430 int ECHCid = G4SDManager::GetSDMpointer()->GetCollectionID(names[idty]);
0431 CaloG4HitCollection *theECHC = (CaloG4HitCollection *)allHC->GetHC(ECHCid);
0432 LogDebug("ValidHcal") << "SimG4HcalValidation:: Hit Collection for " << names[idty] << " of ID " << ECHCid
0433 << " is obtained at " << theECHC;
0434 int theechc_entries = theECHC->entries();
0435 if (ECHCid >= 0 && theECHC != nullptr) {
0436 for (j = 0; j < theechc_entries; j++) {
0437 CaloG4Hit *aHit = (*theECHC)[j];
0438
0439 double e = aHit->getEnergyDeposit() / GeV;
0440 double time = aHit->getTimeSlice();
0441
0442 math::XYZPoint pos = aHit->getPosition();
0443 double theta = pos.theta();
0444 double eta = -log(tan(theta / 2.));
0445 double phi = pos.phi();
0446 HcalNumberingFromDDD::HcalID id = numberingFromDDD->unitID(eta, phi, 1, 1);
0447 uint32_t unitID = org->getUnitID(id);
0448 int subdet, zside, layer, ieta, iphi, lay;
0449 org->unpackHcalIndex(unitID, subdet, zside, layer, ieta, iphi, lay);
0450 subdet = idty + 9;
0451 layer = 0;
0452 unitID = org->packHcalIndex(subdet, zside, layer, ieta, iphi, lay);
0453
0454
0455 if (time >= timeLowlim && time <= timeUplim && e > ehitThreshold) {
0456 enEcal += e;
0457 CaloHit ahit(subdet, lay, e, eta, phi, time, unitID);
0458 hitcache.push_back(ahit);
0459 ++nec;
0460 }
0461
0462 LogDebug("ValidHcal") << "SimG4HcalValidation::debugFill Ecal " << det << " layer " << std::setw(2) << layer
0463 << " lay " << std::setw(2) << lay << " time " << std::setw(6) << time << " theta "
0464 << std::setw(8) << theta << " eta " << std::setw(8) << eta << " phi " << std::setw(8)
0465 << phi << " e " << std::setw(8) << e << " ID 0x" << std::hex << unitID << " ID dec "
0466 << std::dec << (int)unitID;
0467 }
0468 }
0469
0470 LogDebug("ValidHcal") << "SimG4HcalValidation:: " << det << " hits : " << nec;
0471 }
0472 }
0473 }
0474
0475 void SimG4HcalValidation::layerAnalysis(PHcalValidInfoLayer &product) {
0476 int i = 0;
0477 LogDebug("ValidHcal") << "\n ===>>> SimG4HcalValidation: Energy deposit "
0478 << "in MeV\n at EB : " << std::setw(6) << edepEB / MeV << "\n at EE : " << std::setw(6)
0479 << edepEE / MeV << "\n at HB : " << std::setw(6) << edepHB / MeV
0480 << "\n at HE : " << std::setw(6) << edepHE / MeV << "\n at HO : " << std::setw(6)
0481 << edepHO / MeV << "\n ---- SimG4HcalValidation: Energy deposit in";
0482 for (i = 0; i < 5; i++)
0483 LogDebug("ValidHcal") << " Depth " << std::setw(2) << i << " E " << std::setw(8) << edepd[i] / MeV << " MeV";
0484 LogDebug("ValidHcal") << " ---- SimG4HcalValidation: Energy deposit in"
0485 << "layers";
0486 for (i = 0; i < 20; i++)
0487 LogDebug("ValidHcal") << " Layer " << std::setw(2) << i << " E " << std::setw(8) << edepl[i] / MeV << " MeV";
0488
0489 product.fillLayers(edepl, edepd, edepHO, edepHB + edepHE, edepEB + edepEE);
0490
0491
0492 product.fillHF(vhitec, vhithc, enEcal, enHcal);
0493 LogDebug("ValidHcal") << "SimG4HcalValidation::HF hits " << vhitec << " in EC and " << vhithc << " in HC\n"
0494 << " HB/HE " << enEcal << " in EC and " << enHcal << " in HC";
0495
0496
0497 fetchHits(product);
0498
0499 LogDebug("ValidHcal") << " layerAnalysis ===> after fetchHits";
0500 }
0501
0502
0503 void SimG4HcalValidation::nxNAnalysis(PHcalValidInfoNxN &product) {
0504 std::vector<CaloHit> *hits = &hitcache;
0505 std::vector<CaloHit>::iterator hit_itr;
0506
0507 LogDebug("ValidHcal") << "SimG4HcalValidation::NxNAnalysis : entrance ";
0508
0509 collectEnergyRdir(eta0, phi0);
0510
0511
0512 LogDebug("ValidHcal") << " NxNAnalysis : coolectEnergyRdir - Ecal " << een << " Hcal " << hen;
0513
0514 double etot = 0.;
0515 double ee = 0.;
0516 double he = 0.;
0517 double hoe = 0.;
0518
0519 int max = dEta.size();
0520
0521 for (hit_itr = hits->begin(); hit_itr < hits->end(); hit_itr++) {
0522 double e = hit_itr->e();
0523 double t = hit_itr->t();
0524 double eta = hit_itr->eta();
0525 double phi = hit_itr->phi();
0526 int type = hit_itr->det();
0527 int layer = hit_itr->layer();
0528
0529
0530
0531 if (fabs(eta0 - eta) <= dEta[max - 1] && fabs(phi0 - phi) <= dPhi[max - 1]) {
0532 etot += e;
0533 if (type == 10 || type == 11 || type == 12)
0534 ee += e;
0535 if (type == static_cast<int>(HcalBarrel) || type == static_cast<int>(HcalEndcap) ||
0536 type == static_cast<int>(HcalForward)) {
0537 he += e;
0538 if (type == static_cast<int>(HcalBarrel) && layer > 17)
0539 hoe += e;
0540
0541
0542 for (int i = 0; i < max; i++) {
0543 if ((eta0 - eta) <= dEta[i] && (eta0 - eta) >= -dEta[i] && (phi0 - phi) <= dPhi[i] &&
0544 (phi0 - phi) >= -dPhi[i]) {
0545 LogDebug("ValidHcal") << "SimG4HcalValidation:: nxNAnalysis eta0,"
0546 << " phi0 = " << eta0 << " " << phi0 << " type, layer = " << type << "," << layer
0547 << " eta, phi = " << eta << " " << phi;
0548
0549 product.fillTProfileNxN(e, i, t);
0550 break;
0551 }
0552 }
0553
0554 }
0555 }
0556 }
0557
0558 product.fillEcollectNxN(ee, he, hoe, etot);
0559 product.fillHvsE(een, hen, hoen, een + hen);
0560
0561 LogDebug("ValidHcal") << " nxNAnalysis ===> after fillHvsE";
0562 }
0563
0564
0565 void SimG4HcalValidation::jetAnalysis(PHcalValidInfoJets &product) {
0566 std::vector<CaloHit> *hhit = &hitcache;
0567
0568 jetf->setInput(hhit);
0569
0570
0571 std::vector<SimG4HcalHitCluster> *result = jetf->getClusters(hcalOnly);
0572
0573 std::vector<SimG4HcalHitCluster>::iterator clus_itr;
0574
0575 LogDebug("ValidHcal") << "\n ---------- Final list of " << (*result).size() << " clusters ---------------";
0576 for (clus_itr = result->begin(); clus_itr < result->end(); clus_itr++)
0577 LogDebug("ValidHcal") << (*clus_itr);
0578
0579 std::vector<double> enevec, etavec, phivec;
0580
0581 if (!(*result).empty()) {
0582 sort((*result).begin(), (*result).end());
0583
0584 clus_itr = result->begin();
0585 double etac = clus_itr->eta();
0586 double phic = clus_itr->phi();
0587
0588 double ecal_collect = 0.;
0589 if (!hcalOnly) {
0590 ecal_collect = clus_itr->collectEcalEnergyR();
0591 } else {
0592 collectEnergyRdir(etac, phic);
0593 ecal_collect = een;
0594 }
0595 LogDebug("ValidHcal") << " JetAnalysis ===> ecal_collect = " << ecal_collect;
0596
0597
0598 double dist = jetf->rDist(eta0, phi0, etac, phic);
0599 LogDebug("ValidHcal") << " JetAnalysis ===> eta phi deviation = " << dist;
0600 product.fillEtaPhiProfileJet(eta0, phi0, etac, phic, dist);
0601
0602 LogDebug("ValidHcal") << " JetAnalysis ===> after fillEtaPhiProfileJet";
0603
0604 std::vector<CaloHit> *hits = clus_itr->getHits();
0605 std::vector<CaloHit>::iterator hit_itr;
0606
0607 double ee = 0., he = 0., hoe = 0., etot = 0.;
0608
0609
0610 for (hit_itr = hits->begin(); hit_itr < hits->end(); hit_itr++) {
0611 double e = hit_itr->e();
0612 double t = hit_itr->t();
0613 double r = jetf->rDist(&(*clus_itr), &(*hit_itr));
0614
0615
0616 etot += e;
0617 if (hit_itr->det() == 10 || hit_itr->det() == 11 || hit_itr->det() == 12)
0618 ee += e;
0619 if (hit_itr->det() == static_cast<int>(HcalBarrel) || hit_itr->det() == static_cast<int>(HcalEndcap) ||
0620 hit_itr->det() == static_cast<int>(HcalForward)) {
0621 he += e;
0622 if (hit_itr->det() == static_cast<int>(HcalBarrel) && hit_itr->layer() > 17)
0623 hoe += e;
0624 }
0625
0626 if (hit_itr->det() == static_cast<int>(HcalBarrel) || hit_itr->det() == static_cast<int>(HcalEndcap) ||
0627 hit_itr->det() == static_cast<int>(HcalForward)) {
0628 product.fillTProfileJet(he, r, t);
0629 }
0630 }
0631
0632 product.fillEcollectJet(ee, he, hoe, etot);
0633
0634 LogDebug("ValidHcal") << " JetAnalysis ===> after fillEcollectJet: "
0635 << "ee/he/hoe/etot " << ee << "/" << he << "/" << hoe << "/" << etot;
0636
0637
0638 for (clus_itr = result->begin(); clus_itr < result->end(); clus_itr++) {
0639 if (clus_itr->e() > jetThreshold) {
0640 enevec.push_back(clus_itr->e());
0641 etavec.push_back(clus_itr->eta());
0642 phivec.push_back(clus_itr->phi());
0643 }
0644 }
0645 product.fillJets(enevec, etavec, phivec);
0646
0647 LogDebug("ValidHcal") << " JetAnalysis ===> after fillJets\n"
0648 << " JetAnalysis ===> (*result).size() " << (*result).size();
0649
0650
0651 if (etavec.size() > 1) {
0652 if (etavec[0] > -2.5 && etavec[0] < 2.5 && etavec[1] > -2.5 && etavec[1] < 2.5) {
0653 LogDebug("ValidHcal") << " JetAnalysis ===> Di-jet mass enter\n"
0654 << " JetAnalysis ===> Di-jet vectors ";
0655 for (unsigned int i = 0; i < enevec.size(); i++)
0656 LogDebug("ValidHcal") << " e, eta, phi = " << enevec[i] << " " << etavec[i] << " " << phivec[i];
0657
0658 double et0 = enevec[0] / cosh(etavec[0]);
0659 double px0 = et0 * cos(phivec[0]);
0660 double py0 = et0 * sin(phivec[0]);
0661 double pz0 = et0 * sinh(etavec[0]);
0662 double et1 = enevec[1] / cosh(etavec[1]);
0663 double px1 = et1 * cos(phivec[1]);
0664 double py1 = et1 * sin(phivec[1]);
0665 double pz1 = et1 * sinh(etavec[1]);
0666
0667 double dijetmass2 = (enevec[0] + enevec[1]) * (enevec[0] + enevec[1]) - (px1 + px0) * (px1 + px0) -
0668 (py1 + py0) * (py1 + py0) - (pz1 + pz0) * (pz1 + pz0);
0669
0670 LogDebug("ValidHcal") << " JetAnalysis ===> Di-jet massSQ " << dijetmass2;
0671
0672 double dijetmass;
0673 if (dijetmass2 >= 0.)
0674 dijetmass = sqrt(dijetmass2);
0675 else
0676 dijetmass = -sqrt(-1. * dijetmass2);
0677
0678 product.fillDiJets(dijetmass);
0679
0680 LogDebug("ValidHcal") << " JetAnalysis ===> after fillDiJets";
0681 }
0682 }
0683 }
0684 }
0685
0686
0687 void SimG4HcalValidation::fetchHits(PHcalValidInfoLayer &product) {
0688 LogDebug("ValidHcal") << "Enter SimG4HcalValidation::fetchHits with " << hitcache.size() << " hits";
0689 int nHit = hitcache.size();
0690 int hit = 0;
0691 int i;
0692 std::vector<CaloHit>::iterator itr;
0693 std::vector<CaloHit *> lhits(nHit);
0694 for (i = 0, itr = hitcache.begin(); itr != hitcache.end(); i++, itr++) {
0695 uint32_t unitID = itr->id();
0696 int subdet, zside, group, ieta, iphi, lay;
0697 HcalTestNumbering::unpackHcalIndex(unitID, subdet, zside, group, ieta, iphi, lay);
0698 subdet = itr->det();
0699 lay = itr->layer();
0700 group = (subdet & 15) << 20;
0701 group += ((lay - 1) & 31) << 15;
0702 group += (zside & 1) << 14;
0703 group += (ieta & 127) << 7;
0704 group += (iphi & 127);
0705 itr->setId(group);
0706 lhits[i] = &hitcache[i];
0707 LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits:Original " << i << " " << hitcache[i] << "\n"
0708 << "SimG4HcalValidation::fetchHits:Copied " << i << " " << *lhits[i];
0709 }
0710 sort(lhits.begin(), lhits.end(), CaloHitIdMore());
0711 std::vector<CaloHit *>::iterator k1, k2;
0712 for (i = 0, k1 = lhits.begin(); k1 != lhits.end(); i++, k1++)
0713 LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits:Sorted " << i << " " << **k1;
0714 int nHits = 0;
0715 for (i = 0, k1 = lhits.begin(); k1 != lhits.end(); i++, k1++) {
0716 double ehit = (**k1).e();
0717 double t = (**k1).t();
0718 uint32_t unitID = (**k1).id();
0719 int jump = 0;
0720 LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits:Start " << i << " U/T/E"
0721 << " 0x" << std::hex << unitID << std::dec << " " << t << " " << ehit;
0722 for (k2 = k1 + 1; k2 != lhits.end() && (t - (**k2).t()) < 1 && (t - (**k2).t()) > -1 && unitID == (**k2).id();
0723 k2++) {
0724 ehit += (**k2).e();
0725 LogDebug("ValidHcal") << "\t + " << (**k2).e();
0726 jump++;
0727 }
0728 LogDebug("ValidHcal") << "\t = " << ehit << " in " << jump;
0729
0730 double eta = (*k1)->eta();
0731 double phi = (*k1)->phi();
0732 int lay = ((unitID >> 15) & 31) + 1;
0733 int subdet = (unitID >> 20) & 15;
0734 int zside = (unitID >> 14) & 1;
0735 int ieta = (unitID >> 7) & 127;
0736 int iphi = (unitID)&127;
0737
0738
0739 product.fillHits(nHits, lay, subdet, eta, phi, ehit, t);
0740 nHits++;
0741
0742 LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits:Hit " << nHits << " " << i << " ID 0x" << std::hex
0743 << unitID << " det " << std::dec << subdet << " " << lay << " " << zside << " " << ieta
0744 << " " << iphi << " Time " << t << " E " << ehit;
0745
0746 i += jump;
0747 k1 += jump;
0748 }
0749
0750 LogDebug("ValidHcal") << "SimG4HcalValidation::fetchHits called with " << nHit << " hits"
0751 << " and writes out " << nHits << '(' << hit << ") hits";
0752 }
0753
0754 void SimG4HcalValidation::clear() { hitcache.erase(hitcache.begin(), hitcache.end()); }
0755
0756
0757 void SimG4HcalValidation::collectEnergyRdir(const double eta0, const double phi0) {
0758 std::vector<CaloHit> *hits = &hitcache;
0759 std::vector<CaloHit>::iterator hit_itr;
0760
0761 double sume = 0., sumh = 0., sumho = 0.;
0762
0763 for (hit_itr = hits->begin(); hit_itr < hits->end(); hit_itr++) {
0764 double e = hit_itr->e();
0765 double eta = hit_itr->eta();
0766 double phi = hit_itr->phi();
0767
0768 int type = hit_itr->det();
0769
0770 double r = jetf->rDist(eta0, phi0, eta, phi);
0771 if (r < coneSize) {
0772 if (type == 10 || type == 11 || type == 12) {
0773 sume += e;
0774 } else {
0775 sumh += e;
0776 if (type == static_cast<int>(HcalBarrel) && hit_itr->layer() > 17)
0777 sumho += e;
0778 }
0779 }
0780 }
0781
0782 een = sume;
0783 hen = sumh;
0784 hoen = sumho;
0785 }
0786
0787
0788 double SimG4HcalValidation::getHcalScale(std::string det, int layer) const {
0789 double tmp = 0.;
0790
0791 if (det == "HB") {
0792 tmp = scaleHB[layer];
0793 } else if (det == "HES" || det == "HED") {
0794 tmp = scaleHE[layer];
0795 } else if (det == "HF") {
0796 tmp = scaleHF[layer];
0797 }
0798
0799 return tmp;
0800 }
0801
0802 #include "FWCore/Framework/interface/MakerMacros.h"
0803 #include "FWCore/PluginManager/interface/ModuleDef.h"
0804 #include "SimG4Core/Watcher/interface/SimWatcherFactory.h"
0805
0806 DEFINE_SIMWATCHER(SimG4HcalValidation);