CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFEcalRecHitCreator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFEcalRecHitCreator_h
2 #define RecoParticleFlow_PFClusterProducer_PFEcalRecHitCreator_h
3 
8 
11 
12 
18 
25 
26 template <typename Geometry,PFLayer::Layer Layer,int Detector>
28 
29  public:
31  PFRecHitCreatorBase(iConfig,iC)
32  {
34  }
35 
36  void importRecHits(std::auto_ptr<reco::PFRecHitCollection>&out,std::auto_ptr<reco::PFRecHitCollection>& cleaned ,const edm::Event& iEvent,const edm::EventSetup& iSetup) {
37 
38  beginEvent(iEvent,iSetup);
39 
41 
43  iSetup.get<CaloGeometryRecord>().get(geoHandle);
44 
45  // get the ecal geometry
46  const CaloSubdetectorGeometry *gTmp =
47  geoHandle->getSubdetectorGeometry(DetId::Ecal, Detector);
48 
49  const Geometry *ecalGeo =dynamic_cast< const Geometry* > (gTmp);
50 
51  iEvent.getByToken(recHitToken_,recHitHandle);
52  for(const auto& erh : *recHitHandle ) {
53  const DetId& detid = erh.detid();
54  double energy = erh.energy();
55  double time = erh.time();
56 
58  math::XYZVector axis;
59 
60  const CaloCellGeometry *thisCell;
61  thisCell= ecalGeo->getGeometry(detid);
62 
63  // find rechit geometry
64  if(!thisCell) {
65  edm::LogError("PFEcalRecHitCreator")
66  <<"warning detid "<<detid.rawId()
67  <<" not found in geometry"<<std::endl;
68  continue;
69  }
70 
71 
72  const auto point = thisCell->getPosition();
73  position.SetCoordinates ( point.x(),
74  point.y(),
75  point.z() );
76 
77  // the axis vector is the difference
78  const TruncatedPyramid* pyr
79  = dynamic_cast< const TruncatedPyramid* > (thisCell);
80 
81 
82  if( pyr ) {
83  auto const point1 = pyr->getPosition(1);
84  axis.SetCoordinates( point1.x(),
85  point1.y(),
86  point1.z() );
87 
88  auto const point0 = pyr->getPosition(0);
89 
90  math::XYZVector axis0( point0.x(),
91  point0.y(),
92  point0.z() );
93 
94  axis -= axis0;
95  }
96  else continue;
97 
98  reco::PFRecHit rh( detid.rawId(),Layer,
99  energy,
100  position.x(), position.y(), position.z(),
101  axis.x(), axis.y(), axis.z() );
102 
103 
104  //ECAL has no segmentation so put 1
105 
106  const CaloCellGeometry::CornersVec& corners = thisCell->getCorners();
107  assert( corners.size() == 8 );
108 
109  rh.setNECorner( corners[0].x(), corners[0].y(), corners[0].z() );
110  rh.setSECorner( corners[1].x(), corners[1].y(), corners[1].z() );
111  rh.setSWCorner( corners[2].x(), corners[2].y(), corners[2].z() );
112  rh.setNWCorner( corners[3].x(), corners[3].y(), corners[3].z() );
113 
114 
115  bool rcleaned = false;
116  bool keep=true;
117 
118  //Apply Q tests
119  for( const auto& qtest : qualityTests_ ) {
120  if (!qtest->test(rh,erh,rcleaned)) {
121  keep = false;
122  }
123  }
124 
125  if(keep) {
126  rh.setTime(time);
127  rh.setDepth(1);
128  out->push_back(rh);
129  }
130  else if (rcleaned)
131  cleaned->push_back(rh);
132  }
133  }
134 
135 
136 
137  protected:
139 
140 
141 };
142 
143 
146 
147 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
void importRecHits(std::auto_ptr< reco::PFRecHitCollection > &out, std::auto_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
assert(m_qm.get())
SeedingLayerSetsHits::SeedingLayer Layer
Definition: LayerTriplets.h:14
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const int keep
PFEcalRecHitCreator< EcalEndcapGeometry, PFLayer::ECAL_ENDCAP, EcalEndcap > PFEERecHitCreator
int iEvent
Definition: GenABIO.cc:230
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
edm::EDGetTokenT< EcalRecHitCollection > recHitToken_
void beginEvent(const edm::Event &event, const edm::EventSetup &setup)
tuple out
Definition: dbtoconf.py:99
Definition: DetId.h:18
const GlobalPoint getPosition(CCGFloat depth) const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const T & get() const
Definition: EventSetup.h:55
A base class to handle the particular shape of Ecal Xtals. Taken from ORCA Calorimetry Code...
static int position[264][3]
Definition: ReadPGInfo.cc:509
PFEcalRecHitCreator< EcalBarrelGeometry, PFLayer::ECAL_BARREL, EcalBarrel > PFEBRecHitCreator
PFEcalRecHitCreator(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