CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
EcalDetailedTimeRecHitProducer.cc
Go to the documentation of this file.
1 
15 
18 
20 
21 #include <cmath>
22 #include <iostream>
23 #include <memory>
24 
25 #include <vector>
26 
38 
39 #include "CLHEP/Units/GlobalPhysicalConstants.h"
40 #include "CLHEP/Units/GlobalSystemOfUnits.h"
41 
43  EBRecHitCollection_ = consumes<EBRecHitCollection>(ps.getParameter<edm::InputTag>("EBRecHitCollection"));
44  EERecHitCollection_ = consumes<EERecHitCollection>(ps.getParameter<edm::InputTag>("EERecHitCollection"));
45 
46  ebTimeDigiCollection_ = consumes<EcalTimeDigiCollection>(ps.getParameter<edm::InputTag>("EBTimeDigiCollection"));
47  eeTimeDigiCollection_ = consumes<EcalTimeDigiCollection>(ps.getParameter<edm::InputTag>("EETimeDigiCollection"));
48 
49  caloGeometry_ = esConsumes<CaloGeometry, CaloGeometryRecord>();
50 
51  EBDetailedTimeRecHitCollection_ = ps.getParameter<std::string>("EBDetailedTimeRecHitCollection");
52  EEDetailedTimeRecHitCollection_ = ps.getParameter<std::string>("EEDetailedTimeRecHitCollection");
53 
54  correctForVertexZPosition_ = ps.getParameter<bool>("correctForVertexZPosition");
55  useMCTruthVertex_ = ps.getParameter<bool>("useMCTruthVertex");
56  if (correctForVertexZPosition_) {
57  if (not useMCTruthVertex_) {
58  recoVertex_ = consumes<reco::VertexCollection>(ps.getParameter<edm::InputTag>("recoVertex"));
59  } else {
60  simVertex_ = consumes<edm::SimVertexContainer>(ps.getParameter<edm::InputTag>("simVertex"));
61  }
62  }
63 
64  ebTimeLayer_ = ps.getParameter<int>("EBTimeLayer");
65  eeTimeLayer_ = ps.getParameter<int>("EETimeLayer");
66 
67  produces<EBRecHitCollection>(EBDetailedTimeRecHitCollection_);
68  produces<EERecHitCollection>(EEDetailedTimeRecHitCollection_);
69 }
70 
72 
74  using namespace edm;
75  using namespace reco;
76 
78 
79  m_geometry = hGeometry.product();
80 
81  Handle<EBRecHitCollection> pEBRecHits;
82  Handle<EERecHitCollection> pEERecHits;
83 
84  const EBRecHitCollection* EBRecHits = nullptr;
85  const EERecHitCollection* EERecHits = nullptr;
86 
87  evt.getByToken(EBRecHitCollection_, pEBRecHits);
88  if (pEBRecHits.isValid()) {
89  EBRecHits = pEBRecHits.product(); // get a ptr to the product
90 #ifdef DEBUG
91  LogDebug("EcalRecHitDebug") << "total # EB rechits to be re-calibrated: " << EBRecHits->size();
92 #endif
93  }
94 
95  evt.getByToken(EERecHitCollection_, pEERecHits);
96  if (pEERecHits.isValid()) {
97  EERecHits = pEERecHits.product(); // get a ptr to the product
98 #ifdef DEBUG
99  LogDebug("EcalRecHitDebug") << "total # EE uncalibrated rechits to be re-calibrated: " << EERecHits->size();
100 #endif
101  }
102 
103  Handle<EcalTimeDigiCollection> pEBTimeDigis;
104  Handle<EcalTimeDigiCollection> pEETimeDigis;
105 
106  const EcalTimeDigiCollection* ebTimeDigis = nullptr;
107  const EcalTimeDigiCollection* eeTimeDigis = nullptr;
108 
109  evt.getByToken(ebTimeDigiCollection_, pEBTimeDigis);
110  //evt.getByToken( digiProducer_, pEBTimeDigis);
111  if (pEBTimeDigis.isValid()) {
112  ebTimeDigis = pEBTimeDigis.product(); // get a ptr to the produc
113  edm::LogInfo("EcalDetailedTimeRecHitInfo") << "total # ebTimeDigis: " << ebTimeDigis->size();
114  }
115 
116  evt.getByToken(eeTimeDigiCollection_, pEETimeDigis);
117  //evt.getByToken( digiProducer_, pEETimeDigis);
118  if (pEETimeDigis.isValid()) {
119  eeTimeDigis = pEETimeDigis.product(); // get a ptr to the product
120  edm::LogInfo("EcalDetailedTimeRecHitInfo") << "total # eeTimeDigis: " << eeTimeDigis->size();
121  }
122  // collection of rechits to put in the event
123  std::unique_ptr<EBRecHitCollection> EBDetailedTimeRecHits(new EBRecHitCollection);
124  std::unique_ptr<EERecHitCollection> EEDetailedTimeRecHits(new EERecHitCollection);
125 
126  std::unique_ptr<GlobalPoint> vertex;
127 
129  if (!useMCTruthVertex_) {
130  //Get the first reco vertex
131  // get primary vertices
132 
133  edm::Handle<VertexCollection> VertexHandle;
134  evt.getByToken(recoVertex_, VertexHandle);
135 
136  if (VertexHandle.isValid()) {
137  if (!(*VertexHandle).empty()) //at least 1 vertex
138  {
139  const reco::Vertex* myVertex = &(*VertexHandle)[0];
140  vertex = std::make_unique<GlobalPoint>(myVertex->x(), myVertex->y(), myVertex->z());
141  }
142  }
143 
144  } else {
145  edm::Handle<SimVertexContainer> VertexHandle;
146  evt.getByToken(simVertex_, VertexHandle);
147 
148  if (VertexHandle.isValid()) {
149  if (!(*VertexHandle).empty()) //at least 1 vertex
150  {
151  assert((*VertexHandle)[0].vertexId() == 0);
152  const SimVertex* myVertex = &(*VertexHandle)[0];
153  vertex = std::make_unique<GlobalPoint>(
154  myVertex->position().x(), myVertex->position().y(), myVertex->position().z());
155  }
156  }
157  }
158  }
159 
160  if (EBRecHits && ebTimeDigis) {
161  // loop over uncalibrated rechits to make calibrated ones
162  for (EBRecHitCollection::const_iterator it = EBRecHits->begin(); it != EBRecHits->end(); ++it) {
163  EcalRecHit aHit((*it));
164  EcalTimeDigiCollection::const_iterator timeDigi = ebTimeDigis->find((*it).id());
165  if (timeDigi != ebTimeDigis->end()) {
166  if (timeDigi->sampleOfInterest() >= 0) {
167  float myTime = (*timeDigi)[timeDigi->sampleOfInterest()];
168  //Vertex corrected ToF
169  if (vertex) {
170  aHit.setTime(myTime + deltaTimeOfFlight(*vertex, (*it).id(), ebTimeLayer_));
171  } else
172  //Uncorrected ToF
173  aHit.setTime(myTime);
174  }
175  }
176  // leave standard time if no timeDigi is associated (e.g. noise recHits)
177  EBDetailedTimeRecHits->push_back(aHit);
178  }
179  }
180 
181  if (EERecHits && eeTimeDigis) {
182  // loop over uncalibrated rechits to make calibrated ones
183  for (EERecHitCollection::const_iterator it = EERecHits->begin(); it != EERecHits->end(); ++it) {
184  EcalRecHit aHit(*it);
185  EcalTimeDigiCollection::const_iterator timeDigi = eeTimeDigis->find((*it).id());
186  if (timeDigi != eeTimeDigis->end()) {
187  if (timeDigi->sampleOfInterest() >= 0) {
188  float myTime = (*timeDigi)[timeDigi->sampleOfInterest()];
189  //Vertex corrected ToF
190  if (vertex) {
191  aHit.setTime(myTime + deltaTimeOfFlight(*vertex, (*it).id(), eeTimeLayer_));
192  } else
193  //Uncorrected ToF
194  aHit.setTime(myTime);
195  }
196  }
197  EEDetailedTimeRecHits->push_back(aHit);
198  }
199  }
200  // put the collection of recunstructed hits in the event
201  LogInfo("EcalDetailedTimeRecHitInfo") << "total # EB rechits: " << EBDetailedTimeRecHits->size();
202  LogInfo("EcalDetailedTimeRecHitInfo") << "total # EE rechits: " << EEDetailedTimeRecHits->size();
203 
204  evt.put(std::move(EBDetailedTimeRecHits), EBDetailedTimeRecHitCollection_);
205  evt.put(std::move(EEDetailedTimeRecHits), EEDetailedTimeRecHitCollection_);
206 }
207 
209  auto cellGeometry(m_geometry->getGeometry(detId));
210  assert(nullptr != cellGeometry);
211  GlobalPoint layerPos =
212  cellGeometry->getPosition(double(layer) + 0.5); //depth in mm in the middle of the layer position
213  GlobalVector tofVector = layerPos - vertex;
214  return (layerPos.mag() * cm - tofVector.mag() * cm) / (float)c_light;
215 }
216 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double y() const
y coordinate
Definition: Vertex.h:131
edm::EDGetTokenT< EcalTimeDigiCollection > eeTimeDigiCollection_
std::vector< EcalRecHit >::const_iterator const_iterator
assert(be >=bs)
edm::EDGetTokenT< EcalTimeDigiCollection > ebTimeDigiCollection_
constexpr std::array< uint8_t, layerIndexSize > layer
void setTime(float time)
Definition: EcalRecHit.h:71
T mag() const
Definition: PV3DBase.h:64
edm::EDGetTokenT< reco::VertexCollection > recoVertex_
def move
Definition: eostools.py:511
double z() const
z coordinate
Definition: Vertex.h:133
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:21
bool isValid() const
Definition: HandleBase.h:70
Log< level::Info, false > LogInfo
EcalDetailedTimeRecHitProducer(const edm::ParameterSet &ps)
Definition: DetId.h:17
double x() const
x coordinate
Definition: Vertex.h:129
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< EBRecHitCollection > EBRecHitCollection_
double deltaTimeOfFlight(GlobalPoint &vertex, const DetId &detId, int layer) const
edm::EDGetTokenT< edm::SimVertexContainer > simVertex_
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:60
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
edm::EDGetTokenT< EERecHitCollection > EERecHitCollection_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometry_
void produce(edm::Event &evt, const edm::EventSetup &es) override
#define LogDebug(id)