Skip to content

Commit e9b58be

Browse files
authored
Merge pull request #46886 from bsunanda/Phase2-hgx360T
Phase2-hgx360T Make a test code to find the area of HGCal cells through the method "cellArea" of HGCalDDConst
2 parents fbe851c + 03f691f commit e9b58be

File tree

2 files changed

+238
-0
lines changed

2 files changed

+238
-0
lines changed
Lines changed: 162 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,162 @@
1+
// -*- C++ -*-
2+
//
3+
// Package: HGCalTestCellArea
4+
// Class: HGCalTestCellArea
5+
//
6+
/**\class HGCalTestCellArea HGCalTestCellArea.cc
7+
test/HGCalTestCellArea.cc
8+
9+
Description: <one line class summary>
10+
11+
Implementation:
12+
<Notes on implementation>
13+
*/
14+
//
15+
// Original Author: Sunanda Banerjee
16+
// Created: Mon 2024/11/29
17+
//
18+
//
19+
20+
// system include files
21+
#include <fstream>
22+
#include <iostream>
23+
#include <sstream>
24+
#include <string>
25+
#include <vector>
26+
27+
// user include files
28+
#include "FWCore/Framework/interface/Frameworkfwd.h"
29+
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
30+
#include "FWCore/Framework/interface/Event.h"
31+
#include "FWCore/Framework/interface/EventSetup.h"
32+
#include "FWCore/Framework/interface/MakerMacros.h"
33+
#include "FWCore/MessageLogger/interface/MessageLogger.h"
34+
#include "FWCore/ParameterSet/interface/FileInPath.h"
35+
#include "FWCore/ParameterSet/interface/ParameterSet.h"
36+
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
37+
#include "FWCore/Utilities/interface/transform.h"
38+
39+
#include "DataFormats/DetId/interface/DetId.h"
40+
#include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
41+
#include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
42+
#include "Geometry/HGCalCommonData/interface/HGCalGeomUtils.h"
43+
#include "Geometry/Records/interface/IdealGeometryRecord.h"
44+
45+
class HGCalTestCellArea : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
46+
public:
47+
explicit HGCalTestCellArea(const edm::ParameterSet &);
48+
~HGCalTestCellArea() override = default;
49+
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
50+
51+
void beginJob() override {}
52+
void beginRun(edm::Run const &, edm::EventSetup const &) override;
53+
void analyze(edm::Event const &iEvent, edm::EventSetup const &) override {}
54+
void endRun(edm::Run const &, edm::EventSetup const &) override {}
55+
void endJob() override {}
56+
57+
private:
58+
const std::vector<std::string> nameDetectors_;
59+
const std::string fileName_;
60+
const std::vector<edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord>> tok_hgcal_;
61+
std::vector<const HGCalDDDConstants *> hgcCons_;
62+
std::vector<std::pair<DetId, uint32_t>> detIds_;
63+
};
64+
65+
HGCalTestCellArea::HGCalTestCellArea(const edm::ParameterSet &iC)
66+
: nameDetectors_(iC.getParameter<std::vector<std::string>>("nameDetectors")),
67+
fileName_(iC.getParameter<std::string>("fileName")),
68+
tok_hgcal_{edm::vector_transform(nameDetectors_, [this](const std::string &name) {
69+
return esConsumes<HGCalDDDConstants, IdealGeometryRecord, edm::Transition::BeginRun>(edm::ESInputTag{"", name});
70+
})} {
71+
std::ostringstream st1;
72+
for (const auto &name : nameDetectors_)
73+
st1 << " : " << name;
74+
edm::LogVerbatim("HGCGeom") << "Test validity of cells for " << nameDetectors_.size() << " detectors" << st1.str()
75+
<< " with inputs from " << fileName_;
76+
if (!fileName_.empty()) {
77+
edm::FileInPath filetmp("Geometry/HGCalCommonData/data/" + fileName_);
78+
std::string fileName = filetmp.fullPath();
79+
std::ifstream fInput(fileName.c_str());
80+
if (!fInput.good()) {
81+
edm::LogVerbatim("HGCGeom") << "Cannot open file " << fileName;
82+
} else {
83+
char buffer[80];
84+
const std::vector<DetId::Detector> dets = {DetId::HGCalEE, DetId::HGCalHSi, DetId::HGCalHSc};
85+
while (fInput.getline(buffer, 80)) {
86+
std::vector<std::string> items = HGCalGeomUtils::splitString(std::string(buffer));
87+
if (items.size() == 8) {
88+
DetId::Detector det = static_cast<DetId::Detector>(std::atoi(items[0].c_str()));
89+
auto itr = std::find(dets.begin(), dets.end(), det);
90+
if (itr != dets.end()) {
91+
uint32_t pos = static_cast<uint32_t>(itr - dets.begin());
92+
DetId id(0);
93+
if ((det == DetId::HGCalEE) || (det == DetId::HGCalHSi)) {
94+
int type = std::atoi(items[1].c_str());
95+
int zside = std::atoi(items[2].c_str());
96+
int layer = std::atoi(items[3].c_str());
97+
int waferU = std::atoi(items[4].c_str());
98+
int waferV = std::atoi(items[5].c_str());
99+
int cellU = std::atoi(items[6].c_str());
100+
int cellV = std::atoi(items[7].c_str());
101+
id = static_cast<DetId>(HGCSiliconDetId(det, zside, type, layer, waferU, waferV, cellU, cellV));
102+
detIds_.emplace_back(id, pos);
103+
}
104+
}
105+
}
106+
}
107+
fInput.close();
108+
}
109+
}
110+
edm::LogVerbatim("HGCGeom") << "Reads " << detIds_.size() << " ID's from " << fileName_;
111+
for (unsigned int k = 0; k < detIds_.size(); ++k) {
112+
edm::LogVerbatim("HGCGeom") << "[" << k << "] " << HGCSiliconDetId(detIds_[k].first) << " from DDConstant "
113+
<< (detIds_[k].second);
114+
}
115+
}
116+
117+
void HGCalTestCellArea::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
118+
std::vector<std::string> names = {"HGCalEESensitive", "HGCalHESiliconSensitive"};
119+
edm::ParameterSetDescription desc;
120+
desc.add<std::vector<std::string>>("nameDetectors", names);
121+
desc.add<std::string>("fileName", "missD88.txt");
122+
descriptions.add("hgcalTestCellArea", desc);
123+
}
124+
125+
// ------------ method called to produce the data ------------
126+
void HGCalTestCellArea::beginRun(edm::Run const &iRun, edm::EventSetup const &iSetup) {
127+
//initiating hgc Geometry
128+
std::vector<std::string> names = {"HGCalEESensitive", "HGCalHESiliconSensitive"};
129+
std::vector<DetId::Detector> dets = {DetId::HGCalEE, DetId::HGCalHSi};
130+
std::map<DetId::Detector, uint32_t> detMap;
131+
for (uint32_t i = 0; i < nameDetectors_.size(); i++) {
132+
edm::LogVerbatim("HGCGeom") << "Tries to initialize HGCalGeometry and HGCalDDDConstants for " << i << ":"
133+
<< nameDetectors_[i];
134+
const edm::ESHandle<HGCalDDDConstants> &hgcCons = iSetup.getHandle(tok_hgcal_[i]);
135+
if (hgcCons.isValid()) {
136+
hgcCons_.push_back(hgcCons.product());
137+
} else {
138+
edm::LogWarning("HGCGeom") << "Cannot initiate HGCalDDDConstants for " << nameDetectors_[i] << std::endl;
139+
}
140+
auto ii = std::find(names.begin(), names.end(), nameDetectors_[i]);
141+
if (ii != names.end()) {
142+
uint32_t k = static_cast<uint32_t>(ii - names.begin());
143+
detMap[dets[k]] = i;
144+
}
145+
}
146+
edm::LogVerbatim("HGCGeom") << "Loaded HGCalDDConstants for " << detMap.size() << " detectors";
147+
148+
for (auto itr = detMap.begin(); itr != detMap.end(); ++itr)
149+
edm::LogVerbatim("HGCGeom") << "[" << itr->second << "]: " << nameDetectors_[itr->second] << " for Detector "
150+
<< itr->first;
151+
152+
for (unsigned int k = 0; k < detIds_.size(); ++k) {
153+
const HGCalDDDConstants *cons = hgcCons_[detMap[(detIds_[k].first).det()]];
154+
HGCSiliconDetId id(detIds_[k].first);
155+
edm::LogVerbatim("HGCGeom") << "Hit[" << k << "] " << id << " Area " << cons->cellArea(id, false) << " Valid "
156+
<< cons->isValidHex8(
157+
id.layer(), id.waferU(), id.waferV(), id.cellU(), id.cellV(), true);
158+
}
159+
}
160+
161+
// define this as a plug-in
162+
DEFINE_FWK_MODULE(HGCalTestCellArea);
Lines changed: 76 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,76 @@
1+
###############################################################################
2+
# Way to use this:
3+
# cmsRun testHGCalNumbering_cfg.py type=V19
4+
#
5+
# Options for type V16, V17, V17n, V18, V19
6+
#
7+
###############################################################################
8+
import FWCore.ParameterSet.Config as cms
9+
import os, sys, importlib, re
10+
import FWCore.ParameterSet.VarParsing as VarParsing
11+
12+
####################################################################
13+
### SETUP OPTIONS
14+
options = VarParsing.VarParsing('standard')
15+
options.register('type',
16+
"V19",
17+
VarParsing.VarParsing.multiplicity.singleton,
18+
VarParsing.VarParsing.varType.string,
19+
"type of operations: V16, V17, V17n, V18, V19")
20+
21+
### get and parse the command line arguments
22+
options.parseArguments()
23+
print(options)
24+
25+
from Configuration.Eras.Era_Phase2C17I13M9_cff import Phase2C17I13M9
26+
process = cms.Process("HGCalCellArea",Phase2C17I13M9)
27+
28+
geomFile = "Geometry.HGCalCommonData.testHGCal" + options.type + "XML_cfi"
29+
print("Geometry file: ", geomFile)
30+
31+
32+
process.load(geomFile)
33+
process.load("Geometry.HGCalCommonData.hgcalParametersInitialization_cfi")
34+
process.load("Geometry.HGCalCommonData.hgcalNumberingInitialization_cfi")
35+
process.load("Geometry.EcalCommonData.ecalSimulationParameters_cff")
36+
process.load("Geometry.HcalCommonData.hcalDDDSimConstants_cff")
37+
process.load("SimGeneral.HepPDTESSource.pdt_cfi")
38+
process.load('FWCore.MessageService.MessageLogger_cfi')
39+
40+
if hasattr(process,'MessageLogger'):
41+
process.MessageLogger.HGCalGeom=dict()
42+
process.MessageLogger.HGCGeom=dict()
43+
44+
process.load("IOMC.RandomEngine.IOMC_cff")
45+
process.RandomNumberGeneratorService.generator.initialSeed = 456789
46+
47+
process.source = cms.Source("EmptySource")
48+
49+
process.generator = cms.EDProducer("FlatRandomEGunProducer",
50+
PGunParameters = cms.PSet(
51+
PartID = cms.vint32(14),
52+
MinEta = cms.double(-3.5),
53+
MaxEta = cms.double(3.5),
54+
MinPhi = cms.double(-3.14159265359),
55+
MaxPhi = cms.double(3.14159265359),
56+
MinE = cms.double(9.99),
57+
MaxE = cms.double(10.01)
58+
),
59+
AddAntiParticle = cms.bool(False),
60+
Verbosity = cms.untracked.int32(0),
61+
firstRun = cms.untracked.uint32(1)
62+
)
63+
64+
process.maxEvents = cms.untracked.PSet(
65+
input = cms.untracked.int32(1)
66+
)
67+
68+
process.SimpleMemoryCheck = cms.Service("SimpleMemoryCheck",
69+
ignoreTotal = cms.untracked.int32(1),
70+
moduleMemorySummary = cms.untracked.bool(True)
71+
)
72+
73+
process.load("Geometry.HGCalCommonData.hgcalTestCellArea_cfi")
74+
75+
76+
process.p1 = cms.Path(process.generator*process.hgcalTestCellArea)

0 commit comments

Comments
 (0)