CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFHBHERecHitCreatorMaxSample.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFHBHeRecHitCreatorMaxSample_h
2 #define RecoParticleFlow_PFClusterProducer_PFHBHeRecHitCreatorMaxSample_h
3 
5 
11 
17 
22 
24 
25 
27 
28 
30 
31 
33 
34  public:
36  PFRecHitCreatorBase(iConfig,iC)
37  {
39  }
40 
41 
43  }
44 
45 
46  void importRecHits(std::auto_ptr<reco::PFRecHitCollection>&out,std::auto_ptr<reco::PFRecHitCollection>& cleaned ,const edm::Event& iEvent,const edm::EventSetup& iSetup) {
47 
48  beginEvent(iEvent,iSetup);
49 
51 
53  iSetup.get<CaloGeometryRecord>().get(geoHandle);
54 
55 
57  iSetup.get<HcalDbRecord > ().get(conditions);
58 
59  // get the ecal geometry
60  const CaloSubdetectorGeometry *hcalBarrelGeo =
61  geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalBarrel);
62  const CaloSubdetectorGeometry *hcalEndcapGeo =
63  geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalEndcap);
64 
65  iEvent.getByToken(recHitToken_,recHitHandle);
66  for( const auto& erh : *recHitHandle ) {
67  const HcalDetId& detid = (HcalDetId)erh.detid();
69 
70 
71  //CUSTOM ENERGY AND TIME RECO
72 
73  CaloSamples tool;
74  const HcalCalibrations& calibrations = conditions->getHcalCalibrations(detid);
75  const HcalQIECoder* channelCoder = conditions->getHcalCoder(detid);
76  const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
77  HcalCoderDb coder(*channelCoder, *shape);
78  int auxwd1 = erh.auxHBHE(); // TS = 0,1,2,3 info
79  int auxwd2 = erh.aux(); // TS = 4,5,6,7 info
80 
81  int adc[8];
82  int capid[8];
83 
84  adc[0] = (auxwd1) & 0x7F;
85  adc[1] = (auxwd1 >> 7) & 0x7F;
86  adc[2] = (auxwd1 >> 14) & 0x7F;
87  adc[3] = (auxwd1 >> 21) & 0x7F;
88 
89  capid[0] = (auxwd1 >> 28) & 0x3; // rotating through 4 values: 0,1,2,3
90  capid[1] = (capid[0] + 1 <= 3) ? capid[0] + 1 : 0;
91  capid[2] = (capid[1] + 1 <= 3) ? capid[1] + 1 : 0;
92  capid[3] = (capid[2] + 1 <= 3) ? capid[2] + 1 : 0;
93 
94  adc[4] = (auxwd2) & 0x7F;
95  adc[5] = (auxwd2 >> 7) & 0x7F;
96  adc[6] = (auxwd2 >> 14) & 0x7F;
97  adc[7] = (auxwd2 >> 21) & 0x7F;
98 
99  capid[4] = (auxwd2 >> 28) & 0x3;
100  capid[5] = (capid[4] + 1 <= 3) ? capid[4] + 1 : 0;
101  capid[6] = (capid[5] + 1 <= 3) ? capid[5] + 1 : 0;
102  capid[7] = (capid[6] + 1 <= 3) ? capid[6] + 1 : 0;
103 
104  // Pack the above into HBHEDataFrame for subsequent linearization to fC
105  HBHEDataFrame digi(detid);
106  digi.setSize(10);
107  digi.setPresamples(4);
108  for (int ii = 0; ii < 8; ii++) {
109  HcalQIESample s (adc[ii], capid[ii], 0, 0);
110  digi.setSample(ii,s);
111  }
112  coder.adc2fC(digi, tool);
113  HcalGenericDetId hcalGenDetId(erh.id());
114  std::array<double,8> samples;
115 
116 
117  //store the samples over thresholds
118  for (int ii = 0; ii < 8; ii++) {
119  double sampleE = calibrations.respcorrgain(capid[ii]) *
120  (tool[ii] - calibrations.pedestal(capid[ii]));
121  if (sampleE>sampleCut_)
122  samples[ii] = sampleE;
123  else
124  samples[ii] = 0.0;
125 
126  }
127 
128 
129  //run the algorithm
130  double energy=0.0;
131  double energy2=0.0;
132  double time=0.0;
133  double s2=0.0;
134  std::vector<double> hitEnergies;
135  std::vector<double> hitTimes;
136 
137  for (int ii = 0; ii < 8; ii++) {
138  energy=energy+samples[ii];
139  s2=samples[ii]*samples[ii];
140  time = time+(-100. + ii*25.0)*s2;
141  energy2=energy2+s2;
142 
143  if (ii>0 && ii<7) {
144  if (samples[ii]<=samples[ii-1] && samples[ii]<samples[ii+1] && energy>0) {
145  hitEnergies.push_back(energy);
146  hitTimes.push_back(time/energy2);
147  energy=0.0;
148  energy2=0.0;
149  time=0.0;
150  }
151 
152  }
153  else if (ii==7 && energy>0) {
154  hitEnergies.push_back(energy);
155  hitTimes.push_back(time/energy2);
156  energy=0.0;
157  energy2=0.0;
158  time=0.0;
159  }
160 
161  }
162 
163  if (hitEnergies.size()==0)
164  continue;
165 
166  int depth = detid.depth();
168  math::XYZVector axis;
169 
170  const CaloCellGeometry *thisCell=0;
172  switch(esd) {
173  case HcalBarrel:
174  thisCell =hcalBarrelGeo->getGeometry(detid);
175  layer =PFLayer::HCAL_BARREL1;
176  break;
177 
178  case HcalEndcap:
179  thisCell =hcalEndcapGeo->getGeometry(detid);
180  layer =PFLayer::HCAL_ENDCAP;
181  break;
182  default:
183  break;
184  }
185 
186  // find rechit geometry
187  if(!thisCell) {
188  edm::LogError("PFHBHERecHitCreatorMaxSample")
189  <<"warning detid "<<detid.rawId()
190  <<" not found in geometry"<<std::endl;
191  continue;
192  }
193 
194 
195  auto const point = thisCell->getPosition();
196  position.SetCoordinates ( point.x(),
197  point.y(),
198  point.z() );
199 
200 
201  reco::PFRecHit rh( detid.rawId(),layer,
202  energy,
203  position.x(), position.y(), position.z(),
204  0,0,0);
205 
206  rh.setDepth(depth);
207 
208 
209  const CaloCellGeometry::CornersVec& corners = thisCell->getCorners();
210  assert( corners.size() == 8 );
211 
212  rh.setNECorner( corners[0].x(), corners[0].y(), corners[0].z());
213  rh.setSECorner( corners[1].x(), corners[1].y(), corners[1].z());
214  rh.setSWCorner( corners[2].x(), corners[2].y(), corners[2].z());
215  rh.setNWCorner( corners[3].x(), corners[3].y(), corners[3].z());
216 
217  // for (unsigned int i=0;i<hitEnergies.size();++i)
218  // printf(" %f / %f ,",hitEnergies[i],hitTimes[i]);
219 
220  //now find the best hit
221  auto best_hit = std::min_element(hitTimes.begin(),hitTimes.end(),
222  [](const double& a,
223  const double& b){
224  return fabs(a) < fabs(b);
225  });
226 
227 
228  rh.setTime(*best_hit);
229  rh.setEnergy(hitEnergies[std::distance(hitTimes.begin(),best_hit)]);
230  // printf("Best = %f %f\n",rh.energy(),rh.time());
231 
232  //Apply Q tests
233  bool rcleaned = false;
234  bool keep=true;
235 
236  for( const auto& qtest : qualityTests_ ) {
237  if (!qtest->test(rh,erh,rcleaned)) {
238  keep = false;
239  }
240  }
241 
242  if(keep) {
243  out->push_back(rh);
244  }
245  else if (rcleaned)
246  cleaned->push_back(rh);
247  }
248  }
249 
250 
251  protected:
253  const double sampleCut_ = 0.6; // minimalistic threshold just to reduce the iterations
254 
255 };
256 
257 
258 
259 
260 
261 #endif
int adc(sample_type sample)
get the ADC sample (12 bits)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
void setDepth(int depth)
Definition: PFRecHit.h:79
double respcorrgain(int fCapId) const
get response corrected gain for capid=0..3
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
assert(m_qm.get())
double pedestal(int fCapId) const
get pedestal for capid=0..3
int ii
Definition: cuy.py:588
void setSize(int size)
float float float z
virtual const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
tuple s2
Definition: indexGen.py:106
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const int keep
int depth() const
get the tower depth
Definition: HcalDetId.h:40
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::SortedCollection< HBHERecHit > > recHitToken_
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
void setPresamples(int ps)
void setSample(int i, const HcalQIESample &sam)
Definition: HBHEDataFrame.h:50
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const
Definition: HcalCoderDb.cc:44
HcalSubdetector
Definition: HcalAssistant.h:31
void beginEvent(const edm::Event &event, const edm::EventSetup &setup)
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
tuple out
Definition: dbtoconf.py:99
Layer
layer definition
Definition: PFLayer.h:31
size_type size() const
Definition: EZArrayFL.h:81
PFHBHERecHitCreatorMaxSample(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const T & get() const
Definition: EventSetup.h:55
double b
Definition: hdecay.h:120
void importRecHits(std::auto_ptr< reco::PFRecHitCollection > &out, std::auto_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup)
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
Definition: DDAxes.h:10
const CornersVec & getCorners() const
Returns the corner points of this cell&#39;s volume.
const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5