CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFHcalRecHitCreator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFHcalRecHitCreator_h
2 #define RecoParticleFlow_PFClusterProducer_PFHcalRecHitCreator_h
3 
5 
11 
17 
19 
20 template <typename Digi, typename Geometry,PFLayer::Layer Layer,int Detector>
22 
23  public:
25  PFRecHitCreatorBase(iConfig,iC)
26  {
28  }
29 
30  void importRecHits(std::auto_ptr<reco::PFRecHitCollection>&out,std::auto_ptr<reco::PFRecHitCollection>& cleaned ,const edm::Event& iEvent,const edm::EventSetup& iSetup) {
31 
32 
33  beginEvent(iEvent,iSetup);
34 
36 
38  iSetup.get<CaloGeometryRecord>().get(geoHandle);
39 
40  // get the ecal geometry
41  const CaloSubdetectorGeometry *gTmp =
42  geoHandle->getSubdetectorGeometry(DetId::Hcal, Detector);
43 
44  const Geometry *hcalGeo =dynamic_cast< const Geometry* > (gTmp);
45 
46  iEvent.getByToken(recHitToken_,recHitHandle);
47  for( const auto& erh : *recHitHandle ) {
48  const HcalDetId& detid = (HcalDetId)erh.detid();
50 
51  //since hbhe are together kill other detector
52  if (esd !=Detector && Detector != HcalOther )
53  continue;
54 
55 
56  double energy = erh.energy();
57  double time = erh.time();
58  int depth =detid.depth();
59 
61  math::XYZVector axis;
62 
63  const CaloCellGeometry *thisCell;
64  thisCell= hcalGeo->getGeometry(detid);
65 
66  // find rechit geometry
67  if(!thisCell) {
68  edm::LogError("PFHcalRecHitCreator")
69  <<"warning detid "<<detid.rawId()
70  <<" not found in geometry"<<std::endl;
71  continue;
72  }
73 
74  auto const point = thisCell->getPosition();
75  position.SetCoordinates ( point.x(),
76  point.y(),
77  point.z() );
78 
79 
80 
81 
82  reco::PFRecHit rh( detid.rawId(),Layer,
83  energy,
84  position.x(), position.y(), position.z(),
85  0,0,0);
86  rh.setTime(time); //Mike: This we will use later
87  rh.setDepth(depth);
88 
89  const CaloCellGeometry::CornersVec& corners = thisCell->getCorners();
90  assert( corners.size() == 8 );
91 
92  rh.setNECorner( corners[0].x(), corners[0].y(), corners[0].z());
93  rh.setSECorner( corners[1].x(), corners[1].y(), corners[1].z());
94  rh.setSWCorner( corners[2].x(), corners[2].y(), corners[2].z());
95  rh.setNWCorner( corners[3].x(), corners[3].y(), corners[3].z());
96 
97 
98  bool rcleaned = false;
99  bool keep=true;
100 
101  //Apply Q tests
102  for( const auto& qtest : qualityTests_ ) {
103  if (!qtest->test(rh,erh,rcleaned)) {
104  keep = false;
105 
106  }
107  }
108 
109  if(keep) {
110  out->push_back(rh);
111  }
112  else if (rcleaned)
113  cleaned->push_back(rh);
114  }
115  }
116 
117 
118 
119  protected:
121  int hoDepth_;
122 
123 };
124 
130 
131 
132 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
PFHcalRecHitCreator< HORecHit, CaloSubdetectorGeometry, PFLayer::HCAL_BARREL2, HcalOuter > PFHORecHitCreator
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
void importRecHits(std::auto_ptr< reco::PFRecHitCollection > &out, std::auto_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup)
assert(m_qm.get())
PFHcalRecHitCreator< HBHERecHit, CaloSubdetectorGeometry, PFLayer::HCAL_BARREL1, HcalBarrel > PFHBRecHitCreator
float float float z
SeedingLayerSetsHits::SeedingLayer Layer
Definition: LayerTriplets.h:14
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
PFHcalRecHitCreator< HBHERecHit, CaloSubdetectorGeometry, PFLayer::HCAL_ENDCAP, HcalEndcap > PFHERecHitCreator
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
HcalSubdetector
Definition: HcalAssistant.h:31
edm::EDGetTokenT< edm::SortedCollection< Digi > > recHitToken_
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
size_type size() const
Definition: EZArrayFL.h:81
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const T & get() const
Definition: EventSetup.h:55
PFHcalRecHitCreator< HFRecHit, CaloSubdetectorGeometry, PFLayer::HF_EM, HcalForward > PFHFEMRecHitCreator
static int position[264][3]
Definition: ReadPGInfo.cc:509
PFHcalRecHitCreator< HFRecHit, CaloSubdetectorGeometry, PFLayer::HF_HAD, HcalForward > PFHFHADRecHitCreator
Definition: DDAxes.h:10
PFHcalRecHitCreator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
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