CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFHBHERecHitCreator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFHBHeRecHitCreator_h
2 #define RecoParticleFlow_PFClusterProducer_PFHBHeRecHitCreator_h
3 
5 
11 
17 
20 
21  public:
23  PFRecHitCreatorBase(iConfig,iC)
24  {
26  }
27 
28  void importRecHits(std::auto_ptr<reco::PFRecHitCollection>&out,std::auto_ptr<reco::PFRecHitCollection>& cleaned ,const edm::Event& iEvent,const edm::EventSetup& iSetup) {
29 
30  beginEvent(iEvent,iSetup);
31 
33 
35  iSetup.get<CaloGeometryRecord>().get(geoHandle);
36 
37  // get the ecal geometry
38  const CaloSubdetectorGeometry *hcalBarrelGeo =
39  geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalBarrel);
40  const CaloSubdetectorGeometry *hcalEndcapGeo =
41  geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalEndcap);
42 
43  iEvent.getByToken(recHitToken_,recHitHandle);
44  for( const auto& erh : *recHitHandle ) {
45  const HcalDetId& detid = (HcalDetId)erh.detid();
47 
48  double energy = erh.energy();
49  double time = erh.time();
50  int depth = detid.depth();
51 
53  math::XYZVector axis;
54 
55  const CaloCellGeometry *thisCell=0;
57  switch(esd) {
58  case HcalBarrel:
59  thisCell =hcalBarrelGeo->getGeometry(detid);
60  layer =PFLayer::HCAL_BARREL1;
61  break;
62 
63  case HcalEndcap:
64  thisCell =hcalEndcapGeo->getGeometry(detid);
65  layer =PFLayer::HCAL_ENDCAP;
66  break;
67  default:
68  break;
69  }
70 
71  // find rechit geometry
72  if(!thisCell) {
73  edm::LogError("PFHBHERecHitCreator")
74  <<"warning detid "<<detid.rawId()
75  <<" not found in geometry"<<std::endl;
76  continue;
77  }
78 
79  auto const point = thisCell->getPosition();
80  position.SetCoordinates ( point.x(),
81  point.y(),
82  point.z() );
83 
84  reco::PFRecHit rh( detid.rawId(),layer,
85  energy,
86  position.x(), position.y(), position.z(),
87  0,0,0);
88  rh.setTime(time); //Mike: This we will use later
89  rh.setDepth(depth);
90  const CaloCellGeometry::CornersVec& corners = thisCell->getCorners();
91  assert( corners.size() == 8 );
92 
93  rh.setNECorner( corners[0].x(), corners[0].y(), corners[0].z());
94  rh.setSECorner( corners[1].x(), corners[1].y(), corners[1].z());
95  rh.setSWCorner( corners[2].x(), corners[2].y(), corners[2].z());
96  rh.setNWCorner( corners[3].x(), corners[3].y(), corners[3].z());
97 
98 
99  bool rcleaned = false;
100  bool keep=true;
101 
102  //Apply Q tests
103  for( const auto& qtest : qualityTests_ ) {
104  if (!qtest->test(rh,erh,rcleaned)) {
105  keep = false;
106 
107  }
108  }
109 
110  if(keep) {
111  out->push_back(rh);
112  }
113  else if (rcleaned)
114  cleaned->push_back(rh);
115  }
116  }
117 
118 
119 
120  protected:
122 
123 
124 };
125 
126 
127 
128 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
assert(m_qm.get())
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.
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
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
HcalSubdetector
Definition: HcalAssistant.h:31
void beginEvent(const edm::Event &event, const edm::EventSetup &setup)
void setTime(double time)
Definition: PFRecHit.h:78
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
PFHBHERecHitCreator(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
static int position[264][3]
Definition: ReadPGInfo.cc:509
void importRecHits(std::auto_ptr< reco::PFRecHitCollection > &out, std::auto_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup)
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.
edm::EDGetTokenT< edm::SortedCollection< HBHERecHit > > recHitToken_
*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