CMS 3D CMS Logo

HGCalTestPartialWaferRecHits.cc
Go to the documentation of this file.
3 
7 
14 
18 
25 
26 #include <fstream>
27 #include <string>
28 #include <vector>
29 
31 public:
33  ~HGCalTestPartialWaferRecHits() override = default;
34  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
35 
36 protected:
37  void analyze(edm::Event const&, edm::EventSetup const&) override;
38  void beginJob() override {}
39  void endJob() override {}
40 
41 private:
46  std::vector<int> wafers_;
47  std::vector<int> dumpDets_;
48 };
49 
51  : nameDetector_(ps.getParameter<std::string>("detectorName")),
52  missingFile_(ps.getParameter<std::string>("missingFile")),
53  source_(ps.getParameter<edm::InputTag>("source")),
54  recHitSource_(consumes<HGCRecHitCollection>(source_)),
56  edm::LogVerbatim("HGCalSim") << "Test Hit ID using RecHits for " << nameDetector_ << " with module Label: " << source_
57  << " Missing Wafer file " << missingFile_;
58  if (!missingFile_.empty()) {
59  edm::FileInPath filetmp("SimG4CMS/Calo/data/" + missingFile_);
60  std::string fileName = filetmp.fullPath();
61  std::ifstream fInput(fileName.c_str());
62  if (!fInput.good()) {
63  edm::LogVerbatim("HGCalSim") << "Cannot open file " << fileName;
64  } else {
65  char buffer[80];
66  while (fInput.getline(buffer, 80)) {
67  std::vector<std::string> items = CaloSimUtils::splitString(std::string(buffer));
68  if (items.size() > 2) {
69  int layer = std::atoi(items[0].c_str());
70  int waferU = std::atoi(items[1].c_str());
71  int waferV = std::atoi(items[2].c_str());
72  wafers_.emplace_back(HGCalWaferIndex::waferIndex(layer, waferU, waferV, false));
73  } else if (items.size() == 1) {
74  int dumpdet = std::atoi(items[0].c_str());
75  dumpDets_.emplace_back(dumpdet);
76  edm::LogVerbatim("HGCalSim") << nameDetector_ << " Dump detector " << dumpdet;
77  }
78  }
79  edm::LogVerbatim("HGCalSim") << "HGCalTestPartialWaferRecHits::Reads in " << wafers_.size() << ":"
80  << dumpDets_.size() << " wafer|detector information from " << fileName;
81  fInput.close();
82  }
83  }
84 }
85 
88  desc.add<std::string>("detectorName", "HGCalEESensitive");
89  desc.add<std::string>("missingFile", "");
90  desc.add<edm::InputTag>("source", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
91  descriptions.add("hgcalRecHitPartialEE", desc);
92 }
93 
95  // get hgcalGeometry
96  const HGCalGeometry* geom = &iS.getData(tok_hgcGeom_);
97  const HGCalDDDConstants& hgc = geom->topology().dddConstants();
98  int firstLayer = hgc.getLayerOffset();
99  // get the hit collection
100  const edm::Handle<HGCRecHitCollection>& theRecHitContainers = e.getHandle(recHitSource_);
101  bool getHits = (theRecHitContainers.isValid());
102  uint32_t nhits = (getHits) ? theRecHitContainers->size() : 0;
103  uint32_t good(0), allSi(0), all(0);
104  edm::LogVerbatim("HGCalSim") << "HGCalTestPartialWaferRecHits: Input flags Hits " << getHits << " with " << nhits
105  << " hits: Layer Offset " << firstLayer;
106 
107  if (getHits) {
108  // Loop over all hits
109  for (const auto& it : *(theRecHitContainers.product())) {
110  ++all;
111  DetId id(it.id());
112  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
113  ++allSi;
114  HGCSiliconDetId hid(id);
115  const auto& info = hgc.waferInfo(hid.layer(), hid.waferU(), hid.waferV());
116  bool toCheck(false);
117  if (!wafers_.empty()) {
118  int indx = HGCalWaferIndex::waferIndex(firstLayer + hid.layer(), hid.waferU(), hid.waferV(), false);
119  if (std::find(wafers_.begin(), wafers_.end(), indx) != wafers_.end())
120  toCheck = true;
121  } else if (!dumpDets_.empty()) {
122  if ((std::find(dumpDets_.begin(), dumpDets_.end(), static_cast<int>(id.det())) != dumpDets_.end()) &&
123  (info.part != HGCalTypes::WaferFull))
124  toCheck = true;
125  } else {
126  // Only partial wafers
127  toCheck = (info.part != HGCalTypes::WaferFull);
128  }
129  if (toCheck) {
130  ++good;
131  GlobalPoint pos = geom->getPosition(id);
132  bool valid1 = geom->topology().valid(id);
133  bool valid2 = hgc.isValidHex8(hid.layer(), hid.waferU(), hid.waferV(), hid.cellU(), hid.cellV(), false);
134  edm::LogVerbatim("HGCalSim") << "Hit[" << all << ":" << allSi << ":" << good << "]" << hid
135  << " Wafer Type:Part:Orient:Cassette " << info.type << ":" << info.part << ":"
136  << info.orient << ":" << info.cassette << " at (" << pos.x() << ", " << pos.y()
137  << ", " << pos.z() << ") Validity " << valid1 << ":" << valid2;
138  }
139  }
140  }
141  }
142  edm::LogVerbatim("HGCalSim") << "Total hits = " << all << ":" << nhits << " Good DetIds = " << allSi << ":" << good;
143 }
144 
145 //define this as a plug-in
Log< level::Info, true > LogVerbatim
const edm::EDGetTokenT< HGCRecHitCollection > recHitSource_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
static const TGPicture * info(bool iBackgroundIsBlack)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
int32_t waferU(const int32_t index)
size_type size() const
int cellU() const
get the cell #&#39;s in u,v or in x,y
T const * product() const
Definition: Handle.h:70
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
~HGCalTestPartialWaferRecHits() override=default
int layer() const
get the layer #
static constexpr int32_t WaferFull
Definition: HGCalTypes.h:35
std::vector< std::string > splitString(const std::string &)
Definition: CaloSimUtils.cc:3
HGCalTestPartialWaferRecHits(const edm::ParameterSet &ps)
int waferU() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void analyze(edm::Event const &, edm::EventSetup const &) override
int32_t waferIndex(int32_t layer, int32_t waferU, int32_t waferV, bool old=false)
Definition: DetId.h:17
int waferV() const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int cellV() const
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
int32_t waferV(const int32_t index)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::ESGetToken< HGCalGeometry, IdealGeometryRecord > tok_hgcGeom_