CMS 3D CMS Logo

TBHodoActiveVolumeRawInfoProducer.cc
Go to the documentation of this file.
1 /*
2  * \file TBHodoActiveVolumeRawInfoProducer.cc
3  *
4  *
5  */
6 
15 
23 
24 #include <iostream>
25 #include <memory>
26 
28 public:
31 
33  ~TBHodoActiveVolumeRawInfoProducer() override = default;
34 
36  void produce(edm::Event &event, const edm::EventSetup &eventSetup) override;
37 
38 private:
39  double myThreshold;
41 
42  std::unique_ptr<EcalTBHodoscopeGeometry> theTBHodoGeom_;
43 };
44 
45 using namespace cms;
46 using namespace std;
47 
49  : myThreshold(0.05E-3),
50  m_EcalToken(consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", "EcalTBH4BeamHits"))) {
51  produces<EcalTBHodoscopeRawInfo>();
52 
53  theTBHodoGeom_ = std::make_unique<EcalTBHodoscopeGeometry>();
54 }
55 
57  std::unique_ptr<EcalTBHodoscopeRawInfo> product(new EcalTBHodoscopeRawInfo());
58 
59  // caloHit container
60  const edm::Handle<edm::PCaloHitContainer> &pCaloHit = event.getHandle(m_EcalToken);
61  const edm::PCaloHitContainer *caloHits = nullptr;
62  if (pCaloHit.isValid()) {
63  caloHits = pCaloHit.product();
64  LogDebug("EcalTBHodo") << "total # caloHits: " << caloHits->size();
65  } else {
66  edm::LogError("EcalTBHodo") << "Error! can't get the caloHitContainer ";
67  }
68  if (!caloHits) {
69  return;
70  }
71 
72  // detid - energy_sum map
73  std::map<unsigned int, double> energyMap;
74 
75  for (auto &&aHit : *caloHits) {
76  double thisHitEne = aHit.energy();
77 
78  std::map<unsigned int, double>::iterator itmap = energyMap.find(aHit.id());
79  if (itmap == energyMap.end())
80  energyMap.insert(pair<unsigned int, double>(aHit.id(), thisHitEne));
81  else {
82  (*itmap).second += thisHitEne;
83  }
84  }
85 
86  // planes and fibers
87  int nPlanes = theTBHodoGeom_->getNPlanes();
88  int nFibers = theTBHodoGeom_->getNFibres();
89  product->setPlanes(nPlanes);
90 
91  bool firedChannels[4][64];
92  for (int iPlane = 0; iPlane < nPlanes; ++iPlane) {
93  for (int iFiber = 0; iFiber < nFibers; ++iFiber) {
94  firedChannels[iPlane][iFiber] = 0.;
95  }
96  }
97  for (std::map<unsigned int, double>::const_iterator itmap = energyMap.begin(); itmap != energyMap.end(); ++itmap) {
98  if ((*itmap).second > myThreshold) {
99  HodoscopeDetId myHodoDetId = HodoscopeDetId((*itmap).first);
100  firedChannels[myHodoDetId.planeId()][myHodoDetId.fibrId()] = true;
101  }
102  }
103  for (int iPlane = 0; iPlane < nPlanes; ++iPlane) {
104  EcalTBHodoscopePlaneRawHits planeHit(nFibers);
105 
106  for (int iFiber = 0; iFiber < nFibers; ++iFiber) {
107  planeHit.setHit(iFiber, firedChannels[iPlane][iFiber]);
108  }
109  product->setPlane((unsigned int)iPlane, planeHit);
110  }
111 
112  LogDebug("EcalTBHodo") << (*product);
113 
114  event.put(std::move(product));
115 }
116 
std::vector< PCaloHit > PCaloHitContainer
void setHit(unsigned int i, bool status)
T const * product() const
Definition: Handle.h:70
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
Produce digis out of raw data.
Log< level::Error, false > LogError
TBHodoActiveVolumeRawInfoProducer(const edm::ParameterSet &ps)
Constructor.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Namespace of DDCMS conversion namespace.
~TBHodoActiveVolumeRawInfoProducer() override=default
Destructor.
bool isValid() const
Definition: HandleBase.h:70
int fibrId() const
HLT enums.
int planeId() const
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
std::unique_ptr< EcalTBHodoscopeGeometry > theTBHodoGeom_
#define LogDebug(id)
edm::EDGetTokenT< edm::PCaloHitContainer > m_EcalToken